1 |
C $Header: /u/gcmpack/MITgcm/model/src/ini_p_ground.F,v 1.10 2016/03/10 20:52:43 jmc Exp $ |
2 |
C $Name: $ |
3 |
|
4 |
#include "CPP_OPTIONS.h" |
5 |
#undef CHECK_ANALYLIC_THETA |
6 |
|
7 |
CBOP |
8 |
C !ROUTINE: INI_P_GROUND |
9 |
C !INTERFACE: |
10 |
SUBROUTINE INI_P_GROUND(selectMode, |
11 |
& Hfld, Pfld, |
12 |
I myThid ) |
13 |
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 |
IMPLICIT NONE |
25 |
C == Global variables == |
26 |
#include "SIZE.h" |
27 |
#include "GRID.h" |
28 |
#include "EEPARAMS.h" |
29 |
#include "PARAMS.h" |
30 |
#include "SURFACE.h" |
31 |
|
32 |
C !INPUT/OUTPUT PARAMETERS: |
33 |
C == Routine arguments == |
34 |
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 |
C Hfld (input/outp) :: Ground elevation [m] |
38 |
C Pfld (outp/input) :: reference Pressure at the ground [Pa] |
39 |
INTEGER selectMode |
40 |
_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 |
C !LOCAL VARIABLES: |
45 |
C == Local variables == |
46 |
C msgBuf :: Informational/error message 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 |
INTEGER bi,bj,i,j,k, ks |
63 |
_RL Po_surf |
64 |
_RL hRef(2*Nr+1), rHalf(2*Nr+1) |
65 |
LOGICAL findPoSurf |
66 |
|
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 |
_RL psNorm, rMidKp1 |
73 |
_RL ratioRm(Nr), ratioRp(Nr) |
74 |
INTEGER kLev |
75 |
#ifdef CHECK_ANALYLIC_THETA |
76 |
_RL tmpVar(nLevHvR,61) |
77 |
#endif |
78 |
CEOP |
79 |
|
80 |
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 |
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 |
C- Reference Geopotential at Half levels : |
95 |
C Tracer level: hRef(2k) ; Interface_W level: hRef(2k+1) |
96 |
C- Convert phiRef to Z unit : |
97 |
DO k=1,2*Nr+1 |
98 |
hRef(k) = phiRef(k)*recip_gravity |
99 |
ENDDO |
100 |
|
101 |
IF (selectFindRoSurf.EQ.0 .AND. selectMode .GT. 0 ) THEN |
102 |
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 |
IF (Hfld(i,j,bi,bj).GE.hRef(k)) ks = k |
110 |
ENDDO |
111 |
Po_surf = rHalf(ks) + (rHalf(ks+1)-rHalf(ks))* |
112 |
& (Hfld(i,j,bi,bj)-hRef(ks))/(hRef(ks+1)-hRef(ks)) |
113 |
|
114 |
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 |
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 |
C-- endif selectFindRoSurf=0 & selectMode > 0 |
124 |
ENDIF |
125 |
|
126 |
IF ( selectFindRoSurf.EQ.1 ) THEN |
127 |
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)= rF(1)/atm_Po |
134 |
PiHvR(1) = atm_Cp*(pLevHvR(1)**atm_kappa) |
135 |
DO k=1,nLevHvR |
136 |
pLevHvR(k+1)= pLevHvR(1) - 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 |
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 |
#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 |
C-- endif selectFindRoSurf=1 |
167 |
ENDIF |
168 |
|
169 |
IF ( selectFindRoSurf*selectMode .GT. 0) THEN |
170 |
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 |
|
177 |
DO j=1,sNy |
178 |
DO i=1,sNx |
179 |
phiLoc = Hfld(i,j,bi,bj) - hRef(1) |
180 |
IF ( phiLoc .LE. 0. _d 0 ) THEN |
181 |
Pfld(i,j,bi,bj) = rF(1) |
182 |
ELSE |
183 |
yLatLoc = yC(i,j,bi,bj) |
184 |
CALL ANALYLIC_THETA( yLatLoc, pMidHvR, |
185 |
& thetaHvR, nLevHvR, myThid ) |
186 |
zLoc = 0. |
187 |
DO k=1,nLevHvR |
188 |
IF (zLoc.GE.0.) THEN |
189 |
C- compute phi/g corresponding to next p_level: |
190 |
dzLoc = dPiHvR(k)*thetaHvR(k)*recip_gravity |
191 |
IF ( phiLoc .LE. zLoc+dzLoc ) THEN |
192 |
C- compute the normalized surf. Pressure psNorm |
193 |
PiLoc = PiHvR(k) |
194 |
& - gravity*( phiLoc - zLoc )/thetaHvR(k) |
195 |
psNorm = (PiLoc/atm_Cp)**recip_kappa |
196 |
C- use linear interpolation: |
197 |
c psNorm = pLevHvR(k) |
198 |
c & - dpHvR*( phiLoc - zLoc )/dzLoc |
199 |
zLoc = -1. |
200 |
ELSE |
201 |
zLoc = zLoc + dzLoc |
202 |
ENDIF |
203 |
ENDIF |
204 |
ENDDO |
205 |
IF (zLoc.GE.0.) THEN |
206 |
WRITE(msgBuf,'(2A)') |
207 |
& 'INI_P_GROUND: FAIL in trying to find Pfld:', |
208 |
& ' selectMode,i,j,bi,bj=',selectMode,i,j,bi,bj |
209 |
CALL PRINT_ERROR( msgBuf, myThid ) |
210 |
WRITE(msgBuf,'(A,F7.1,2A,F6.4,A,F8.0)') |
211 |
& 'INI_P_GROUND: Hfld=', Hfld(i,j,bi,bj), ' exceeds', |
212 |
& ' Zloc(lowest P=', pLevHvR(1+nLevHvR),' )=',zLoc |
213 |
CALL PRINT_ERROR( msgBuf, myThid ) |
214 |
STOP 'ABNORMAL END: S/R INI_P_GROUND' |
215 |
ELSE |
216 |
Pfld(i,j,bi,bj) = psNorm*atm_Po |
217 |
ENDIF |
218 |
ENDIF |
219 |
ENDDO |
220 |
ENDDO |
221 |
|
222 |
IF (selectMode.EQ.2 .AND. integr_GeoPot.NE.1) THEN |
223 |
C--------- |
224 |
C Modify pressure to account for non fully linear relation between |
225 |
C Phi & p (due to numerical truncation of the Finite Difference scheme) |
226 |
C--------- |
227 |
DO j=1,sNy |
228 |
DO i=1,sNx |
229 |
Po_surf = Pfld(i,j,bi,bj) |
230 |
IF ( Po_surf.LT.rC(1) .AND. Po_surf.GT.rC(Nr) ) THEN |
231 |
findPoSurf = .TRUE. |
232 |
DO k=1,Nr |
233 |
IF ( findPoSurf .AND. Po_surf.GE.rC(k) ) THEN |
234 |
Po_surf = rC(k) + (Po_surf-rC(k))/ratioRm(k) |
235 |
findPoSurf = .FALSE. |
236 |
ENDIF |
237 |
rMidKp1 = rF(k+1) |
238 |
IF (k.LT.Nr) rMidKp1 = (rC(k)+rC(k+1))*0.5 _d 0 |
239 |
IF ( findPoSurf .AND. Po_surf.GE.rMidKp1 ) THEN |
240 |
Po_surf = rC(k) + (Po_surf-rC(k))/ratioRp(k) |
241 |
findPoSurf = .FALSE. |
242 |
ENDIF |
243 |
ENDDO |
244 |
IF ( findPoSurf ) |
245 |
& STOP 'S/R INI_P_GROUND: Pb with selectMode=2' |
246 |
ENDIF |
247 |
Pfld(i,j,bi,bj) = Po_surf |
248 |
ENDDO |
249 |
ENDDO |
250 |
C--------- |
251 |
ENDIF |
252 |
|
253 |
C- end bi,bj loop. |
254 |
ENDDO |
255 |
ENDDO |
256 |
|
257 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
258 |
C-- endif selectFindRoSurf*selectMode > 0 |
259 |
ENDIF |
260 |
|
261 |
IF (selectMode .LT. 0) THEN |
262 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
263 |
C--- Compute Hfld = Phi(Pfld)/g. |
264 |
|
265 |
DO bj = myByLo(myThid), myByHi(myThid) |
266 |
DO bi = myBxLo(myThid), myBxHi(myThid) |
267 |
C- start bi,bj loop: |
268 |
|
269 |
C-- Compute Hfld from g*Hfld = hRef(Po_surf) |
270 |
DO j=1,sNy |
271 |
DO i=1,sNx |
272 |
C- compute phiLoc = hRef(Ro_surf): |
273 |
ks = kSurfC(i,j,bi,bj) |
274 |
IF (ks.LE.Nr) THEN |
275 |
C- sigma coord. case (=> uniform kSurfC), find true ks: |
276 |
IF ( selectSigmaCoord.NE.0 ) THEN |
277 |
DO k=2,Nr |
278 |
IF ( Pfld(i,j,bi,bj).LT.rF(k) ) ks = k |
279 |
ENDDO |
280 |
ENDIF |
281 |
IF ( Pfld(i,j,bi,bj).GE.rC(ks) ) THEN |
282 |
phiLoc = hRef(2*ks) |
283 |
& + (hRef(2*ks-1)-hRef(2*ks)) |
284 |
& *(Pfld(i,j,bi,bj)-rC(ks))/(rHalf(2*ks-1)-rHalf(2*ks)) |
285 |
ELSE |
286 |
phiLoc = hRef(2*ks) |
287 |
& + (hRef(2*ks+1)-hRef(2*ks)) |
288 |
& *(Pfld(i,j,bi,bj)-rC(ks))/(rHalf(2*ks+1)-rHalf(2*ks)) |
289 |
ENDIF |
290 |
Hfld(i,j,bi,bj) = phiLoc |
291 |
ELSE |
292 |
Hfld(i,j,bi,bj) = 0. |
293 |
ENDIF |
294 |
ENDDO |
295 |
ENDDO |
296 |
|
297 |
IF (selectFindRoSurf.EQ.1) THEN |
298 |
C----- |
299 |
C goal: evaluate phi0surf (used for computing geopotential_prime = Phi - PhiRef): |
300 |
C phi0surf = g*Ztopo-PhiRef(Ro_surf) if no truncation was applied to Ro_surf; |
301 |
C but because of hFacMin, surf.ref.pressure (=Ro_surf) is truncated and |
302 |
C phi0surf = Phi(Theta-Analytic,P=Ro_surf) - PhiRef(P=Ro_surf) |
303 |
C----- |
304 |
C-- Compute Hfld from g*Hfld = Phi[Po_surf,theta(yLat,p)]: |
305 |
DO j=1,sNy |
306 |
DO i=1,sNx |
307 |
zLoc = hRef(1) |
308 |
IF ( Pfld(i,j,bi,bj) .LT. rF(1) ) THEN |
309 |
Po_surf = Pfld(i,j,bi,bj) |
310 |
C--------- |
311 |
C Modify pressure to account for non fully linear relation between |
312 |
C Phi & p (due to numerical truncation of the Finite Difference scheme) |
313 |
IF (selectMode.EQ.-2 .AND. integr_GeoPot.NE.1) THEN |
314 |
IF ( Po_surf.LT.rC(1) .AND. Po_surf.GT.rC(Nr) ) THEN |
315 |
findPoSurf = .TRUE. |
316 |
DO k=1,Nr |
317 |
IF ( findPoSurf .AND. Po_surf.GE.rC(k) ) THEN |
318 |
Po_surf = rC(k) + (Po_surf-rC(k))*ratioRm(k) |
319 |
findPoSurf = .FALSE. |
320 |
ENDIF |
321 |
IF ( findPoSurf .AND. Po_surf.GE.rF(k+1) ) THEN |
322 |
Po_surf = rC(k) + (Po_surf-rC(k))*ratioRp(k) |
323 |
findPoSurf = .FALSE. |
324 |
ENDIF |
325 |
ENDDO |
326 |
ENDIF |
327 |
ENDIF |
328 |
C--------- |
329 |
psNorm = Po_surf/atm_Po |
330 |
kLev = 1 + INT( (pLevHvR(1)-psNorm)/dpHvR ) |
331 |
yLatLoc = yC(i,j,bi,bj) |
332 |
CALL ANALYLIC_THETA( yLatLoc, pMidHvR, |
333 |
& thetaHvR, kLev, myThid ) |
334 |
C- compute height at level pLev(kLev) just below Pfld: |
335 |
DO k=1,kLev-1 |
336 |
dzLoc = dPiHvR(k)*thetaHvR(k)*recip_gravity |
337 |
zLoc = zLoc + dzLoc |
338 |
ENDDO |
339 |
dzLoc = ( PiHvR(kLev)-atm_Cp*(psNorm**atm_kappa) ) |
340 |
& * thetaHvR(kLev)*recip_gravity |
341 |
zLoc = zLoc + dzLoc |
342 |
ENDIF |
343 |
C- compute phi0surf = Phi[Po_surf,theta(yLat,p)] - PhiRef(Po_surf) |
344 |
phi0surf(i,j,bi,bj) = gravity*(zLoc - Hfld(i,j,bi,bj)) |
345 |
C- save Phi[Po_surf,theta(yLat,p)] in Hfld (replacing PhiRef(Po_surf)): |
346 |
Hfld(i,j,bi,bj) = zLoc |
347 |
ENDDO |
348 |
ENDDO |
349 |
C- endif selectFindRoSurf=1 |
350 |
ENDIF |
351 |
|
352 |
C- end bi,bj loop. |
353 |
ENDDO |
354 |
ENDDO |
355 |
|
356 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
357 |
C-- endif selectMode < 0 |
358 |
ENDIF |
359 |
|
360 |
RETURN |
361 |
END |
362 |
|
363 |
CBOP |
364 |
C !SUBROUTINE: ANALYLIC_THETA |
365 |
C !INTERFACE: |
366 |
SUBROUTINE ANALYLIC_THETA( yLat, pNlev, |
367 |
O thetaLev, |
368 |
I kSize, myThid ) |
369 |
|
370 |
C !DESCRIPTION: \bv |
371 |
C *==========================================================* |
372 |
C | SUBROUTINE ANALYLIC_THETA |
373 |
C | o Conpute analyticaly the potential temperature Theta |
374 |
C | as a function of Latitude (yLat) and normalized |
375 |
C | pressure pNlev. |
376 |
C | approximatively match the N-S symetric, zonal-mean and |
377 |
C | annual-mean NCEP climatology in the troposphere. |
378 |
C *==========================================================* |
379 |
C \ev |
380 |
|
381 |
C !USES: |
382 |
IMPLICIT NONE |
383 |
|
384 |
C == Global variables == |
385 |
#include "SIZE.h" |
386 |
#include "EEPARAMS.h" |
387 |
#include "PARAMS.h" |
388 |
|
389 |
C !INPUT/OUTPUT PARAMETERS: |
390 |
C == Routine arguments == |
391 |
C yLat :: latitude (degre) |
392 |
C pNlev :: normalized pressure levels |
393 |
C kSize :: dimension of pNlev & ANALYLIC_THETA |
394 |
C myThid :: Thread number for this instance of the routine |
395 |
INTEGER kSize |
396 |
_RL yLat |
397 |
_RL pNlev (kSize) |
398 |
_RL thetaLev(kSize) |
399 |
INTEGER myThid |
400 |
|
401 |
C !LOCAL VARIABLES: |
402 |
C == Local variables == |
403 |
INTEGER k |
404 |
_RL yyA, yyB, yyC, yyAd, yyBd, yyCd |
405 |
_RL cAtmp, cBtmp, ttdC |
406 |
_RL ppN0, ppN1, ppN2, ppN3a, ppN3b, ppN4 |
407 |
_RL ttp1, ttp2, ttp3, ttp4, ttp5 |
408 |
_RL yAtmp, yBtmp, yCtmp, yDtmp |
409 |
_RL ttp2y, ttp4y, a1tmp |
410 |
_RL ppl, ppm, pph, ppr |
411 |
CEOP |
412 |
|
413 |
DATA yyA , yyB , yyC , yyAd , yyBd , yyCd |
414 |
& / 45. _d 0, 65. _d 0, 65. _d 0, .9 _d 0, .9 _d 0, 10. _d 0 / |
415 |
DATA cAtmp , cBtmp , ttdC |
416 |
& / 2.6 _d 0, 1.5 _d 0, 3.3 _d 0 / |
417 |
DATA ppN0 , ppN1 , ppN2 , ppN3a , ppN3b , ppN4 |
418 |
& / .1 _d 0, .19 _d 0, .3 _d 0, .9 _d 0, .7 _d 0, .925 _d 0 / |
419 |
DATA ttp1 , ttp2 , ttp3 , ttp4 , ttp5 |
420 |
& / 350. _d 0, 342. _d 0, 307. _d 0, 301. _d 0, 257. _d 0 / |
421 |
|
422 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
423 |
|
424 |
yAtmp = ABS(yLat) - yyA |
425 |
yAtmp = yyA + MIN(0. _d 0,yAtmp/yyAd) + MAX(yAtmp, 0. _d 0) |
426 |
yAtmp = COS( deg2rad*MAX(yAtmp, 0. _d 0) ) |
427 |
yBtmp = ABS(yLat) - yyB |
428 |
yBtmp = yyB + yBtmp/yyBd |
429 |
yBtmp = COS( deg2rad*MAX( 0. _d 0, MIN(yBtmp,90. _d 0) ) ) |
430 |
yCtmp = ABS(yLat) - yyC |
431 |
yCtmp = MAX( 0. _d 0, 1. _d 0 - (yCtmp/yyCd)**2 ) |
432 |
yDtmp = ppN3a +(ppN3b - ppN3a)*yCtmp |
433 |
ttp2y = ttp3 + (ttp2-ttp3)*yAtmp**cAtmp |
434 |
ttp4y = ttp5 + (ttp4-ttp5)*yBtmp**cBtmp |
435 |
a1tmp = (ttp1-ttp2y)*ppN1*ppN2/(ppN2-ppN1) |
436 |
DO k=1,kSize |
437 |
ppl = MIN(pNlev(k),ppN1) |
438 |
ppm = MIN(MAX(pNlev(k),ppN1),ppN2) |
439 |
pph = MAX(pNlev(k),ppN2) |
440 |
ppr =( ppN0 + ABS(ppl-ppN0) - ppN1 )/(ppN2-ppN1) |
441 |
thetaLev(k) = |
442 |
& ( (1. _d 0 -ppr)*ttp1*ppN1**atm_kappa |
443 |
& + ppr*ttp2y*ppN2**atm_kappa |
444 |
& )*ppl**(-atm_kappa) |
445 |
& + a1tmp*(1. _d 0 /ppm - 1. _d 0/ppN1) |
446 |
& + (ttp4y-ttp2y)*(pph-ppN2)/(ppN4-ppN2) |
447 |
& + (ttdC+yCtmp)*MAX(0. _d 0,pNlev(k)-yDtmp)/(1-yDtmp) |
448 |
ENDDO |
449 |
|
450 |
C--------------------------------------------------- |
451 |
C matlab script, input: pN, yp=abs(yLat) |
452 |
C pN0=.1; pN1=.19 ; pN2=.3; pN4=.925; |
453 |
C pm=min(max(pN,pN1),pN2); pp=max(pN,pN2); |
454 |
C pl=min(pN,pN1); pr=(pN0+abs(pl-pN0)-pN1)/(pN2-pN1); |
455 |
C |
456 |
C yA=yp-45; yA=45+min(0,yA/.9)+max(0,yA); yA=max(0,yA); cosyA=cos(yA*rad) ; |
457 |
C yB=yp-65; yB=65+yB/.9; yB=min(max(0,yB),90); cosyB=cos(yB*rad) ; |
458 |
C tp1=350*ones(nyA,1); |
459 |
C tp2=307+(342-307)*cosyA.^2.6; |
460 |
C tp4=257+(301-257)*cosyB.^1.5; |
461 |
C a1=(tp1-tp2)*pN1*pN2/(pN2-pN1); |
462 |
C pF=max(0,1.-((yp-65)/10).^2); pT=.9+(0.7-.9)*pF; |
463 |
C |
464 |
C tA0=((1-pr(k))*tp1(j)*pN1^kappa+pr(k)*tp2(j)*pN2^kappa)*pl(k)^-kappa; |
465 |
C tA1=a1(j)*(1./pm(k)-1./pN1); |
466 |
C tA2=(tp4(j)-tp2(j))*(pp(k)-pN2)/(pN4-pN2); |
467 |
C tA3=(3.3+pF(j))*max(0,pN(k)-pT(j))/(1-pT(j)); |
468 |
C tAn(j,k)=tA0+tA1+tA2+tA3; |
469 |
C--------------------------------------------------- |
470 |
|
471 |
RETURN |
472 |
END |