/[MITgcm]/MITgcm_contrib/jscott/igsm/src/daily_new.F
ViewVC logotype

Annotation of /MITgcm_contrib/jscott/igsm/src/daily_new.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.3 - (hide annotations) (download)
Mon Apr 23 21:20:18 2007 UTC (18 years, 3 months ago) by jscott
Branch: MAIN
Changes since 1.2: +40 -6 lines
bring igsm atmos code up to date

1 jscott 1.3 c source sokolov users 16842 Mar 29 17:41 daily_new.F
2 jscott 1.1 #include "ctrparam.h"
3    
4     ! ==========================================================
5     !
6     ! MD2G04.F: Lots of utility functions.
7     !
8     ! ----------------------------------------------------------
9     !
10     ! Revision History:
11     !
12     ! When Who What
13     ! ---- ---------- -------
14     ! 073100 Chien Wang repack based on CliChem3 & M24x11,
15     ! and add cpp.
16     !
17     ! ==========================================================
18    
19    
20     SUBROUTINE DAILY_NEW 1001.
21     C**** 1002.
22     C**** THIS SUBROUTINE PERFORMS THOSE FUNCTIONS OF THE PROGRAM WHICH 1003.
23     C**** TAKE PLACE AT THE BEGINNING OF A NEW DAY. 1004.
24     C**** 1005.
25    
26     #include "BD2G04.COM" 1006.
27    
28     COMMON/SPEC2/KM,KINC,COEK,C3LAND(IO0,JM0),C3OICE(IO0,JM0) 1006.1
29     * ,C3LICE(IO0,JM0),WMGE(IO0,JM0),TSSFC(IM0,JM0,4) 1006.2
30     COMMON U,V,T,P,Q 1007.
31     COMMON/WORK2/Z1OOLD(IO0,JM0),XO(IO0,JM0,3),XZO(IO0,JM0) 1008.
32     COMMON/OLDZO/ZMLOLD(IO0,JM0)
33     DIMENSION AMONTH(12),JDOFM(13) 1009.
34     CHARACTER*4 AMONTH 1009.1
35     DIMENSION XA(1,JM0),XB(1,JM0),OI(IO0,JM0),XOI(IO0,JM0) 1009.5
36     dimension sst1(JM0,3),sst2(JM0,3),dsst(JM0,3),intem(3),
37     & sstmin(12,2)
38     & ,miceo(JM0)
39     common/qfl/QFLUX(JM0,0:13),ZOAV(JM0),QFLUXT(JM0)
40 jscott 1.3 ! common/TSUR/TSURFC(JM0,0:13),TSURFT(JM0),TSURFD(JM0),DTSURF(JM0)
41     #include "TSRF.COM"
42 jscott 1.1 common/fixcld/cldssm(JM0,LM0,0:13),cldmcm(JM0,LM0,0:13),
43     & CLDSST(JM0,LM0),
44     & CLDMCT(JM0,LM0)
45     common/surps/srps(JM0+3),nsrps
46     #if ( defined OCEAN_3D || defined ML_2D)
47 jscott 1.2 #include "AGRID.h"
48 jscott 1.1 #endif
49 jscott 1.3 #if ( defined CLM )
50     #include "CLM.h"
51     #endif
52 jscott 1.1 LOGICAL HPRNT
53     common/conprn/HPRNT,JPR,LPR
54     data ifirst /1/
55     data intem /1,4,5/
56     data sstmin /-1.56,-1.56,-0.75,6*0.0,2*-0.75,-1.56,
57     * 3*0.0,2*-0.75,3*-1.56,-0.75,3*0.0/
58     DATA AMONTH/'JAN','FEB','MAR','APR','MAY','JUNE','JULY','AUG', 1010.
59     * 'SEP','OCT','NOV','DEC'/ 1011.
60     DATA JDOFM/0,31,59,90,120,151,181,212,243,273,304,334,365/ 1012.
61     DATA JDPERY/365/,JMPERY/12/,EDPERY/365./,Z1I/.1/,RHOI/916.6/ 1013.
62     C**** ORBITAL PARAMETERS FOR EARTH FOR YEAR 2000 A.D. 1014.
63     DATA SOLS/173./,APHEL/186./,OBLIQ/23.44/,ECCN/.0167/ 1015.
64     DATA OMEGT/282.9/
65     C**** 1016.
66     C**** THE GLOBAL MEAN PRESSURE IS KEPT CONSTANT AT PSF MILLIBARS 1017.
67     C**** 1018.
68     C**** CALCULATE THE CURRENT GLOBAL MEAN PRESSURE 1019.
69     100 SMASS=0. 1020.
70     c print *,' from Daily KOCEAN=',KOCEAN
71     nsrps=nsrps+1
72     DO 120 J=1,JM 1021.
73     SPRESS=0. 1022.
74     DO 110 I=1,IM 1023.
75     110 SPRESS=SPRESS+P(I,J) 1024.
76     srps(J)=srps(J)+P(1,J)
77     SMASS=SMASS+SPRESS*DXYP(J) 1025.
78     if(J.EQ.JM/2)PBARSH=SMASS
79     120 continue
80     PBAR=SMASS/AREAG+PTOP 1026.
81     PBARNH=2.*(SMASS-PBARSH)/AREAG
82     PBARSH=2.*PBARSH/AREAG
83     srps(JM+1)=srps(JM+1)+PBARSH
84     srps(JM+2)=srps(JM+2)+PBARNH
85     srps(JM+3)=srps(JM+3)+PBAR-PTOP
86     #if ( defined OCEAN_3D)
87     Cjrs do j=1,jm
88     Cjrs surfpr(j)=surfpr(j)+P(1,J)
89     Cjrs enddo
90     #endif
91     C**** CORRECT PRESSURE FIELD FOR ANY LOSS OF MASS BY TRUNCATION ERROR 1027.
92     130 DELTAP=PSF-PBAR 1028.
93     if(DELTAP.gt.1.)then
94     print *,' from Daily DELTAP=',DELTAP
95     print *,' PBAR=',PBAR,' PBARNH=',PBARNH,' PBARSH=',PBARSH
96     endif
97     c GO TO 1140
98     DO 140 J=1,JM 1029.
99     DO 140 I=1,IM 1030.
100     140 P(I,J)=P(I,J)+DELTAP 1031.
101     DOPK=1. 1032.
102     1140 continue
103     C WRITE (6,901) DELTAP 1033.
104     C**** 1034.
105     C**** CALCULATE THE DAILY CALENDAR 1035.
106     C**** 1036.
107     200 JYEAR=IYEAR+(IDAY-1)/JDPERY 1037.
108     JDAY=IDAY-(JYEAR-IYEAR)*JDPERY 1038.
109     DO 210 MONTH=1,JMPERY 1039.
110     IF(JDAY.LE.JDOFM(MONTH+1)) GO TO 220 1040.
111     210 CONTINUE 1041.
112     220 JDATE=JDAY-JDOFM(MONTH) 1042.
113     JMONTH=AMONTH(MONTH) 1043.
114     C**** CALCULATE SOLAR ANGLES AND ORBIT POSITION 1044.
115     if(ifirst.eq.1.or.HPRNT)then
116     print *,' DAILY_ATM IDAY=',IDAY,' IYEAR=',IYEAR
117     print *,' JYEAR=',JYEAR,' JDAY=',JDAY
118     print *,' JDATE=',JDATE,' JMONTH=',JMONTH
119     print *,'OBLIQ=',OBLIQ
120     ifirst=0
121     endif
122     JDSAVE=JDAY 1044.5
123     JDATES=JDATE 1044.51
124     MONSAV=MONTH 1044.52
125     c JDAY=197 1044.53
126     c JDATE=16 1044.54
127     c MONTH=7 1044.55
128     ! RSDIST=(1.+ECCN*COS(TWOPI*(JDAY-APHEL)/EDPERY))**2 1045.
129     ! DEC=COS(TWOPI*(JDAY-SOLS)/EDPERY)*OBLIQ*TWOPI/360. 1046.
130     ! SIND=SIN(DEC) 1047.
131     ! COSD=COS(DEC) 1048.
132     ! 03/03/06
133     ! Fixed calculation of incoming solar radiation
134     CALL ORBIT (OBLIQ,ECCN,OMEGT,JDAY-0.5,RSDIST,SIND,COSD,LAMBDA)
135     if(JDATE.le.16)then
136     do 7231 j=1,JM
137     do 7231 L=1,LM
138     CLDSST(j,L)=((16-JDATE)*cldssm(j,L,MONTH-1)+
139     * (JDATE+15)*cldssm(j,L,MONTH))/31.
140     CLDMCT(j,L)=((16-JDATE)*cldmcm(j,L,MONTH-1)+
141     * (JDATE+15)*cldmcm(j,L,MONTH))/31.
142     7231 continue
143     else
144     do 7241 j=1,JM
145     do 7241 L=1,LM
146     CLDSST(j,L)=((JDATE-16)*cldssm(j,L,MONTH+1)+
147     * (31-JDATE+16)*cldssm(j,L,MONTH))/31.
148     CLDMCT(j,L)=((JDATE-16)*cldmcm(j,L,MONTH+1)+
149     * (31-JDATE+16)*cldmcm(j,L,MONTH))/31.
150     7241 continue
151     endif
152     #if (defined OCEAN_3D || defined ML_2D)
153     if(JDATE.le.16)then
154     do 723 j=1,JM
155     TSURFT(j)=((16-JDATE)*TSURFC(j,MONTH-1)+
156     * (JDATE+15)*TSURFC(j,MONTH))/31.
157 jscott 1.3 TLANDT(j)=((16-JDATE)*TLANDC(j,MONTH-1)+
158     * (JDATE+15)*TLANDC(j,MONTH))/31.
159 jscott 1.1 723 continue
160     else
161     do 724 j=1,JM
162     TSURFT(j)=((JDATE-16)*TSURFC(j,MONTH+1)+
163     * (31-JDATE+16)*TSURFC(j,MONTH))/31.
164 jscott 1.3 TLANDT(j)=((JDATE-16)*TLANDC(j,MONTH+1)+
165     * (31-JDATE+16)*TLANDC(j,MONTH))/31.
166 jscott 1.1 724 continue
167     endif
168 jscott 1.3 ! print *,'From daily_new TSURFD'
169     ! print *,TSURFD
170     ! print *,'TSURFT'
171     ! print *,TSURFT
172     ! print *,'From daily_new TLANDD'
173     ! print *,TLANDD
174     ! print *,'TLANDT'
175     ! print *,TLANDT
176     do 725 j=1,JM
177     DT2MGL(j)=TSURFD(j)-TSURFT(j)
178     DT2MLD(j)=TLANDD(j)-TLANDT(j)
179     TSURFD(j)=0.
180     TLANDD(j)=0.
181     725 continue
182     #if ( defined CLM )
183     DT2MLAND=0.
184     AREAL=0.
185     do j=1,jm
186     DT2MLAND=DT2MLAND+DT2MLD(J)*DXYP(j)*FDATA(1,j,2)
187     AREAL=AREAL+DXYP(j)*FDATA(1,j,2)
188     end do !j
189     DT2MLAND=DT2MLAND/AREAL
190     ! print *,'DT2MLD'
191     ! print *,DT2MLD
192     if(JDATE.eq.1)then
193     print *,'JDATE=',JDATE,' DT2MLAND=',DT2MLAND
194     endif
195     ! print *,'AREAL=',AREAL
196     #endif
197    
198 jscott 1.1 #endif
199    
200     RETURN 1108.5
201     C**** 1109.
202     ENTRY DAILY_NEW0 1110.
203     c IF(TAU.GT.TAUI+DT/7200.) GO TO 200 1111.
204     c GO TO 100 1112.
205     go to 200
206     C***** 1113.
207     901 FORMAT ('0PRESSURE ADDED IN GMP IS',F10.6/) 1114.
208     902 FORMAT ('0MEAN SURFACE PRESSURE OF THE ATMOSPHERE IS',F10.4) 1115.
209     910 FORMAT('1',33A4/) 1116.
210     915 FORMAT (47X,'DAY',I5,', HR',I3,' (',I2,A5,I5,')',F8.1) 1117.
211     920 FORMAT('1') 1118.
212     END 1119.
213     SUBROUTINE ORBIT (OBLIQ,ECCN,OMEGT,DAY,SDIST,SIND,COSD,LAMBDA) 8201.
214     C**** 8202.
215     C**** ORBIT receives the orbital parameters and time of year, and 8203.
216     C**** returns the distance from the sun and its declination angle. 8204.
217     C**** The reference for the following caculations is: V.M.Blanco 8205.
218     C**** and S.W.McCuskey, 1961, "Basic Physics of the Solar System", 8206.
219     C**** pages 135 - 151. 8207.
220     C**** 8208.
221     C**** Program authors: Gary L. Russell and Robert J. Suozzo, 12/13/85 8209.
222     C**** 8210.
223     C**** All computations are in double-precision; 8211.
224     C**** but the arguments are single-precision. 8212.
225     C**** Input: OBLIQ = latitude of tropics in degrees 8213.
226     C**** ECCEN = eccentricity of the orbital ellipse 8214.
227     C**** OMEGT = angle from vernal equinox to perihelion in degrees 8215.
228     C**** DAY = day of the year in days; 0 = Jan 1, hour 0 8216.
229     C**** 8217.
230     C**** Constants: EDAYPY = Earth days per year = 365 8218.
231     C**** VERQNX = occurence of vernal equinox = day 79 = Mar 21 8219.
232     C**** 8220.
233     C**** Intermediate quantities: 8221.
234     C**** PERIHE = perihelion during the year in temporal radians 8222.
235     C**** MA = mean anomaly in temporal radians = 2J DAY/365 - PERIHE8223.
236     C**** EA = eccentric anomaly in radians 8224.
237     C**** TA = true anomaly in radians 8225.
238     C**** BSEMI = semi minor axis in units of the semi major axis 8226.
239     C**** GREENW = longitude of Greenwich in the Earth's reference frame 8227.
240     C**** 8228.
241     C**** Output: DIST = distance to the sun in units of the semi major axis8229.
242     C**** SDIST = square of DIST 8229.5
243     C**** SIND = sine of the declination angle 8230.
244     C**** COSD = cosine of the declination angle 8231.
245     C**** LAMBDA = sun longitude in Earth's rotating reference frame 8232.
246     C**** 8233.
247     IMPLICIT REAL*8 (A-H,O-Z) 8234.
248     REAL*8 MA 8235.
249     C REAL*4 SIND,COSD,SDIST,LAMBDA,OBLIQ,ECCN,OMEGT,DAY 8236.
250     C**** 8237.
251     PI = 3.14159265358979D0 8238.
252     EDAYPY = 365. 8239.
253     VERQNX = 79. 8240.
254     OMEGA=OMEGT*(PI/180.D0) 8241.
255     DOBLIQ=OBLIQ*(PI/180.D0) 8242.
256     ECCEN=ECCN 8243.
257     C**** 8244.
258     C**** Determine time of perihelion using Kepler's equation: 8245.
259     C**** PERIHE-VERQNX = OMEGA - ECCEN sin(OMEGA) 8246.
260     C**** 8247.
261     PERIHE = OMEGA-ECCEN*SIN(OMEGA)+VERQNX*2.*PI/365. 8248.
262     C PERIHE = DMOD(PERIHE,2.*PI) 8249.
263     MA = 2.*PI*DAY/365.-PERIHE 8250.
264     MA = DMOD(MA,2.*PI) 8251.
265     C**** 8252.
266     C**** Numerically solve Kepler's equation: MA = EA - ECCEN sin(EA) 8253.
267     C**** 8254.
268     EA = MA+ECCEN*(SIN(MA)+ECCEN*SIN(2.*MA)/2.) 8255.
269     110 DEA = (MA-EA+ECCEN*SIN(MA))/(1.-ECCEN*COS(EA)) 8256.
270     EA = EA+DEA 8257.
271     IF (DABS(DEA).GT.1.D-8) GO TO 110 8258.
272     C**** 8259.
273     C**** Calculate the distance to the sun and the true anomaly 8260.
274     C**** 8261.
275     BSEMI = DSQRT(1.-ECCEN*ECCEN) 8262.
276     COSEA = COS(EA) 8263.
277     SINEA = SIN(EA) 8264.
278     SDIST = (1.-ECCEN*COSEA)*(1.-ECCEN*COSEA) 8265.
279     TA = DATAN2(SINEA*BSEMI,COSEA-ECCEN) 8266.
280     C**** 8267.
281     C**** Change the reference frame to be the Earth's equatorial plane 8268.
282     C**** with the Earth at the center and the positive x axis parallel to 8269.
283     C**** the ray from the sun to the Earth were it at vernal equinox. 8270.
284     C**** The distance from the current Earth to that ray (or x axis) is: 8271.
285     C**** DIST sin(TA+OMEGA). The sun is located at: 8272.
286     C**** 8273.
287     C**** SUN = (-DIST cos(TA+OMEGA), 8274.
288     C**** -DIST sin(TA+OMEGA) cos(OBLIQ), 8275.
289     C**** DIST sin(TA+OMEGA) sin(OBLIQ)) 8276.
290     C**** SIND = sin(TA+OMEGA) sin(OBLIQ) 8277.
291     C**** COSD = sqrt(1-SIND**2) 8278.
292     C**** LAMBDA = atan[tan(TA+OMEGA) cos(OBLIQ)] - GREENW 8279.
293     C**** GREENW = 2*3.14159 DAY (EDAYPY-1)/EDAYPY 8280.
294     C**** 8281.
295     SINDD = SIN(TA+OMEGA)*SIN(DOBLIQ) 8282.
296     COSD = DSQRT(1.-SINDD*SINDD) 8283.
297     SIND = SINDD 8284.
298     C GREENW = 2.*PI*(DAY-VERQNX)*(EDAYPY+1.)/EDAYPY 8285.
299     C SUNX = -COS(TA+OMEGA) 8286.
300     C SUNY = -SIN(TA+OMEGA)*COS(DOBLIQ) 8287.
301     C LAMBDA = DATAN2(SUNY,SUNX)-GREENW 8288.
302     C LAMBDA = DMOD(LAMBDA,2.*PI) 8289.
303     C**** 8290.
304     RETURN 8291.
305     END 8292.

  ViewVC Help
Powered by ViewVC 1.1.22