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