1 |
C $Header: /u/gcmpack/MITgcm/pkg/exf/exf_bulkformulae.F,v 1.23 2010/05/21 10:08:44 mlosch Exp $ |
2 |
C $Name: $ |
3 |
|
4 |
#include "EXF_OPTIONS.h" |
5 |
|
6 |
CBOP |
7 |
C !ROUTINE: EXF_BULKFORMULAE |
8 |
C !INTERFACE: |
9 |
SUBROUTINE EXF_BULKFORMULAE( myTime, myIter, myThid ) |
10 |
|
11 |
C !DESCRIPTION: \bv |
12 |
C *==========================================================* |
13 |
C | SUBROUTINE EXF_BULKFORMULAE |
14 |
C | o Calculate bulk formula fluxes over open ocean |
15 |
C | either following |
16 |
C | Large and Pond, 1981 & 1982 |
17 |
C | or (if defined ALLOW_BULK_LARGEYEAGER04) |
18 |
C | Large and Yeager, 2004, NCAR/TN-460+STR. |
19 |
C *==========================================================* |
20 |
C \ev |
21 |
C |
22 |
C |
23 |
C NOTES: |
24 |
C ====== |
25 |
C |
26 |
C See EXF_OPTIONS.h for a description of the various possible |
27 |
C ocean-model forcing configurations. |
28 |
C |
29 |
C The bulk formulae of pkg/exf are not valid for sea-ice covered |
30 |
C oceans but they can be used in combination with a sea-ice model, |
31 |
C for example, pkg/seaice, to specify open water flux contributions. |
32 |
C |
33 |
C ================================================================== |
34 |
C |
35 |
C for Large and Pond, 1981 & 1982 |
36 |
C |
37 |
C The calculation of the bulk surface fluxes has been adapted from |
38 |
C the NCOM model which uses the formulae given in Large and Pond |
39 |
C (1981 & 1982 ) |
40 |
C |
41 |
C |
42 |
C Header taken from NCOM version: ncom1.4.1 |
43 |
C ----------------------------------------- |
44 |
C |
45 |
C Following procedures and coefficients in Large and Pond |
46 |
C (1981 ; 1982) |
47 |
C |
48 |
C Output: Bulk estimates of the turbulent surface fluxes. |
49 |
C ------- |
50 |
C |
51 |
C hs - sensible heat flux (W/m^2), into ocean |
52 |
C hl - latent heat flux (W/m^2), into ocean |
53 |
C |
54 |
C Input: |
55 |
C ------ |
56 |
C |
57 |
C us - mean wind speed (m/s) at height hu (m) |
58 |
C th - mean air temperature (K) at height ht (m) |
59 |
C qh - mean air humidity (kg/kg) at height hq (m) |
60 |
C sst - sea surface temperature (K) |
61 |
C tk0 - Kelvin temperature at 0 Celsius (K) |
62 |
C |
63 |
C Assume 1) a neutral 10m drag coefficient = |
64 |
C |
65 |
C cdn = .0027/u10 + .000142 + .0000764 u10 |
66 |
C |
67 |
C 2) a neutral 10m stanton number = |
68 |
C |
69 |
C ctn = .0327 sqrt(cdn), unstable |
70 |
C ctn = .0180 sqrt(cdn), stable |
71 |
C |
72 |
C 3) a neutral 10m dalton number = |
73 |
C |
74 |
C cen = .0346 sqrt(cdn) |
75 |
C |
76 |
C 4) the saturation humidity of air at |
77 |
C |
78 |
C t(k) = exf_BulkqSat(t) (kg/m^3) |
79 |
C |
80 |
C Note: 1) here, tstar = <wt>/u*, and qstar = <wq>/u*. |
81 |
C 2) wind speeds should all be above a minimum speed, |
82 |
C say 0.5 m/s. |
83 |
C 3) with optional iteration loop, niter=3, should suffice. |
84 |
C 4) this version is for analyses inputs with hu = 10m and |
85 |
C ht = hq. |
86 |
C 5) sst enters in Celsius. |
87 |
C |
88 |
C ================================================================== |
89 |
C |
90 |
C started: Christian Eckert eckert@mit.edu 27-Aug-1999 |
91 |
C |
92 |
C changed: Christian Eckert eckert@mit.edu 14-Jan-2000 |
93 |
C - restructured the original version in order to have a |
94 |
C better interface to the MITgcmUV. |
95 |
C |
96 |
C Christian Eckert eckert@mit.edu 12-Feb-2000 |
97 |
C - Changed Routine names (package prefix: exf_) |
98 |
C |
99 |
C Patrick Heimbach, heimbach@mit.edu 04-May-2000 |
100 |
C - changed the handling of precip and sflux with respect |
101 |
C to CPP options ALLOW_BULKFORMULAE and ALLOW_ATM_TEMP |
102 |
C - included some CPP flags ALLOW_BULKFORMULAE to make |
103 |
C sure ALLOW_ATM_TEMP, ALLOW_ATM_WIND are used only in |
104 |
C conjunction with defined ALLOW_BULKFORMULAE |
105 |
C - statement functions discarded |
106 |
C |
107 |
C Ralf.Giering@FastOpt.de 25-Mai-2000 |
108 |
C - total rewrite using new subroutines |
109 |
C |
110 |
C Detlef Stammer: include river run-off. Nov. 21, 2001 |
111 |
C |
112 |
C heimbach@mit.edu, 10-Jan-2002 |
113 |
C - changes to enable field swapping |
114 |
C |
115 |
C mods for pkg/seaice: menemenlis@jpl.nasa.gov 20-Dec-2002 |
116 |
C |
117 |
C martin.losch@awi.de: merged with exf_bulk_largeyeager04, 21-May-2010 |
118 |
C |
119 |
C ================================================================== |
120 |
C |
121 |
C for Large and Yeager, 2004 |
122 |
C |
123 |
C === Turbulent Fluxes : |
124 |
C * use the approach "B": shift coeff to height & stability of the |
125 |
C atmosphere state (instead of "C": shift temp & humid to the height |
126 |
C of wind, then shift the coeff to this height & stability of the atmos). |
127 |
C * similar to EXF (except over sea-ice) ; default parameter values |
128 |
C taken from Large & Yeager. |
129 |
C * assume that Qair & Tair inputs are from the same height (zq=zt) |
130 |
C * formulae in short: |
131 |
C wind stress = (ust,vst) = rhoA * Cd * Ws * (del.u,del.v) |
132 |
C Sensib Heat flux = fsha = rhoA * Ch * Ws * del.T * CpAir |
133 |
C Latent Heat flux = flha = rhoA * Ce * Ws * del.Q * Lvap |
134 |
C = -Evap * Lvap |
135 |
C with Ws = wind speed = sqrt(del.u^2 +del.v^2) ; |
136 |
C del.T = Tair - Tsurf ; del.Q = Qair - Qsurf ; |
137 |
C Cd,Ch,Ce = drag coefficient, Stanton number and Dalton number |
138 |
C respectively [no-units], function of height & stability |
139 |
|
140 |
C !USES: |
141 |
IMPLICIT NONE |
142 |
C === Global variables === |
143 |
#include "EEPARAMS.h" |
144 |
#include "SIZE.h" |
145 |
#include "PARAMS.h" |
146 |
#include "DYNVARS.h" |
147 |
#include "GRID.h" |
148 |
|
149 |
#include "EXF_PARAM.h" |
150 |
#include "EXF_FIELDS.h" |
151 |
#include "EXF_CONSTANTS.h" |
152 |
|
153 |
#ifdef ALLOW_AUTODIFF_TAMC |
154 |
#include "tamc.h" |
155 |
#endif |
156 |
|
157 |
C !INPUT/OUTPUT PARAMETERS: |
158 |
C input: |
159 |
C myTime :: Current time in simulation |
160 |
C myIter :: Current iteration number in simulation |
161 |
C myThid :: My Thread Id number |
162 |
_RL myTime |
163 |
INTEGER myIter |
164 |
INTEGER myThid |
165 |
C output: |
166 |
CEOP |
167 |
|
168 |
#ifdef ALLOW_BULKFORMULAE |
169 |
C == external Functions |
170 |
|
171 |
C == Local variables == |
172 |
C i,j :: grid point indices |
173 |
C bi,bj :: tile indices |
174 |
C ssq :: saturation specific humidity [kg/kg] in eq. with Sea-Surface water |
175 |
INTEGER i,j,bi,bj |
176 |
|
177 |
_RL czol |
178 |
_RL Tsf ! surface temperature [K] |
179 |
_RL wsm ! limited wind speed [m/s] (> umin) |
180 |
_RL t0 ! virtual temperature [K] |
181 |
C these need to be 2D-arrays for vectorizing code |
182 |
_RL tstar (1:sNx,1:sNy) ! turbulent temperature scale [K] |
183 |
_RL qstar (1:sNx,1:sNy) ! turbulent humidity scale [kg/kg] |
184 |
_RL ustar (1:sNx,1:sNy) ! friction velocity [m/s] |
185 |
_RL tau (1:sNx,1:sNy) ! surface stress coef = rhoA * Ws * sqrt(Cd) |
186 |
_RL rdn (1:sNx,1:sNy) ! neutral, zref (=10m) values of rd |
187 |
_RL rd (1:sNx,1:sNy) ! = sqrt(Cd) [-] |
188 |
_RL delq (1:sNx,1:sNy) ! specific humidity difference [kg/kg] |
189 |
_RL deltap(1:sNx,1:sNy) |
190 |
C |
191 |
_RL dzTmp |
192 |
_RL SSTtmp |
193 |
_RL ssq |
194 |
_RL re ! = Ce / sqrt(Cd) [-] |
195 |
_RL rh ! = Ch / sqrt(Cd) [-] |
196 |
_RL ren, rhn ! neutral, zref (=10m) values of re, rh |
197 |
_RL usn, usm |
198 |
_RL stable ! = 1 if stable ; = 0 if unstable |
199 |
_RL huol ! stability parameter at zwd [-] (=z/Monin-Obuklov length) |
200 |
_RL htol ! stability parameter at zth [-] |
201 |
_RL hqol |
202 |
_RL x ! stability function [-] |
203 |
_RL xsq ! = x^2 [-] |
204 |
_RL psimh ! momentum stability function |
205 |
_RL psixh ! latent & sensib. stability function |
206 |
_RL zwln ! = log(zwd/zref) |
207 |
_RL ztln ! = log(zth/zref) |
208 |
_RL windStress ! surface wind-stress (@ grid cell center) |
209 |
_RL tmpbulk |
210 |
_RL recip_rhoConstFresh |
211 |
INTEGER ksrf, ksrfp1 |
212 |
INTEGER iter |
213 |
C solve4Stress :: if F, by-pass momentum turb. flux (wind-stress) calculation |
214 |
LOGICAL solve4Stress |
215 |
#ifdef ALLOW_ATM_WIND |
216 |
PARAMETER ( solve4Stress = .TRUE. ) |
217 |
#endif |
218 |
|
219 |
#ifdef ALLOW_AUTODIFF_TAMC |
220 |
integer ikey_1 |
221 |
integer ikey_2 |
222 |
#endif |
223 |
|
224 |
#ifndef ALLOW_ATM_WIND |
225 |
#ifdef ALLOW_BULK_LARGEYEAGER04 |
226 |
solve4Stress = wspeedfile .NE. ' ' |
227 |
#else |
228 |
solve4Stress = .FALSE. |
229 |
#endif |
230 |
#endif |
231 |
|
232 |
C-- Set surface parameters : |
233 |
zwln = LOG(hu/zref) |
234 |
ztln = LOG(ht/zref) |
235 |
czol = hu*karman*gravity_mks |
236 |
C-- abbreviation |
237 |
recip_rhoConstFresh = 1. _d 0/rhoConstFresh |
238 |
|
239 |
c Loop over tiles. |
240 |
#ifdef ALLOW_AUTODIFF_TAMC |
241 |
C-- HPF directive to help TAMC |
242 |
CHPF$ INDEPENDENT |
243 |
#endif |
244 |
DO bj = myByLo(myThid),myByHi(myThid) |
245 |
#ifdef ALLOW_AUTODIFF_TAMC |
246 |
C-- HPF directive to help TAMC |
247 |
CHPF$ INDEPENDENT |
248 |
#endif |
249 |
DO bi = myBxLo(myThid),myBxHi(myThid) |
250 |
ksrf = 1 |
251 |
ksrfp1 = 2 |
252 |
DO j = 1,sNy |
253 |
DO i = 1,sNx |
254 |
|
255 |
#ifdef ALLOW_AUTODIFF_TAMC |
256 |
act1 = bi - myBxLo(myThid) |
257 |
max1 = myBxHi(myThid) - myBxLo(myThid) + 1 |
258 |
act2 = bj - myByLo(myThid) |
259 |
max2 = myByHi(myThid) - myByLo(myThid) + 1 |
260 |
act3 = myThid - 1 |
261 |
max3 = nTx*nTy |
262 |
act4 = ikey_dynamics - 1 |
263 |
|
264 |
ikey_1 = i |
265 |
& + sNx*(j-1) |
266 |
& + sNx*sNy*act1 |
267 |
& + sNx*sNy*max1*act2 |
268 |
& + sNx*sNy*max1*max2*act3 |
269 |
& + sNx*sNy*max1*max2*max3*act4 |
270 |
#endif |
271 |
|
272 |
#ifdef ALLOW_ATM_TEMP |
273 |
|
274 |
#ifdef ALLOW_AUTODIFF_TAMC |
275 |
deltap(i,j) = 0. _d 0 |
276 |
delq(i,j) = 0. _d 0 |
277 |
#endif |
278 |
C--- Compute turbulent surface fluxes |
279 |
C- Pot. Temp and saturated specific humidity |
280 |
IF ( atemp(i,j,bi,bj) .NE. 0. _d 0 ) THEN |
281 |
C- Surface Temp. |
282 |
Tsf = theta(i,j,ksrf,bi,bj) + cen2kel |
283 |
IF ( sstExtrapol.GT.0. _d 0 ) THEN |
284 |
SSTtmp = sstExtrapol |
285 |
& *( theta(i,j,ksrf,bi,bj)-theta(i,j,ksrfp1,bi,bj) ) |
286 |
& * maskC(i,j,ksrfp1,bi,bj) |
287 |
Tsf = Tsf + MAX( SSTtmp, 0. _d 0 ) |
288 |
ENDIF |
289 |
c |
290 |
tmpbulk = cvapor_fac*exp(-cvapor_exp/Tsf) |
291 |
ssq = saltsat*tmpbulk/atmrho |
292 |
deltap(i,j) = atemp(i,j,bi,bj) + gamma_blk*ht - Tsf |
293 |
delq(i,j) = aqh(i,j,bi,bj) - ssq |
294 |
C-- no negative evaporation over ocean: |
295 |
IF ( noNegativeEvap ) delq(i,j) = MIN( 0. _d 0, delq(i,j) ) |
296 |
|
297 |
C-- initial guess for exchange coefficients: |
298 |
C take U_N = del.U ; stability from del.Theta ; |
299 |
stable = exf_half + SIGN(exf_half, deltap(i,j)) |
300 |
|
301 |
#ifndef ALLOW_ATM_WIND |
302 |
IF ( solve4Stress ) THEN |
303 |
#endif |
304 |
#ifdef ALLOW_AUTODIFF_TAMC |
305 |
CADJ STORE sh(i,j,bi,bj) = comlev1_exf_1, key = ikey_1 |
306 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
307 |
C-- Wind speed |
308 |
wsm = sh(i,j,bi,bj) |
309 |
tmpbulk = exf_scal_BulkCdn * |
310 |
& ( cdrag_1/wsm + cdrag_2 + cdrag_3*wsm ) |
311 |
rdn(i,j) = SQRT(tmpbulk) |
312 |
ustar(i,j) = rdn(i,j)*wsm |
313 |
#ifndef ALLOW_ATM_WIND |
314 |
ELSE |
315 |
rdn(i,j) = 0. _d 0 |
316 |
windStress = wStress(i,j,bi,bj) |
317 |
ustar(i,j) = sqrt(windStress/atmrho) |
318 |
c tau(i,j) = windStress/ustar(i,j) |
319 |
tau(i,j) = sqrt(windStress*atmrho) |
320 |
ENDIF |
321 |
#endif /* ALLOW_ATM_WIND */ |
322 |
|
323 |
C-- initial guess for exchange other coefficients: |
324 |
rhn = (exf_one-stable)*cstanton_1 + stable*cstanton_2 |
325 |
ren = cDalton |
326 |
C-- calculate turbulent scales |
327 |
tstar(i,j)=rhn*deltap(i,j) |
328 |
qstar(i,j)=ren*delq(i,j) |
329 |
|
330 |
ELSE |
331 |
C atemp(i,j,bi,bj) .EQ. 0. |
332 |
tstar (i,j) = 0. _d 0 |
333 |
qstar (i,j) = 0. _d 0 |
334 |
ustar (i,j) = 0. _d 0 |
335 |
tau (i,j) = 0. _d 0 |
336 |
rdn (i,j) = 0. _d 0 |
337 |
ENDIF |
338 |
ENDDO |
339 |
ENDDO |
340 |
DO iter=1,niter_bulk |
341 |
DO j = 1,sNy |
342 |
DO i = 1,sNx |
343 |
IF ( atemp(i,j,bi,bj) .NE. 0. _d 0 ) THEN |
344 |
C--- iterate with psi-functions to find transfer coefficients |
345 |
|
346 |
#ifdef ALLOW_AUTODIFF_TAMC |
347 |
ikey_2 = i |
348 |
& + sNx*(j-1) |
349 |
& + sNx*sNy*(iter-1) |
350 |
& + sNx*sNy*niter_bulk*act1 |
351 |
& + sNx*sNy*niter_bulk*max1*act2 |
352 |
& + sNx*sNy*niter_bulk*max1*max2*act3 |
353 |
& + sNx*sNy*niter_bulk*max1*max2*max3*act4 |
354 |
CADJ STORE rdn (i,j) = comlev1_exf_2, key = ikey_2 |
355 |
CADJ STORE ustar (i,j) = comlev1_exf_2, key = ikey_2 |
356 |
CADJ STORE qstar (i,j) = comlev1_exf_2, key = ikey_2 |
357 |
CADJ STORE tstar (i,j) = comlev1_exf_2, key = ikey_2 |
358 |
CADJ STORE sh (i,j,bi,bj) = comlev1_exf_2, key = ikey_2 |
359 |
CADJ STORE wspeed(i,j,bi,bj) = comlev1_exf_2, key = ikey_2 |
360 |
#endif |
361 |
|
362 |
t0 = atemp(i,j,bi,bj)* |
363 |
& (exf_one + humid_fac*aqh(i,j,bi,bj)) |
364 |
huol = ( tstar(i,j)/t0 + |
365 |
& qstar(i,j)/(exf_one/humid_fac+aqh(i,j,bi,bj)) |
366 |
& )*czol/(ustar(i,j)*ustar(i,j)) |
367 |
#ifdef ALLOW_BULK_LARGEYEAGER04 |
368 |
tmpbulk = MIN(abs(huol),10. _d 0) |
369 |
huol = SIGN(tmpbulk , huol) |
370 |
#else |
371 |
C-- Large&Pond1981: |
372 |
huol = max(huol,zolmin) |
373 |
#endif /* ALLOW_BULK_LARGEYEAGER04 */ |
374 |
htol = huol*ht/hu |
375 |
hqol = huol*hq/hu |
376 |
stable = exf_half + sign(exf_half, huol) |
377 |
|
378 |
c Evaluate all stability functions assuming hq = ht. |
379 |
#ifndef ALLOW_ATM_WIND |
380 |
IF ( solve4Stress ) THEN |
381 |
#endif |
382 |
#ifdef ALLOW_BULK_LARGEYEAGER04 |
383 |
C-- Large&Yeager04: |
384 |
xsq = SQRT( ABS(exf_one - huol*16. _d 0) ) |
385 |
#else |
386 |
C-- Large&Pond1981: |
387 |
xsq = max(sqrt(abs(exf_one - 16.*huol)),exf_one) |
388 |
#endif /* ALLOW_BULK_LARGEYEAGER04 */ |
389 |
x = SQRT(xsq) |
390 |
psimh = -psim_fac*huol*stable |
391 |
& +(exf_one-stable)* |
392 |
& ( LOG( (exf_one + exf_two*x + xsq)* |
393 |
& (exf_one+xsq)*.125 _d 0 ) |
394 |
& -exf_two*ATAN(x) + exf_half*pi ) |
395 |
#ifndef ALLOW_ATM_WIND |
396 |
ELSE |
397 |
psimh = 0. _d 0 |
398 |
ENDIF |
399 |
#endif |
400 |
#ifdef ALLOW_BULK_LARGEYEAGER04 |
401 |
C-- Large&Yeager04: |
402 |
xsq = SQRT( ABS(exf_one - htol*16. _d 0) ) |
403 |
#else |
404 |
C-- Large&Pond1981: |
405 |
xsq = max(sqrt(abs(exf_one - 16.*htol)),exf_one) |
406 |
#endif /* ALLOW_BULK_LARGEYEAGER04 */ |
407 |
psixh = -psim_fac*htol*stable + (exf_one-stable)* |
408 |
& ( exf_two*LOG(exf_half*(exf_one+xsq)) ) |
409 |
|
410 |
#ifndef ALLOW_ATM_WIND |
411 |
IF ( solve4Stress ) THEN |
412 |
#endif |
413 |
C- Shift wind speed using old coefficient |
414 |
#ifdef ALLOW_BULK_LARGEYEAGER04 |
415 |
C-- Large&Yeager04: |
416 |
dzTmp = (zwln-psimh)/karman |
417 |
usn = wspeed(i,j,bi,bj)/(exf_one + rdn(i,j)*dzTmp ) |
418 |
#else |
419 |
C-- Large&Pond1981: |
420 |
c rd(i,j)= rdn(i,j)/(exf_one - rdn(i,j)/karman*psimh ) |
421 |
c usn = sh(i,j,bi,bj)*rd(i,j)/rdn(i,j) |
422 |
c ML: the original formulation above is replaced below to be |
423 |
c similar to largeyeager04, but does not give the same results, strange |
424 |
usn = sh(i,j,bi,bj)/(exf_one - rdn(i,j)/karman*psimh) |
425 |
#endif /* ALLOW_BULK_LARGEYEAGER04 */ |
426 |
usm = MAX(usn, umin) |
427 |
|
428 |
C- Update the 10m, neutral stability transfer coefficients (momentum) |
429 |
tmpbulk = exf_scal_BulkCdn * |
430 |
& ( cdrag_1/usm + cdrag_2 + cdrag_3*usm ) |
431 |
rdn(i,j) = SQRT(tmpbulk) |
432 |
#ifdef ALLOW_BULK_LARGEYEAGER04 |
433 |
rd(i,j) = rdn(i,j)/(exf_one + rdn(i,j)*dzTmp) |
434 |
#else |
435 |
rd(i,j) = rdn(i,j)/(exf_one - rdn(i,j)/karman*psimh) |
436 |
#endif /* ALLOW_BULK_LARGEYEAGER04 */ |
437 |
ustar(i,j) = rd(i,j)*sh(i,j,bi,bj) |
438 |
C- Coeff: |
439 |
tau(i,j) = atmrho*rd(i,j)*wspeed(i,j,bi,bj) |
440 |
#ifndef ALLOW_ATM_WIND |
441 |
ENDIF |
442 |
#endif |
443 |
|
444 |
C- Update the 10m, neutral stability transfer coefficients (sens&evap) |
445 |
rhn = (exf_one-stable)*cstanton_1 + stable*cstanton_2 |
446 |
ren = cDalton |
447 |
|
448 |
C- Shift all coefficients to the measurement height and stability. |
449 |
rh = rhn/(exf_one + rhn*(ztln-psixh)/karman) |
450 |
re = ren/(exf_one + ren*(ztln-psixh)/karman) |
451 |
|
452 |
C-- Update ustar, tstar, qstar using updated, shifted coefficients. |
453 |
qstar(i,j) = re*delq(i,j) |
454 |
tstar(i,j) = rh*deltap(i,j) |
455 |
|
456 |
ENDIF |
457 |
ENDDO |
458 |
ENDDO |
459 |
c end of iteration loop |
460 |
ENDDO |
461 |
DO j = 1,sNy |
462 |
DO i = 1,sNx |
463 |
IF ( atemp(i,j,bi,bj) .NE. 0. _d 0 ) THEN |
464 |
|
465 |
#ifdef ALLOW_AUTODIFF_TAMC |
466 |
ikey_1 = i |
467 |
& + sNx*(j-1) |
468 |
& + sNx*sNy*act1 |
469 |
& + sNx*sNy*max1*act2 |
470 |
& + sNx*sNy*max1*max2*act3 |
471 |
& + sNx*sNy*max1*max2*max3*act4 |
472 |
CADJ STORE qstar(i,j) = comlev1_exf_1, key = ikey_1 |
473 |
CADJ STORE tstar(i,j) = comlev1_exf_1, key = ikey_1 |
474 |
CADJ STORE tau (i,j) = comlev1_exf_1, key = ikey_1 |
475 |
CADJ STORE rd (i,j) = comlev1_exf_1, key = ikey_1 |
476 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
477 |
|
478 |
C- Turbulent Fluxes |
479 |
hs(i,j,bi,bj) = atmcp*tau(i,j)*tstar(i,j) |
480 |
hl(i,j,bi,bj) = flamb*tau(i,j)*qstar(i,j) |
481 |
#ifndef EXF_READ_EVAP |
482 |
C change sign and convert from kg/m^2/s to m/s via rhoConstFresh |
483 |
c evap(i,j,bi,bj) = -recip_rhonil*tau(i,j)*qstar(i,j) |
484 |
evap(i,j,bi,bj) = -recip_rhoConstFresh*tau(i,j)*qstar(i,j) |
485 |
C but older version was using rhonil instead: |
486 |
c evap(i,j,bi,bj) = -recip_rhonil*tau(i,j)*qstar(i,j) |
487 |
#endif |
488 |
#ifdef ALLOW_ATM_WIND |
489 |
c ustress(i,j,bi,bj) = tau(i,j)*rd(i,j)*ws*cw(i,j,bi,bj) |
490 |
c vstress(i,j,bi,bj) = tau(i,j)*rd(i,j)*ws*sw(i,j,bi,bj) |
491 |
C- jmc: below is how it should be written ; different from above when |
492 |
C both wind-speed & 2 compon. of the wind are specified, and in |
493 |
C- this case, formula below is better. |
494 |
tmpbulk = tau(i,j)*rd(i,j) |
495 |
ustress(i,j,bi,bj) = tmpbulk*uwind(i,j,bi,bj) |
496 |
vstress(i,j,bi,bj) = tmpbulk*vwind(i,j,bi,bj) |
497 |
#endif |
498 |
|
499 |
c IF ( ATEMP(i,j,bi,bj) .NE. 0. ) |
500 |
else |
501 |
#ifdef ALLOW_ATM_WIND |
502 |
ustress(i,j,bi,bj) = 0. _d 0 |
503 |
vstress(i,j,bi,bj) = 0. _d 0 |
504 |
#endif /* ALLOW_ATM_WIND */ |
505 |
hflux (i,j,bi,bj) = 0. _d 0 |
506 |
evap (i,j,bi,bj) = 0. _d 0 |
507 |
hs (i,j,bi,bj) = 0. _d 0 |
508 |
hl (i,j,bi,bj) = 0. _d 0 |
509 |
|
510 |
c IF ( ATEMP(i,j,bi,bj) .NE. 0. ) |
511 |
endif |
512 |
|
513 |
#else /* ifndef ALLOW_ATM_TEMP */ |
514 |
#ifdef ALLOW_ATM_WIND |
515 |
wsm = sh(i,j,bi,bj) |
516 |
tmpbulk = exf_scal_BulkCdn * |
517 |
& ( cdrag_1/wsm + cdrag_2 + cdrag_3*wsm ) |
518 |
ustress(i,j,bi,bj) = atmrho*tmpbulk*wspeed(i,j,bi,bj)* |
519 |
& uwind(i,j,bi,bj) |
520 |
vstress(i,j,bi,bj) = atmrho*tmpbulk*wspeed(i,j,bi,bj)* |
521 |
& vwind(i,j,bi,bj) |
522 |
#endif |
523 |
#endif /* ifndef ALLOW_ATM_TEMP */ |
524 |
ENDDO |
525 |
ENDDO |
526 |
ENDDO |
527 |
ENDDO |
528 |
|
529 |
#endif /* ALLOW_BULKFORMULAE */ |
530 |
|
531 |
RETURN |
532 |
END |