/[PAMELA software]/DarthVader/OrbitalInfo/src/igrf_sub.for
ViewVC logotype

Annotation of /DarthVader/OrbitalInfo/src/igrf_sub.for

Parent Directory Parent Directory | Revision Log Revision Log


Revision 1.4 - (hide annotations) (download)
Tue Aug 11 12:58:22 2009 UTC (15 years, 3 months ago) by mocchiut
Branch: MAIN
CVS Tags: v9r00, v9r01
Changes since 1.3: +17 -17 lines
Compilation warnings with gcc >= 4.3 fixed

1 mocchiut 1.1 subroutine igrf_sub(xlat,xlong,year,height,
2     & xl,icode,dip,dec)
3     c----------------------------------------------------------------
4     c INPUT:
5 mocchiut 1.4 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 mocchiut 1.1 c OUTPUT:
11 mocchiut 1.4 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 mocchiut 1.1 c----------------------------------------------------------------
17    
18     REAL LATI,LONGI
19     COMMON/GENER/ UMR,ERA,AQUAD,BQUAD
20     SAVE /GENER/
21     C
22     CALL INITIZE
23 mocchiut 1.4 ibbb=0
24 mocchiut 1.1 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 mocchiut 1.4 DEC=ASIN(BEAST/SQRT(BEAST*BEAST+BNORTH*BNORTH))/UMR
36 mocchiut 1.1 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 mocchiut 1.3 IF((R-RMAX).le.0.) goto 44
312     IF((R-RMAX).gt.0.) goto 45
313 mocchiut 1.1 45 ICODE=2
314     RADIK=RADIK-12.*(R-RMAX)**2
315     44 IF(RADIK+RADIK.LE.ORADIK) GOTO 10
316     TERM=SQRT(RADIK)*FF*((E2*T+E1)*T+E0)/(RQ+ZQ)
317     FI=FI+STP*(OTERM+TERM)
318     ORADIK=RADIK
319     OTERM=TERM
320     STP=R*STEQ
321     Z=Z+STP
322     GOTO4
323     C*****PREDICTOR (FIELD LINE TRACING)
324     5 P(1,N+1)=P(1,N)+STEP12*(23.*P(4,N)-16.*P(4,N-1)+5.*P(4,N-2))
325     P(2,N+1)=P(2,N)+STEP12*(23.*P(5,N)-16.*P(5,N-1)+5.*P(5,N-2))
326     P(3,N+1)=P(3,N)+STEP
327     CALL STOER(P(1,N+1),BQ3,R3)
328     C*****SEARCH FOR LOWEST MAGNETIC FIELD STRENGTH
329     IF(BQ3.LT.BEQU) THEN
330     IEQU=N+1
331     BEQU=BQ3
332     ENDIF
333     3 CONTINUE
334     10 IF(IEQU.lt.2) IEQU=2
335     SP(1)=P(1,IEQU-1)
336     SP(2)=P(2,IEQU-1)
337     SP(3)=P(3,IEQU-1)
338     IF(ORADIK.LT.1E-15)GOTO11
339     FI=FI+STP/0.75*OTERM*ORADIK/(ORADIK-RADIK)
340     C
341     C-- The minimal allowable value of FI was changed from 1E-15 to 1E-12,
342     C-- because 1E-38 is the minimal allowable arg. for ALOG in our envir.
343     C-- D. Bilitza, Nov 87.
344     C
345 mocchiut 1.4 11 FI=0.5*ABS(FI)/SQRT(B0)+1E-12
346 mocchiut 1.1 C
347     C*****COMPUTE L FROM B AND I. SAME AS CARMEL IN INVAR.
348     C
349     C-- Correct dipole moment is used here. D. Bilitza, Nov 87.
350     C
351     DIMOB0=DIMO/B0
352 mocchiut 1.4 arg1=alog(FI)
353     arg2=alog(DIMOB0)
354     c arg = FI*FI*FI/DIMOB0
355     c if(abs(arg).gt.88.0) arg=88.0
356 mocchiut 1.1 XX=3*arg1-arg2
357     IF(XX.GT.23.0) GOTO 776
358     IF(XX.GT.11.7) GOTO 775
359     IF(XX.GT.+3.0) GOTO 774
360     IF(XX.GT.-3.0) GOTO 773
361     IF(XX.GT.-22.) GOTO 772
362 mocchiut 1.3 c 771 GG=3.33338E-1*XX+3.0062102E-1
363     GG=3.33338E-1*XX+3.0062102E-1
364 mocchiut 1.1 GOTO777
365     772 GG=((((((((-8.1537735E-14*XX+8.3232531E-13)*XX+1.0066362E-9)*XX+
366     18.1048663E-8)*XX+3.2916354E-6)*XX+8.2711096E-5)*XX+1.3714667E-3)*
367     2XX+1.5017245E-2)*XX+4.3432642E-1)*XX+6.2337691E-1
368     GOTO777
369     773 GG=((((((((2.6047023E-10*XX+2.3028767E-9)*XX-2.1997983E-8)*XX-
370     15.3977642E-7)*XX-3.3408822E-6)*XX+3.8379917E-5)*XX+1.1784234E-3)*
371     2XX+1.4492441E-2)*XX+4.3352788E-1)*XX+6.228644E-1
372     GOTO777
373     774 GG=((((((((6.3271665E-10*XX-3.958306E-8)*XX+9.9766148E-07)*XX-
374     11.2531932E-5)*XX+7.9451313E-5)*XX-3.2077032E-4)*XX+2.1680398E-3)*
375     2XX+1.2817956E-2)*XX+4.3510529E-1)*XX+6.222355E-1
376     GOTO777
377     775 GG=(((((2.8212095E-8*XX-3.8049276E-6)*XX+2.170224E-4)*XX-6.7310339
378     1E-3)*XX+1.2038224E-1)*XX-1.8461796E-1)*XX+2.0007187E0
379     GOTO777
380     776 GG=XX-3.0460681E0
381     777 FL=EXP(ALOG((1.+EXP(GG))*DIMOB0)/3.0)
382     RETURN
383     C*****APPROXIMATION FOR HIGH VALUES OF L.
384     30 ICODE=3
385     T=-P(3,N-1)/STEP
386     FL=1./(ABS(((C3*T+C2)*T+C1)*T+C0)+1E-15)
387     RETURN
388     END
389     C
390     C
391     SUBROUTINE STOER(P,BQ,R)
392     C*******************************************************************
393     C* SUBROUTINE USED FOR FIELD LINE TRACING IN SHELLG *
394     C* CALLS ENTRY POINT FELDI IN GEOMAGNETIC FIELD SUBROUTINE FELDG *
395     C*******************************************************************
396     DIMENSION P(7),U(3,3)
397     COMMON XI(3),H(144)
398     C*****XM,YM,ZM ARE GEOMAGNETIC CARTESIAN INVERSE CO-ORDINATES
399     ZM=P(3)
400     FLI=P(1)*P(1)+P(2)*P(2)+1E-15
401     R=0.5*(FLI+SQRT(FLI*FLI+(ZM+ZM)**2))
402     RQ=R*R
403     WR=SQRT(R)
404     XM=P(1)*WR
405     YM=P(2)*WR
406     C*****TRANSFORM TO GEOGRAPHIC CO-ORDINATE SYSTEM
407     DATA U/ +0.3511737,-0.9148385,-0.1993679,
408     A +0.9335804,+0.3583680,+0.0000000,
409     B +0.0714471,-0.1861260,+0.9799247/
410     XI(1)=XM*U(1,1)+YM*U(1,2)+ZM*U(1,3)
411     XI(2)=XM*U(2,1)+YM*U(2,2)+ZM*U(2,3)
412     XI(3)=XM*U(3,1) +ZM*U(3,3)
413     C*****COMPUTE DERIVATIVES
414     C CALL FELDI(XI,H)
415     CALL FELDI
416     Q=H(1)/RQ
417     DX=H(3)+H(3)+Q*XI(1)
418     DY=H(4)+H(4)+Q*XI(2)
419     DZ=H(2)+H(2)+Q*XI(3)
420     C*****TRANSFORM BACK TO GEOMAGNETIC CO-ORDINATE SYSTEM
421     DXM=U(1,1)*DX+U(2,1)*DY+U(3,1)*DZ
422     DYM=U(1,2)*DX+U(2,2)*DY
423     DZM=U(1,3)*DX+U(2,3)*DY+U(3,3)*DZ
424     DR=(XM*DXM+YM*DYM+ZM*DZM)/R
425     C*****FORM SLOWLY VARYING EXPRESSIONS
426     P(4)=(WR*DXM-0.5*P(1)*DR)/(R*DZM)
427     P(5)=(WR*DYM-0.5*P(2)*DR)/(R*DZM)
428     DSQ=RQ*(DXM*DXM+DYM*DYM+DZM*DZM)
429     BQ=DSQ*RQ*RQ
430     P(6)=SQRT(DSQ/(RQ+3.*ZM*ZM))
431     P(7)=P(6)*(RQ+ZM*ZM)/(RQ*DZM)
432     RETURN
433     END
434     C
435     C
436     SUBROUTINE FELDG(GLAT,GLON,ALT,BNORTH,BEAST,BDOWN,BABS)
437     C-------------------------------------------------------------------
438     C CALCULATES EARTH MAGNETIC FIELD FROM SPHERICAL HARMONICS MODEL
439     C REF: G. KLUGE, EUROPEAN SPACE OPERATIONS CENTRE, INTERNAL NOTE 61,
440     C 1970.
441     C--------------------------------------------------------------------
442     C CHANGES (D. BILITZA, NOV 87):
443     C - FIELD COEFFICIENTS IN BINARY DATA FILES INSTEAD OF BLOCK DATA
444     C - CALCULATES DIPOL MOMENT
445     C--------------------------------------------------------------------
446     C INPUT: ENTRY POINT FELDG
447     C GLAT GEODETIC LATITUDE IN DEGREES (NORTH)
448     C GLON GEODETIC LONGITUDE IN DEGREES (EAST)
449     C ALT ALTITUDE IN KM ABOVE SEA LEVEL
450     C
451     C ENTRY POINT FELDC
452     C V(3) CARTESIAN COORDINATES IN EARTH RADII (6371.2 KM)
453     C X-AXIS POINTING TO EQUATOR AT 0 LONGITUDE
454     C Y-AXIS POINTING TO EQUATOR AT 90 LONG.
455     C Z-AXIS POINTING TO NORTH POLE
456     C
457     C COMMON BLANK AND ENTRY POINT FELDI ARE NEEDED WHEN USED
458     C IN CONNECTION WITH L-CALCULATION PROGRAM SHELLG.
459     C
460     C COMMON /MODEL/ AND /GENER/
461     C UMR = ATAN(1.0)*4./180. <DEGREE>*UMR=<RADIANT>
462     C ERA EARTH RADIUS FOR NORMALIZATION OF CARTESIAN
463     C COORDINATES (6371.2 KM)
464     C AQUAD, BQUAD SQUARE OF MAJOR AND MINOR HALF AXIS FOR
465     C EARTH ELLIPSOID AS RECOMMENDED BY INTERNATIONAL
466     C ASTRONOMICAL UNION (6378.160, 6356.775 KM).
467     C NMAX MAXIMUM ORDER OF SPHERICAL HARMONICS
468     C TIME YEAR (DECIMAL: 1973.5) FOR WHICH MAGNETIC
469     C FIELD IS TO BE CALCULATED
470     C G(M) NORMALIZED FIELD COEFFICIENTS (SEE FELDCOF)
471     C M=NMAX*(NMAX+2)
472     C------------------------------------------------------------------------
473     C OUTPUT: BABS MAGNETIC FIELD STRENGTH IN GAUSS
474     C BNORTH, BEAST, BDOWN COMPONENTS OF THE FIELD WITH RESPECT
475     C TO THE LOCAL GEODETIC COORDINATE SYSTEM, WITH AXIS
476     C POINTING IN THE TANGENTIAL PLANE TO THE NORTH, EAST
477     C AND DOWNWARD.
478     C-----------------------------------------------------------------------
479     DIMENSION V(3),B(3),G(144)
480     CHARACTER*258 NAME
481     INTEGER NMAX
482     REAL TIME
483     COMMON XI(3),H(144)
484     COMMON/MODEL/ G,NMAX,TIME,NAME
485     SAVE/MODEL/
486     COMMON/GENER/ UMR,ERA,AQUAD,BQUAD
487     SAVE/GENER/
488     C
489     C-- IS RECORDS ENTRY POINT
490     C
491     C*****ENTRY POINT FELDG TO BE USED WITH GEODETIC CO-ORDINATES
492     IS=1
493     RLAT=GLAT*UMR
494     CT=SIN(RLAT)
495     ST=COS(RLAT)
496     D=SQRT(AQUAD-(AQUAD-BQUAD)*CT*CT)
497     RLON=GLON*UMR
498     CP=COS(RLON)
499     SP=SIN(RLON)
500     ZZZ=(ALT+BQUAD/D)*CT/ERA
501     RHO=(ALT+AQUAD/D)*ST/ERA
502     XXX=RHO*CP
503     YYY=RHO*SP
504     GOTO10
505     ENTRY FELDC(V,B)
506     C*****ENTRY POINT FELDC TO BE USED WITH CARTESIAN CO-ORDINATES
507     IS=2
508     XXX=V(1)
509     YYY=V(2)
510     ZZZ=V(3)
511     10 RQ=1./(XXX*XXX+YYY*YYY+ZZZ*ZZZ)
512     XI(1)=XXX*RQ
513     XI(2)=YYY*RQ
514     XI(3)=ZZZ*RQ
515     GOTO20
516     ENTRY FELDI
517     C*****ENTRY POINT FELDI USED FOR L COMPUTATION
518     IS=3
519     20 IHMAX=NMAX*NMAX+1
520     LAST=IHMAX+NMAX+NMAX
521     IMAX=NMAX+NMAX-1
522     DO 8 I=IHMAX,LAST
523     8 H(I)=G(I)
524     DO 6 K=1,3,2
525     I=IMAX
526     IH=IHMAX
527     1 IL=IH-I
528     F=2./FLOAT(I-K+2)
529     X=XI(1)*F
530     Y=XI(2)*F
531     Z=XI(3)*(F+F)
532     I=I-2
533 mocchiut 1.3 c print *,' I ',I
534     IF((I-1).lt.0) goto 5
535     IF((I-1).eq.0) goto 4
536     IF((I-1).gt.0) goto 2
537 mocchiut 1.1 2 DO 3 M=3,I,2
538     H(IL+M+1)=G(IL+M+1)+Z*H(IH+M+1)+X*(H(IH+M+3)-H(IH+M-1))
539     A -Y*(H(IH+M+2)+H(IH+M-2))
540     3 H(IL+M)=G(IL+M)+Z*H(IH+M)+X*(H(IH+M+2)-H(IH+M-2))
541     A +Y*(H(IH+M+3)+H(IH+M-1))
542     4 H(IL+2)=G(IL+2)+Z*H(IH+2)+X*H(IH+4)-Y*(H(IH+3)+H(IH))
543     H(IL+1)=G(IL+1)+Z*H(IH+1)+Y*H(IH+4)+X*(H(IH+3)-H(IH))
544     5 H(IL)=G(IL)+Z*H(IH)+2.*(X*H(IH+1)+Y*H(IH+2))
545     IH=IL
546     IF(I.GE.K)GOTO1
547     6 CONTINUE
548     IF(IS.EQ.3)RETURN
549     S=.5*H(1)+2.*(H(2)*XI(3)+H(3)*XI(1)+H(4)*XI(2))
550     T=(RQ+RQ)*SQRT(RQ)
551     BXXX=T*(H(3)-S*XXX)
552     BYYY=T*(H(4)-S*YYY)
553     BZZZ=T*(H(2)-S*ZZZ)
554     IF(IS.EQ.2)GOTO7
555     BABS=SQRT(BXXX*BXXX+BYYY*BYYY+BZZZ*BZZZ)
556     BEAST=BYYY*CP-BXXX*SP
557     BRHO=BYYY*SP+BXXX*CP
558     BNORTH=BZZZ*ST-BRHO*CT
559     BDOWN=-BZZZ*CT-BRHO*ST
560     RETURN
561     7 B(1)=BXXX
562     B(2)=BYYY
563     B(3)=BZZZ
564     RETURN
565     END
566     C
567     C
568 mocchiut 1.3 SUBROUTINE FELDCOF(YEAR,DIMO)
569 mocchiut 1.1 C------------------------------------------------------------------------
570 mocchiut 1.3 C DETERMINES COEFFICIENTS AND DIPOL MOMENT FROM IGRF MODELS
571     C
572     C INPUT: YEAR DECIMAL YEAR FOR WHICH GEOMAGNETIC FIELD IS TO
573     C BE CALCULATED
574     C OUTPUT: DIMO GEOMAGNETIC DIPOL MOMENT IN GAUSS (NORMALIZED
575     C TO EARTH'S RADIUS) AT THE TIME (YEAR)
576     C D. BILITZA, NSSDC, GSFC, CODE 633, GREENBELT, MD 20771,
577     C (301)286-9536 NOV 1987.
578     C ### updated to IGRF-2000 version -dkb- 5/31/2000
579     C ### updated to IGRF-2005 version -dkb- 3/24/2000
580 mocchiut 1.1 C-----------------------------------------------------------------------
581 mocchiut 1.3 CHARACTER*258 FIL1, FIL2
582     CHARACTER*258 FILMOD
583     C ### FILMOD, DTEMOD arrays +1
584     c DIMENSION GH1(144),GH2(120),GHA(144),FILMOD(14),DTEMOD(14)
585     DIMENSION GH1(144),GH2(120),GHA(144),FILMOD(3),DTEMOD(3)
586     DOUBLE PRECISION X,F0,F
587     INTEGER L1,L2,L3
588     INTEGER NMAX
589     REAL TIME
590     CHARACTER *258 P1,P2,P3
591     COMMON/PPATH/ L1,L2,L3,P1, P2, P3
592     SAVE/PPATH/
593     COMMON/MODEL/ GH1,NMAX,TIME,FIL1
594     SAVE/MODEL/
595     COMMON/GENER/ UMR,ERAD,AQUAD,BQUAD
596     SAVE/GENER/
597     C ### updated to 2005
598     C CHARACTER COEFPATH*80, COEF1*80, COEF2*80, COEF3*80
599    
600     c COEFPATH = 'OrbitalInfo/src/'
601     c COEF1 = 'dgrf00.dat'
602     c COEF2 = 'igrf05.dat'
603     c COEF3 = 'igrf05s.dat'
604     c COEF1 = COEFPATH(1:16)//COEF1
605     c COEF2 = COEFPATH(1:16)//COEF2
606     c COEF3 = COEFPATH(1:16)//COEF3
607     c FILMOD(1) = COEF1
608     c FILMOD(2) = COEF2
609     c FILMOD(3) = COEF3
610     c print *, "qui"
611 mocchiut 1.1 FILMOD(1) = P1(1:L1)
612     FILMOD(2) = P2(1:L2)
613     FILMOD(3) = P3(1:L3)
614 mocchiut 1.3 c print *, "qua"
615     c FILMOD(1) = 'OrbitalInfo/src/dgrf00.dat'
616     c FILMOD(2) = 'OrbitalInfo/src/igrf05.dat'
617     c FILMOD(3) = 'OrbitalInfo/src/igrf05s.dat'
618     c WRITE(*,*) FILMOD(1)
619     c WRITE(*,*) FILMOD(2)
620     c WRITE(*,*) FILMOD(3)
621 mocchiut 1.1 c DATA FILMOD / 'dgrf00.dat', 'igrf05.dat', 'igrf05s.dat'/
622     DATA DTEMOD / 2000., 2005., 2010./
623 mocchiut 1.3 c
624     c DATA FILMOD /'dgrf45.dat', 'dgrf50.dat',
625 mocchiut 1.1 c 1 'dgrf55.dat', 'dgrf60.dat', 'dgrf65.dat',
626     c 2 'dgrf70.dat', 'dgrf75.dat', 'dgrf80.dat',
627     c 3 'dgrf85.dat', 'dgrf90.dat', 'dgrf95.dat',
628     c 4 'dgrf00.dat','igrf05.dat','igrf05s.dat'/
629 mocchiut 1.3 c DATA DTEMOD / 1945., 1950., 1955., 1960., 1965., 1970.,
630 mocchiut 1.1 c 1 1975., 1980., 1985., 1990., 1995., 2000.,2005.,2010./
631 mocchiut 1.3 C
632     C ### numye = numye + 1 ; is number of years represented by IGRF
633     C
634     c NUMYE=13
635     NUMYE=2
636     c print *, "quo"
637    
638     C
639     C IS=0 FOR SCHMIDT NORMALIZATION IS=1 GAUSS NORMALIZATION
640     C IU IS INPUT UNIT NUMBER FOR IGRF COEFFICIENT SETS
641     C
642     IU = 10
643     IS = 0
644     C-- DETERMINE IGRF-YEARS FOR INPUT-YEAR
645     TIME = YEAR
646     IYEA = INT(YEAR/5.)*5
647     c L = (IYEA - 1945)/5 + 1
648     L = (IYEA - 2000)/5 + 1
649     IF(L.LT.1) L=1
650     IF(L.GT.NUMYE) L=NUMYE
651     DTE1 = DTEMOD(L)
652     FIL1 = FILMOD(L)
653     DTE2 = DTEMOD(L+1)
654     FIL2 = FILMOD(L+1)
655     c print *, "que"
656     C-- GET IGRF COEFFICIENTS FOR THE BOUNDARY YEARS
657     CALL GETSHC (IU, FIL1, NMAX1, ERAD, GH1, IER)
658     IF (IER .NE. 0) STOP
659     CALL GETSHC (IU, FIL2, NMAX2, ERAD, GH2, IER)
660     IF (IER .NE. 0) STOP
661     c print *, "quj"
662     C-- DETERMINE IGRF COEFFICIENTS FOR YEAR
663     IF (L .LE. NUMYE-1) THEN
664     CALL INTERSHC (YEAR, DTE1, NMAX1, GH1, DTE2,
665     1 NMAX2, GH2, NMAX, GHA)
666     ELSE
667     CALL EXTRASHC (YEAR, DTE1, NMAX1, GH1, NMAX2,
668     1 GH2, NMAX, GHA)
669     ENDIF
670     c print *, "quw"
671     C-- DETERMINE MAGNETIC DIPOL MOMENT AND COEFFIECIENTS G
672     F0=0.D0
673     DO 1234 J=1,3
674     F = GHA(J) * 1.D-5
675     F0 = F0 + F * F
676     1234 CONTINUE
677     DIMO = DSQRT(F0)
678    
679     GH1(1) = 0.0
680     I=2
681     F0=1.D-5
682     IF(IS.EQ.0) F0=-F0
683     SQRT2=SQRT(2.)
684 mocchiut 1.1
685 mocchiut 1.3 c print *, "quq"
686    
687 mocchiut 1.1 DO 9 N=1,NMAX
688 mocchiut 1.3 X = N
689     F0 = F0 * X * X / (4.D0 * X - 2.D0)
690     IF(IS.EQ.0) F0 = F0 * (2.D0 * X - 1.D0) / X
691     F = F0 * 0.5D0
692     IF(IS.EQ.0) F = F * SQRT2
693     GH1(I) = GHA(I-1) * F0
694     I = I+1
695     DO 9 M=1,N
696     F = F * (X + M) / (X - M + 1.D0)
697     IF(IS.EQ.0) F = F * DSQRT((X - M + 1.D0) / (X + M))
698     GH1(I) = GHA(I-1) * F
699     GH1(I+1) = GHA(I) * F
700     I=I+2
701     9 CONTINUE
702     RETURN
703     END
704     C
705     C
706     SUBROUTINE GETSHC (IU, FSPEC, NMAX, ERAD, GH, IER)
707 mocchiut 1.1
708     C ===============================================================
709     C
710 mocchiut 1.3 C Version 1.01
711     C
712     C Reads spherical harmonic coefficients from the specified
713     C file into an array.
714     C
715     C Input:
716     C IU - Logical unit number
717     C FSPEC - File specification
718     C
719     C Output:
720     C NMAX - Maximum degree and order of model
721     C ERAD - Earth's radius associated with the spherical
722     C harmonic coefficients, in the same units as
723     C elevation
724     C GH - Schmidt quasi-normal internal spherical
725     C harmonic coefficients
726     C IER - Error number: = 0, no error
727     C = -2, records out of order
728     C = FORTRAN run-time error number
729     C
730     C A. Zunde
731     C USGS, MS 964, Box 25046 Federal Center, Denver, CO 80225
732     C
733     C ===============================================================
734    
735     CHARACTER FSPEC*(*), FOUT*258
736     DIMENSION GH(*)
737     C ---------------------------------------------------------------
738     C Open coefficient file. Read past first header record.
739     C Read degree and order of model and Earth's radius.
740     C ---------------------------------------------------------------
741 mocchiut 1.1 WRITE(FOUT,667) FSPEC
742 mocchiut 1.3 c 667 FORMAT('/usr/local/etc/httpd/cgi-bin/natasha/IRI/',A12)
743     667 FORMAT(A258)
744     c print *," gui"
745 mocchiut 1.1 OPEN (IU, FILE=FOUT, STATUS='OLD', IOSTAT=IER, ERR=999)
746 mocchiut 1.3 c print *," gua"
747 mocchiut 1.1 READ (IU, *, IOSTAT=IER, ERR=999)
748 mocchiut 1.3 c print *," gue"
749 mocchiut 1.1 READ (IU, *, IOSTAT=IER, ERR=999) NMAX, ERAD
750 mocchiut 1.3 c print *," guo"
751 mocchiut 1.1 C ---------------------------------------------------------------
752     C Read the coefficient file, arranged as follows:
753     C
754     C N M G H
755     C ----------------------
756     C / 1 0 GH(1) -
757     C / 1 1 GH(2) GH(3)
758     C / 2 0 GH(4) -
759     C / 2 1 GH(5) GH(6)
760     C NMAX*(NMAX+3)/2 / 2 2 GH(7) GH(8)
761     C records \ 3 0 GH(9) -
762     C \ . . . .
763     C \ . . . .
764     C NMAX*(NMAX+2) \ . . . .
765     C elements in GH \ NMAX NMAX . .
766     C
767     C N and M are, respectively, the degree and order of the
768     C coefficient.
769     C ---------------------------------------------------------------
770    
771     I = 0
772 mocchiut 1.3 DO 2211 NN = 1, NMAX
773     DO 2233 MM = 0, NN
774 mocchiut 1.1 READ (IU, *, IOSTAT=IER, ERR=999) N, M, G, H
775     IF (NN .NE. N .OR. MM .NE. M) THEN
776     IER = -2
777     GOTO 999
778     ENDIF
779     I = I + 1
780     GH(I) = G
781     IF (M .NE. 0) THEN
782     I = I + 1
783     GH(I) = H
784     ENDIF
785     2233 CONTINUE
786     2211 CONTINUE
787 mocchiut 1.3 c print *," guj"
788 mocchiut 1.1
789     999 CLOSE (IU)
790 mocchiut 1.3 c print *," guw IER",IER
791     if ( IER .eq. -1 ) IER = 0 ! gfortran 4.1.2 bug workaround... hoping not to create problems with other versions
792 mocchiut 1.1
793     RETURN
794     END
795     C
796     C
797     SUBROUTINE INTERSHC (DATE, DTE1, NMAX1, GH1, DTE2,
798     1 NMAX2, GH2, NMAX, GH)
799    
800     C ===============================================================
801     C
802     C Version 1.01
803     C
804     C Interpolates linearly, in time, between two spherical
805     C harmonic models.
806     C
807     C Input:
808     C DATE - Date of resulting model (in decimal year)
809     C DTE1 - Date of earlier model
810     C NMAX1 - Maximum degree and order of earlier model
811     C GH1 - Schmidt quasi-normal internal spherical
812     C harmonic coefficients of earlier model
813     C DTE2 - Date of later model
814     C NMAX2 - Maximum degree and order of later model
815     C GH2 - Schmidt quasi-normal internal spherical
816     C harmonic coefficients of later model
817     C
818     C Output:
819     C GH - Coefficients of resulting model
820     C NMAX - Maximum degree and order of resulting model
821     C
822     C A. Zunde
823     C USGS, MS 964, Box 25046 Federal Center, Denver, CO 80225
824     C
825     C ===============================================================
826    
827     DIMENSION GH1(*), GH2(*), GH(*)
828    
829     C ---------------------------------------------------------------
830     C The coefficients (GH) of the resulting model, at date
831     C DATE, are computed by linearly interpolating between the
832     C coefficients of the earlier model (GH1), at date DTE1,
833     C and those of the later model (GH2), at date DTE2. If one
834     C model is smaller than the other, the interpolation is
835     C performed with the missing coefficients assumed to be 0.
836     C ---------------------------------------------------------------
837    
838     FACTOR = (DATE - DTE1) / (DTE2 - DTE1)
839    
840     IF (NMAX1 .EQ. NMAX2) THEN
841     K = NMAX1 * (NMAX1 + 2)
842     NMAX = NMAX1
843     ELSE IF (NMAX1 .GT. NMAX2) THEN
844     K = NMAX2 * (NMAX2 + 2)
845     L = NMAX1 * (NMAX1 + 2)
846 mocchiut 1.3 DO 1122 I = K + 1, L
847 mocchiut 1.1 1122 GH(I) = GH1(I) + FACTOR * (-GH1(I))
848     NMAX = NMAX1
849     ELSE
850     K = NMAX1 * (NMAX1 + 2)
851     L = NMAX2 * (NMAX2 + 2)
852 mocchiut 1.3 DO 1133 I = K + 1, L
853 mocchiut 1.1 1133 GH(I) = FACTOR * GH2(I)
854     NMAX = NMAX2
855     ENDIF
856    
857 mocchiut 1.3 DO 1144 I = 1, K
858 mocchiut 1.1 1144 GH(I) = GH1(I) + FACTOR * (GH2(I) - GH1(I))
859    
860     RETURN
861     END
862     C
863     C
864     SUBROUTINE EXTRASHC (DATE, DTE1, NMAX1, GH1, NMAX2,
865     1 GH2, NMAX, GH)
866    
867     C ===============================================================
868     C
869     C Version 1.01
870     C
871     C Extrapolates linearly a spherical harmonic model with a
872     C rate-of-change model.
873     C
874     C Input:
875     C DATE - Date of resulting model (in decimal year)
876     C DTE1 - Date of base model
877     C NMAX1 - Maximum degree and order of base model
878     C GH1 - Schmidt quasi-normal internal spherical
879     C harmonic coefficients of base model
880     C NMAX2 - Maximum degree and order of rate-of-change
881     C model
882     C GH2 - Schmidt quasi-normal internal spherical
883     C harmonic coefficients of rate-of-change model
884     C
885     C Output:
886     C GH - Coefficients of resulting model
887     C NMAX - Maximum degree and order of resulting model
888     C
889     C A. Zunde
890     C USGS, MS 964, Box 25046 Federal Center, Denver, CO 80225
891     C
892     C ===============================================================
893    
894     DIMENSION GH1(*), GH2(*), GH(*)
895    
896     C ---------------------------------------------------------------
897     C The coefficients (GH) of the resulting model, at date
898     C DATE, are computed by linearly extrapolating the coef-
899     C ficients of the base model (GH1), at date DTE1, using
900     C those of the rate-of-change model (GH2), at date DTE2. If
901     C one model is smaller than the other, the extrapolation is
902     C performed with the missing coefficients assumed to be 0.
903     C ---------------------------------------------------------------
904    
905     FACTOR = (DATE - DTE1)
906    
907     IF (NMAX1 .EQ. NMAX2) THEN
908     K = NMAX1 * (NMAX1 + 2)
909     NMAX = NMAX1
910     ELSE IF (NMAX1 .GT. NMAX2) THEN
911     K = NMAX2 * (NMAX2 + 2)
912     L = NMAX1 * (NMAX1 + 2)
913 mocchiut 1.3 DO 1155 I = K + 1, L
914 mocchiut 1.1 1155 GH(I) = GH1(I)
915     NMAX = NMAX1
916     ELSE
917     K = NMAX1 * (NMAX1 + 2)
918     L = NMAX2 * (NMAX2 + 2)
919 mocchiut 1.3 DO 1166 I = K + 1, L
920 mocchiut 1.1 1166 GH(I) = FACTOR * GH2(I)
921     NMAX = NMAX2
922     ENDIF
923    
924 mocchiut 1.3 DO 1177 I = 1, K
925 mocchiut 1.1 1177 GH(I) = GH1(I) + FACTOR * GH2(I)
926    
927     RETURN
928     END
929     C
930     C
931     SUBROUTINE INITIZE(TP1,TL1,TP2,TL2,TP3,TL3)
932     C----------------------------------------------------------------
933     C Initializes the parameters in COMMON/GENER/
934     C
935     C UMR = ATAN(1.0)*4./180. <DEGREE>*UMR=<RADIANT>
936     C ERA EARTH RADIUS FOR NORMALIZATION OF CARTESIAN
937     C COORDINATES (6371.2 KM)
938     C EREQU MAJOR HALF AXIS FOR EARTH ELLIPSOID (6378.160 KM)
939     C ERPOL MINOR HALF AXIS FOR EARTH ELLIPSOID (6356.775 KM)
940     C AQUAD SQUARE OF MAJOR HALF AXIS FOR EARTH ELLIPSOID
941     C BQUAD SQUARE OF MINOR HALF AXIS FOR EARTH ELLIPSOID
942     C
943     C ERA, EREQU and ERPOL as recommended by the INTERNATIONAL
944     C ASTRONOMICAL UNION .
945     C-----------------------------------------------------------------
946     INTEGER TL1,TL2,TL3
947     CHARACTER *258 TP1,TP2,TP3
948     INTEGER L1,L2,L3
949     CHARACTER *258 P1,P2,P3
950     COMMON/PPATH/ L1,L2,L3,P1, P2, P3
951     SAVE/PPATH/
952    
953 mocchiut 1.2 COMMON/GENER/UMR,ERA,AQUAD,BQUAD
954 mocchiut 1.1 SAVE/GENER/
955    
956     L1 = TL1
957     L2 = TL2
958     L3 = TL3
959 mocchiut 1.2 P1 = TP1(1:L1)
960     P2 = TP2(1:L2)
961     P3 = TP3(1:L3)
962 mocchiut 1.1
963     ERA=6371.2
964     EREQU=6378.16
965     ERPOL=6356.775
966     AQUAD=EREQU*EREQU
967     BQUAD=ERPOL*ERPOL
968     UMR=ATAN(1.0)*4./180.
969     RETURN
970     END

  ViewVC Help
Powered by ViewVC 1.1.23