1 |
subroutine igrf_sub(xlat,xlong,year,height, |
2 |
& xl,icode,dip,dec) |
3 |
c---------------------------------------------------------------- |
4 |
c INPUT: |
5 |
c xlat geodatic latitude in degrees |
6 |
c xlong geodatic longitude in degrees |
7 |
c year decimal year (year+month/12.0-0.5 or year+day-of-year/365 |
8 |
c or 366 if leap year) |
9 |
c height height in km |
10 |
c OUTPUT: |
11 |
c xl L value |
12 |
c icode =1 L is correct; =2 L is not correct; |
13 |
c =3 an approximation is used |
14 |
c dip geomagnetic inclination in degrees |
15 |
c dec geomagnetic declination in degress |
16 |
c---------------------------------------------------------------- |
17 |
|
18 |
REAL LATI,LONGI |
19 |
COMMON/GENER/ UMR,ERA,AQUAD,BQUAD |
20 |
SAVE /GENER/ |
21 |
C |
22 |
CALL INITIZE |
23 |
ibbb=0 |
24 |
ALOG2=ALOG(2.) |
25 |
ISTART=1 |
26 |
lati=xlat |
27 |
longi=xlong |
28 |
c |
29 |
C----------------CALCULATE PROFILES----------------------------------- |
30 |
c |
31 |
CALL FELDCOF(YEAR,DIMO) |
32 |
CALL FELDG(LATI,LONGI,HEIGHT,BNORTH,BEAST,BDOWN,BABS) |
33 |
CALL SHELLG(LATI,LONGI,HEIGHT,DIMO,XL,ICODE,BAB1) |
34 |
DIP=ASIN(BDOWN/BABS)/UMR |
35 |
DEC=ASIN(BEAST/SQRT(BEAST*BEAST+BNORTH*BNORTH))/UMR |
36 |
RETURN |
37 |
END |
38 |
c |
39 |
c |
40 |
C SHELLIG.FOR, Version 2.0, January 1992 |
41 |
C |
42 |
C 11/01/91-DKB- SHELLG: lowest starting point for B0 search is 2 |
43 |
C 1/27/92-DKB- Adopted to IGRF-91 coeffcients model |
44 |
C 2/05/92-DKB- Reduce variable-names: INTER(P)SHC,EXTRA(P)SHC,INITI(ALI)ZE |
45 |
C 8/08/95-DKB- Updated to IGRF-45-95; new coeff. DGRF90, IGRF95, IGRF95S |
46 |
C 5/31/00-DKB- Updated to IGRF-45-00; new coeff.: IGRF00, IGRF00s |
47 |
C 3/24/05-DKB- Updated to IGRF-45-10; new coeff.: IGRF05, IGRF05s |
48 |
C 4/25/05-DKB- CALL FELDI and DO 1111 I=1,7 [Alexey Petrov] |
49 |
C |
50 |
C********************************************************************* |
51 |
C SUBROUTINES FINDB0, SHELLG, STOER, FELDG, FELDCOF, GETSHC, * |
52 |
C INTERSHC, EXTRASHC, INITIZE * |
53 |
C********************************************************************* |
54 |
C********************************************************************* |
55 |
C |
56 |
C |
57 |
SUBROUTINE FINDB0(STPS,BDEL,VALUE,BEQU,RR0) |
58 |
C-------------------------------------------------------------------- |
59 |
C FINDS SMALLEST MAGNETIC FIELD STRENGTH ON FIELD LINE |
60 |
C |
61 |
C INPUT: STPS STEP SIZE FOR FIELD LINE TRACING |
62 |
C COMMON/FIDB0/ |
63 |
C SP DIPOLE ORIENTED COORDINATES FORM SHELLG; P(1,*), |
64 |
C P(2,*), P(3,*) CLOSEST TO MAGNETIC EQUATOR |
65 |
C BDEL REQUIRED ACCURACY = [ B(LAST) - BEQU ] / BEQU |
66 |
C B(LAST) IS FIELD STRENGTH BEFORE BEQU |
67 |
C |
68 |
C OUTPUT: VALUE =.FALSE., IF BEQU IS NOT MINIMAL VALUE ON FIELD LINE |
69 |
C BEQU MAGNETIC FIELD STRENGTH AT MAGNETIC EQUATOR |
70 |
C RR0 EQUATORIAL RADIUS NORMALIZED TO EARTH RADIUS |
71 |
C BDEL FINAL ACHIEVED ACCURACY |
72 |
C-------------------------------------------------------------------- |
73 |
DIMENSION P(8,4),SP(3) |
74 |
LOGICAL VALUE |
75 |
COMMON/FIDB0/ SP |
76 |
SAVE /FIDB0/ |
77 |
C |
78 |
STEP=STPS |
79 |
IRUN=0 |
80 |
7777 IRUN=IRUN+1 |
81 |
IF(IRUN.GT.5) THEN |
82 |
VALUE=.FALSE. |
83 |
GOTO 8888 |
84 |
ENDIF |
85 |
C*********************FIRST THREE POINTS |
86 |
P(1,2)=SP(1) |
87 |
P(2,2)=SP(2) |
88 |
P(3,2)=SP(3) |
89 |
STEP=-SIGN(STEP,P(3,2)) |
90 |
CALL STOER(P(1,2),BQ2,R2) |
91 |
P(1,3)=P(1,2)+0.5*STEP*P(4,2) |
92 |
P(2,3)=P(2,2)+0.5*STEP*P(5,2) |
93 |
P(3,3)=P(3,2)+0.5*STEP |
94 |
CALL STOER(P(1,3),BQ3,R3) |
95 |
P(1,1)=P(1,2)-STEP*(2.*P(4,2)-P(4,3)) |
96 |
P(2,1)=P(2,2)-STEP*(2.*P(5,2)-P(5,3)) |
97 |
P(3,1)=P(3,2)-STEP |
98 |
CALL STOER(P(1,1),BQ1,R1) |
99 |
P(1,3)=P(1,2)+STEP*(20.*P(4,3)-3.*P(4,2)+P(4,1))/18. |
100 |
P(2,3)=P(2,2)+STEP*(20.*P(5,3)-3.*P(5,2)+P(5,1))/18. |
101 |
P(3,3)=P(3,2)+STEP |
102 |
CALL STOER(P(1,3),BQ3,R3) |
103 |
C******************INVERT SENSE IF REQUIRED |
104 |
IF(BQ3.LE.BQ1) GOTO 2 |
105 |
STEP=-STEP |
106 |
R3=R1 |
107 |
BQ3=BQ1 |
108 |
DO 1 I=1,5 |
109 |
ZZ=P(I,1) |
110 |
P(I,1)=P(I,3) |
111 |
1 P(I,3)=ZZ |
112 |
C******************INITIALIZATION |
113 |
2 STEP12=STEP/12. |
114 |
VALUE=.TRUE. |
115 |
BMIN=1.E4 |
116 |
BOLD=1.E4 |
117 |
C******************CORRECTOR (FIELD LINE TRACING) |
118 |
N=0 |
119 |
5555 P(1,3)=P(1,2)+STEP12*(5.*P(4,3)+8.*P(4,2)-P(4,1)) |
120 |
N=N+1 |
121 |
P(2,3)=P(2,2)+STEP12*(5.*P(5,3)+8.*P(5,2)-P(5,1)) |
122 |
C******************PREDICTOR (FIELD LINE TRACING) |
123 |
P(1,4)=P(1,3)+STEP12*(23.*P(4,3)-16.*P(4,2)+5.*P(4,1)) |
124 |
P(2,4)=P(2,3)+STEP12*(23.*P(5,3)-16.*P(5,2)+5.*P(5,1)) |
125 |
P(3,4)=P(3,3)+STEP |
126 |
CALL STOER(P(1,4),BQ3,R3) |
127 |
DO 1111 J=1,3 |
128 |
C DO 1111 I=1,8 |
129 |
DO 1111 I=1,7 |
130 |
1111 P(I,J)=P(I,J+1) |
131 |
B=SQRT(BQ3) |
132 |
IF(B.LT.BMIN) BMIN=B |
133 |
IF(B.LE.BOLD) THEN |
134 |
BOLD=B |
135 |
ROLD=1./R3 |
136 |
SP(1)=P(1,4) |
137 |
SP(2)=P(2,4) |
138 |
SP(3)=P(3,4) |
139 |
GOTO 5555 |
140 |
ENDIF |
141 |
IF(BOLD.NE.BMIN) THEN |
142 |
VALUE=.FALSE. |
143 |
ENDIF |
144 |
BDELTA=(B-BOLD)/BOLD |
145 |
IF(BDELTA.GT.BDEL) THEN |
146 |
STEP=STEP/10. |
147 |
GOTO 7777 |
148 |
ENDIF |
149 |
8888 RR0=ROLD |
150 |
BEQU=BOLD |
151 |
BDEL=BDELTA |
152 |
RETURN |
153 |
END |
154 |
C |
155 |
C |
156 |
SUBROUTINE SHELLG(GLAT,GLON,ALT,DIMO,FL,ICODE,B0) |
157 |
C-------------------------------------------------------------------- |
158 |
C CALCULATES L-VALUE FOR SPECIFIED GEODAETIC COORDINATES, ALTITUDE |
159 |
C AND GEMAGNETIC FIELD MODEL. |
160 |
C REF: G. KLUGE, EUROPEAN SPACE OPERATIONS CENTER, INTERNAL NOTE |
161 |
C NO. 67, 1970. |
162 |
C G. KLUGE, COMPUTER PHYSICS COMMUNICATIONS 3, 31-35, 1972 |
163 |
C-------------------------------------------------------------------- |
164 |
C CHANGES (D. BILITZA, NOV 87): |
165 |
C - USING CORRECT DIPOL MOMENT I.E.,DIFFERENT COMMON/MODEL/ |
166 |
C - USING IGRF EARTH MAGNETIC FIELD MODELS FROM 1945 TO 1990 |
167 |
C-------------------------------------------------------------------- |
168 |
C INPUT: ENTRY POINT SHELLG |
169 |
C GLAT GEODETIC LATITUDE IN DEGREES (NORTH) |
170 |
C GLON GEODETIC LONGITUDE IN DEGREES (EAST) |
171 |
C ALT ALTITUDE IN KM ABOVE SEA LEVEL |
172 |
C |
173 |
C ENTRY POINT SHELLC |
174 |
C V(3) CARTESIAN COORDINATES IN EARTH RADII (6371.2 KM) |
175 |
C X-AXIS POINTING TO EQUATOR AT 0 LONGITUDE |
176 |
C Y-AXIS POINTING TO EQUATOR AT 90 LONG. |
177 |
C Z-AXIS POINTING TO NORTH POLE |
178 |
C |
179 |
C DIMO DIPOL MOMENT IN GAUSS (NORMALIZED TO EARTH RADIUS) |
180 |
C |
181 |
C COMMON |
182 |
C X(3) NOT USED |
183 |
C H(144) FIELD MODEL COEFFICIENTS ADJUSTED FOR SHELLG |
184 |
C----------------------------------------------------------------------- |
185 |
C OUTPUT: FL L-VALUE |
186 |
C ICODE =1 NORMAL COMPLETION |
187 |
C =2 UNPHYSICAL CONJUGATE POINT (FL MEANINGLESS) |
188 |
C =3 SHELL PARAMETER GREATER THAN LIMIT UP TO |
189 |
C WHICH ACCURATE CALCULATION IS REQUIRED; |
190 |
C APPROXIMATION IS USED. |
191 |
C B0 MAGNETIC FIELD STRENGTH IN GAUSS |
192 |
C----------------------------------------------------------------------- |
193 |
DIMENSION V(3),U(3,3),P(8,100),SP(3) |
194 |
COMMON X(3),H(144) |
195 |
COMMON/FIDB0/ SP |
196 |
SAVE /FIDB0/ |
197 |
COMMON/GENER/ UMR,ERA,AQUAD,BQUAD |
198 |
SAVE /GENER/ |
199 |
C |
200 |
C-- RMIN, RMAX ARE BOUNDARIES FOR IDENTIFICATION OF ICODE=2 AND 3 |
201 |
C-- STEP IS STEP SIZE FOR FIELD LINE TRACING |
202 |
C-- STEQ IS STEP SIZE FOR INTEGRATION |
203 |
C |
204 |
DATA RMIN,RMAX /0.05,1.01/ |
205 |
DATA STEP,STEQ /0.20,0.03/ |
206 |
BEQU=1.E10 |
207 |
C*****ENTRY POINT SHELLG TO BE USED WITH GEODETIC CO-ORDINATES |
208 |
RLAT=GLAT*UMR |
209 |
CT=SIN(RLAT) |
210 |
ST=COS(RLAT) |
211 |
D=SQRT(AQUAD-(AQUAD-BQUAD)*CT*CT) |
212 |
X(1)=(ALT+AQUAD/D)*ST/ERA |
213 |
X(3)=(ALT+BQUAD/D)*CT/ERA |
214 |
RLON=GLON*UMR |
215 |
X(2)=X(1)*SIN(RLON) |
216 |
X(1)=X(1)*COS(RLON) |
217 |
GOTO9 |
218 |
ENTRY SHELLC(V,FL,B0) |
219 |
C*****ENTRY POINT SHELLC TO BE USED WITH CARTESIAN CO-ORDINATES |
220 |
X(1)=V(1) |
221 |
X(2)=V(2) |
222 |
X(3)=V(3) |
223 |
C*****CONVERT TO DIPOL-ORIENTED CO-ORDINATES |
224 |
DATA U/ +0.3511737,-0.9148385,-0.1993679, |
225 |
A +0.9335804,+0.3583680,+0.0000000, |
226 |
B +0.0714471,-0.1861260,+0.9799247/ |
227 |
9 RQ=1./(X(1)*X(1)+X(2)*X(2)+X(3)*X(3)) |
228 |
R3H=SQRT(RQ*SQRT(RQ)) |
229 |
P(1,2)=(X(1)*U(1,1)+X(2)*U(2,1)+X(3)*U(3,1))*R3H |
230 |
P(2,2)=(X(1)*U(1,2)+X(2)*U(2,2) )*R3H |
231 |
P(3,2)=(X(1)*U(1,3)+X(2)*U(2,3)+X(3)*U(3,3))*RQ |
232 |
C*****FIRST THREE POINTS OF FIELD LINE |
233 |
STEP=-SIGN(STEP,P(3,2)) |
234 |
CALL STOER(P(1,2),BQ2,R2) |
235 |
B0=SQRT(BQ2) |
236 |
P(1,3)=P(1,2)+0.5*STEP*P(4,2) |
237 |
P(2,3)=P(2,2)+0.5*STEP*P(5,2) |
238 |
P(3,3)=P(3,2)+0.5*STEP |
239 |
CALL STOER(P(1,3),BQ3,R3) |
240 |
P(1,1)=P(1,2)-STEP*(2.*P(4,2)-P(4,3)) |
241 |
P(2,1)=P(2,2)-STEP*(2.*P(5,2)-P(5,3)) |
242 |
P(3,1)=P(3,2)-STEP |
243 |
CALL STOER(P(1,1),BQ1,R1) |
244 |
P(1,3)=P(1,2)+STEP*(20.*P(4,3)-3.*P(4,2)+P(4,1))/18. |
245 |
P(2,3)=P(2,2)+STEP*(20.*P(5,3)-3.*P(5,2)+P(5,1))/18. |
246 |
P(3,3)=P(3,2)+STEP |
247 |
CALL STOER(P(1,3),BQ3,R3) |
248 |
C*****INVERT SENSE IF REQUIRED |
249 |
IF(BQ3.LE.BQ1)GOTO2 |
250 |
STEP=-STEP |
251 |
R3=R1 |
252 |
BQ3=BQ1 |
253 |
DO 1 I=1,7 |
254 |
ZZ=P(I,1) |
255 |
P(I,1)=P(I,3) |
256 |
1 P(I,3)=ZZ |
257 |
C*****SEARCH FOR LOWEST MAGNETIC FIELD STRENGTH |
258 |
2 IF(BQ1.LT.BEQU) THEN |
259 |
BEQU=BQ1 |
260 |
IEQU=1 |
261 |
ENDIF |
262 |
IF(BQ2.LT.BEQU) THEN |
263 |
BEQU=BQ2 |
264 |
IEQU=2 |
265 |
ENDIF |
266 |
IF(BQ3.LT.BEQU) THEN |
267 |
BEQU=BQ3 |
268 |
IEQU=3 |
269 |
ENDIF |
270 |
C*****INITIALIZATION OF INTEGRATION LOOPS |
271 |
STEP12=STEP/12. |
272 |
STEP2=STEP+STEP |
273 |
STEQ=SIGN(STEQ,STEP) |
274 |
FI=0. |
275 |
ICODE=1 |
276 |
ORADIK=0. |
277 |
OTERM=0. |
278 |
STP=R2*STEQ |
279 |
Z=P(3,2)+STP |
280 |
STP=STP/0.75 |
281 |
P(8,1)=STEP2*(P(1,1)*P(4,1)+P(2,1)*P(5,1)) |
282 |
P(8,2)=STEP2*(P(1,2)*P(4,2)+P(2,2)*P(5,2)) |
283 |
C*****MAIN LOOP (FIELD LINE TRACING) |
284 |
DO 3 N=3,3333 |
285 |
C*****CORRECTOR (FIELD LINE TRACING) |
286 |
P(1,N)=P(1,N-1)+STEP12*(5.*P(4,N)+8.*P(4,N-1)-P(4,N-2)) |
287 |
P(2,N)=P(2,N-1)+STEP12*(5.*P(5,N)+8.*P(5,N-1)-P(5,N-2)) |
288 |
C*****PREPARE EXPANSION COEFFICIENTS FOR INTERPOLATION |
289 |
C*****OF SLOWLY VARYING QUANTITIES |
290 |
P(8,N)=STEP2*(P(1,N)*P(4,N)+P(2,N)*P(5,N)) |
291 |
C0=P(1,N-1)**2+P(2,N-1)**2 |
292 |
C1=P(8,N-1) |
293 |
C2=(P(8,N)-P(8,N-2))*0.25 |
294 |
C3=(P(8,N)+P(8,N-2)-C1-C1)/6.0 |
295 |
D0=P(6,N-1) |
296 |
D1=(P(6,N)-P(6,N-2))*0.5 |
297 |
D2=(P(6,N)+P(6,N-2)-D0-D0)*0.5 |
298 |
E0=P(7,N-1) |
299 |
E1=(P(7,N)-P(7,N-2))*0.5 |
300 |
E2=(P(7,N)+P(7,N-2)-E0-E0)*0.5 |
301 |
C*****INNER LOOP (FOR QUADRATURE) |
302 |
4 T=(Z-P(3,N-1))/STEP |
303 |
IF(T.GT.1.)GOTO5 |
304 |
HLI=0.5*(((C3*T+C2)*T+C1)*T+C0) |
305 |
ZQ=Z*Z |
306 |
R=HLI+SQRT(HLI*HLI+ZQ) |
307 |
IF(R.LE.RMIN)GOTO30 |
308 |
RQ=R*R |
309 |
FF=SQRT(1.+3.*ZQ/RQ) |
310 |
RADIK=B0-((D2*T+D1)*T+D0)*R*RQ*FF |
311 |
IF(R-RMAX)44,44,45 |
312 |
45 ICODE=2 |
313 |
RADIK=RADIK-12.*(R-RMAX)**2 |
314 |
44 IF(RADIK+RADIK.LE.ORADIK) GOTO 10 |
315 |
TERM=SQRT(RADIK)*FF*((E2*T+E1)*T+E0)/(RQ+ZQ) |
316 |
FI=FI+STP*(OTERM+TERM) |
317 |
ORADIK=RADIK |
318 |
OTERM=TERM |
319 |
STP=R*STEQ |
320 |
Z=Z+STP |
321 |
GOTO4 |
322 |
C*****PREDICTOR (FIELD LINE TRACING) |
323 |
5 P(1,N+1)=P(1,N)+STEP12*(23.*P(4,N)-16.*P(4,N-1)+5.*P(4,N-2)) |
324 |
P(2,N+1)=P(2,N)+STEP12*(23.*P(5,N)-16.*P(5,N-1)+5.*P(5,N-2)) |
325 |
P(3,N+1)=P(3,N)+STEP |
326 |
CALL STOER(P(1,N+1),BQ3,R3) |
327 |
C*****SEARCH FOR LOWEST MAGNETIC FIELD STRENGTH |
328 |
IF(BQ3.LT.BEQU) THEN |
329 |
IEQU=N+1 |
330 |
BEQU=BQ3 |
331 |
ENDIF |
332 |
3 CONTINUE |
333 |
10 IF(IEQU.lt.2) IEQU=2 |
334 |
SP(1)=P(1,IEQU-1) |
335 |
SP(2)=P(2,IEQU-1) |
336 |
SP(3)=P(3,IEQU-1) |
337 |
IF(ORADIK.LT.1E-15)GOTO11 |
338 |
FI=FI+STP/0.75*OTERM*ORADIK/(ORADIK-RADIK) |
339 |
C |
340 |
C-- The minimal allowable value of FI was changed from 1E-15 to 1E-12, |
341 |
C-- because 1E-38 is the minimal allowable arg. for ALOG in our envir. |
342 |
C-- D. Bilitza, Nov 87. |
343 |
C |
344 |
11 FI=0.5*ABS(FI)/SQRT(B0)+1E-12 |
345 |
C |
346 |
C*****COMPUTE L FROM B AND I. SAME AS CARMEL IN INVAR. |
347 |
C |
348 |
C-- Correct dipole moment is used here. D. Bilitza, Nov 87. |
349 |
C |
350 |
DIMOB0=DIMO/B0 |
351 |
arg1=alog(FI) |
352 |
arg2=alog(DIMOB0) |
353 |
c arg = FI*FI*FI/DIMOB0 |
354 |
c if(abs(arg).gt.88.0) arg=88.0 |
355 |
XX=3*arg1-arg2 |
356 |
IF(XX.GT.23.0) GOTO 776 |
357 |
IF(XX.GT.11.7) GOTO 775 |
358 |
IF(XX.GT.+3.0) GOTO 774 |
359 |
IF(XX.GT.-3.0) GOTO 773 |
360 |
IF(XX.GT.-22.) GOTO 772 |
361 |
771 GG=3.33338E-1*XX+3.0062102E-1 |
362 |
GOTO777 |
363 |
772 GG=((((((((-8.1537735E-14*XX+8.3232531E-13)*XX+1.0066362E-9)*XX+ |
364 |
18.1048663E-8)*XX+3.2916354E-6)*XX+8.2711096E-5)*XX+1.3714667E-3)* |
365 |
2XX+1.5017245E-2)*XX+4.3432642E-1)*XX+6.2337691E-1 |
366 |
GOTO777 |
367 |
773 GG=((((((((2.6047023E-10*XX+2.3028767E-9)*XX-2.1997983E-8)*XX- |
368 |
15.3977642E-7)*XX-3.3408822E-6)*XX+3.8379917E-5)*XX+1.1784234E-3)* |
369 |
2XX+1.4492441E-2)*XX+4.3352788E-1)*XX+6.228644E-1 |
370 |
GOTO777 |
371 |
774 GG=((((((((6.3271665E-10*XX-3.958306E-8)*XX+9.9766148E-07)*XX- |
372 |
11.2531932E-5)*XX+7.9451313E-5)*XX-3.2077032E-4)*XX+2.1680398E-3)* |
373 |
2XX+1.2817956E-2)*XX+4.3510529E-1)*XX+6.222355E-1 |
374 |
GOTO777 |
375 |
775 GG=(((((2.8212095E-8*XX-3.8049276E-6)*XX+2.170224E-4)*XX-6.7310339 |
376 |
1E-3)*XX+1.2038224E-1)*XX-1.8461796E-1)*XX+2.0007187E0 |
377 |
GOTO777 |
378 |
776 GG=XX-3.0460681E0 |
379 |
777 FL=EXP(ALOG((1.+EXP(GG))*DIMOB0)/3.0) |
380 |
RETURN |
381 |
C*****APPROXIMATION FOR HIGH VALUES OF L. |
382 |
30 ICODE=3 |
383 |
T=-P(3,N-1)/STEP |
384 |
FL=1./(ABS(((C3*T+C2)*T+C1)*T+C0)+1E-15) |
385 |
RETURN |
386 |
END |
387 |
C |
388 |
C |
389 |
SUBROUTINE STOER(P,BQ,R) |
390 |
C******************************************************************* |
391 |
C* SUBROUTINE USED FOR FIELD LINE TRACING IN SHELLG * |
392 |
C* CALLS ENTRY POINT FELDI IN GEOMAGNETIC FIELD SUBROUTINE FELDG * |
393 |
C******************************************************************* |
394 |
DIMENSION P(7),U(3,3) |
395 |
COMMON XI(3),H(144) |
396 |
C*****XM,YM,ZM ARE GEOMAGNETIC CARTESIAN INVERSE CO-ORDINATES |
397 |
ZM=P(3) |
398 |
FLI=P(1)*P(1)+P(2)*P(2)+1E-15 |
399 |
R=0.5*(FLI+SQRT(FLI*FLI+(ZM+ZM)**2)) |
400 |
RQ=R*R |
401 |
WR=SQRT(R) |
402 |
XM=P(1)*WR |
403 |
YM=P(2)*WR |
404 |
C*****TRANSFORM TO GEOGRAPHIC CO-ORDINATE SYSTEM |
405 |
DATA U/ +0.3511737,-0.9148385,-0.1993679, |
406 |
A +0.9335804,+0.3583680,+0.0000000, |
407 |
B +0.0714471,-0.1861260,+0.9799247/ |
408 |
XI(1)=XM*U(1,1)+YM*U(1,2)+ZM*U(1,3) |
409 |
XI(2)=XM*U(2,1)+YM*U(2,2)+ZM*U(2,3) |
410 |
XI(3)=XM*U(3,1) +ZM*U(3,3) |
411 |
C*****COMPUTE DERIVATIVES |
412 |
C CALL FELDI(XI,H) |
413 |
CALL FELDI |
414 |
Q=H(1)/RQ |
415 |
DX=H(3)+H(3)+Q*XI(1) |
416 |
DY=H(4)+H(4)+Q*XI(2) |
417 |
DZ=H(2)+H(2)+Q*XI(3) |
418 |
C*****TRANSFORM BACK TO GEOMAGNETIC CO-ORDINATE SYSTEM |
419 |
DXM=U(1,1)*DX+U(2,1)*DY+U(3,1)*DZ |
420 |
DYM=U(1,2)*DX+U(2,2)*DY |
421 |
DZM=U(1,3)*DX+U(2,3)*DY+U(3,3)*DZ |
422 |
DR=(XM*DXM+YM*DYM+ZM*DZM)/R |
423 |
C*****FORM SLOWLY VARYING EXPRESSIONS |
424 |
P(4)=(WR*DXM-0.5*P(1)*DR)/(R*DZM) |
425 |
P(5)=(WR*DYM-0.5*P(2)*DR)/(R*DZM) |
426 |
DSQ=RQ*(DXM*DXM+DYM*DYM+DZM*DZM) |
427 |
BQ=DSQ*RQ*RQ |
428 |
P(6)=SQRT(DSQ/(RQ+3.*ZM*ZM)) |
429 |
P(7)=P(6)*(RQ+ZM*ZM)/(RQ*DZM) |
430 |
RETURN |
431 |
END |
432 |
C |
433 |
C |
434 |
SUBROUTINE FELDG(GLAT,GLON,ALT,BNORTH,BEAST,BDOWN,BABS) |
435 |
C------------------------------------------------------------------- |
436 |
C CALCULATES EARTH MAGNETIC FIELD FROM SPHERICAL HARMONICS MODEL |
437 |
C REF: G. KLUGE, EUROPEAN SPACE OPERATIONS CENTRE, INTERNAL NOTE 61, |
438 |
C 1970. |
439 |
C-------------------------------------------------------------------- |
440 |
C CHANGES (D. BILITZA, NOV 87): |
441 |
C - FIELD COEFFICIENTS IN BINARY DATA FILES INSTEAD OF BLOCK DATA |
442 |
C - CALCULATES DIPOL MOMENT |
443 |
C-------------------------------------------------------------------- |
444 |
C INPUT: ENTRY POINT FELDG |
445 |
C GLAT GEODETIC LATITUDE IN DEGREES (NORTH) |
446 |
C GLON GEODETIC LONGITUDE IN DEGREES (EAST) |
447 |
C ALT ALTITUDE IN KM ABOVE SEA LEVEL |
448 |
C |
449 |
C ENTRY POINT FELDC |
450 |
C V(3) CARTESIAN COORDINATES IN EARTH RADII (6371.2 KM) |
451 |
C X-AXIS POINTING TO EQUATOR AT 0 LONGITUDE |
452 |
C Y-AXIS POINTING TO EQUATOR AT 90 LONG. |
453 |
C Z-AXIS POINTING TO NORTH POLE |
454 |
C |
455 |
C COMMON BLANK AND ENTRY POINT FELDI ARE NEEDED WHEN USED |
456 |
C IN CONNECTION WITH L-CALCULATION PROGRAM SHELLG. |
457 |
C |
458 |
C COMMON /MODEL/ AND /GENER/ |
459 |
C UMR = ATAN(1.0)*4./180. <DEGREE>*UMR=<RADIANT> |
460 |
C ERA EARTH RADIUS FOR NORMALIZATION OF CARTESIAN |
461 |
C COORDINATES (6371.2 KM) |
462 |
C AQUAD, BQUAD SQUARE OF MAJOR AND MINOR HALF AXIS FOR |
463 |
C EARTH ELLIPSOID AS RECOMMENDED BY INTERNATIONAL |
464 |
C ASTRONOMICAL UNION (6378.160, 6356.775 KM). |
465 |
C NMAX MAXIMUM ORDER OF SPHERICAL HARMONICS |
466 |
C TIME YEAR (DECIMAL: 1973.5) FOR WHICH MAGNETIC |
467 |
C FIELD IS TO BE CALCULATED |
468 |
C G(M) NORMALIZED FIELD COEFFICIENTS (SEE FELDCOF) |
469 |
C M=NMAX*(NMAX+2) |
470 |
C------------------------------------------------------------------------ |
471 |
C OUTPUT: BABS MAGNETIC FIELD STRENGTH IN GAUSS |
472 |
C BNORTH, BEAST, BDOWN COMPONENTS OF THE FIELD WITH RESPECT |
473 |
C TO THE LOCAL GEODETIC COORDINATE SYSTEM, WITH AXIS |
474 |
C POINTING IN THE TANGENTIAL PLANE TO THE NORTH, EAST |
475 |
C AND DOWNWARD. |
476 |
C----------------------------------------------------------------------- |
477 |
DIMENSION V(3),B(3),G(144) |
478 |
CHARACTER*258 NAME |
479 |
INTEGER NMAX |
480 |
REAL TIME |
481 |
COMMON XI(3),H(144) |
482 |
COMMON/MODEL/ G,NMAX,TIME,NAME |
483 |
SAVE/MODEL/ |
484 |
COMMON/GENER/ UMR,ERA,AQUAD,BQUAD |
485 |
SAVE/GENER/ |
486 |
C |
487 |
C-- IS RECORDS ENTRY POINT |
488 |
C |
489 |
C*****ENTRY POINT FELDG TO BE USED WITH GEODETIC CO-ORDINATES |
490 |
IS=1 |
491 |
RLAT=GLAT*UMR |
492 |
CT=SIN(RLAT) |
493 |
ST=COS(RLAT) |
494 |
D=SQRT(AQUAD-(AQUAD-BQUAD)*CT*CT) |
495 |
RLON=GLON*UMR |
496 |
CP=COS(RLON) |
497 |
SP=SIN(RLON) |
498 |
ZZZ=(ALT+BQUAD/D)*CT/ERA |
499 |
RHO=(ALT+AQUAD/D)*ST/ERA |
500 |
XXX=RHO*CP |
501 |
YYY=RHO*SP |
502 |
GOTO10 |
503 |
ENTRY FELDC(V,B) |
504 |
C*****ENTRY POINT FELDC TO BE USED WITH CARTESIAN CO-ORDINATES |
505 |
IS=2 |
506 |
XXX=V(1) |
507 |
YYY=V(2) |
508 |
ZZZ=V(3) |
509 |
10 RQ=1./(XXX*XXX+YYY*YYY+ZZZ*ZZZ) |
510 |
XI(1)=XXX*RQ |
511 |
XI(2)=YYY*RQ |
512 |
XI(3)=ZZZ*RQ |
513 |
GOTO20 |
514 |
ENTRY FELDI |
515 |
C*****ENTRY POINT FELDI USED FOR L COMPUTATION |
516 |
IS=3 |
517 |
20 IHMAX=NMAX*NMAX+1 |
518 |
LAST=IHMAX+NMAX+NMAX |
519 |
IMAX=NMAX+NMAX-1 |
520 |
DO 8 I=IHMAX,LAST |
521 |
8 H(I)=G(I) |
522 |
DO 6 K=1,3,2 |
523 |
I=IMAX |
524 |
IH=IHMAX |
525 |
1 IL=IH-I |
526 |
F=2./FLOAT(I-K+2) |
527 |
X=XI(1)*F |
528 |
Y=XI(2)*F |
529 |
Z=XI(3)*(F+F) |
530 |
I=I-2 |
531 |
IF(I-1)5,4,2 |
532 |
2 DO 3 M=3,I,2 |
533 |
H(IL+M+1)=G(IL+M+1)+Z*H(IH+M+1)+X*(H(IH+M+3)-H(IH+M-1)) |
534 |
A -Y*(H(IH+M+2)+H(IH+M-2)) |
535 |
3 H(IL+M)=G(IL+M)+Z*H(IH+M)+X*(H(IH+M+2)-H(IH+M-2)) |
536 |
A +Y*(H(IH+M+3)+H(IH+M-1)) |
537 |
4 H(IL+2)=G(IL+2)+Z*H(IH+2)+X*H(IH+4)-Y*(H(IH+3)+H(IH)) |
538 |
H(IL+1)=G(IL+1)+Z*H(IH+1)+Y*H(IH+4)+X*(H(IH+3)-H(IH)) |
539 |
5 H(IL)=G(IL)+Z*H(IH)+2.*(X*H(IH+1)+Y*H(IH+2)) |
540 |
IH=IL |
541 |
IF(I.GE.K)GOTO1 |
542 |
6 CONTINUE |
543 |
IF(IS.EQ.3)RETURN |
544 |
S=.5*H(1)+2.*(H(2)*XI(3)+H(3)*XI(1)+H(4)*XI(2)) |
545 |
T=(RQ+RQ)*SQRT(RQ) |
546 |
BXXX=T*(H(3)-S*XXX) |
547 |
BYYY=T*(H(4)-S*YYY) |
548 |
BZZZ=T*(H(2)-S*ZZZ) |
549 |
IF(IS.EQ.2)GOTO7 |
550 |
BABS=SQRT(BXXX*BXXX+BYYY*BYYY+BZZZ*BZZZ) |
551 |
BEAST=BYYY*CP-BXXX*SP |
552 |
BRHO=BYYY*SP+BXXX*CP |
553 |
BNORTH=BZZZ*ST-BRHO*CT |
554 |
BDOWN=-BZZZ*CT-BRHO*ST |
555 |
RETURN |
556 |
7 B(1)=BXXX |
557 |
B(2)=BYYY |
558 |
B(3)=BZZZ |
559 |
RETURN |
560 |
END |
561 |
C |
562 |
C |
563 |
SUBROUTINE FELDCOF(YEAR,DIMO) |
564 |
C------------------------------------------------------------------------ |
565 |
C DETERMINES COEFFICIENTS AND DIPOL MOMENT FROM IGRF MODELS |
566 |
C |
567 |
C INPUT: YEAR DECIMAL YEAR FOR WHICH GEOMAGNETIC FIELD IS TO |
568 |
C BE CALCULATED |
569 |
C OUTPUT: DIMO GEOMAGNETIC DIPOL MOMENT IN GAUSS (NORMALIZED |
570 |
C TO EARTH'S RADIUS) AT THE TIME (YEAR) |
571 |
C D. BILITZA, NSSDC, GSFC, CODE 633, GREENBELT, MD 20771, |
572 |
C (301)286-9536 NOV 1987. |
573 |
C ### updated to IGRF-2000 version -dkb- 5/31/2000 |
574 |
C ### updated to IGRF-2005 version -dkb- 3/24/2000 |
575 |
C----------------------------------------------------------------------- |
576 |
CHARACTER*258 FIL1, FIL2 |
577 |
CHARACTER*258 FILMOD |
578 |
C ### FILMOD, DTEMOD arrays +1 |
579 |
c DIMENSION GH1(144),GH2(120),GHA(144),FILMOD(14),DTEMOD(14) |
580 |
DIMENSION GH1(144),GH2(120),GHA(144),FILMOD(3),DTEMOD(3) |
581 |
DOUBLE PRECISION X,F0,F |
582 |
INTEGER L1,L2,L3 |
583 |
INTEGER NMAX |
584 |
REAL TIME |
585 |
CHARACTER *258 P1,P2,P3 |
586 |
COMMON/PPATH/ L1,L2,L3,P1, P2, P3 |
587 |
SAVE/PPATH/ |
588 |
COMMON/MODEL/ GH1,NMAX,TIME,FIL1 |
589 |
SAVE/MODEL/ |
590 |
COMMON/GENER/ UMR,ERAD,AQUAD,BQUAD |
591 |
SAVE/GENER/ |
592 |
C ### updated to 2005 |
593 |
C CHARACTER COEFPATH*80, COEF1*80, COEF2*80, COEF3*80 |
594 |
|
595 |
c COEFPATH = 'OrbitalInfo/src/' |
596 |
c COEF1 = 'dgrf00.dat' |
597 |
c COEF2 = 'igrf05.dat' |
598 |
c COEF3 = 'igrf05s.dat' |
599 |
c COEF1 = COEFPATH(1:16)//COEF1 |
600 |
c COEF2 = COEFPATH(1:16)//COEF2 |
601 |
c COEF3 = COEFPATH(1:16)//COEF3 |
602 |
c FILMOD(1) = COEF1 |
603 |
c FILMOD(2) = COEF2 |
604 |
c FILMOD(3) = COEF3 |
605 |
FILMOD(1) = P1(1:L1) |
606 |
FILMOD(2) = P2(1:L2) |
607 |
FILMOD(3) = P3(1:L3) |
608 |
c FILMOD(1) = 'OrbitalInfo/src/dgrf00.dat' |
609 |
c FILMOD(2) = 'OrbitalInfo/src/igrf05.dat' |
610 |
c FILMOD(3) = 'OrbitalInfo/src/igrf05s.dat' |
611 |
c WRITE(*,*) FILMOD(1) |
612 |
c WRITE(*,*) FILMOD(2) |
613 |
c WRITE(*,*) FILMOD(3) |
614 |
c DATA FILMOD / 'dgrf00.dat', 'igrf05.dat', 'igrf05s.dat'/ |
615 |
DATA DTEMOD / 2000., 2005., 2010./ |
616 |
c |
617 |
c DATA FILMOD /'dgrf45.dat', 'dgrf50.dat', |
618 |
c 1 'dgrf55.dat', 'dgrf60.dat', 'dgrf65.dat', |
619 |
c 2 'dgrf70.dat', 'dgrf75.dat', 'dgrf80.dat', |
620 |
c 3 'dgrf85.dat', 'dgrf90.dat', 'dgrf95.dat', |
621 |
c 4 'dgrf00.dat','igrf05.dat','igrf05s.dat'/ |
622 |
c DATA DTEMOD / 1945., 1950., 1955., 1960., 1965., 1970., |
623 |
c 1 1975., 1980., 1985., 1990., 1995., 2000.,2005.,2010./ |
624 |
C |
625 |
C ### numye = numye + 1 ; is number of years represented by IGRF |
626 |
C |
627 |
c NUMYE=13 |
628 |
NUMYE=2 |
629 |
|
630 |
C |
631 |
C IS=0 FOR SCHMIDT NORMALIZATION IS=1 GAUSS NORMALIZATION |
632 |
C IU IS INPUT UNIT NUMBER FOR IGRF COEFFICIENT SETS |
633 |
C |
634 |
IU = 10 |
635 |
IS = 0 |
636 |
C-- DETERMINE IGRF-YEARS FOR INPUT-YEAR |
637 |
TIME = YEAR |
638 |
IYEA = INT(YEAR/5.)*5 |
639 |
c L = (IYEA - 1945)/5 + 1 |
640 |
L = (IYEA - 2000)/5 + 1 |
641 |
IF(L.LT.1) L=1 |
642 |
IF(L.GT.NUMYE) L=NUMYE |
643 |
DTE1 = DTEMOD(L) |
644 |
FIL1 = FILMOD(L) |
645 |
DTE2 = DTEMOD(L+1) |
646 |
FIL2 = FILMOD(L+1) |
647 |
C-- GET IGRF COEFFICIENTS FOR THE BOUNDARY YEARS |
648 |
CALL GETSHC (IU, FIL1, NMAX1, ERAD, GH1, IER) |
649 |
IF (IER .NE. 0) STOP |
650 |
CALL GETSHC (IU, FIL2, NMAX2, ERAD, GH2, IER) |
651 |
IF (IER .NE. 0) STOP |
652 |
C-- DETERMINE IGRF COEFFICIENTS FOR YEAR |
653 |
IF (L .LE. NUMYE-1) THEN |
654 |
CALL INTERSHC (YEAR, DTE1, NMAX1, GH1, DTE2, |
655 |
1 NMAX2, GH2, NMAX, GHA) |
656 |
ELSE |
657 |
CALL EXTRASHC (YEAR, DTE1, NMAX1, GH1, NMAX2, |
658 |
1 GH2, NMAX, GHA) |
659 |
ENDIF |
660 |
C-- DETERMINE MAGNETIC DIPOL MOMENT AND COEFFIECIENTS G |
661 |
F0=0.D0 |
662 |
DO 1234 J=1,3 |
663 |
F = GHA(J) * 1.D-5 |
664 |
F0 = F0 + F * F |
665 |
1234 CONTINUE |
666 |
DIMO = DSQRT(F0) |
667 |
|
668 |
GH1(1) = 0.0 |
669 |
I=2 |
670 |
F0=1.D-5 |
671 |
IF(IS.EQ.0) F0=-F0 |
672 |
SQRT2=SQRT(2.) |
673 |
|
674 |
DO 9 N=1,NMAX |
675 |
X = N |
676 |
F0 = F0 * X * X / (4.D0 * X - 2.D0) |
677 |
IF(IS.EQ.0) F0 = F0 * (2.D0 * X - 1.D0) / X |
678 |
F = F0 * 0.5D0 |
679 |
IF(IS.EQ.0) F = F * SQRT2 |
680 |
GH1(I) = GHA(I-1) * F0 |
681 |
I = I+1 |
682 |
DO 9 M=1,N |
683 |
F = F * (X + M) / (X - M + 1.D0) |
684 |
IF(IS.EQ.0) F = F * DSQRT((X - M + 1.D0) / (X + M)) |
685 |
GH1(I) = GHA(I-1) * F |
686 |
GH1(I+1) = GHA(I) * F |
687 |
I=I+2 |
688 |
9 CONTINUE |
689 |
RETURN |
690 |
END |
691 |
C |
692 |
C |
693 |
SUBROUTINE GETSHC (IU, FSPEC, NMAX, ERAD, GH, IER) |
694 |
|
695 |
C =============================================================== |
696 |
C |
697 |
C Version 1.01 |
698 |
C |
699 |
C Reads spherical harmonic coefficients from the specified |
700 |
C file into an array. |
701 |
C |
702 |
C Input: |
703 |
C IU - Logical unit number |
704 |
C FSPEC - File specification |
705 |
C |
706 |
C Output: |
707 |
C NMAX - Maximum degree and order of model |
708 |
C ERAD - Earth's radius associated with the spherical |
709 |
C harmonic coefficients, in the same units as |
710 |
C elevation |
711 |
C GH - Schmidt quasi-normal internal spherical |
712 |
C harmonic coefficients |
713 |
C IER - Error number: = 0, no error |
714 |
C = -2, records out of order |
715 |
C = FORTRAN run-time error number |
716 |
C |
717 |
C A. Zunde |
718 |
C USGS, MS 964, Box 25046 Federal Center, Denver, CO 80225 |
719 |
C |
720 |
C =============================================================== |
721 |
|
722 |
CHARACTER FSPEC*(*), FOUT*258 |
723 |
DIMENSION GH(*) |
724 |
C --------------------------------------------------------------- |
725 |
C Open coefficient file. Read past first header record. |
726 |
C Read degree and order of model and Earth's radius. |
727 |
C --------------------------------------------------------------- |
728 |
WRITE(FOUT,667) FSPEC |
729 |
c 667 FORMAT('/usr/local/etc/httpd/cgi-bin/natasha/IRI/',A12) |
730 |
667 FORMAT(A258) |
731 |
OPEN (IU, FILE=FOUT, STATUS='OLD', IOSTAT=IER, ERR=999) |
732 |
READ (IU, *, IOSTAT=IER, ERR=999) |
733 |
READ (IU, *, IOSTAT=IER, ERR=999) NMAX, ERAD |
734 |
C --------------------------------------------------------------- |
735 |
C Read the coefficient file, arranged as follows: |
736 |
C |
737 |
C N M G H |
738 |
C ---------------------- |
739 |
C / 1 0 GH(1) - |
740 |
C / 1 1 GH(2) GH(3) |
741 |
C / 2 0 GH(4) - |
742 |
C / 2 1 GH(5) GH(6) |
743 |
C NMAX*(NMAX+3)/2 / 2 2 GH(7) GH(8) |
744 |
C records \ 3 0 GH(9) - |
745 |
C \ . . . . |
746 |
C \ . . . . |
747 |
C NMAX*(NMAX+2) \ . . . . |
748 |
C elements in GH \ NMAX NMAX . . |
749 |
C |
750 |
C N and M are, respectively, the degree and order of the |
751 |
C coefficient. |
752 |
C --------------------------------------------------------------- |
753 |
|
754 |
I = 0 |
755 |
DO 2211 NN = 1, NMAX |
756 |
DO 2233 MM = 0, NN |
757 |
READ (IU, *, IOSTAT=IER, ERR=999) N, M, G, H |
758 |
IF (NN .NE. N .OR. MM .NE. M) THEN |
759 |
IER = -2 |
760 |
GOTO 999 |
761 |
ENDIF |
762 |
I = I + 1 |
763 |
GH(I) = G |
764 |
IF (M .NE. 0) THEN |
765 |
I = I + 1 |
766 |
GH(I) = H |
767 |
ENDIF |
768 |
2233 CONTINUE |
769 |
2211 CONTINUE |
770 |
|
771 |
999 CLOSE (IU) |
772 |
|
773 |
|
774 |
RETURN |
775 |
END |
776 |
C |
777 |
C |
778 |
SUBROUTINE INTERSHC (DATE, DTE1, NMAX1, GH1, DTE2, |
779 |
1 NMAX2, GH2, NMAX, GH) |
780 |
|
781 |
C =============================================================== |
782 |
C |
783 |
C Version 1.01 |
784 |
C |
785 |
C Interpolates linearly, in time, between two spherical |
786 |
C harmonic models. |
787 |
C |
788 |
C Input: |
789 |
C DATE - Date of resulting model (in decimal year) |
790 |
C DTE1 - Date of earlier model |
791 |
C NMAX1 - Maximum degree and order of earlier model |
792 |
C GH1 - Schmidt quasi-normal internal spherical |
793 |
C harmonic coefficients of earlier model |
794 |
C DTE2 - Date of later model |
795 |
C NMAX2 - Maximum degree and order of later model |
796 |
C GH2 - Schmidt quasi-normal internal spherical |
797 |
C harmonic coefficients of later model |
798 |
C |
799 |
C Output: |
800 |
C GH - Coefficients of resulting model |
801 |
C NMAX - Maximum degree and order of resulting model |
802 |
C |
803 |
C A. Zunde |
804 |
C USGS, MS 964, Box 25046 Federal Center, Denver, CO 80225 |
805 |
C |
806 |
C =============================================================== |
807 |
|
808 |
DIMENSION GH1(*), GH2(*), GH(*) |
809 |
|
810 |
C --------------------------------------------------------------- |
811 |
C The coefficients (GH) of the resulting model, at date |
812 |
C DATE, are computed by linearly interpolating between the |
813 |
C coefficients of the earlier model (GH1), at date DTE1, |
814 |
C and those of the later model (GH2), at date DTE2. If one |
815 |
C model is smaller than the other, the interpolation is |
816 |
C performed with the missing coefficients assumed to be 0. |
817 |
C --------------------------------------------------------------- |
818 |
|
819 |
FACTOR = (DATE - DTE1) / (DTE2 - DTE1) |
820 |
|
821 |
IF (NMAX1 .EQ. NMAX2) THEN |
822 |
K = NMAX1 * (NMAX1 + 2) |
823 |
NMAX = NMAX1 |
824 |
ELSE IF (NMAX1 .GT. NMAX2) THEN |
825 |
K = NMAX2 * (NMAX2 + 2) |
826 |
L = NMAX1 * (NMAX1 + 2) |
827 |
DO 1122 I = K + 1, L |
828 |
1122 GH(I) = GH1(I) + FACTOR * (-GH1(I)) |
829 |
NMAX = NMAX1 |
830 |
ELSE |
831 |
K = NMAX1 * (NMAX1 + 2) |
832 |
L = NMAX2 * (NMAX2 + 2) |
833 |
DO 1133 I = K + 1, L |
834 |
1133 GH(I) = FACTOR * GH2(I) |
835 |
NMAX = NMAX2 |
836 |
ENDIF |
837 |
|
838 |
DO 1144 I = 1, K |
839 |
1144 GH(I) = GH1(I) + FACTOR * (GH2(I) - GH1(I)) |
840 |
|
841 |
RETURN |
842 |
END |
843 |
C |
844 |
C |
845 |
SUBROUTINE EXTRASHC (DATE, DTE1, NMAX1, GH1, NMAX2, |
846 |
1 GH2, NMAX, GH) |
847 |
|
848 |
C =============================================================== |
849 |
C |
850 |
C Version 1.01 |
851 |
C |
852 |
C Extrapolates linearly a spherical harmonic model with a |
853 |
C rate-of-change model. |
854 |
C |
855 |
C Input: |
856 |
C DATE - Date of resulting model (in decimal year) |
857 |
C DTE1 - Date of base model |
858 |
C NMAX1 - Maximum degree and order of base model |
859 |
C GH1 - Schmidt quasi-normal internal spherical |
860 |
C harmonic coefficients of base model |
861 |
C NMAX2 - Maximum degree and order of rate-of-change |
862 |
C model |
863 |
C GH2 - Schmidt quasi-normal internal spherical |
864 |
C harmonic coefficients of rate-of-change model |
865 |
C |
866 |
C Output: |
867 |
C GH - Coefficients of resulting model |
868 |
C NMAX - Maximum degree and order of resulting model |
869 |
C |
870 |
C A. Zunde |
871 |
C USGS, MS 964, Box 25046 Federal Center, Denver, CO 80225 |
872 |
C |
873 |
C =============================================================== |
874 |
|
875 |
DIMENSION GH1(*), GH2(*), GH(*) |
876 |
|
877 |
C --------------------------------------------------------------- |
878 |
C The coefficients (GH) of the resulting model, at date |
879 |
C DATE, are computed by linearly extrapolating the coef- |
880 |
C ficients of the base model (GH1), at date DTE1, using |
881 |
C those of the rate-of-change model (GH2), at date DTE2. If |
882 |
C one model is smaller than the other, the extrapolation is |
883 |
C performed with the missing coefficients assumed to be 0. |
884 |
C --------------------------------------------------------------- |
885 |
|
886 |
FACTOR = (DATE - DTE1) |
887 |
|
888 |
IF (NMAX1 .EQ. NMAX2) THEN |
889 |
K = NMAX1 * (NMAX1 + 2) |
890 |
NMAX = NMAX1 |
891 |
ELSE IF (NMAX1 .GT. NMAX2) THEN |
892 |
K = NMAX2 * (NMAX2 + 2) |
893 |
L = NMAX1 * (NMAX1 + 2) |
894 |
DO 1155 I = K + 1, L |
895 |
1155 GH(I) = GH1(I) |
896 |
NMAX = NMAX1 |
897 |
ELSE |
898 |
K = NMAX1 * (NMAX1 + 2) |
899 |
L = NMAX2 * (NMAX2 + 2) |
900 |
DO 1166 I = K + 1, L |
901 |
1166 GH(I) = FACTOR * GH2(I) |
902 |
NMAX = NMAX2 |
903 |
ENDIF |
904 |
|
905 |
DO 1177 I = 1, K |
906 |
1177 GH(I) = GH1(I) + FACTOR * GH2(I) |
907 |
|
908 |
RETURN |
909 |
END |
910 |
C |
911 |
C |
912 |
SUBROUTINE INITIZE(TP1,TL1,TP2,TL2,TP3,TL3) |
913 |
C---------------------------------------------------------------- |
914 |
C Initializes the parameters in COMMON/GENER/ |
915 |
C |
916 |
C UMR = ATAN(1.0)*4./180. <DEGREE>*UMR=<RADIANT> |
917 |
C ERA EARTH RADIUS FOR NORMALIZATION OF CARTESIAN |
918 |
C COORDINATES (6371.2 KM) |
919 |
C EREQU MAJOR HALF AXIS FOR EARTH ELLIPSOID (6378.160 KM) |
920 |
C ERPOL MINOR HALF AXIS FOR EARTH ELLIPSOID (6356.775 KM) |
921 |
C AQUAD SQUARE OF MAJOR HALF AXIS FOR EARTH ELLIPSOID |
922 |
C BQUAD SQUARE OF MINOR HALF AXIS FOR EARTH ELLIPSOID |
923 |
C |
924 |
C ERA, EREQU and ERPOL as recommended by the INTERNATIONAL |
925 |
C ASTRONOMICAL UNION . |
926 |
C----------------------------------------------------------------- |
927 |
INTEGER TL1,TL2,TL3 |
928 |
CHARACTER *258 TP1,TP2,TP3 |
929 |
INTEGER L1,L2,L3 |
930 |
CHARACTER *258 P1,P2,P3 |
931 |
COMMON/PPATH/ L1,L2,L3,P1, P2, P3 |
932 |
SAVE/PPATH/ |
933 |
|
934 |
COMMON/GENER/UMR,ERA,AQUAD,BQUAD |
935 |
SAVE/GENER/ |
936 |
|
937 |
L1 = TL1 |
938 |
L2 = TL2 |
939 |
L3 = TL3 |
940 |
P1 = TP1(1:L1) |
941 |
P2 = TP2(1:L2) |
942 |
P3 = TP3(1:L3) |
943 |
|
944 |
ERA=6371.2 |
945 |
EREQU=6378.16 |
946 |
ERPOL=6356.775 |
947 |
AQUAD=EREQU*EREQU |
948 |
BQUAD=ERPOL*ERPOL |
949 |
UMR=ATAN(1.0)*4./180. |
950 |
RETURN |
951 |
END |