1 |
jmc |
1.5 |
C $Header: /u/gcmpack/MITgcm/model/src/ini_p_ground.F,v 1.4 2002/12/10 02:55:47 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 |
|
|
_RL ddPI, Po_surf |
64 |
|
|
_RL phiRef(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.5 |
IF (selectMode*selectFindRoSurf .LE. 0) THEN |
95 |
jmc |
1.1 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
96 |
|
|
C- Compute Reference Geopotential at Half levels : |
97 |
|
|
C Tracer level: phiRef(2K) ; Interface_W level: phiRef(2K+1) |
98 |
|
|
|
99 |
|
|
phiRef(1) = 0. |
100 |
|
|
|
101 |
jmc |
1.4 |
IF (integr_GeoPot.EQ.1) THEN |
102 |
jmc |
1.1 |
C- Finite Volume Form, linear by half level : |
103 |
|
|
DO K=1,2*Nr |
104 |
|
|
Ks = (K+1)/2 |
105 |
jmc |
1.4 |
ddPI=atm_Cp*( ((rHalf( K )/atm_Po)**atm_kappa) |
106 |
|
|
& -((rHalf(K+1)/atm_Po)**atm_kappa) ) |
107 |
jmc |
1.1 |
phiRef(K+1) = phiRef(K)+ddPI*tRef(Ks) |
108 |
|
|
ENDDO |
109 |
|
|
C------ |
110 |
|
|
ELSE |
111 |
|
|
C- Finite Difference Form, linear between Tracer level : |
112 |
jmc |
1.4 |
C works with integr_GeoPot = 0, 2 or 3 |
113 |
jmc |
1.1 |
K = 1 |
114 |
jmc |
1.4 |
ddPI=atm_Cp*( ((rF(K)/atm_Po)**atm_kappa) |
115 |
|
|
& -((rC(K)/atm_Po)**atm_kappa) ) |
116 |
jmc |
1.1 |
phiRef(2*K) = phiRef(1) + ddPI*tRef(K) |
117 |
|
|
DO K=1,Nr-1 |
118 |
jmc |
1.4 |
ddPI=atm_Cp*( ((rC( K )/atm_Po)**atm_kappa) |
119 |
|
|
& -((rC(K+1)/atm_Po)**atm_kappa) ) |
120 |
jmc |
1.1 |
phiRef(2*K+1) = phiRef(2*K) + ddPI*0.5*tRef(K) |
121 |
|
|
phiRef(2*K+2) = phiRef(2*K) |
122 |
|
|
& + ddPI*0.5*(tRef(K)+tRef(K+1)) |
123 |
|
|
ENDDO |
124 |
|
|
K = Nr |
125 |
jmc |
1.4 |
ddPI=atm_Cp*( ((rC( K )/atm_Po)**atm_kappa) |
126 |
|
|
& -((rF(K+1)/atm_Po)**atm_kappa) ) |
127 |
jmc |
1.1 |
phiRef(2*K+1) = phiRef(2*K) + ddPI*tRef(K) |
128 |
|
|
C------ |
129 |
|
|
ENDIF |
130 |
|
|
|
131 |
|
|
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
132 |
|
|
|
133 |
|
|
C- Convert phiRef to Z unit : |
134 |
|
|
DO K=1,2*Nr+1 |
135 |
|
|
phiRef(K) = phiRef(K)*recip_gravity |
136 |
|
|
ENDDO |
137 |
|
|
|
138 |
jmc |
1.3 |
_BEGIN_MASTER( myThid ) |
139 |
jmc |
1.1 |
C- Write to check : |
140 |
jmc |
1.3 |
WRITE(standardMessageUnit,'(2A)') |
141 |
|
|
& 'INI_P_GROUND: PhiRef/g [m] at Center (integer) ', |
142 |
|
|
& 'and Interface (half-int.) levels:' |
143 |
jmc |
1.1 |
DO K=1,2*Nr+1 |
144 |
jmc |
1.3 |
WRITE(standardMessageUnit,'(A,F5.1,A,F10.1,A,F12.3)') |
145 |
|
|
& ' K=',K*0.5,' ; r=',rHalf(K),' ; phiRef/g=', phiRef(K) |
146 |
jmc |
1.1 |
ENDDO |
147 |
jmc |
1.3 |
_END_MASTER( myThid ) |
148 |
jmc |
1.1 |
|
149 |
jmc |
1.5 |
C-- endif selectMode*selectFindRoSurf =< 0 |
150 |
jmc |
1.4 |
ENDIF |
151 |
|
|
|
152 |
jmc |
1.5 |
IF (selectFindRoSurf.EQ.0 .AND. selectMode .GT. 0 ) THEN |
153 |
jmc |
1.1 |
C- Find Po_surf : Linear between 2 half levels : |
154 |
|
|
DO bj = myByLo(myThid), myByHi(myThid) |
155 |
|
|
DO bi = myBxLo(myThid), myBxHi(myThid) |
156 |
|
|
DO j=1,sNy |
157 |
|
|
DO i=1,sNx |
158 |
|
|
Ks = 1 |
159 |
|
|
DO K=1,2*Nr |
160 |
|
|
IF (Hfld(i,j,bi,bj).GE.phiRef(K)) Ks = K |
161 |
|
|
ENDDO |
162 |
|
|
Po_surf = rHalf(Ks) + (rHalf(Ks+1)-rHalf(Ks))* |
163 |
|
|
& (Hfld(i,j,bi,bj)-phiRef(Ks))/(phiRef(Ks+1)-phiRef(Ks)) |
164 |
|
|
|
165 |
|
|
c IF (Hfld(i,j,bi,bj).LT.phiRef(1)) Po_surf= rHalf(1) |
166 |
jmc |
1.4 |
c IF (Hfld(i,j,bi,bj).GT.phiRef(2*Nr+1)) Po_surf=rHalf(2*Nr+1) |
167 |
jmc |
1.1 |
Pfld(i,j,bi,bj) = Po_surf |
168 |
|
|
ENDDO |
169 |
|
|
ENDDO |
170 |
|
|
ENDDO |
171 |
|
|
ENDDO |
172 |
|
|
|
173 |
|
|
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
174 |
jmc |
1.5 |
C-- endif selectFindRoSurf=0 & selectMode > 0 |
175 |
jmc |
1.4 |
ENDIF |
176 |
|
|
|
177 |
jmc |
1.5 |
IF ( selectFindRoSurf.EQ.1 ) THEN |
178 |
jmc |
1.4 |
C-- define high resolution Pressure discretization: |
179 |
|
|
|
180 |
|
|
recip_kappa = 1. _d 0 / atm_kappa |
181 |
|
|
plowHvR = 0.4 _d 0 |
182 |
|
|
dpHvR = nLevHvR |
183 |
|
|
dpHvR = (1. - plowHvR) / dpHvR |
184 |
|
|
pLevHvR(1)= Ro_SeaLevel/atm_Po |
185 |
|
|
PiHvR(1) = atm_Cp*(pLevHvR(1)**atm_kappa) |
186 |
|
|
DO k=1,nLevHvR |
187 |
|
|
pLevHvR(k+1)= pLevHvR(1) - float(k)*dpHvR |
188 |
|
|
PiHvR(k+1) = atm_Cp*(pLevHvR(k+1)**atm_kappa) |
189 |
|
|
pMidHvR(k)= (pLevHvR(k)+pLevHvR(k+1))*0.5 _d 0 |
190 |
|
|
dPiHvR(k) = PiHvR(k) - PiHvR(k+1) |
191 |
|
|
ENDDO |
192 |
|
|
|
193 |
jmc |
1.5 |
C-- to modify pressure when using non fully linear relation between Phi & p |
194 |
|
|
C (Integr_GeoPot=2 & Tracer Point at middle between 2 interfaces) |
195 |
|
|
DO k=1,Nr |
196 |
|
|
ratioRm(k) = 1. _d 0 |
197 |
|
|
ratioRp(k) = 1. _d 0 |
198 |
|
|
IF (k.GT.1 ) ratioRm(k) = 0.5 _d 0*drC(k)/(rF(k)-rC(k)) |
199 |
|
|
IF (k.LT.Nr) ratioRp(k) = 0.5 _d 0*drC(k+1)/(rC(k)-rF(k+1)) |
200 |
|
|
ENDDO |
201 |
|
|
|
202 |
jmc |
1.4 |
#ifdef CHECK_ANALYLIC_THETA |
203 |
|
|
_BEGIN_MASTER( mythid ) |
204 |
|
|
DO j=1,61 |
205 |
|
|
yLatLoc =-90.+(j-1)*3. |
206 |
|
|
CALL ANALYLIC_THETA( yLatLoc , pMidHvR, |
207 |
|
|
& tmpVar(1,j), nLevHvR, mythid) |
208 |
|
|
ENDDO |
209 |
|
|
OPEN(88,FILE='analytic_theta', |
210 |
|
|
& STATUS='unknown',FORM='unformatted') |
211 |
|
|
WRITE(88) tmpVar |
212 |
|
|
CLOSE(88) |
213 |
|
|
_END_MASTER( mythid ) |
214 |
|
|
STOP 'CHECK_ANALYLIC_THETA' |
215 |
|
|
#endif /* CHECK_ANALYLIC_THETA */ |
216 |
|
|
|
217 |
jmc |
1.5 |
C-- endif selectFindRoSurf=1 |
218 |
jmc |
1.4 |
ENDIF |
219 |
|
|
|
220 |
jmc |
1.5 |
IF ( selectFindRoSurf*selectMode .GT. 0) THEN |
221 |
jmc |
1.4 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
222 |
|
|
C- Find Po_surf such as g*Hfld = Phi[Po_surf,theta(yLat,p)]: |
223 |
|
|
|
224 |
|
|
DO bj = myByLo(myThid), myByHi(myThid) |
225 |
|
|
DO bi = myBxLo(myThid), myBxHi(myThid) |
226 |
|
|
C- start bi,bj loop: |
227 |
jmc |
1.5 |
|
228 |
jmc |
1.4 |
DO j=1,sNy |
229 |
|
|
DO i=1,sNx |
230 |
|
|
IF ( Hfld(i,j,bi,bj) .LE. 0. _d 0) THEN |
231 |
|
|
Pfld(i,j,bi,bj) = Ro_SeaLevel |
232 |
|
|
ELSE |
233 |
|
|
yLatLoc = yC(i,j,bi,bj) |
234 |
|
|
CALL ANALYLIC_THETA( yLatLoc , pMidHvR, |
235 |
|
|
& thetaHvR, nLevHvR, mythid) |
236 |
|
|
zLoc = 0. |
237 |
|
|
DO k=1,nLevHvR |
238 |
|
|
IF (zLoc.GE.0.) THEN |
239 |
|
|
C- compute phi/g corresponding to next p_level: |
240 |
|
|
dzLoc = dPiHvR(k)*thetaHvR(k)*recip_gravity |
241 |
|
|
IF ( Hfld(i,j,bi,bj) .LE. zLoc+dzLoc ) THEN |
242 |
jmc |
1.5 |
C- compute the normalized surf. Pressure psNorm |
243 |
jmc |
1.4 |
PiLoc = PiHvR(k) |
244 |
|
|
& - gravity*(Hfld(i,j,bi,bj)-zLoc)/thetaHvR(k) |
245 |
jmc |
1.5 |
psNorm = (PiLoc/atm_Cp)**recip_kappa |
246 |
jmc |
1.4 |
C- use linear interpolation: |
247 |
jmc |
1.5 |
c psNorm = pLevHvR(k) |
248 |
|
|
c & - dpHvR*(Hfld(i,j,bi,bj)-zLoc)/dzLoc |
249 |
jmc |
1.4 |
zLoc = -1. |
250 |
|
|
ELSE |
251 |
|
|
zLoc = zLoc + dzLoc |
252 |
|
|
ENDIF |
253 |
|
|
ENDIF |
254 |
|
|
ENDDO |
255 |
|
|
IF (zLoc.GE.0.) THEN |
256 |
|
|
WRITE(msgBuf,'(2A)') |
257 |
|
|
& 'INI_P_GROUND: FAIL in trying to find Pfld:', |
258 |
|
|
& ' selectMode,i,j,bi,bj=',selectMode,i,j,bi,bj |
259 |
|
|
CALL PRINT_ERROR( msgBuf , myThid) |
260 |
|
|
WRITE(msgBuf,'(A,F7.1,2A,F6.4,A,F8.0)') |
261 |
|
|
& 'INI_P_GROUND: Hfld=', Hfld(i,j,bi,bj), ' exceeds', |
262 |
|
|
& ' Zloc(lowest P=', pLevHvR(1+nLevHvR),' )=',zLoc |
263 |
|
|
CALL PRINT_ERROR( msgBuf , myThid) |
264 |
|
|
STOP 'ABNORMAL END: S/R INI_P_GROUND' |
265 |
|
|
ELSE |
266 |
jmc |
1.5 |
Pfld(i,j,bi,bj) = psNorm*atm_Po |
267 |
jmc |
1.4 |
ENDIF |
268 |
|
|
ENDIF |
269 |
|
|
ENDDO |
270 |
|
|
ENDDO |
271 |
jmc |
1.5 |
|
272 |
|
|
IF (selectMode.EQ.2 .AND. integr_GeoPot.NE.1) THEN |
273 |
|
|
C--------- |
274 |
|
|
C Modify pressure to account for non fully linear relation between |
275 |
|
|
C Phi & p (due to numerical truncation of the Finite Difference scheme) |
276 |
|
|
C--------- |
277 |
|
|
DO j=1,sNy |
278 |
|
|
DO i=1,sNx |
279 |
|
|
Po_surf = Pfld(i,j,bi,bj) |
280 |
|
|
IF ( Po_surf.LT.rC(1) .AND. Po_surf.GT.rC(Nr) ) THEN |
281 |
|
|
findPoSurf = .TRUE. |
282 |
|
|
DO k=1,Nr |
283 |
|
|
IF ( findPoSurf .AND. Po_surf.GE.rC(k) ) THEN |
284 |
|
|
Po_surf = rC(k) + (Po_surf-rC(k))/ratioRm(k) |
285 |
|
|
findPoSurf = .FALSE. |
286 |
|
|
ENDIF |
287 |
|
|
rMidKp1 = rF(k+1) |
288 |
|
|
IF (k.LT.Nr) rMidKp1 = (rC(k)+rC(k+1))*0.5 _d 0 |
289 |
|
|
IF ( findPoSurf .AND. Po_surf.GE.rMidKp1 ) THEN |
290 |
|
|
Po_surf = rC(k) + (Po_surf-rC(k))/ratioRp(k) |
291 |
|
|
findPoSurf = .FALSE. |
292 |
|
|
ENDIF |
293 |
|
|
ENDDO |
294 |
|
|
IF ( findPoSurf ) |
295 |
|
|
& STOP 'S/R INI_P_GROUND: Pb with selectMode=2' |
296 |
|
|
ENDIF |
297 |
|
|
Pfld(i,j,bi,bj) = Po_surf |
298 |
|
|
ENDDO |
299 |
|
|
ENDDO |
300 |
|
|
C--------- |
301 |
|
|
ENDIF |
302 |
|
|
|
303 |
jmc |
1.4 |
C- end bi,bj loop. |
304 |
|
|
ENDDO |
305 |
|
|
ENDDO |
306 |
|
|
|
307 |
|
|
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
308 |
jmc |
1.5 |
C-- endif selectFindRoSurf*selectMode > 0 |
309 |
jmc |
1.4 |
ENDIF |
310 |
|
|
|
311 |
jmc |
1.5 |
IF (selectMode .LT. 0) THEN |
312 |
jmc |
1.4 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
313 |
jmc |
1.5 |
C--- Compute Hfld = Phi(Pfld)/g. |
314 |
jmc |
1.4 |
|
315 |
|
|
DO bj = myByLo(myThid), myByHi(myThid) |
316 |
|
|
DO bi = myBxLo(myThid), myBxHi(myThid) |
317 |
|
|
C- start bi,bj loop: |
318 |
|
|
|
319 |
jmc |
1.5 |
C-- Compute Hfld from g*Hfld = PhiRef(Po_surf) |
320 |
jmc |
1.4 |
DO j=1,sNy |
321 |
|
|
DO i=1,sNx |
322 |
|
|
C- compute phiLoc = PhiRef(Ro_surf): |
323 |
|
|
ks = ksurfC(i,j,bi,bj) |
324 |
|
|
IF (ks.LE.Nr) THEN |
325 |
|
|
IF ( Pfld(i,j,bi,bj).GE.rC(ks) ) THEN |
326 |
|
|
phiLoc = phiRef(2*ks) |
327 |
|
|
& + (phiRef(2*ks-1)-phiRef(2*ks)) |
328 |
|
|
& *(Pfld(i,j,bi,bj)-rC(ks))/(rHalf(2*ks-1)-rHalf(2*ks)) |
329 |
|
|
ELSE |
330 |
|
|
phiLoc = phiRef(2*ks) |
331 |
|
|
& + (phiRef(2*ks+1)-phiRef(2*ks)) |
332 |
|
|
& *(Pfld(i,j,bi,bj)-rC(ks))/(rHalf(2*ks+1)-rHalf(2*ks)) |
333 |
|
|
ENDIF |
334 |
jmc |
1.5 |
Hfld(i,j,bi,bj) = phiLoc |
335 |
jmc |
1.4 |
ELSE |
336 |
|
|
Hfld(i,j,bi,bj) = 0. |
337 |
|
|
ENDIF |
338 |
|
|
ENDDO |
339 |
|
|
ENDDO |
340 |
jmc |
1.5 |
|
341 |
|
|
IF (selectFindRoSurf.EQ.1) THEN |
342 |
|
|
C----- |
343 |
|
|
C goal: evaluate phi0surf (used for computing geopotential'= Phi - PhiRef): |
344 |
|
|
C phi0surf = g*Ztopo-PhiRef(Ro_surf) if no truncation was applied to Ro_surf; |
345 |
|
|
C but because of hFacMin, surf.ref.pressure (=Ro_surf) is truncated and |
346 |
|
|
C phi0surf = Phi(Theta-Analytic,P=Ro_surf) - PhiRef(P=Ro_surf) |
347 |
|
|
C----- |
348 |
|
|
C-- Compute Hfld from g*Hfld = Phi[Po_surf,theta(yLat,p)]: |
349 |
|
|
DO j=1,sNy |
350 |
|
|
DO i=1,sNx |
351 |
|
|
zLoc = 0. |
352 |
|
|
IF ( Pfld(i,j,bi,bj) .LT. Ro_SeaLevel) THEN |
353 |
|
|
Po_surf = Pfld(i,j,bi,bj) |
354 |
|
|
C--------- |
355 |
|
|
C Modify pressure to account for non fully linear relation between |
356 |
|
|
C Phi & p (due to numerical truncation of the Finite Difference scheme) |
357 |
|
|
IF (selectMode.EQ.-2 .AND. integr_GeoPot.NE.1) THEN |
358 |
|
|
IF ( Po_surf.LT.rC(1) .AND. Po_surf.GT.rC(Nr) ) THEN |
359 |
|
|
findPoSurf = .TRUE. |
360 |
|
|
DO k=1,Nr |
361 |
|
|
IF ( findPoSurf .AND. Po_surf.GE.rC(k) ) THEN |
362 |
|
|
Po_surf = rC(k) + (Po_surf-rC(k))*ratioRm(k) |
363 |
|
|
findPoSurf = .FALSE. |
364 |
|
|
ENDIF |
365 |
|
|
IF ( findPoSurf .AND. Po_surf.GE.rF(k+1) ) THEN |
366 |
|
|
Po_surf = rC(k) + (Po_surf-rC(k))*ratioRp(k) |
367 |
|
|
findPoSurf = .FALSE. |
368 |
|
|
ENDIF |
369 |
|
|
ENDDO |
370 |
|
|
ENDIF |
371 |
|
|
ENDIF |
372 |
|
|
C--------- |
373 |
|
|
psNorm = Po_surf/atm_Po |
374 |
|
|
kLev = 1 + INT( (pLevHvR(1)-psNorm)/dpHvR ) |
375 |
|
|
yLatLoc = yC(i,j,bi,bj) |
376 |
|
|
CALL ANALYLIC_THETA( yLatLoc , pMidHvR, |
377 |
|
|
& thetaHvR, kLev, mythid) |
378 |
|
|
C- compute height at level pLev(kLev) just below Pfld: |
379 |
|
|
DO k=1,kLev-1 |
380 |
|
|
dzLoc = dPiHvR(k)*thetaHvR(k)*recip_gravity |
381 |
|
|
zLoc = zLoc + dzLoc |
382 |
|
|
ENDDO |
383 |
|
|
dzLoc = ( PiHvR(kLev)-atm_Cp*(psNorm**atm_kappa) ) |
384 |
|
|
& * thetaHvR(kLev)*recip_gravity |
385 |
|
|
zLoc = zLoc + dzLoc |
386 |
|
|
ENDIF |
387 |
|
|
C- compute phi0surf = Phi[Po_surf,theta(yLat,p)] - PhiRef(Po_surf) |
388 |
|
|
phi0surf(i,j,bi,bj) = gravity*(zLoc - Hfld(i,j,bi,bj)) |
389 |
|
|
C- save Phi[Po_surf,theta(yLat,p)] in Hfld (replacing PhiRef(Po_surf)): |
390 |
|
|
Hfld(i,j,bi,bj) = zLoc |
391 |
|
|
ENDDO |
392 |
|
|
ENDDO |
393 |
|
|
C- endif selectFindRoSurf=1 |
394 |
|
|
ENDIF |
395 |
|
|
|
396 |
jmc |
1.4 |
C- end bi,bj loop. |
397 |
|
|
ENDDO |
398 |
|
|
ENDDO |
399 |
|
|
|
400 |
|
|
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
401 |
jmc |
1.5 |
C-- endif selectMode < 0 |
402 |
jmc |
1.4 |
ENDIF |
403 |
|
|
|
404 |
|
|
RETURN |
405 |
|
|
END |
406 |
|
|
|
407 |
|
|
CBOP |
408 |
|
|
C !SUBROUTINE: ANALYLIC_THETA |
409 |
|
|
C !INTERFACE: |
410 |
|
|
SUBROUTINE ANALYLIC_THETA( yLat, pNlev, |
411 |
|
|
O thetaLev, |
412 |
|
|
I kSize,myThid) |
413 |
|
|
|
414 |
|
|
C !DESCRIPTION: \bv |
415 |
|
|
C *==========================================================* |
416 |
|
|
C | SUBROUTINE ANALYLIC_THETA |
417 |
|
|
C | o Conpute analyticaly the potential temperature Theta |
418 |
|
|
C | as a function of Latitude (yLat) and normalized |
419 |
|
|
C | pressure pNlev. |
420 |
|
|
C | approximatively match the N-S symetric, zonal-mean and |
421 |
|
|
C | annual-mean NCEP climatology in the troposphere. |
422 |
|
|
C *==========================================================* |
423 |
|
|
C \ev |
424 |
|
|
|
425 |
|
|
C !USES: |
426 |
|
|
IMPLICIT NONE |
427 |
|
|
|
428 |
|
|
C == Global variables == |
429 |
|
|
#include "SIZE.h" |
430 |
|
|
#include "EEPARAMS.h" |
431 |
|
|
#include "PARAMS.h" |
432 |
|
|
|
433 |
|
|
C !INPUT/OUTPUT PARAMETERS: |
434 |
|
|
C == Routine arguments == |
435 |
|
|
C yLat :: latitude (degre) |
436 |
|
|
C pNlev :: normalized pressure levels |
437 |
|
|
C kSize :: dimension of pNlev & ANALYLIC_THETA |
438 |
|
|
C myThid :: Thread number for this instance of the routine |
439 |
|
|
INTEGER kSize |
440 |
|
|
_RL yLat |
441 |
|
|
_RL pNlev (kSize) |
442 |
|
|
_RL thetaLev(kSize) |
443 |
|
|
INTEGER myThid |
444 |
|
|
|
445 |
|
|
C !LOCAL VARIABLES: |
446 |
|
|
C == Local variables == |
447 |
|
|
INTEGER k |
448 |
|
|
_RL yyA, yyB, yyC, yyAd, yyBd, yyCd |
449 |
|
|
_RL cAtmp, cBtmp, ttdC |
450 |
|
|
_RL ppN0, ppN1, ppN2, ppN3a, ppN3b, ppN4 |
451 |
|
|
_RL ttp1, ttp2, ttp3, ttp4, ttp5 |
452 |
|
|
_RL yAtmp, yBtmp, yCtmp, yDtmp |
453 |
|
|
_RL ttp2y, ttp4y, a1tmp |
454 |
|
|
_RL ppl, ppm, pph, ppr |
455 |
|
|
|
456 |
|
|
CEOP |
457 |
|
|
|
458 |
|
|
DATA yyA , yyB , yyC , yyAd , yyBd , yyCd |
459 |
|
|
& / 45. _d 0, 65. _d 0, 65. _d 0, .9 _d 0, .9 _d 0, 10. _d 0 / |
460 |
|
|
DATA cAtmp , cBtmp , ttdC |
461 |
|
|
& / 2.6 _d 0, 1.5 _d 0, 3.3 _d 0 / |
462 |
|
|
DATA ppN0 , ppN1 , ppN2 , ppN3a , ppN3b , ppN4 |
463 |
|
|
& / .1 _d 0, .19 _d 0, .3 _d 0, .9 _d 0, .7 _d 0, .925 _d 0 / |
464 |
|
|
DATA ttp1 , ttp2 , ttp3 , ttp4 , ttp5 |
465 |
|
|
& / 350. _d 0, 342. _d 0, 307. _d 0, 301. _d 0, 257. _d 0 / |
466 |
|
|
|
467 |
|
|
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
468 |
|
|
|
469 |
|
|
yAtmp = ABS(yLat) - yyA |
470 |
|
|
yAtmp = yyA + MIN(0. _d 0,yAtmp/yyAd) + MAX(yAtmp, 0. _d 0) |
471 |
|
|
yAtmp = COS( deg2rad*MAX(yAtmp, 0. _d 0) ) |
472 |
|
|
yBtmp = ABS(yLat) - yyB |
473 |
|
|
yBtmp = yyB + yBtmp/yyBd |
474 |
|
|
yBtmp = COS( deg2rad*MAX( 0. _d 0, MIN(yBtmp,90. _d 0) ) ) |
475 |
|
|
yCtmp = ABS(yLat) - yyC |
476 |
|
|
yCtmp = MAX( 0. _d 0, 1. _d 0 - (yCtmp/yyCd)**2 ) |
477 |
|
|
yDtmp = ppN3a +(ppN3b - ppN3a)*yCtmp |
478 |
|
|
ttp2y = ttp3 + (ttp2-ttp3)*yAtmp**cAtmp |
479 |
|
|
ttp4y = ttp5 + (ttp4-ttp5)*yBtmp**cBtmp |
480 |
|
|
a1tmp = (ttp1-ttp2y)*ppN1*ppN2/(ppN2-ppN1) |
481 |
|
|
DO k=1,kSize |
482 |
|
|
ppl = MIN(pNlev(k),ppN1) |
483 |
|
|
ppm = MIN(MAX(pNlev(k),ppN1),ppN2) |
484 |
|
|
pph = MAX(pNlev(k),ppN2) |
485 |
|
|
ppr =( ppN0 + ABS(ppl-ppN0) - ppN1 )/(ppN2-ppN1) |
486 |
|
|
thetaLev(k) = |
487 |
|
|
& ( (1. _d 0 -ppr)*ttp1*ppN1**atm_kappa |
488 |
|
|
& + ppr*ttp2y*ppN2**atm_kappa |
489 |
|
|
& )*ppl**(-atm_kappa) |
490 |
|
|
& + a1tmp*(1. _d 0 /ppm - 1. _d 0/ppN1) |
491 |
|
|
& + (ttp4y-ttp2y)*(pph-ppN2)/(ppN4-ppN2) |
492 |
|
|
& + (ttdC+yCtmp)*MAX(0. _d 0,pNlev(k)-yDtmp)/(1-yDtmp) |
493 |
|
|
ENDDO |
494 |
jmc |
1.1 |
|
495 |
jmc |
1.4 |
C--------------------------------------------------- |
496 |
|
|
C matlab script, input: pN, yp=abs(yLat) |
497 |
|
|
C pN0=.1; pN1=.19 ; pN2=.3; pN4=.925; |
498 |
|
|
C pm=min(max(pN,pN1),pN2); pp=max(pN,pN2); |
499 |
|
|
C pl=min(pN,pN1); pr=(pN0+abs(pl-pN0)-pN1)/(pN2-pN1); |
500 |
|
|
C |
501 |
|
|
C yA=yp-45; yA=45+min(0,yA/.9)+max(0,yA); yA=max(0,yA); cosyA=cos(yA*rad) ; |
502 |
|
|
C yB=yp-65; yB=65+yB/.9; yB=min(max(0,yB),90); cosyB=cos(yB*rad) ; |
503 |
|
|
C tp1=350*ones(nyA,1); |
504 |
|
|
C tp2=307+(342-307)*cosyA.^2.6; |
505 |
|
|
C tp4=257+(301-257)*cosyB.^1.5; |
506 |
|
|
C a1=(tp1-tp2)*pN1*pN2/(pN2-pN1); |
507 |
|
|
C pF=max(0,1.-((yp-65)/10).^2); pT=.9+(0.7-.9)*pF; |
508 |
|
|
C |
509 |
|
|
C tA0=((1-pr(k))*tp1(j)*pN1^kappa+pr(k)*tp2(j)*pN2^kappa)*pl(k)^-kappa; |
510 |
|
|
C tA1=a1(j)*(1./pm(k)-1./pN1); |
511 |
|
|
C tA2=(tp4(j)-tp2(j))*(pp(k)-pN2)/(pN4-pN2); |
512 |
|
|
C tA3=(3.3+pF(j))*max(0,pN(k)-pT(j))/(1-pT(j)); |
513 |
|
|
C tAn(j,k)=tA0+tA1+tA2+tA3; |
514 |
|
|
C--------------------------------------------------- |
515 |
jmc |
1.1 |
|
516 |
|
|
RETURN |
517 |
|
|
END |