forked from NCAR/apex_fortran
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathapex.f90
2968 lines (2803 loc) · 112 KB
/
apex.f90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
module apex
!
! This software is part of the NCAR TIE-GCM. Use is governed by the
! Open Source Academic Research License Agreement contained in the file
! tiegcmlicense.txt.
!
! 2018 May: A.D. Richmond modifications to add full base vectors f3,g1,g2,g3;
! to update comments, to unify specification of constants, and to add
! subroutine apex_setup that simplifies initialization for many purposes.
! 2013 April: B. Foster (NCAR/HAO)
! This is a refactored version of the legacy apex code, originally written
! by Art Richmond and Roy Barnes in the 1993-2000 timeframe, using
! field-line tracing routines developed in 1973 by Wally Clark at NOAA.
! This version is written in free-format fortran90. Subroutines and
! module data may be use-associated from this module.
!
! Reference: Richmond, A. D., Ionospheric Electrodynamics Using Magnetic Apex
! Coordinates, J. Geomag. Geoelectr., 47, 191-212, 1995,
! https://doi.org/10.5636/jgg.47.191
!
! A typical calling sequence for a code calling this module is as follows:
!
! subroutine apex_setup(date,altmax):
! Set up 3D arrays over the Earth for geomagnetic-field epoch date,
! up to an altitude of at least altmax,
! for interpolation use when calling apex_mall, apex_q2g, or apex_m2g.
! (This must be called before apex_mall, apex_q2g, or apex_m2g)
! The resolution parameter nvert is set to 40 here, which results in an array
! of lat,lon dimensions 121,201, adequate for most purposes.
! The nalt dimension depends on nvert and altmax; for nvert=40 and altmax
! between 910.2 and 1124.3, nalt=7 . For very large altmax the maximum
! nalt is nvert+1.
!
! subroutine apex_mall:
! Calculate modified Apex coordinates, basis vectors described in
! Richmond [1995], and other magnetic field parameters
! (often called from lat,lon,alt nested loop)
!
! subroutine apex_q2g or apex_m2g:
! Convert from quasi-dipole or modified magnetic-apex coordinates
! to geodetic coordinates.
!
! Subroutines called by those above can be invoked directly to do coordinate
! conversions; to get geomagnetic-field quantities and grid-point locations;
! and to calculate the subsolar location, cosine of the solar zenith angle,
! and magnetic local time.
!
implicit none
real,parameter :: re = 6371.2, & ! Reference radius of IGRF
eps = 1.e-5, & ! Small number
fltnvrs = 298.257223563 ! Inverse flatness of geoid
real,allocatable,save :: &
xarray(:,:,:), & ! cos(quasi-dipole latitude)*cos(apex longitude)
! minus dipole value
yarray(:,:,:), & ! cos(quasi-dipole latitude)*sin(apex longitude)
! minus dipole value
zarray(:,:,:), & ! sin(quasi-dipole latitude)
! minus dipole value
varray(:,:,:) ! magnetic potential, normalized by radially varying
! dipole value at S pole,
! minus dipole value
!
integer :: nglat,nglon,ngalt ! numbers of lat,lon,alt grid points,
! set in apex_mka to nlat,nlon,nalt as determined by ggrid.
!
! This grid (geolat,geolon,geoalt is equivalent to gdlat,gdlon,gdalt,
! as passed to apex_mka.
real,allocatable,save :: geolat(:), geolon(:), geoalt(:)
integer,parameter :: nmax=13 ! maximum degree of IGRF coefficients
integer,parameter :: ncoef = nmax*nmax + 2*nmax + 1 ! 196
real,dimension(ncoef) :: &
gb, & ! Coefficients for magnetic field calculation
gv ! Coefficients for magnetic potential calculation
!
real, parameter :: &
dtr=0.0174532925199432957692369076847, & ! degrees to radians
rtd=57.2957795130823208767981548147 ! radians to degrees
real :: &
pola ! Limiting latitude magnitude (deg); when the geographic
! latitude is poleward of pola, then xarray,yarray,zarray,varray
! are forced to be constant for all longitudes at each altitude.
real,parameter :: & ! Formerly common /APXCON/
req = 6378.137, & ! Equatorial earth radius
precise = 7.6e-11, & ! Precision factor
glatlim = 89.9, & ! Limit above which gradients are recalculated
xmiss = -32767., & ! Set missing values to this
req2 = req*req, &
req2e2 = req2*(2. - 1./fltnvrs)/fltnvrs
!
real :: & ! Formerly /APXDIPL/ and /DIPOLE/
colat, & ! Geocentric colatitude of geomagnetic dipole north pole (deg)
elon, & ! East longitude of geomagnetic dipole north pole (deg)
vp, & ! Magnitude, in T.m, of dipole component of magnetic
! potential at geomagnetic pole and geocentric radius re
ctp,stp ! cosine and sine of colat
!
real :: & ! Formerly /FLDCOMD/
bx, & ! X comp. of field vector at the current tracing point (Gauss)
by, & ! Y comp. of field vector at the current tracing point (Gauss)
bz, & ! Z comp. of field vector at the current tracing point (Gauss)
bb ! Magnitude of field vector at the current tracing point (Gauss)
real :: & ! Formerly /APXIN/
yapx(3,3) ! Matrix of cartesian coordinates (loaded columnwise)
!
! /ITRA/ was only in subs linapx and itrace, so it could be removed from module
! data if those routines are properly restructured to pass the /ITRA/
! variables.
!
integer :: & ! Formerly /ITRA/
nstp ! Step count. Incremented in sub linapx.
real :: &
y(3), & ! Array containing current tracing point cartesian coordinates.
yp(3), & ! Array containing previous tracing point cartesian coordinates.
sgn, & ! Determines direction of trace. Set in subprogram linapx
ds ! Step size (Km) Computed in subprogram linapx.
real :: & ! limits beyond which east-west gradients are computed
glatmn,glatmx ! differently to avoid potential underflow (apex_mka)
contains
!-----------------------------------------------------------------------
subroutine apex_setup(date,altmax)
! Set up 3D arrays over the Earth for geomagnetic-field epoch date,
! up to an altitude of at least altmax,
! for use when calling apex_mall, apex_q2g, or apex_m2g.
real, intent(in) :: &
date, & ! Date in year and fraction (e.g., 2018.5)
altmax ! Maximum altitude that will be used (km).
!
integer,parameter :: nvert = 40 ! Resolution parameter, corresponding
! to the maximum number of vertical grid increments when altmax=infinity.
! Points are spaced uniformly in 1/r between the Earth's surface and an
! altitude that is at least altmax (or a very large value if altmax=infinity).
integer,parameter :: &
mxlat = 3*nvert+1, & ! number of latitude grid points
mxlon = 5*nvert+1, & ! number of longitude grid points
mxalt = nvert+1 ! maximum number of height grid points
! Local
integer :: nlat,nlon,nalt ! Numbers of lat, lon, alt grid points
real :: gplat(mxlat),gplon(mxlon),gpalt(mxalt) ! grid latitudes (deg),
! longitudes (deg), and altitudes (km), respectively
integer :: ier ! Error flag, non-zero if there is a problem.
nlat = 0
nlon = 0
nalt = 0
call ggrid(nvert,-90.,90.,-180.,180.,0.,altmax, &
gplat,gplon,gpalt,mxlat,mxlon,mxalt,nlat,nlon,nalt)
call apex_mka(date,gplat,gplon,gpalt,nlat,nlon,nalt,ier)
if (ier /= 0) then
write(6,"('>>> apex_setup error from apex_mka')")
stop 'apex_setup apex_mka'
endif
!
!
end subroutine apex_setup
!-----------------------------------------------------------------------
subroutine ggrid(nvert,glatmin,glatmax,glonmin,glonmax,altmin,altmax, &
gplat,gplon,gpalt,mxlat,mxlon,mxalt,nlat,nlon,nalt)
!
! Given desired range of geographic latitude, longitude and altitude,
! choose an appropriate grid that can be used in subsequent calls to
! subs apex_mka, apex_mall, apex_q2g.
!
! Input args:
integer,intent(in) :: nvert,mxlat,mxlon,mxalt
real,intent(in) :: glatmin,glatmax,glonmin,glonmax,altmin,altmax
!
! Output args:
integer,intent(out) :: nlat,nlon,nalt
real,intent(out) :: gplat(mxlat),gplon(mxlon),gpalt(mxalt)
!
! Local:
real :: dlon,dlat,diht,dnv,glonmaxx,x
integer :: nlatmin,nlatmax,nlonmin,nlonmax,naltmin,naltmax
integer :: i,j,k,kk
!
! Check inputs:
if (glatmin > glatmax) then
write(6,"('>>> ggrid: glatmin=',f9.2,' must be <= glatmax=',f9.2)") &
glatmin,glatmax
stop 'ggrid'
endif
if (glonmin > glonmax) then
write(6,"('>>> ggrid: glonmin=',f9.2,' must be <= glonmax=',f9.2)") &
glonmin,glonmax
stop 'ggrid'
endif
if (altmin > altmax) then
write(6,"('>>> ggrid: altmin=',f9.2,' must be <= altmax=',f9.2)") &
altmin,altmax
stop 'ggrid'
endif
!
! Init outputs:
nlat = 0 ; nlon = 0 ; nalt = 0
gplat = 0. ; gplon = 0. ; gpalt = 0.
!
dnv = float(nvert)
dlon = 360. / (5.*dnv)
dlat = 180. / (3.*dnv)
diht = 1. / dnv
nlatmin = max(int((glatmin+90.)/dlat),0)
nlatmax = min(int((glatmax+90.)/dlat+1.),3*nvert)
nlonmin = max(int((glonmin+180.)/dlon),0)
glonmaxx = min(glonmax,glonmin+360.)
nlonmax = min(int((glonmaxx+180.)/dlon+1.),10*nvert)
x = re/(re+altmax)/diht-eps
naltmin = max(x,1.)
naltmin = min(naltmin,nvert-1)
x = re/(re+altmin)/diht+eps
i = x + 1.
naltmax = min(i,nvert)
nlat = nlatmax - nlatmin + 1
nlon = nlonmax - nlonmin + 1
nlon = min(nlon,5*nvert+1)
nalt = naltmax - naltmin + 1
do j=1,nlat
gplat(j) = dlat*float(nlatmin+j-1) - 90.
enddo
gplat(nlat) = amin1(gplat(nlat),90.)
do i=1,nlon
gplon(i) = dlon*float(nlonmin+i-1) - 180.
enddo
do k=1,nalt
kk = naltmax - k +1
gpalt(k) = re*(float(nvert-kk) - eps) / (float(kk)+eps)
enddo
if (gplon(nlon-1) >= glonmax) nlon = nlon-1
gpalt(1) = max(gpalt(1),0.)
! write(6,"('ggrid: nlat=',i4,' gplat=',/,(6f9.2))") nlat,gplat
! write(6,"('ggrid: nlon=',i4,' gplon=',/,(6f9.2))") nlon,gplon
! write(6,"('ggrid: nalt=',i4,' gpalt=',/,(6f9.2))") nalt,gpalt
end subroutine ggrid
!-----------------------------------------------------------------------
subroutine apex_mka(date,gplat,gplon,gpalt,nlat,nlon,nalt,ier)
!
! Given a 3d lat,lon,altitude grid, calculate xarray,yarray,zarray,varray in
! module data above. These arrays are used later for calculating quantities
! involving gradients of Apex coordinates, such as base vectors in the
! Modified-Apex and Quasi-Dipole systems.
!
! This defines module 3d data xarray,yarray,zarray,varray
!
! Input args:
real,intent(in) :: date ! year and fraction
integer,intent(in) :: nlat,nlon,nalt ! dimensions of 3d grid
real,intent(inout) :: gplat(nlat),gplon(nlon),gpalt(nalt)
!
! Output args:
integer,intent(out) :: ier
!
! Local:
integer :: i,j,k,kpol,istat
real :: reqore,rqorm1,cp,ct,st,sp,stmcpm,stmspm,ctm
real :: aht,alat,phia,bmag,xmag,ymag,zdown,vmp ! apex_sub output
real :: vnor,rp,reqam1,a,slp,clp,phiar
!
! pola:
! Limiting latitude magnitude (deg); when the geographic latitude is
! poleward of pola, then xarray,yarray,zarray,varray are forced to
! be constant for all longitudes at each altitude.
!
! This makes pola = 89.9995 for precise = 7.6e-11:
pola = 90.-sqrt(precise)*rtd ! Pole angle (deg)
ier = 0
!
! Allocate 3d xarray,yarray,zarray,varray arrays:
! Deallocate if necessary in case this routine is called more than once
! (e.g., called by each of coupled models with different grids)
!
if (allocated(xarray)) deallocate(xarray)
allocate(xarray(nlat,nlon,nalt),stat=istat)
if (istat /= 0) stop 'allocate xarray'
xarray = 0.
if (allocated(yarray)) deallocate(yarray)
allocate(yarray(nlat,nlon,nalt),stat=istat)
if (istat /= 0) stop 'allocate yarray'
yarray = 0.
if (allocated(zarray)) deallocate(zarray)
allocate(zarray(nlat,nlon,nalt),stat=istat)
if (istat /= 0) stop 'allocate zarray'
zarray = 0.
if (allocated(varray)) deallocate(varray)
allocate(varray(nlat,nlon,nalt),stat=istat)
if (istat /= 0) stop 'allocate varray'
varray = 0.
!
! Set geographic grid in module data for later reference:
! (these also are not deallocated by this module)
!
nglon=nlon ; nglat=nlat ; ngalt=nalt
if (allocated(geolat)) deallocate(geolat)
if (allocated(geolon)) deallocate(geolon)
if (allocated(geoalt)) deallocate(geoalt)
allocate(geolat(nglat),stat=istat)
allocate(geolon(nglon),stat=istat)
allocate(geoalt(ngalt),stat=istat)
geolat(:) = gplat(:)
geolon(:) = gplon(:)
geoalt(:) = gpalt(:)
!
! Set coefficients gb,gv (module data) for requested year:
!
call cofrm(date)
! 20180531 ADR: cofrm also calculates dipole quantities colat,elon,vp,ctp,stp
! write(6,"('apex_mka after cofrm: ncoef=',i4,' gb=',/,(6f12.3))") ncoef,gb
! write(6,"('apex_mka after cofrm: ncoef=',i4,' gv=',/,(6f12.3))") ncoef,gv
reqore = req/re
rqorm1 = reqore-1.
do j=1,nlat
ct = sin(gplat(j)*dtr)
st = cos(gplat(j)*dtr)
kpol = 0
if (abs(gplat(j)) > pola) kpol = 1
do i=1,nlon
cp = cos((gplon(i)-elon)*dtr)
sp = sin((gplon(i)-elon)*dtr)
!
! ctm is pseudodipole component of z
! -ctm is pseudodipole component of v
! stmcpm is pseudodipole component of x
! stmspm is pseudodipole component of y
!
ctm = ctp*ct + stp*st*cp
stmcpm = st*ctp*cp - ct*stp
stmspm = st*sp
do k=1,nalt
call apex_sub(gplat(j),gplon(i),gpalt(k),&
aht,alat,phia,bmag,xmag,ymag,zdown,vmp)
vnor = vmp/vp
rp = 1. + gpalt(k)/re
varray(j,i,k) = vnor*rp*rp + ctm
reqam1 = req*(aht-1.)
slp = sqrt(max(reqam1-gpalt(k),0.)/(reqam1+re))
!
! Reverse sign of slp in southern magnetic hemisphere
!
if (zdown.lt.0.) slp = -slp
clp = sqrt (rp/(reqore*aht-rqorm1))
phiar = phia*dtr
xarray(j,i,k) = clp*cos (phiar) - stmcpm
yarray(j,i,k) = clp*sin (phiar) - stmspm
zarray(j,i,k) = slp - ctm
enddo ! k=1,nalt
if (kpol==1.and.i > 1) then
xarray(j,i,:) = xarray(j,1,:)
yarray(j,i,:) = yarray(j,1,:)
zarray(j,i,:) = zarray(j,1,:)
varray(j,i,:) = varray(j,1,:)
cycle
endif
enddo ! i=1,nlon
enddo ! j=1,nlat
!
! Establish for this grid polar latitude limits beyond which east-west
! gradients are computed differently to avoid potential underflow
! (glatmx,glatmn are in module data, glatlim is parameter constant)
!
glatmx = max( glatlim,gplat(nlat-2))
glatmn = min(-glatlim,gplat(2))
end subroutine apex_mka
!-----------------------------------------------------------------------
subroutine apex_mall(glat,glon,alt,hr, b,bhat,bmag,si,alon,xlatm,vmp,W,&
D,be3,sim,d1,d2,d3,e1,e2,e3,xlatqd,F,f1,f2,f3,g1,g2,g3,ier)
!
! Compute Modified Apex coordinates, quasi-dipole coordinates,
! base vectors and other parameters by interpolation from
! precalculated arrays. Subroutine apex_setup or apex_mka must be called
! before calling this subroutine.
!
! Args:
real,intent(in) :: & ! Input
glat ,& ! Geographic (geodetic) latitude (deg)
glon ,& ! Geographic (geodetic) longitude (deg)
alt ,& ! Altitude (km)
hr ! Reference altitude (km)
real,intent(out) :: & ! Output
b(3) ,& ! Magnetic field components (east, north, up), in nT
bhat(3) ,& ! components (east, north, up) of unit vector along
! geomagnetic field direction
bmag ,& ! Magnitude of magnetic field (nT)
si ,& ! sine of magnetic inclination
alon ,& ! Apex longitude = modified apex longitude =
! quasi-dipole longitude (deg)
xlatm ,& ! Modified Apex latitude (deg)
vmp ,& ! Magnetic potential (T.m)
W ,& ! W of Richmond [1995] reference above, in km**2 /nT
! (i.e., 10**15 m**2 /T)
D ,& ! D of Richmond [1995]
be3 ,& ! B_e3 of Richmond [1995] (= Bmag/D), in nT
sim ,& ! sin(I_m) described in Richmond [1995]
xlatqd ,& ! Quasi-dipole latitude (deg)
F ! F described in Richmond [1995] for quasi-dipole
! coordinates
!
real,dimension(3),intent(out) :: d1,d2,d3,e1,e2,e3,f1,f2,f3,g1,g2,g3 ! Components of base vectors
integer,intent(out) :: ier ! error return
!
! Local:
real :: glonloc,cth,sth,glatx,clm,r3_2
real :: fx,fy,fz,fv
real :: dfxdth,dfydth,dfzdth,dfvdth, &
dfxdln,dfydln,dfzdln,dfvdln, &
dfxdh ,dfydh ,dfzdh ,dfvdh
real,dimension(3) :: gradx,grady,gradz,gradv, grclm,clmgrp,rgrlp
real :: & ! dummies for polar calls to intrp
fxdum,fydum,fzdum,fvdum, &
dmxdth,dmydth,dmzdth,dmvdth, &
dmxdh,dmydh,dmzdh,dmvdh
!
! Init:
!
ier = 0
glonloc = glon
call intrp (glat,glonloc,alt, geolat,geolon,geoalt,nglat,nglon,ngalt, &
fx,fy,fz,fv, &
dfxdth,dfydth,dfzdth,dfvdth, &
dfxdln,dfydln,dfzdln,dfvdln, &
dfxdh ,dfydh ,dfzdh ,dfvdh, ier)
if (ier /= 0) then
call setmiss(xmiss,xlatm,alon,vmp,b,bmag,be3,sim,si,F,D,W, &
bhat,d1,d2,d3,e1,e2,e3,f1,f2,f3,g1,g2,g3)
write(6,"('apex_mall called setmiss: glat,glon,alt=',3f12.3)") &
glat,glon,alt
return
endif
call adpl(glat,glonloc,cth,sth,fx,fy,fz,fv, &
dfxdth,dfydth,dfzdth,dfvdth,dfxdln,dfydln,dfzdln,dfvdln)
call gradxyzv(alt,cth,sth, &
dfxdth,dfydth,dfzdth,dfvdth,dfxdln,dfydln,dfzdln,dfvdln, &
dfxdh,dfydh,dfzdh,dfvdh,gradx,grady,gradz,gradv)
!
! If the point is very close to either the North or South
! geographic pole, recompute the east-west gradients after
! stepping a small distance from the pole.
!
if (glat > glatmx .or. glat < glatmn) then
glatx = glatmx
if (glat < 0.) glatx = glatmn
call intrp (glatx,glonloc,alt, geolat,geolon,geoalt,nglat,nglon,ngalt, &
fxdum,fydum,fzdum,fvdum, &
dmxdth,dmydth,dmzdth,dmvdth,dfxdln,dfydln,dfzdln, &
dfvdln,dmxdh,dmydh,dmzdh,dmvdh, ier)
call adpl(glatx,glonloc,cth,sth,fxdum,fydum,fzdum,fvdum, &
dmxdth,dmydth,dmzdth,dmvdth,dfxdln,dfydln,dfzdln,dfvdln)
call grapxyzv(alt,cth,sth,dfxdln,dfydln,dfzdln,dfvdln, &
gradx,grady,gradz,gradv)
endif
call gradlpv(hr,alt,fx,fy,fz,fv,gradx,grady,gradz,gradv, &
xlatm,alon,vmp,grclm,clmgrp,xlatqd,rgrlp,b,clm,r3_2)
call basevec(hr,xlatm,grclm,clmgrp,rgrlp,b,clm,r3_2, &
bmag,sim,si,F,D,W,bhat,d1,d2,d3,e1,e2,e3,f1,f2,f3,g1,g2,g3)
be3 = bmag/D
ier = 0
end subroutine apex_mall
!-----------------------------------------------------------------------
subroutine apex_q2g(qdlat_in,qdlon,alt,gdlat,gdlon,ier)
!
! Convert from quasi-dipole to geodetic coordinates. This subroutine
! (input magnetic, output geodetic) is the functional inverse of
! subroutine apex_mall (input geodetic, output magnetic). Sub apex_setup
! or apex_mka must be called before this routine.
!
! 20160811 ADR If qdlat_in is +/- 90, move it a little off the pole so
! that when gdlat, gdlon are used in apxmall the directions of the base
! vectors vary appropriately with longitude.
!
! Args:
real,intent(in) :: & ! inputs
qdlat_in, & ! quasi-dipole latitude (deg)
qdlon, & ! quasi-dipole longitude (deg)
alt ! altitude (km)
real,intent(out) :: & ! outputs
gdlat, & ! geodetic latitude (deg)
gdlon ! geodetic longitude (deg)
integer,intent(out) :: ier ! error return
!
! Local:
real :: x0,y0,z0,xnorm,xdif,ydif,zdif,dist2,hgrd2e,hgrd2n,hgrd2,&
angdist,distlon,glatx,cal,sal,coslm,slm,cad,sad,slp,clm2,slm2,&
sad2,cal2,clp2,clp,dylon,qdlat
real :: ylat,ylon ! first guess output by gm2gc, input to intrp
integer :: iter
integer,parameter :: niter=20
real :: & ! output of sub intrp
fx,fy,fz,fv, & ! interpolated values of x,y,z,v
dfxdth,dfydth,dfzdth,dfvdth, & ! derivatives of x,y,z,v wrt colatitude
dfxdln,dfydln,dfzdln,dfvdln, & ! derivatives of x,y,z,v wrt longitude
dfxdh ,dfydh ,dfzdh ,dfvdh ! derivatives of x,y,z,v wrt altitude
real :: & ! dummies for polar calls to intrp
fxdum,fydum,fzdum,fvdum, &
dmxdth,dmydth,dmzdth,dmvdth, &
dmxdh,dmydh,dmzdh,dmvdh
real :: cth,sth ! output of adpl
character(len=5) :: edge
ier = 0 ; gdlat = 0. ; gdlon = 0.
! Keep qdlat away from poles
qdlat = amax1(qdlat_in,-90.+sqrt(precise)*rtd)
qdlat = amin1(qdlat , 90.-sqrt(precise)*rtd)
!
! Determine quasi-cartesian coordinates on a unit sphere of the
! desired magnetic lat,lon in quasi-dipole coordinates.
!
x0 = cos (qdlat*dtr) * cos (qdlon*dtr)
y0 = cos (qdlat*dtr) * sin (qdlon*dtr)
z0 = sin (qdlat*dtr)
!
! Initial guess: use centered dipole, convert to geocentric coords
!
call gm2gc (qdlat,qdlon,ylat,ylon)
!
! Iterate until (angular distance)**2 (units: radians) is within
! precise of location (qdlat,qdlon) on a unit sphere.
! (precise is a parameter in module data)
!
do iter=1,niter
!
! geolat,lon,alt and nglat,lon,alt are in module data (set by apex_mka)
!
call intrp (ylat,ylon,alt, geolat,geolon,geoalt,nglat,nglon,ngalt, &
fx,fy,fz,fv, &
dfxdth,dfydth,dfzdth,dfvdth, &
dfxdln,dfydln,dfzdln,dfvdln, &
dfxdh ,dfydh ,dfzdh ,dfvdh, ier)
if (ier /= 0) then
write(6,"('>>> apex_q2g error from intrp')")
stop 'qpex_q2g intrp'
endif
!
! Add-back of pseudodipole component to x,y,z,v and their derivatives.
!
call adpl(ylat,ylon,cth,sth,fx,fy,fz,fv, &
dfxdth,dfydth,dfzdth,dfvdth,dfxdln,dfydln,dfzdln,dfvdln)
distlon = cos(ylat*dtr)
if (ylat > glatmx .or. ylat < glatmn) then ! glatmx,glatmn are module data
glatx = glatmx
if (ylat.lt.0.) glatx = glatmn
distlon = cos (glatx*dtr)
call intrp (glatx,ylon,alt, geolat,geolon,geoalt,nglat,nglon,ngalt, &
fxdum,fydum,fzdum,fvdum, &
dmxdth,dmydth,dmzdth,dmvdth,dfxdln,dfydln,dfzdln, &
dfvdln,dmxdh,dmydh,dmzdh,dmvdh, ier)
if (ier /= 0) then
write(6,"('>>> apex_q2g error from polar intrp')")
stop 'qpex_q2g intrp'
endif
call adpl(glatx,ylon,cth,sth,fxdum,fydum,fzdum,fvdum, &
dmxdth,dmydth,dmzdth,dmvdth,dfxdln,dfydln,dfzdln,dfvdln)
endif
!
! At this point, fx,fy,fz are approximate quasi-cartesian
! coordinates on a unit sphere for the quasi-dipole coordinates
! corresponding to the geodetic coordinates ylat, ylon.
! Normalize the vector length of (fx,fy,fz) to unity using xnorm
! so that the resultant vector can be directly compared with the
! target vector (x0,y0,z0).
!
xnorm = sqrt(fx*fx + fy*fy + fz*fz)
xdif = fx/xnorm - x0
ydif = fy/xnorm - y0
zdif = fz/xnorm - z0
!
! dist2 = square of distance between normalized (fx,fy,fz) and x0,y0,z0.
!
dist2 = xdif*xdif + ydif*ydif + zdif*zdif
if (dist2 <= precise) then
ier = 0
gdlat = ylat
gdlon = ylon
return
endif
!
! hgrd2* = one-half of east or north gradient of dist2 on unit sphere.
!
hgrd2e = (xdif*dfxdln + ydif*dfydln + zdif*dfzdln)/distlon
hgrd2n = -(xdif*dfxdth + ydif*dfydth + zdif*dfzdth)
hgrd2 = sqrt(hgrd2e*hgrd2e + hgrd2n*hgrd2n)
!
! angdist = magnitude of angular distance to be moved for new guess
! of ylat, ylon.
!
angdist = dist2/hgrd2
!
! Following spherical trigonometry moves ylat,ylon to new location,
! in direction of grad(dist2), by amount angdist.
!
cal = -hgrd2n/hgrd2
sal = -hgrd2e/hgrd2
coslm = cos(ylat*dtr)
slm = sin(ylat*dtr)
cad = cos(angdist)
sad = sin(angdist)
slp = slm*cad + coslm*sad*cal
clm2 = coslm*coslm
slm2 = slm*slm
sad2 = sad*sad
cal2 = cal*caL
clp2 = clm2 + slm2*sad2 - 2.*slm*cad*coslm*sad*cal -clm2*sad2*cal2
clp = sqrt (max(0.,clp2))
ylat = atan2(slp,clp)*rtd
!
! Restrict latitude iterations to stay within the interpolation grid
! limits, but let intrp find any longitude exceedence. This is only
! an issue when the interpolation grid does not cover the entire
! magnetic pole region.
!
ylat = min(ylat,geolat(nglat))
ylat = max(ylat,geolat(1))
dylon = atan2 (sad*sal,cad*coslm-sad*slm*cal)*rtd
ylon = ylon + dylon
if (ylon > geolon(nglon)) ylon = ylon - 360.
if (ylon < geolon(1)) ylon = ylon + 360.
enddo ! iter=1,niter
write(6,"('>>> apex_q2g: ',i3,' iterations only reduced the angular')") niter
write(6,"(' difference to ',f10.5,' degrees, where test criterion')") &
sqrt(dist2)*rtd
write(6,"(' is ',f10.5,' degrees.')") sqrt(precise)*rtd
edge = ' '
if (ylat == geolat(nglat)) edge = 'north'
if (ylat == geolat(1)) edge = 'south'
if (edge /= ' ') then
write(6,"('Coordinates are on the ',a,' edge of the interpolation grid ')") edge
write(6,"('and latitude is constrained to stay within grid limits when iterating.')")
endif
ier = 1
end subroutine apex_q2g
!-----------------------------------------------------------------------
subroutine gradxyzv(alt,cth,sth, &
dfxdth,dfydth,dfzdth,dfvdth,dfxdln,dfydln,dfzdln,dfvdln, &
dfxdh,dfydh,dfzdh,dfvdh,gradx,grady,gradz,gradv)
!
! Calculates east,north,up components of gradients of x,y,z,v in
! geodetic coordinates. All gradients are in inverse km.
! 940803 A. D. Richmond; updated 180510
!
! Args:
real,intent(in) :: alt,cth,sth
real,dimension(3),intent(out) :: gradx,grady,gradz,gradv
real,intent(in) :: &
dfxdth,dfydth,dfzdth,dfvdth, &
dfxdln,dfydln,dfzdln,dfvdln, &
dfxdh,dfydh,dfzdh,dfvdh
!
! Local:
real :: d,d2,rho,dddthod,drhodth,dzetdth,ddisdth
!
d2 = req2 - req2e2*cth*cth
d = sqrt(d2)
rho = sth*(alt + req2/d)
dddthod = req2e2*cth*sth/d2
drhodth = alt*cth + (req2/d)*(cth-sth*dddthod)
dzetdth =-alt*sth - ((req2-req2e2)/d)*(sth+cth*dddthod)
ddisdth = sqrt(drhodth*drhodth + dzetdth*dzetdth)
gradx(1) = dfxdln/rho
grady(1) = dfydln/rho
gradz(1) = dfzdln/rho
gradv(1) = dfvdln/rho
gradx(2) = -dfxdth/ddisdth
grady(2) = -dfydth/ddisdth
gradz(2) = -dfzdth/ddisdth
gradv(2) = -dfvdth/ddisdth
gradx(3) = dfxdh
grady(3) = dfydh
gradz(3) = dfzdh
gradv(3) = dfvdh
end subroutine gradxyzv
!-----------------------------------------------------------------------
subroutine grapxyzv(alt,cth,sth, &
dfxdln,dfydln,dfzdln,dfvdln,gradx,grady,gradz,gradv)
!
! Calculates east component of gradient near pole.
!
! Args:
real,intent(in) :: alt,cth,sth
real,intent(in) :: dfxdln,dfydln,dfzdln,dfvdln
real,dimension(3),intent(out) :: gradx,grady,gradz,gradv
!
! Local:
real :: d,d2,rho,dddthod,drhodth,dzetdth,ddisdth
!
d2 = req2 - req2e2*cth*cth
d = sqrt(d2)
rho = sth*(alt + req2/d)
dddthod = req2e2*cth*sth/d2
drhodth = alt*cth + (req2/d)*(cth-sth*dddthod)
dzetdth =-alt*sth - ((req2-req2e2)/d)*(sth+cth*dddthod)
ddisdth = sqrt(drhodth*drhodth + dzetdth*dzetdth)
gradx(1) = dfxdln/rho
grady(1) = dfydln/rho
gradz(1) = dfzdln/rho
gradv(1) = dfvdln/rho
end subroutine grapxyzv
!-----------------------------------------------------------------------
subroutine gradlpv(hr,alt,fx,fy,fz,fv,gradx,grady,gradz,gradv, &
xlatm,xlonm,vmp,grclm,clmgrp,qdlat,rgrlp,b,clm,r3_2)
!
! Uses gradients of x,y,z,v to compute geomagnetic field and
! gradients of apex latitude, longitude.
!
! Args:
real,intent(in) :: & ! scalar inputs
hr, & ! reference altitude (km)
alt, & ! altitude (km)
fx,fy,fz,fv ! interpolated values of x,y,z,v, plus
! pseudodipole component
real,dimension(3),intent(in) :: & ! 3-component inputs
gradx,grady,gradz,gradv ! interpolated gradients of x,y,z,v,
! including pseudodipole components (east,north,up)
!
! Local:
integer :: i
real :: rr,r,rn,sqrror,cpm,spm,bo,rn2,x2py2,xnorm,xlp,slp,clp,grclp
real,intent(out) :: & ! scalar outputs
xlatm, & ! modified apex latitude (lambda_m), degrees
xlonm, & ! apex longitude (phi_a), degrees
vmp, & ! magnetic potential, in T.m.
qdlat, & ! quasi-dipole latitude, degrees
clm, & ! cos(lambda_m)
r3_2 ! ((re + alt)/(re + hr))**(3/2)
real,dimension(3),intent(out) :: & ! 3-component outputs
grclm, & ! grad(cos(lambda_m)), in km-1
clmgrp, & ! cos(lambda_m)*grad(phi_a), in km-1
rgrlp, & ! (re + alt)*grad(lambda')
b ! magnetic field, in nT
xlatm=0. ; xlonm=0. ; vmp=0 ; grclm=0. ; clmgrp=0. ; rgrlp = 0. ; b=0.
clm=0. ; r3_2=0. ; qdlat=0.
rr = re + hr
r = re + alt
rn = r/re
sqrror = sqrt(rr/r)
r3_2 = 1./sqrror/sqrror/sqrror
xlonm = atan2(fy,fx)
cpm = cos(xlonm)
spm = sin(xlonm)
xlonm = rtd*xlonm ! output
bo = vp*1.e6 ! vp is module data; 1.e6 converts T to nT and km-1 to m-1
rn2 = rn*rn
vmp = vp*fv/rn2 ! output
b(1) = -bo*gradv(1)/rn2
b(2) = -bo*gradv(2)/rn2
b(3) = -bo*(gradv(3)-2.*fv/r)/rn2
x2py2 = fx*fx + fy*fy
xnorm = sqrt(x2py2 + fz*fz)
xlp = atan2(fz,sqrt(x2py2))
slp = sin(xlp)
clp = cos(xlp)
qdlat = xlp*rtd ! output
clm = sqrror*clp ! output
if (clm > 1.) then
write(6,"('>>> gradlpv: hr=',f12.3,' alt=',f12.3)") hr,alt
write(6,"(' Point lies below field line that peaks at reference height.')")
stop 'gradlpv'
endif
xlatm = rtd*acos(clm)
!
! If southern magnetic hemisphere, reverse sign of xlatm
!
if (slp < 0.) xlatm = -xlatm
do i=1,3
grclp = cpm*gradx(i) + spm*grady(i)
rgrlp(i) = r*(clp*gradz(i) - slp*grclp)
grclm(i) = sqrror*grclp
clmgrp(i) = sqrror*(cpm*grady(i)-spm*gradx(i))
enddo
grclm(3) = grclm(3) - sqrror*clp/(2.*r)
end subroutine gradlpv
!-----------------------------------------------------------------------
subroutine basevec(hr,xlatm,grclm,clmgrp,rgrlp,b,clm,r3_2, &
bmag,sim,si,F,D,W,bhat,d1,d2,d3,e1,e2,e3,f1,f2,f3,g1,g2,g3)
!
! Computes base vectors and other parameters for apex coordinates.
! Vector components: east, north, up
!
! Args:
real,intent(in) :: & ! scalar inputs
hr, & ! reference altitude
xlatm, & ! modified apex latitude (deg)
clm, & ! cos(lambda_m)
r3_2 ! ((re + altitude)/(re + hr))**(3/2)
real,dimension(3),intent(in) :: & ! 3-component inputs
grclm, & ! grad(cos(lambda_m)), in km-1
clmgrp, & ! cos(lambda_m)*grad(phi_a), in km-1
rgrlp, & ! (re + altitude)*grad(lambda')
b ! magnetic field, in nT
real,intent(out) :: & ! scalar output
bmag, & ! magnitude of magnetic field, in nT
sim, & ! sin(I_m) of Richmond reference
si, & ! sin(I)
F, & ! F of Richmond reference
D, & ! D of Richmond reference
W ! W of Richmond reference
real,dimension(3),intent(out) :: & ! 3-component outputs
bhat, & ! unit vector along geomagnetic field direction
d1,d2,d3,e1,e2,e3,f1,f2,f3,g1,g2,g3 ! base vectors of Richmond reference
!
! Local:
integer :: i
real :: rr,simoslm,d1db,d2db
rr = re + hr
simoslm = 2./sqrt(4. - 3.*clm*clm)
sim = simoslm*sin(xlatm*dtr)
bmag = sqrt(b(1)*b(1) + b(2)*b(2) + b(3)*b(3))
d1db = 0.
d2db = 0.
do i=1,3
bhat(i) = b(i)/bmag
d1(i) = rr*clmgrp(i)
d1db = d1db + d1(i)*bhat(i)
d2(i) = rr*simoslm*grclm(i)
d2db = d2db + d2(i)*bhat(i)
enddo
!
! Ensure that d1,d2 are exactly perpendicular to B:
!
do i=1,3
d1(i) = d1(i) - d1db*bhat(i)
d2(i) = d2(i) - d2db*bhat(i)
enddo
e3(1) = d1(2)*d2(3) - d1(3)*d2(2)
e3(2) = d1(3)*d2(1) - d1(1)*d2(3)
e3(3) = d1(1)*d2(2) - d1(2)*d2(1)
D = bhat(1)*e3(1) + bhat(2)*e3(2) + bhat(3)*e3(3)
do i=1,3
d3(i) = bhat(i)/D
e3(i) = bhat(i)*D ! Ensure that e3 lies along bhat.
enddo
e1(1) = d2(2)*d3(3) - d2(3)*d3(2)
e1(2) = d2(3)*d3(1) - d2(1)*d3(3)
e1(3) = d2(1)*d3(2) - d2(2)*d3(1)
e2(1) = d3(2)*d1(3) - d3(3)*d1(2)
e2(2) = d3(3)*d1(1) - d3(1)*d1(3)
e2(3) = d3(1)*d1(2) - d3(2)*d1(1)
W = rr*rr*clm*abs(sim)/(bmag*D)
si = -bhat(3)
f1(1) = rgrlp(2)
f1(2) = -rgrlp(1)
f2(1) = -d1(2)*r3_2
f2(2) = d1(1)*r3_2
F = f1(1)*f2(2) - f1(2)*f2(1)
! Added output 2015 February 17
f1(3) = 0.
f2(3) = 0.
g1(1) = r3_2*d1(1)/F
g1(2) = r3_2*d1(2)/F
g1(3) = r3_2*d1(3)/F
g2(1) = rgrlp(1)/F
g2(2) = rgrlp(2)/F
g2(3) = rgrlp(3)/F
g3(1) = 0.
g3(2) = 0.
g3(3) = F
f3(1) = g1(2)*g2(3) - g1(3)*g2(2)
f3(2) = g1(3)*g2(1) - g1(1)*g2(3)
f3(3) = g1(1)*g2(2) - g1(2)*g2(1)
end subroutine basevec
!-----------------------------------------------------------------------
subroutine apex_sub(dlat,dlon,alt,aht,alat,alon,bmag,xmag,ymag,zmag,vmp)
!
! Calculate apex radius, latitude, longitude; and magnetic field and
! scalar magnetic potential.
!
! INPUTS:
! dlat = Geodetic latitude in degrees
! dlon = Geodetic longitude in degrees
! alt = Altitude in km
!
! RETURNS:
! aht = (Apex height + req)/req, where req = equatorial Earth radius
! aht is analogous to the L value in invariant coordinates.
! alat = Apex latitude in degrees (negative in S magnetic hemisphere)
! alon = Apex longitude (geomagnetic longitude of apex) in degrees
! bmag = geomagnetic field magnitude (nT)
! xmag = geomagnetic field component (nT): north
! ymag = geomagnetic field component (nT): east
! zmag = geomagnetic field component (nT): downward
! vmp = geomagnetic potential (T-m)
!
! Args:
real,intent(inout) :: dlat,dlon,alt
real,intent(out) :: aht,alat,alon,bmag,xmag,ymag,zmag,vmp
!
! Local:
real :: x,y,z,xre,yre,zre
integer :: iflag
vmp = 0.
!
! Last 7 args of linapx are output:
!
call linapx(dlat,dlon,alt, aht,alat,alon,xmag,ymag,zmag,bmag)
xmag = xmag*1.e5
ymag = ymag*1.e5
zmag = zmag*1.e5
bmag = bmag*1.e5
call gd2cart (dlat,dlon,alt,x,y,z)
iflag = 3
xre = x/re ; yre = y/re ; zre = z/re
call feldg(iflag,xre,yre,zre,bx,by,bz,vmp)
end subroutine apex_sub
!-----------------------------------------------------------------------
subroutine linapx(gdlat,glon,alt,aht,alat,alon,xmag,ymag,zmag,fmag)
!
! Transform geographic coordinates to Apex coordinates.
! Based on code originally written in 1973 by Wally Clark at NOAA
! Modified by A.D. Richmond, 1994, to return meaningful values when the
! maximum number of steps has been reached, by treating the field line
! above that high altitude as dipolar.
! Rewritten in Fortran 90, eliminating common blocks, by B. Foster, 2013
! Reference: Stassinopoulos E.G., Mead Gilbert D., X-841-72-17 (1971),
! NASA GSFC, Greenbelt, Maryland
!
! Input Args:
!
real,intent(inout) :: & ! These may be changed by convrt, depending on iflag
gdlat, & ! latitude of starting point (deg)
glon, & ! longitude of starting point (deg)
alt ! height of starting point (km)
!
! Output Args:
!
real,intent(out) :: &
aht, & ! (Apex height+req)/req, where req is equatorial earth radius
alat, & ! Apex latitude (deg)
alon, & ! Apex longitude (deg)
xmag, & ! North component of magnetic field at starting point
ymag, & ! East component of magnetic field at starting point
zmag, & ! Down component of magnetic field at starting point
fmag ! Magnetic field magnitude at starting point
!
! Trace the geomagnetic field line from the given location to find the
! apex of the field line. Before starting iterations to trace along
! the field line: (1) Establish a step size (ds, arc length in km)
! based on the geomagnetic dipole latitude; (2) determine the step
! direction from the sign of the vertical component of the geomagnetic
! field; and (3) convert to geocentric cartesian coordinates. Each
! iteration increments a step count (nstp) and calls itrace to move
! along the the field line until reaching the iteration count limit