1 |
jmc |
1.10 |
C $Header: /u/gcmpack/MITgcm/pkg/aim_v23/phy_driver.F,v 1.9 2010/10/26 20:59:53 dfer Exp $ |
2 |
jmc |
1.1 |
C $Name: $ |
3 |
|
|
|
4 |
|
|
#include "AIM_OPTIONS.h" |
5 |
|
|
|
6 |
jmc |
1.8 |
CBOP |
7 |
|
|
C !ROUTINE: PHY_DRIVER |
8 |
|
|
C !INTERFACE: |
9 |
jmc |
1.7 |
SUBROUTINE PHY_DRIVER( tYear, usePkgDiag, |
10 |
|
|
I bi, bj, myTime, myIter, myThid ) |
11 |
jmc |
1.1 |
|
12 |
jmc |
1.8 |
C !DESCRIPTION: \bv |
13 |
jmc |
1.1 |
C------------------------ |
14 |
|
|
C-- SUBROUTINE PHYDRIVER (tYear, myTime, bi, bj, myThid ) |
15 |
|
|
C-- Purpose: stand-alone driver for physical parametrization routines |
16 |
|
|
C-- Input : TYEAR : fraction of year (0 = 1jan.00, 1 = 31dec.24) |
17 |
|
|
C-- grid-point model fields in common block: PHYGR1 |
18 |
|
|
C-- forcing fields in common blocks : LSMASK, FORFIX, FORCIN |
19 |
|
|
C-- Output : Diagnosed upper-air variables in common block: PHYGR2 |
20 |
|
|
C-- Diagnosed surface variables in common block: PHYGR3 |
21 |
|
|
C-- Physical param. tendencies in common block: PHYTEN |
22 |
|
|
C-- Surface and upper boundary fluxes in common block: FLUXES |
23 |
|
|
C------- |
24 |
jmc |
1.8 |
C Note: tendencies are not /dpFac here but later in AIM_AIM2DYN |
25 |
jmc |
1.1 |
C------- |
26 |
jmc |
1.8 |
C from SPEDDY code: (part of original code left with c_FM) |
27 |
|
|
C * S/R PHYPAR : except interp. dynamical Var. from Spectral of grid point |
28 |
|
|
C here dynamical var. are loaded within S/R AIM_DYN2AIM. |
29 |
|
|
C * S/R FORDATE: only the CALL SOL_OZ (done once / day in SPEEDY) |
30 |
|
|
C------------------------ |
31 |
|
|
C \ev |
32 |
jmc |
1.1 |
|
33 |
jmc |
1.8 |
C !USES: |
34 |
jmc |
1.1 |
IMPLICIT NONE |
35 |
|
|
|
36 |
jmc |
1.8 |
C == Global variables === |
37 |
jmc |
1.1 |
|
38 |
|
|
C-- size for MITgcm & Physics package : |
39 |
jmc |
1.7 |
#include "AIM_SIZE.h" |
40 |
jmc |
1.1 |
#include "EEPARAMS.h" |
41 |
jmc |
1.3 |
|
42 |
|
|
C-- Physics package |
43 |
|
|
#include "AIM_PARAMS.h" |
44 |
jmc |
1.1 |
#include "AIM_GRID.h" |
45 |
jmc |
1.10 |
#include "AIM_CO2.h" |
46 |
jmc |
1.1 |
|
47 |
|
|
C Constants + functions of sigma and latitude |
48 |
|
|
#include "com_physcon.h" |
49 |
|
|
|
50 |
|
|
C Model variables, tendencies and fluxes on gaussian grid |
51 |
|
|
#include "com_physvar.h" |
52 |
|
|
|
53 |
|
|
C Surface forcing fields (time-inv. or functions of seasonal cycle) |
54 |
|
|
#include "com_forcing.h" |
55 |
|
|
|
56 |
|
|
C Constants for forcing fields: |
57 |
|
|
#include "com_forcon.h" |
58 |
|
|
|
59 |
jmc |
1.8 |
C Radiation scheme variables |
60 |
jmc |
1.1 |
#include "com_radvar.h" |
61 |
|
|
|
62 |
jmc |
1.3 |
C Radiation constants |
63 |
|
|
#include "com_radcon.h" |
64 |
|
|
|
65 |
jmc |
1.1 |
C Logical flags |
66 |
|
|
c_FM include "com_lflags.h" |
67 |
|
|
|
68 |
jmc |
1.8 |
C !INPUT/OUTPUT PARAMETERS: |
69 |
|
|
C == Routine arguments == |
70 |
|
|
C tYear :: Fraction into year |
71 |
|
|
C usePkgDiag :: logical flag, true if using Diagnostics PKG |
72 |
|
|
C bi, bj :: Tile index |
73 |
|
|
C myTime :: Current time of simulation ( s ) |
74 |
|
|
C myIter :: Current iteration number in simulation |
75 |
|
|
C myThid :: Number of this instance of the routine |
76 |
|
|
_RL tYear |
77 |
jmc |
1.7 |
LOGICAL usePkgDiag |
78 |
|
|
INTEGER bi,bj |
79 |
jmc |
1.8 |
_RL myTime |
80 |
jmc |
1.7 |
INTEGER myIter, myThid |
81 |
jmc |
1.8 |
CEOP |
82 |
jmc |
1.1 |
|
83 |
|
|
#ifdef ALLOW_AIM |
84 |
jmc |
1.8 |
C !FUNCTIONS: |
85 |
|
|
C !LOCAL VARIABLES: |
86 |
|
|
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
87 |
|
|
C-- Local Variables originally (Speedy) in common bloc (com_physvar.h): |
88 |
|
|
C TG1 :: absolute temperature |
89 |
|
|
C QG1 :: specific humidity (g/kg) |
90 |
|
|
C VsurfSq :: Square of surface wind speed (grid position = as T,Q) |
91 |
|
|
C SE :: dry static energy <- replaced by Pot.Temp. |
92 |
|
|
C QSAT :: saturation specific humidity (g/kg) |
93 |
|
|
C PSG :: surface pressure (normalized) |
94 |
|
|
_RL TG1 (NGP,NLEV) |
95 |
|
|
_RL QG1 (NGP,NLEV) |
96 |
|
|
_RL VsurfSq(NGP) |
97 |
|
|
_RL SE (NGP,NLEV) |
98 |
|
|
_RL QSAT (NGP,NLEV) |
99 |
|
|
_RL PSG (NGP) |
100 |
|
|
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
101 |
jmc |
1.1 |
C-- Local variables: |
102 |
jmc |
1.10 |
C absLW_CO2 :: LW absorbtion in CO2 band (uniform value) |
103 |
|
|
C kGround :: Ground level index (2-dim) |
104 |
jmc |
1.3 |
C dpFac :: cell delta_P fraction (3-dim) |
105 |
|
|
C dTskin :: temp. correction for daily-cycle heating [K] |
106 |
jmc |
1.6 |
C T1s :: near-surface air temperature (from Pot.Temp) |
107 |
|
|
C DENVV :: surface flux (sens,lat.) coeff. (=Rho*|V|) [kg/m2/s] |
108 |
|
|
C Shf0 :: sensible heat flux over freezing surf. |
109 |
|
|
C dShf :: sensible heat flux derivative relative to surf. temp |
110 |
jmc |
1.3 |
C Evp0 :: evaporation computed over freezing surface (Ts=0.oC) |
111 |
jmc |
1.8 |
C dEvp :: evaporation derivative relative to surf. temp |
112 |
jmc |
1.3 |
C Slr0 :: upward long wave radiation over freezing surf. |
113 |
|
|
C dSlr :: upward long wave rad. derivative relative to surf. temp |
114 |
|
|
C sFlx :: net surface flux (+=down) function of surf. temp Ts: |
115 |
|
|
C 0: Flux(Ts=0.oC) ; 1: Flux(Ts^n) ; 2: d.Flux/d.Ts(Ts^n) |
116 |
jmc |
1.1 |
LOGICAL LRADSW |
117 |
|
|
INTEGER ICLTOP(NGP) |
118 |
|
|
INTEGER kGround(NGP) |
119 |
jmc |
1.10 |
_RL absLW_CO2 |
120 |
jmc |
1.1 |
_RL dpFac(NGP,NLEV) |
121 |
|
|
c_FM REAL RPS(NGP), ST4S(NGP) |
122 |
|
|
_RL ST4S(NGP) |
123 |
|
|
_RL PSG_1(NGP), RPS_1 |
124 |
jmc |
1.6 |
_RL dTskin(NGP), T1s(NGP), DENVV(NGP) |
125 |
|
|
_RL Shf0(NGP), dShf(NGP), Evp0(NGP), dEvp(NGP) |
126 |
|
|
_RL Slr0(NGP), dSlr(NGP), sFlx(NGP,0:2) |
127 |
dfer |
1.9 |
_RL UPSWG(NGP) |
128 |
jmc |
1.1 |
|
129 |
|
|
INTEGER J, K |
130 |
|
|
|
131 |
jmc |
1.6 |
#ifdef ALLOW_CLR_SKY_DIAG |
132 |
|
|
_RL dummyR(NGP) |
133 |
|
|
INTEGER dummyI(NGP) |
134 |
|
|
#endif |
135 |
jmc |
1.1 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
136 |
|
|
|
137 |
|
|
C-- 1. Compute grid-point fields |
138 |
|
|
|
139 |
|
|
C- 1.1 Convert model spectral variables to grid-point variables |
140 |
|
|
|
141 |
|
|
CALL AIM_DYN2AIM( |
142 |
|
|
O TG1, QG1, SE, VsurfSq, PSG, dpFac, kGround, |
143 |
|
|
I bi, bj, myTime, myIter, myThid ) |
144 |
|
|
|
145 |
|
|
C- 1.2 Compute thermodynamic variables |
146 |
|
|
|
147 |
|
|
C- 1.2.a Surface pressure (ps), 1/ps and surface temperature |
148 |
|
|
RPS_1 = 1. _d 0 |
149 |
|
|
DO J=1,NGP |
150 |
|
|
PSG_1(J)=1. _d 0 |
151 |
|
|
c_FM PSG(J)=EXP(PSLG1(J)) |
152 |
|
|
c_FM RPS(J)=1./PSG(J) |
153 |
|
|
ENDDO |
154 |
|
|
|
155 |
|
|
C 1.2.b Dry static energy |
156 |
|
|
C <= replaced by Pot.Temp in aim_dyn2aim |
157 |
|
|
c DO K=1,NLEV |
158 |
|
|
c DO J=1,NGP |
159 |
jmc |
1.8 |
c_FM SE(J,K)=CP*TG1(J,K)+PHIG1(J,K) |
160 |
jmc |
1.1 |
c ENDDO |
161 |
|
|
c ENDDO |
162 |
|
|
|
163 |
|
|
C 1.2.c Relative humidity and saturation spec. humidity |
164 |
|
|
|
165 |
|
|
DO K=1,NLEV |
166 |
|
|
c_FM CALL SHTORH (1,NGP,TG1(1,K),PSG,SIG(K),QG1(1,K), |
167 |
|
|
c_FM & RH(1,K),QSAT(1,K)) |
168 |
|
|
CALL SHTORH (1,NGP,TG1(1,K),PSG_1,SIG(K),QG1(1,K), |
169 |
|
|
O RH(1,K,myThid),QSAT(1,K), |
170 |
|
|
I myThid) |
171 |
|
|
ENDDO |
172 |
|
|
|
173 |
|
|
C-- 2. Precipitation |
174 |
|
|
|
175 |
|
|
C 2.1 Deep convection |
176 |
|
|
|
177 |
|
|
c_FM CALL CONVMF (PSG,SE,QG1,QSAT, |
178 |
|
|
c_FM & ICLTOP,CBMF,PRECNV,TT_CNV,QT_CNV) |
179 |
|
|
CALL CONVMF (PSG,dpFac,SE,QG1,QSAT, |
180 |
|
|
O ICLTOP,CBMF(1,myThid),PRECNV(1,myThid), |
181 |
|
|
O TT_CNV(1,1,myThid),QT_CNV(1,1,myThid), |
182 |
|
|
I kGround,bi,bj,myThid) |
183 |
|
|
|
184 |
|
|
DO K=2,NLEV |
185 |
|
|
DO J=1,NGP |
186 |
|
|
TT_CNV(J,K,myThid)=TT_CNV(J,K,myThid)*RPS_1*GRDSCP(K) |
187 |
|
|
QT_CNV(J,K,myThid)=QT_CNV(J,K,myThid)*RPS_1*GRDSIG(K) |
188 |
|
|
ENDDO |
189 |
|
|
ENDDO |
190 |
|
|
|
191 |
|
|
C 2.2 Large-scale condensation |
192 |
|
|
|
193 |
|
|
c_FM CALL LSCOND (PSG,QG1,QSAT, |
194 |
|
|
c_FM & PRECLS,TT_LSC,QT_LSC) |
195 |
|
|
CALL LSCOND (PSG,dpFac,QG1,QSAT, |
196 |
|
|
O PRECLS(1,myThid),TT_LSC(1,1,myThid), |
197 |
|
|
O QT_LSC(1,1,myThid), |
198 |
|
|
I kGround,bi,bj,myThid) |
199 |
|
|
|
200 |
jmc |
1.3 |
IF ( aim_energPrecip ) THEN |
201 |
|
|
C 2.3 Snow Precipitation (update TT_CNV & TT_LSC) |
202 |
jmc |
1.8 |
CALL SNOW_PRECIP ( |
203 |
jmc |
1.3 |
I PSG, dpFac, SE, ICLTOP, |
204 |
|
|
I PRECNV(1,myThid), QT_CNV(1,1,myThid), |
205 |
|
|
I PRECLS(1,myThid), QT_LSC(1,1,myThid), |
206 |
|
|
U TT_CNV(1,1,myThid), TT_LSC(1,1,myThid), |
207 |
|
|
O EnPrec(1,myThid), |
208 |
|
|
I kGround,bi,bj,myThid) |
209 |
|
|
ELSE |
210 |
|
|
DO J=1,NGP |
211 |
|
|
EnPrec(J,myThid) = 0. _d 0 |
212 |
|
|
ENDDO |
213 |
|
|
ENDIF |
214 |
|
|
|
215 |
jmc |
1.1 |
C-- 3. Radiation (shortwave and longwave) and surface fluxes |
216 |
|
|
|
217 |
|
|
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
218 |
|
|
C --> from FORDATE (in SPEEDY) : |
219 |
|
|
|
220 |
|
|
C 3.0 Compute Incomming shortwave rad. (from FORDATE in SPEEDY) |
221 |
|
|
|
222 |
|
|
c_FM CALL SOL_OZ (SOLC,TYEAR) |
223 |
|
|
CALL SOL_OZ (SOLC,tYear, snLat(1,myThid), csLat(1,myThid), |
224 |
|
|
O FSOL, OZONE, OZUPP, ZENIT, STRATZ, |
225 |
|
|
I bi,bj,myThid) |
226 |
|
|
|
227 |
|
|
C <-- from FORDATE (in SPEEDY). |
228 |
|
|
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
229 |
|
|
|
230 |
|
|
C 3.1 Compute shortwave tendencies and initialize lw transmissivity |
231 |
|
|
|
232 |
jmc |
1.10 |
C Set LW absorption in CO2 band |
233 |
|
|
IF ( aim_select_pCO2.EQ.1 .OR. aim_select_pCO2.EQ.3 ) THEN |
234 |
|
|
absLW_CO2 = ABLCO2 |
235 |
|
|
& + aim_abs_pCO2*LOG( aim_pCO2/aim_ref_pCO2 ) |
236 |
|
|
ELSE |
237 |
|
|
absLW_CO2 = ABLCO2 |
238 |
|
|
ENDIF |
239 |
|
|
|
240 |
jmc |
1.1 |
C The sw radiation may be called at selected time steps |
241 |
|
|
LRADSW = .TRUE. |
242 |
jmc |
1.8 |
|
243 |
jmc |
1.1 |
IF (LRADSW) THEN |
244 |
jmc |
1.8 |
|
245 |
jmc |
1.1 |
c_FM CALL RADSW (PSG,QG1,RH,ALB1, |
246 |
|
|
c_FM & ICLTOP,CLOUDC,TSR,SSR,TT_RSW) |
247 |
jmc |
1.6 |
ICLTOP(1) = 1 |
248 |
jmc |
1.3 |
CALL RADSW (PSG,dpFac,QG1,RH(1,1,myThid),ALB1(1,0,myThid), |
249 |
jmc |
1.1 |
I FSOL, OZONE, OZUPP, ZENIT, STRATZ, |
250 |
|
|
O TAU2, STRATC, |
251 |
|
|
O ICLTOP,CLOUDC(1,myThid), |
252 |
dfer |
1.9 |
O TSR(1,myThid),SSR(1,0,myThid), |
253 |
|
|
O UPSWG,TT_RSW(1,1,myThid), |
254 |
jmc |
1.10 |
I absLW_CO2, kGround, bi, bj, myThid ) |
255 |
jmc |
1.8 |
|
256 |
jmc |
1.1 |
DO J=1,NGP |
257 |
|
|
CLTOP(J,myThid)=SIGH(ICLTOP(J)-1)*PSG_1(J) |
258 |
|
|
ENDDO |
259 |
jmc |
1.8 |
|
260 |
jmc |
1.1 |
DO K=1,NLEV |
261 |
|
|
DO J=1,NGP |
262 |
|
|
TT_RSW(J,K,myThid)=TT_RSW(J,K,myThid)*RPS_1*GRDSCP(K) |
263 |
|
|
ENDDO |
264 |
|
|
ENDDO |
265 |
jmc |
1.8 |
|
266 |
dfer |
1.9 |
#ifdef ALLOW_DIAGNOSTICS |
267 |
|
|
IF ( usePkgDiag ) THEN |
268 |
|
|
CALL DIAGNOSTICS_FILL( UPSWG, |
269 |
|
|
& 'UPSWG ', 1, 1 , 3,bi,bj, myThid ) |
270 |
|
|
ENDIF |
271 |
|
|
#endif |
272 |
|
|
|
273 |
jmc |
1.1 |
ENDIF |
274 |
|
|
|
275 |
|
|
C 3.2 Compute downward longwave fluxes |
276 |
jmc |
1.8 |
|
277 |
jmc |
1.1 |
c_FM CALL RADLW (-1,TG1,TS,ST4S, |
278 |
|
|
c_FM & OLR,SLR,TT_RLW) |
279 |
|
|
CALL RADLW (-1,TG1,TS(1,myThid),ST4S, |
280 |
|
|
& OZUPP, STRATC, TAU2, FLUX, ST4A, |
281 |
jmc |
1.3 |
O OLR(1,myThid),SLR(1,0,myThid),TT_RLW(1,1,myThid), |
282 |
jmc |
1.1 |
I kGround,bi,bj,myThid) |
283 |
|
|
|
284 |
jmc |
1.3 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
285 |
jmc |
1.1 |
C 3.3. Compute surface fluxes and land skin temperature |
286 |
jmc |
1.8 |
|
287 |
jmc |
1.1 |
c_FM CALL SUFLUX (PSG,UG1,VG1,TG1,QG1,RH,PHIG1, |
288 |
|
|
c_FM & PHIS0,FMASK1,STL1,SST1,SOILW1,SSR,SLR, |
289 |
|
|
c_FM & USTR,VSTR,SHF,EVAP,ST4S, |
290 |
jmc |
1.8 |
c_FM & TS,TSKIN,U0,V0,T0,Q0) |
291 |
jmc |
1.3 |
CALL SUFLUX_PREP( |
292 |
|
|
I PSG, TG1, QG1, RH(1,1,myThid), SE, VsurfSq, |
293 |
jmc |
1.1 |
I WVSurf(1,myThid),csLat(1,myThid),fOrogr(1,myThid), |
294 |
jmc |
1.3 |
I FMASK1(1,1,myThid),STL1(1,myThid),SST1(1,myThid), |
295 |
|
|
I sti1(1,myThid), SSR(1,0,myThid), |
296 |
jmc |
1.6 |
O SPEED0(1,myThid),DRAG(1,0,myThid),DENVV, |
297 |
|
|
O dTskin,T1s,T0(1,myThid),Q0(1,myThid), |
298 |
jmc |
1.1 |
I kGround,bi,bj,myThid) |
299 |
|
|
|
300 |
jmc |
1.3 |
CALL SUFLUX_LAND ( |
301 |
|
|
I PSG, FMASK1(1,1,myThid), EMISFC, |
302 |
|
|
I STL1(1,myThid), dTskin, |
303 |
|
|
I SOILW1(1,myThid), SSR(1,1,myThid), SLR(1,0,myThid), |
304 |
jmc |
1.6 |
I T1s, T0(1,myThid), Q0(1,myThid), DENVV, |
305 |
jmc |
1.3 |
O SHF(1,1,myThid), EVAP(1,1,myThid), SLR(1,1,myThid), |
306 |
jmc |
1.6 |
O Shf0, dShf, Evp0, dEvp, Slr0, dSlr, sFlx, |
307 |
jmc |
1.3 |
O TS(1,myThid), TSKIN(1,myThid), |
308 |
|
|
I bi,bj,myThid) |
309 |
jmc |
1.8 |
#ifdef ALLOW_LAND |
310 |
jmc |
1.3 |
CALL AIM_LAND_IMPL( |
311 |
jmc |
1.5 |
I FMASK1(1,1,myThid), dTskin, |
312 |
jmc |
1.6 |
I Shf0, dShf, Evp0, dEvp, Slr0, dSlr, |
313 |
|
|
U sFlx, STL1(1,myThid), |
314 |
|
|
U SHF(1,1,myThid), EVAP(1,1,myThid), SLR(1,1,myThid), |
315 |
|
|
O dTsurf(1,1,myThid), |
316 |
jmc |
1.3 |
I bi, bj, myTime, myIter, myThid) |
317 |
|
|
#endif /* ALLOW_LAND */ |
318 |
|
|
|
319 |
|
|
CALL SUFLUX_OCEAN( |
320 |
|
|
I PSG, FMASK1(1,2,myThid), |
321 |
|
|
I SST1(1,myThid), |
322 |
|
|
I SSR(1,2,myThid), SLR(1,0,myThid), |
323 |
jmc |
1.6 |
O T1s, T0(1,myThid), Q0(1,myThid), DENVV, |
324 |
jmc |
1.3 |
O SHF(1,2,myThid), EVAP(1,2,myThid), SLR(1,2,myThid), |
325 |
|
|
I bi,bj,myThid) |
326 |
|
|
|
327 |
|
|
IF ( aim_splitSIOsFx ) THEN |
328 |
|
|
CALL SUFLUX_SICE ( |
329 |
|
|
I PSG, FMASK1(1,3,myThid), EMISFC, |
330 |
|
|
I STI1(1,myThid), dTskin, |
331 |
|
|
I SSR(1,3,myThid), SLR(1,0,myThid), |
332 |
jmc |
1.6 |
I T1s, T0(1,myThid), Q0(1,myThid), DENVV, |
333 |
jmc |
1.3 |
O SHF(1,3,myThid), EVAP(1,3,myThid), SLR(1,3,myThid), |
334 |
jmc |
1.6 |
O Shf0, dShf, Evp0, dEvp, Slr0, dSlr, sFlx, |
335 |
jmc |
1.3 |
O TS(1,myThid), TSKIN(1,myThid), |
336 |
|
|
I bi,bj,myThid) |
337 |
jmc |
1.8 |
#ifdef ALLOW_THSICE |
338 |
jmc |
1.4 |
CALL AIM_SICE_IMPL( |
339 |
|
|
I FMASK1(1,3,myThid), SSR(1,3,myThid), sFlx, |
340 |
jmc |
1.6 |
I Shf0, dShf, Evp0, dEvp, Slr0, dSlr, |
341 |
|
|
U STI1(1,myThid), |
342 |
|
|
U SHF(1,3,myThid), EVAP(1,3,myThid), SLR(1,3,myThid), |
343 |
|
|
O dTsurf(1,3,myThid), |
344 |
jmc |
1.4 |
I bi, bj, myTime, myIter, myThid) |
345 |
|
|
#endif /* ALLOW_THSICE */ |
346 |
jmc |
1.3 |
ELSE |
347 |
|
|
DO J=1,NGP |
348 |
jmc |
1.6 |
SHF (J,3,myThid) = 0. _d 0 |
349 |
jmc |
1.3 |
EVAP(J,3,myThid) = 0. _d 0 |
350 |
|
|
SLR (J,3,myThid) = 0. _d 0 |
351 |
|
|
ENDDO |
352 |
|
|
ENDIF |
353 |
|
|
|
354 |
|
|
CALL SUFLUX_POST( |
355 |
jmc |
1.8 |
I FMASK1(1,1,myThid), EMISFC, |
356 |
|
|
I STL1(1,myThid), SST1(1,myThid), sti1(1,myThid), |
357 |
jmc |
1.3 |
I dTskin, SLR(1,0,myThid), |
358 |
jmc |
1.6 |
I T0(1,myThid), Q0(1,myThid), DENVV, |
359 |
jmc |
1.8 |
U DRAG(1,0,myThid), SHF(1,0,myThid), |
360 |
jmc |
1.3 |
U EVAP(1,0,myThid), SLR(1,1,myThid), |
361 |
|
|
O ST4S, TS(1,myThid), TSKIN(1,myThid), |
362 |
|
|
I bi,bj,myThid) |
363 |
jmc |
1.7 |
|
364 |
|
|
#ifdef ALLOW_DIAGNOSTICS |
365 |
|
|
IF ( usePkgDiag ) THEN |
366 |
|
|
CALL DIAGNOSTICS_FILL( SLR(1,0,myThid), |
367 |
|
|
& 'DWNLWG ', 1, 1 , 3,bi,bj, myThid ) |
368 |
|
|
ENDIF |
369 |
|
|
#endif |
370 |
jmc |
1.3 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
371 |
|
|
|
372 |
jmc |
1.1 |
C 3.4 Compute upward longwave fluxes, convert them to tendencies |
373 |
|
|
C and add shortwave tendencies |
374 |
|
|
|
375 |
|
|
c_FM CALL RADLW (1,TG1,TS,ST4S, |
376 |
|
|
c_FM & OLR,SLR,TT_RLW) |
377 |
|
|
CALL RADLW (1,TG1,TS(1,myThid),ST4S, |
378 |
|
|
& OZUPP, STRATC, TAU2, FLUX, ST4A, |
379 |
jmc |
1.3 |
O OLR(1,myThid),SLR(1,0,myThid),TT_RLW(1,1,myThid), |
380 |
jmc |
1.1 |
I kGround,bi,bj,myThid) |
381 |
jmc |
1.8 |
|
382 |
jmc |
1.1 |
DO K=1,NLEV |
383 |
|
|
DO J=1,NGP |
384 |
|
|
TT_RLW(J,K,myThid)=TT_RLW(J,K,myThid)*RPS_1*GRDSCP(K) |
385 |
|
|
c_FM TTEND (J,K)=TTEND(J,K)+TT_RSW(J,K)+TT_RLW(J,K) |
386 |
|
|
ENDDO |
387 |
jmc |
1.8 |
ENDDO |
388 |
jmc |
1.1 |
|
389 |
jmc |
1.6 |
#ifdef ALLOW_CLR_SKY_DIAG |
390 |
|
|
C 3.5 Compute clear-sky radiation (for diagnostics only) |
391 |
|
|
IF ( aim_clrSkyDiag ) THEN |
392 |
jmc |
1.8 |
|
393 |
jmc |
1.6 |
C 3.5.1 Compute shortwave tendencies |
394 |
|
|
dummyI(1) = -1 |
395 |
|
|
CALL RADSW (PSG,dpFac,QG1,RH(1,1,myThid),ALB1(1,0,myThid), |
396 |
|
|
I FSOL, OZONE, OZUPP, ZENIT, STRATZ, |
397 |
|
|
O TAU2, STRATC, |
398 |
|
|
O dummyI, dummyR, |
399 |
dfer |
1.9 |
O TSWclr(1,myThid), SSWclr(1,myThid), UPSWG, TT_SWclr(1,1,myThid), |
400 |
jmc |
1.10 |
I absLW_CO2, kGround, bi, bj, myThid ) |
401 |
jmc |
1.8 |
|
402 |
dfer |
1.9 |
#ifdef ALLOW_DIAGNOSTICS |
403 |
|
|
IF ( usePkgDiag ) THEN |
404 |
|
|
CALL DIAGNOSTICS_FILL( UPSWG, |
405 |
|
|
& 'UPSWGclr', 1, 1 , 3,bi,bj, myThid ) |
406 |
|
|
ENDIF |
407 |
|
|
#endif |
408 |
|
|
|
409 |
jmc |
1.6 |
C 3.5.2 Compute downward longwave fluxes |
410 |
jmc |
1.8 |
|
411 |
jmc |
1.6 |
CALL RADLW (-1,TG1,TS(1,myThid),ST4S, |
412 |
|
|
& OZUPP, STRATC, TAU2, FLUX, ST4A, |
413 |
|
|
O OLWclr(1,myThid), SLWclr(1,myThid), TT_LWclr(1,1,myThid), |
414 |
|
|
I kGround,bi,bj,myThid) |
415 |
|
|
|
416 |
|
|
C 3.5.3 Compute upward longwave fluxes, convert them to tendencies |
417 |
|
|
|
418 |
|
|
CALL RADLW (1,TG1,TS(1,myThid),ST4S, |
419 |
|
|
& OZUPP, STRATC, TAU2, FLUX, ST4A, |
420 |
|
|
O OLWclr(1,myThid), SLWclr(1,myThid), TT_LWclr(1,1,myThid), |
421 |
|
|
I kGround,bi,bj,myThid) |
422 |
jmc |
1.8 |
|
423 |
jmc |
1.6 |
DO K=1,NLEV |
424 |
|
|
DO J=1,NGP |
425 |
|
|
TT_SWclr(J,K,myThid)=TT_SWclr(J,K,myThid)*RPS_1*GRDSCP(K) |
426 |
|
|
TT_LWclr(J,K,myThid)=TT_LWclr(J,K,myThid)*RPS_1*GRDSCP(K) |
427 |
|
|
ENDDO |
428 |
jmc |
1.8 |
ENDDO |
429 |
jmc |
1.6 |
|
430 |
|
|
ENDIF |
431 |
|
|
#endif /* ALLOW_CLR_SKY_DIAG */ |
432 |
|
|
|
433 |
jmc |
1.1 |
C-- 4. PBL interactions with lower troposphere |
434 |
|
|
|
435 |
|
|
C 4.1 Vertical diffusion and shallow convection |
436 |
jmc |
1.8 |
|
437 |
jmc |
1.1 |
c_FM CALL VDIFSC (UG1,VG1,SE,RH,QG1,QSAT,PHIG1, |
438 |
|
|
c_FM & UT_PBL,VT_PBL,TT_PBL,QT_PBL) |
439 |
|
|
CALL VDIFSC (dpFac, SE, RH(1,1,myThid), QG1, QSAT, |
440 |
|
|
O TT_PBL(1,1,myThid),QT_PBL(1,1,myThid), |
441 |
|
|
I kGround,bi,bj,myThid) |
442 |
jmc |
1.8 |
|
443 |
jmc |
1.1 |
C 4.2 Add tendencies due to surface fluxes |
444 |
jmc |
1.8 |
|
445 |
jmc |
1.1 |
DO J=1,NGP |
446 |
|
|
c_FM UT_PBL(J,NLEV)=UT_PBL(J,NLEV)+USTR(J,3)*RPS(J)*GRDSIG(NLEV) |
447 |
|
|
c_FM VT_PBL(J,NLEV)=VT_PBL(J,NLEV)+VSTR(J,3)*RPS(J)*GRDSIG(NLEV) |
448 |
|
|
c_FM TT_PBL(J,NLEV)=TT_PBL(J,NLEV)+ SHF(J,3)*RPS(J)*GRDSCP(NLEV) |
449 |
|
|
c_FM QT_PBL(J,NLEV)=QT_PBL(J,NLEV)+EVAP(J,3)*RPS(J)*GRDSIG(NLEV) |
450 |
|
|
K = kGround(J) |
451 |
|
|
IF ( K.GT.0 ) THEN |
452 |
|
|
TT_PBL(J,K,myThid) = TT_PBL(J,K,myThid) |
453 |
jmc |
1.3 |
& + SHF(J,0,myThid) *RPS_1*GRDSCP(K) |
454 |
jmc |
1.1 |
QT_PBL(J,K,myThid) = QT_PBL(J,K,myThid) |
455 |
jmc |
1.3 |
& + EVAP(J,0,myThid)*RPS_1*GRDSIG(K) |
456 |
jmc |
1.1 |
ENDIF |
457 |
|
|
ENDDO |
458 |
jmc |
1.8 |
|
459 |
jmc |
1.1 |
c_FM DO K=1,NLEV |
460 |
|
|
c_FM DO J=1,NGP |
461 |
|
|
c_FM UTEND(J,K)=UTEND(J,K)+UT_PBL(J,K) |
462 |
|
|
c_FM VTEND(J,K)=VTEND(J,K)+VT_PBL(J,K) |
463 |
|
|
c_FM TTEND(J,K)=TTEND(J,K)+TT_PBL(J,K) |
464 |
|
|
c_FM QTEND(J,K)=QTEND(J,K)+QT_PBL(J,K) |
465 |
|
|
c_FM ENDDO |
466 |
jmc |
1.8 |
c_FM ENDDO |
467 |
jmc |
1.1 |
|
468 |
jmc |
1.8 |
#endif /* ALLOW_AIM */ |
469 |
jmc |
1.1 |
|
470 |
|
|
RETURN |
471 |
|
|
END |