/[MITgcm]/MITgcm/model/src/ini_p_ground.F
ViewVC logotype

Annotation of /MITgcm/model/src/ini_p_ground.F

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


Revision 1.7 - (hide annotations) (download)
Thu Dec 15 01:11:52 2005 UTC (18 years, 6 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint58b_post, checkpoint62c, checkpoint59, checkpoint58, checkpoint58f_post, checkpoint58d_post, checkpoint58a_post, checkpoint57z_post, checkpoint62a, checkpoint58y_post, checkpoint58t_post, checkpoint58m_post, checkpoint60, checkpoint61, checkpoint62, checkpoint58w_post, checkpoint58o_post, checkpoint58p_post, checkpoint58q_post, checkpoint58e_post, mitgcm_mapl_00, checkpoint58r_post, checkpoint58n_post, checkpoint59q, checkpoint59p, checkpoint59r, checkpoint59e, checkpoint59d, checkpoint59g, checkpoint59f, checkpoint59a, checkpoint59c, checkpoint59b, checkpoint59m, checkpoint59l, checkpoint59o, checkpoint59n, checkpoint59i, checkpoint59h, checkpoint59k, checkpoint58k_post, checkpoint62b, checkpoint58v_post, checkpoint58l_post, checkpoint61f, checkpoint58g_post, checkpoint58x_post, checkpoint61n, checkpoint59j, checkpoint58h_post, checkpoint58j_post, checkpoint61q, checkpoint61e, checkpoint58i_post, checkpoint58c_post, checkpoint58u_post, checkpoint58s_post, checkpoint61g, checkpoint61d, checkpoint61b, checkpoint61c, checkpoint61a, checkpoint61o, checkpoint61l, checkpoint61m, checkpoint61j, checkpoint61k, checkpoint61h, checkpoint61i, checkpoint61v, checkpoint61w, checkpoint61t, checkpoint61u, checkpoint61r, checkpoint61s, checkpoint61p, checkpoint61z, checkpoint61x, checkpoint61y
Changes since 1.6: +2 -2 lines
remove unused variable (reduce number of warnings).

1 jmc 1.7 C $Header: /u/gcmpack/MITgcm/model/src/ini_p_ground.F,v 1.6 2005/12/05 14:37:41 jmc Exp $
2 jmc 1.1 C $Name: $
3    
4     #include "CPP_OPTIONS.h"
5 jmc 1.4 #undef CHECK_ANALYLIC_THETA
6 jmc 1.1
7 cnh 1.2 CBOP
8     C !ROUTINE: INI_P_GROUND
9     C !INTERFACE:
10 jmc 1.4 SUBROUTINE INI_P_GROUND(selectMode,
11     & Hfld, Pfld,
12 jmc 1.1 I myThid )
13 cnh 1.2 C !DESCRIPTION: \bv
14     C *==========================================================*
15     C | SUBROUTINE INI_P_GROUND
16     C | o Convert Topography [m] to (reference) Surface Pressure
17     C | according to tRef profile,
18     C | using same discretisation as in calc_phi_hyd
19     C |
20     C *==========================================================*
21     C \ev
22    
23     C !USES:
24 jmc 1.1 IMPLICIT NONE
25     C == Global variables ==
26     #include "SIZE.h"
27     #include "GRID.h"
28     #include "EEPARAMS.h"
29     #include "PARAMS.h"
30 jmc 1.4 #include "SURFACE.h"
31 jmc 1.1
32 cnh 1.2 C !INPUT/OUTPUT PARAMETERS:
33 jmc 1.1 C == Routine arguments ==
34 jmc 1.5 C selectMode :: > 0 = find Pfld from Hfld ; < 0 = compute Hfld from Pfld
35     C selectFindRoSurf = 0 : use Theta_Ref profile
36     C selectFindRoSurf = 1 : use analytic fct Theta(Lat,P)
37 jmc 1.4 C Hfld (input/outp) :: Ground elevation [m]
38     C Pfld (outp/input) :: reference Pressure at the ground [Pa]
39     INTEGER selectMode
40 jmc 1.1 _RS Hfld(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
41     _RS Pfld(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
42     INTEGER myThid
43    
44 cnh 1.2 C !LOCAL VARIABLES:
45 jmc 1.1 C == Local variables ==
46 jmc 1.4 C msgBuf :: Informational/error meesage buffer
47     C-
48     C For an accurate definition of the reference surface pressure,
49     C define a High vertical resolution (in P):
50     C nLevHvR :: Number of P-level used for High vertical Resolution (HvR)
51     C plowHvR :: Lower bound of the High vertical Resolution
52     C dpHvR :: normalized pressure increment (HvR)
53     C pLevHvR :: normalized P-level of the High vertical Resolution
54     C pMidHvR :: normalized mid point level (HvR)
55     C thetaHvR :: potential temperature at mid point level (HvR)
56     C PiHvR :: Exner function at P-level
57     C dPiHvR :: Exner function difference between 2 P-levels
58     C recip_kappa :: 1/kappa = Cp/R
59     C PiLoc, zLoc, dzLoc, yLatLoc, phiLoc :: hold on temporary values
60     C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
61     CHARACTER*(MAX_LEN_MBUF) msgBuf
62 jmc 1.1 INTEGER bi,bj,i,j,K, Ks
63 jmc 1.7 _RL Po_surf
64 jmc 1.6 _RL hRef(2*Nr+1), rHalf(2*Nr+1)
65 jmc 1.5 LOGICAL findPoSurf
66 jmc 1.4
67     INTEGER nLevHvR
68     PARAMETER ( nLevHvR = 60 )
69     _RL plowHvR, dpHvR, pLevHvR(nLevHvR+1), pMidHvR(nLevHvR)
70     _RL thetaHvR(nLevHvR), PiHvR(nLevHvR+1), dPiHvR(nLevHvR)
71     _RL recip_kappa, PiLoc, zLoc, dzLoc, yLatLoc, phiLoc
72 jmc 1.5 _RL psNorm, rMidKp1
73     _RL ratioRm(Nr), ratioRp(Nr)
74 jmc 1.4 INTEGER kLev
75     #ifdef CHECK_ANALYLIC_THETA
76 jmc 1.5 _RL tmpVar(nLevHvR,61)
77 jmc 1.4 #endif
78 cnh 1.2 CEOP
79 jmc 1.1
80 jmc 1.4 IF ( selectFindRoSurf.LT.0 .OR. selectFindRoSurf.GT.1 ) THEN
81     WRITE(msgBuf,'(A,I2,A)')
82     & 'INI_P_GROUND: selectFindRoSurf =', selectFindRoSurf,
83     & ' <== bad value !'
84     CALL PRINT_ERROR( msgBuf , myThid)
85     STOP 'INI_P_GROUND'
86     ENDIF
87    
88 jmc 1.1 DO K=1,Nr
89     rHalf(2*K-1) = rF(K)
90     rHalf(2*K) = rC(K)
91     ENDDO
92     rHalf(2*Nr+1) = rF(Nr+1)
93    
94 jmc 1.6 C- Reference Geopotential at Half levels :
95     C Tracer level: hRef(2K) ; Interface_W level: hRef(2K+1)
96 jmc 1.1 C- Convert phiRef to Z unit :
97     DO K=1,2*Nr+1
98 jmc 1.6 hRef(K) = phiRef(K)*recip_gravity
99 jmc 1.1 ENDDO
100 jmc 1.4
101 jmc 1.5 IF (selectFindRoSurf.EQ.0 .AND. selectMode .GT. 0 ) THEN
102 jmc 1.1 C- Find Po_surf : Linear between 2 half levels :
103     DO bj = myByLo(myThid), myByHi(myThid)
104     DO bi = myBxLo(myThid), myBxHi(myThid)
105     DO j=1,sNy
106     DO i=1,sNx
107     Ks = 1
108     DO K=1,2*Nr
109 jmc 1.6 IF (Hfld(i,j,bi,bj).GE.hRef(K)) Ks = K
110 jmc 1.1 ENDDO
111     Po_surf = rHalf(Ks) + (rHalf(Ks+1)-rHalf(Ks))*
112 jmc 1.6 & (Hfld(i,j,bi,bj)-hRef(Ks))/(hRef(Ks+1)-hRef(Ks))
113 jmc 1.1
114 jmc 1.6 c IF (Hfld(i,j,bi,bj).LT.hRef(1)) Po_surf= rHalf(1)
115     c IF (Hfld(i,j,bi,bj).GT.hRef(2*Nr+1)) Po_surf=rHalf(2*Nr+1)
116 jmc 1.1 Pfld(i,j,bi,bj) = Po_surf
117     ENDDO
118     ENDDO
119     ENDDO
120     ENDDO
121    
122     C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
123 jmc 1.5 C-- endif selectFindRoSurf=0 & selectMode > 0
124 jmc 1.4 ENDIF
125    
126 jmc 1.5 IF ( selectFindRoSurf.EQ.1 ) THEN
127 jmc 1.4 C-- define high resolution Pressure discretization:
128    
129     recip_kappa = 1. _d 0 / atm_kappa
130     plowHvR = 0.4 _d 0
131     dpHvR = nLevHvR
132     dpHvR = (1. - plowHvR) / dpHvR
133     pLevHvR(1)= Ro_SeaLevel/atm_Po
134     PiHvR(1) = atm_Cp*(pLevHvR(1)**atm_kappa)
135     DO k=1,nLevHvR
136     pLevHvR(k+1)= pLevHvR(1) - float(k)*dpHvR
137     PiHvR(k+1) = atm_Cp*(pLevHvR(k+1)**atm_kappa)
138     pMidHvR(k)= (pLevHvR(k)+pLevHvR(k+1))*0.5 _d 0
139     dPiHvR(k) = PiHvR(k) - PiHvR(k+1)
140     ENDDO
141    
142 jmc 1.5 C-- to modify pressure when using non fully linear relation between Phi & p
143     C (Integr_GeoPot=2 & Tracer Point at middle between 2 interfaces)
144     DO k=1,Nr
145     ratioRm(k) = 1. _d 0
146     ratioRp(k) = 1. _d 0
147     IF (k.GT.1 ) ratioRm(k) = 0.5 _d 0*drC(k)/(rF(k)-rC(k))
148     IF (k.LT.Nr) ratioRp(k) = 0.5 _d 0*drC(k+1)/(rC(k)-rF(k+1))
149     ENDDO
150    
151 jmc 1.4 #ifdef CHECK_ANALYLIC_THETA
152     _BEGIN_MASTER( mythid )
153     DO j=1,61
154     yLatLoc =-90.+(j-1)*3.
155     CALL ANALYLIC_THETA( yLatLoc , pMidHvR,
156     & tmpVar(1,j), nLevHvR, mythid)
157     ENDDO
158     OPEN(88,FILE='analytic_theta',
159     & STATUS='unknown',FORM='unformatted')
160     WRITE(88) tmpVar
161     CLOSE(88)
162     _END_MASTER( mythid )
163     STOP 'CHECK_ANALYLIC_THETA'
164     #endif /* CHECK_ANALYLIC_THETA */
165    
166 jmc 1.5 C-- endif selectFindRoSurf=1
167 jmc 1.4 ENDIF
168    
169 jmc 1.5 IF ( selectFindRoSurf*selectMode .GT. 0) THEN
170 jmc 1.4 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
171     C- Find Po_surf such as g*Hfld = Phi[Po_surf,theta(yLat,p)]:
172    
173     DO bj = myByLo(myThid), myByHi(myThid)
174     DO bi = myBxLo(myThid), myBxHi(myThid)
175     C- start bi,bj loop:
176 jmc 1.5
177 jmc 1.4 DO j=1,sNy
178     DO i=1,sNx
179     IF ( Hfld(i,j,bi,bj) .LE. 0. _d 0) THEN
180     Pfld(i,j,bi,bj) = Ro_SeaLevel
181     ELSE
182     yLatLoc = yC(i,j,bi,bj)
183     CALL ANALYLIC_THETA( yLatLoc , pMidHvR,
184     & thetaHvR, nLevHvR, mythid)
185     zLoc = 0.
186     DO k=1,nLevHvR
187     IF (zLoc.GE.0.) THEN
188     C- compute phi/g corresponding to next p_level:
189     dzLoc = dPiHvR(k)*thetaHvR(k)*recip_gravity
190     IF ( Hfld(i,j,bi,bj) .LE. zLoc+dzLoc ) THEN
191 jmc 1.5 C- compute the normalized surf. Pressure psNorm
192 jmc 1.4 PiLoc = PiHvR(k)
193     & - gravity*(Hfld(i,j,bi,bj)-zLoc)/thetaHvR(k)
194 jmc 1.5 psNorm = (PiLoc/atm_Cp)**recip_kappa
195 jmc 1.4 C- use linear interpolation:
196 jmc 1.5 c psNorm = pLevHvR(k)
197     c & - dpHvR*(Hfld(i,j,bi,bj)-zLoc)/dzLoc
198 jmc 1.4 zLoc = -1.
199     ELSE
200     zLoc = zLoc + dzLoc
201     ENDIF
202     ENDIF
203     ENDDO
204     IF (zLoc.GE.0.) THEN
205     WRITE(msgBuf,'(2A)')
206     & 'INI_P_GROUND: FAIL in trying to find Pfld:',
207     & ' selectMode,i,j,bi,bj=',selectMode,i,j,bi,bj
208     CALL PRINT_ERROR( msgBuf , myThid)
209     WRITE(msgBuf,'(A,F7.1,2A,F6.4,A,F8.0)')
210     & 'INI_P_GROUND: Hfld=', Hfld(i,j,bi,bj), ' exceeds',
211     & ' Zloc(lowest P=', pLevHvR(1+nLevHvR),' )=',zLoc
212     CALL PRINT_ERROR( msgBuf , myThid)
213     STOP 'ABNORMAL END: S/R INI_P_GROUND'
214     ELSE
215 jmc 1.5 Pfld(i,j,bi,bj) = psNorm*atm_Po
216 jmc 1.4 ENDIF
217     ENDIF
218     ENDDO
219     ENDDO
220 jmc 1.5
221     IF (selectMode.EQ.2 .AND. integr_GeoPot.NE.1) THEN
222     C---------
223     C Modify pressure to account for non fully linear relation between
224     C Phi & p (due to numerical truncation of the Finite Difference scheme)
225     C---------
226     DO j=1,sNy
227     DO i=1,sNx
228     Po_surf = Pfld(i,j,bi,bj)
229     IF ( Po_surf.LT.rC(1) .AND. Po_surf.GT.rC(Nr) ) THEN
230     findPoSurf = .TRUE.
231     DO k=1,Nr
232     IF ( findPoSurf .AND. Po_surf.GE.rC(k) ) THEN
233     Po_surf = rC(k) + (Po_surf-rC(k))/ratioRm(k)
234     findPoSurf = .FALSE.
235     ENDIF
236     rMidKp1 = rF(k+1)
237     IF (k.LT.Nr) rMidKp1 = (rC(k)+rC(k+1))*0.5 _d 0
238     IF ( findPoSurf .AND. Po_surf.GE.rMidKp1 ) THEN
239     Po_surf = rC(k) + (Po_surf-rC(k))/ratioRp(k)
240     findPoSurf = .FALSE.
241     ENDIF
242     ENDDO
243     IF ( findPoSurf )
244     & STOP 'S/R INI_P_GROUND: Pb with selectMode=2'
245     ENDIF
246     Pfld(i,j,bi,bj) = Po_surf
247     ENDDO
248     ENDDO
249     C---------
250     ENDIF
251    
252 jmc 1.4 C- end bi,bj loop.
253     ENDDO
254     ENDDO
255    
256     C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
257 jmc 1.5 C-- endif selectFindRoSurf*selectMode > 0
258 jmc 1.4 ENDIF
259    
260 jmc 1.5 IF (selectMode .LT. 0) THEN
261 jmc 1.4 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
262 jmc 1.5 C--- Compute Hfld = Phi(Pfld)/g.
263 jmc 1.4
264     DO bj = myByLo(myThid), myByHi(myThid)
265     DO bi = myBxLo(myThid), myBxHi(myThid)
266     C- start bi,bj loop:
267    
268 jmc 1.6 C-- Compute Hfld from g*Hfld = hRef(Po_surf)
269 jmc 1.4 DO j=1,sNy
270     DO i=1,sNx
271 jmc 1.6 C- compute phiLoc = hRef(Ro_surf):
272 jmc 1.4 ks = ksurfC(i,j,bi,bj)
273     IF (ks.LE.Nr) THEN
274     IF ( Pfld(i,j,bi,bj).GE.rC(ks) ) THEN
275 jmc 1.6 phiLoc = hRef(2*ks)
276     & + (hRef(2*ks-1)-hRef(2*ks))
277 jmc 1.4 & *(Pfld(i,j,bi,bj)-rC(ks))/(rHalf(2*ks-1)-rHalf(2*ks))
278     ELSE
279 jmc 1.6 phiLoc = hRef(2*ks)
280     & + (hRef(2*ks+1)-hRef(2*ks))
281 jmc 1.4 & *(Pfld(i,j,bi,bj)-rC(ks))/(rHalf(2*ks+1)-rHalf(2*ks))
282     ENDIF
283 jmc 1.5 Hfld(i,j,bi,bj) = phiLoc
284 jmc 1.4 ELSE
285     Hfld(i,j,bi,bj) = 0.
286     ENDIF
287     ENDDO
288     ENDDO
289 jmc 1.5
290     IF (selectFindRoSurf.EQ.1) THEN
291     C-----
292     C goal: evaluate phi0surf (used for computing geopotential'= Phi - PhiRef):
293     C phi0surf = g*Ztopo-PhiRef(Ro_surf) if no truncation was applied to Ro_surf;
294     C but because of hFacMin, surf.ref.pressure (=Ro_surf) is truncated and
295     C phi0surf = Phi(Theta-Analytic,P=Ro_surf) - PhiRef(P=Ro_surf)
296     C-----
297     C-- Compute Hfld from g*Hfld = Phi[Po_surf,theta(yLat,p)]:
298     DO j=1,sNy
299     DO i=1,sNx
300     zLoc = 0.
301     IF ( Pfld(i,j,bi,bj) .LT. Ro_SeaLevel) THEN
302     Po_surf = Pfld(i,j,bi,bj)
303     C---------
304     C Modify pressure to account for non fully linear relation between
305     C Phi & p (due to numerical truncation of the Finite Difference scheme)
306     IF (selectMode.EQ.-2 .AND. integr_GeoPot.NE.1) THEN
307     IF ( Po_surf.LT.rC(1) .AND. Po_surf.GT.rC(Nr) ) THEN
308     findPoSurf = .TRUE.
309     DO k=1,Nr
310     IF ( findPoSurf .AND. Po_surf.GE.rC(k) ) THEN
311     Po_surf = rC(k) + (Po_surf-rC(k))*ratioRm(k)
312     findPoSurf = .FALSE.
313     ENDIF
314     IF ( findPoSurf .AND. Po_surf.GE.rF(k+1) ) THEN
315     Po_surf = rC(k) + (Po_surf-rC(k))*ratioRp(k)
316     findPoSurf = .FALSE.
317     ENDIF
318     ENDDO
319     ENDIF
320     ENDIF
321     C---------
322     psNorm = Po_surf/atm_Po
323     kLev = 1 + INT( (pLevHvR(1)-psNorm)/dpHvR )
324     yLatLoc = yC(i,j,bi,bj)
325     CALL ANALYLIC_THETA( yLatLoc , pMidHvR,
326     & thetaHvR, kLev, mythid)
327     C- compute height at level pLev(kLev) just below Pfld:
328     DO k=1,kLev-1
329     dzLoc = dPiHvR(k)*thetaHvR(k)*recip_gravity
330     zLoc = zLoc + dzLoc
331     ENDDO
332     dzLoc = ( PiHvR(kLev)-atm_Cp*(psNorm**atm_kappa) )
333     & * thetaHvR(kLev)*recip_gravity
334     zLoc = zLoc + dzLoc
335     ENDIF
336     C- compute phi0surf = Phi[Po_surf,theta(yLat,p)] - PhiRef(Po_surf)
337     phi0surf(i,j,bi,bj) = gravity*(zLoc - Hfld(i,j,bi,bj))
338     C- save Phi[Po_surf,theta(yLat,p)] in Hfld (replacing PhiRef(Po_surf)):
339     Hfld(i,j,bi,bj) = zLoc
340     ENDDO
341     ENDDO
342     C- endif selectFindRoSurf=1
343     ENDIF
344    
345 jmc 1.4 C- end bi,bj loop.
346     ENDDO
347     ENDDO
348    
349     C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
350 jmc 1.5 C-- endif selectMode < 0
351 jmc 1.4 ENDIF
352    
353     RETURN
354     END
355    
356     CBOP
357     C !SUBROUTINE: ANALYLIC_THETA
358     C !INTERFACE:
359     SUBROUTINE ANALYLIC_THETA( yLat, pNlev,
360     O thetaLev,
361     I kSize,myThid)
362    
363     C !DESCRIPTION: \bv
364     C *==========================================================*
365     C | SUBROUTINE ANALYLIC_THETA
366     C | o Conpute analyticaly the potential temperature Theta
367     C | as a function of Latitude (yLat) and normalized
368     C | pressure pNlev.
369     C | approximatively match the N-S symetric, zonal-mean and
370     C | annual-mean NCEP climatology in the troposphere.
371     C *==========================================================*
372     C \ev
373    
374     C !USES:
375     IMPLICIT NONE
376    
377     C == Global variables ==
378     #include "SIZE.h"
379     #include "EEPARAMS.h"
380     #include "PARAMS.h"
381    
382     C !INPUT/OUTPUT PARAMETERS:
383     C == Routine arguments ==
384     C yLat :: latitude (degre)
385     C pNlev :: normalized pressure levels
386     C kSize :: dimension of pNlev & ANALYLIC_THETA
387     C myThid :: Thread number for this instance of the routine
388     INTEGER kSize
389     _RL yLat
390     _RL pNlev (kSize)
391     _RL thetaLev(kSize)
392     INTEGER myThid
393    
394     C !LOCAL VARIABLES:
395     C == Local variables ==
396     INTEGER k
397     _RL yyA, yyB, yyC, yyAd, yyBd, yyCd
398     _RL cAtmp, cBtmp, ttdC
399     _RL ppN0, ppN1, ppN2, ppN3a, ppN3b, ppN4
400     _RL ttp1, ttp2, ttp3, ttp4, ttp5
401     _RL yAtmp, yBtmp, yCtmp, yDtmp
402     _RL ttp2y, ttp4y, a1tmp
403     _RL ppl, ppm, pph, ppr
404    
405     CEOP
406    
407     DATA yyA , yyB , yyC , yyAd , yyBd , yyCd
408     & / 45. _d 0, 65. _d 0, 65. _d 0, .9 _d 0, .9 _d 0, 10. _d 0 /
409     DATA cAtmp , cBtmp , ttdC
410     & / 2.6 _d 0, 1.5 _d 0, 3.3 _d 0 /
411     DATA ppN0 , ppN1 , ppN2 , ppN3a , ppN3b , ppN4
412     & / .1 _d 0, .19 _d 0, .3 _d 0, .9 _d 0, .7 _d 0, .925 _d 0 /
413     DATA ttp1 , ttp2 , ttp3 , ttp4 , ttp5
414     & / 350. _d 0, 342. _d 0, 307. _d 0, 301. _d 0, 257. _d 0 /
415    
416     C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
417    
418     yAtmp = ABS(yLat) - yyA
419     yAtmp = yyA + MIN(0. _d 0,yAtmp/yyAd) + MAX(yAtmp, 0. _d 0)
420     yAtmp = COS( deg2rad*MAX(yAtmp, 0. _d 0) )
421     yBtmp = ABS(yLat) - yyB
422     yBtmp = yyB + yBtmp/yyBd
423     yBtmp = COS( deg2rad*MAX( 0. _d 0, MIN(yBtmp,90. _d 0) ) )
424     yCtmp = ABS(yLat) - yyC
425     yCtmp = MAX( 0. _d 0, 1. _d 0 - (yCtmp/yyCd)**2 )
426     yDtmp = ppN3a +(ppN3b - ppN3a)*yCtmp
427     ttp2y = ttp3 + (ttp2-ttp3)*yAtmp**cAtmp
428     ttp4y = ttp5 + (ttp4-ttp5)*yBtmp**cBtmp
429     a1tmp = (ttp1-ttp2y)*ppN1*ppN2/(ppN2-ppN1)
430     DO k=1,kSize
431     ppl = MIN(pNlev(k),ppN1)
432     ppm = MIN(MAX(pNlev(k),ppN1),ppN2)
433     pph = MAX(pNlev(k),ppN2)
434     ppr =( ppN0 + ABS(ppl-ppN0) - ppN1 )/(ppN2-ppN1)
435     thetaLev(k) =
436     & ( (1. _d 0 -ppr)*ttp1*ppN1**atm_kappa
437     & + ppr*ttp2y*ppN2**atm_kappa
438     & )*ppl**(-atm_kappa)
439     & + a1tmp*(1. _d 0 /ppm - 1. _d 0/ppN1)
440     & + (ttp4y-ttp2y)*(pph-ppN2)/(ppN4-ppN2)
441     & + (ttdC+yCtmp)*MAX(0. _d 0,pNlev(k)-yDtmp)/(1-yDtmp)
442     ENDDO
443 jmc 1.1
444 jmc 1.4 C---------------------------------------------------
445     C matlab script, input: pN, yp=abs(yLat)
446     C pN0=.1; pN1=.19 ; pN2=.3; pN4=.925;
447     C pm=min(max(pN,pN1),pN2); pp=max(pN,pN2);
448     C pl=min(pN,pN1); pr=(pN0+abs(pl-pN0)-pN1)/(pN2-pN1);
449     C
450     C yA=yp-45; yA=45+min(0,yA/.9)+max(0,yA); yA=max(0,yA); cosyA=cos(yA*rad) ;
451     C yB=yp-65; yB=65+yB/.9; yB=min(max(0,yB),90); cosyB=cos(yB*rad) ;
452     C tp1=350*ones(nyA,1);
453     C tp2=307+(342-307)*cosyA.^2.6;
454     C tp4=257+(301-257)*cosyB.^1.5;
455     C a1=(tp1-tp2)*pN1*pN2/(pN2-pN1);
456     C pF=max(0,1.-((yp-65)/10).^2); pT=.9+(0.7-.9)*pF;
457     C
458     C tA0=((1-pr(k))*tp1(j)*pN1^kappa+pr(k)*tp2(j)*pN2^kappa)*pl(k)^-kappa;
459     C tA1=a1(j)*(1./pm(k)-1./pN1);
460     C tA2=(tp4(j)-tp2(j))*(pp(k)-pN2)/(pN4-pN2);
461     C tA3=(3.3+pF(j))*max(0,pN(k)-pT(j))/(1-pT(j));
462     C tAn(j,k)=tA0+tA1+tA2+tA3;
463     C---------------------------------------------------
464 jmc 1.1
465     RETURN
466     END

  ViewVC Help
Powered by ViewVC 1.1.22