1 |
c $Header: /u/gcmpack/MITgcm/pkg/exf/exf_getffields.F,v 1.6 2002/12/18 19:23:29 cheisey Exp $ |
2 |
|
3 |
#include "EXF_CPPOPTIONS.h" |
4 |
|
5 |
subroutine exf_GetFFields( mycurrenttime, mycurrentiter, mythid ) |
6 |
|
7 |
c ================================================================== |
8 |
c SUBROUTINE exf_GetFFields |
9 |
c ================================================================== |
10 |
c |
11 |
c o Get the surface fluxes either from file or as derived from bulk |
12 |
c formulae that use the atmospheric state. |
13 |
c |
14 |
c The calculation of the bulk surface fluxes has been adapted from |
15 |
c the NCOM model which uses the formulae given in Large and Pond |
16 |
c (1981 & 1982 ) |
17 |
c |
18 |
c |
19 |
c Header taken from NCOM version: ncom1.4.1 |
20 |
c ----------------------------------------- |
21 |
c |
22 |
c Following procedures and coefficients in Large and Pond |
23 |
c (1981 ; 1982) |
24 |
c |
25 |
c Output: Bulk estimates of the turbulent surface fluxes. |
26 |
c ------- |
27 |
c |
28 |
c hs - sensible heat flux (W/m^2), into ocean |
29 |
c hl - latent heat flux (W/m^2), into ocean |
30 |
c |
31 |
c Input: |
32 |
c ------ |
33 |
c |
34 |
c us - mean wind speed (m/s) at height hu (m) |
35 |
c th - mean air temperature (K) at height ht (m) |
36 |
c qh - mean air humidity (kg/kg) at height hq (m) |
37 |
c sst - sea surface temperature (K) |
38 |
c tk0 - Kelvin temperature at 0 Celsius (K) |
39 |
c |
40 |
c Assume 1) a neutral 10m drag coefficient = |
41 |
c |
42 |
c cdn = .0027/u10 + .000142 + .0000764 u10 |
43 |
c |
44 |
c 2) a neutral 10m stanton number = |
45 |
c |
46 |
c ctn = .0327 sqrt(cdn), unstable |
47 |
c ctn = .0180 sqrt(cdn), stable |
48 |
c |
49 |
c 3) a neutral 10m dalton number = |
50 |
c |
51 |
c cen = .0346 sqrt(cdn) |
52 |
c |
53 |
c 4) the saturation humidity of air at |
54 |
c |
55 |
c t(k) = exf_BulkqSat(t) (kg/m^3) |
56 |
c |
57 |
c Note: 1) here, tstar = <wt>/u*, and qstar = <wq>/u*. |
58 |
c 2) wind speeds should all be above a minimum speed, |
59 |
c say 0.5 m/s. |
60 |
c 3) with optional iteration loop, niter=3, should suffice. |
61 |
c 4) this version is for analyses inputs with hu = 10m and |
62 |
c ht = hq. |
63 |
c 5) sst enters in Celsius. |
64 |
c |
65 |
c ------------------------------------ |
66 |
c |
67 |
c By setting CPP options in the header file *EXF_CPPOPTIONS.h* it |
68 |
c is possible to combine data sets in four different ways: |
69 |
c |
70 |
c The following options are available: |
71 |
c |
72 |
c ALLOW_ATM_TEMP (UAT) |
73 |
c ALLOW_ATM_WIND (UAW) |
74 |
c |
75 |
c |
76 |
c UAT | UAW | action |
77 |
c ---------------------------------------------------- |
78 |
c undefined | undefined | Use surface fluxes. |
79 |
c undefined | defined | Assume cdn(u) given to |
80 |
c | | infer the wind stress. |
81 |
c defined | undefined | Compute wind field from |
82 |
c | | given stress assuming a |
83 |
c | | linear relation. |
84 |
c defined | defined | Use the bulk formulae. |
85 |
c ---------------------------------------------------- |
86 |
c |
87 |
c Implementations of the bulk formulae exist for the follwing |
88 |
c versions of the MITgcm: |
89 |
c |
90 |
c MITgcm : Patrick Heimbach |
91 |
c MITgcmUV: Christian Eckert |
92 |
c |
93 |
c started: Christian Eckert eckert@mit.edu 27-Aug-1999 |
94 |
c |
95 |
c changed: Christian Eckert eckert@mit.edu 14-Jan-2000 |
96 |
c - restructured the original version in order to have a |
97 |
c better interface to the MITgcmUV. |
98 |
c |
99 |
c Christian Eckert eckert@mit.edu 12-Feb-2000 |
100 |
c - Changed Routine names (package prefix: exf_) |
101 |
c |
102 |
c Patrick Heimbach, heimbach@mit.edu 04-May-2000 |
103 |
c - changed the handling of precip and sflux with respect |
104 |
c to CPP options ALLOW_BULKFORMULAE and ALLOW_ATM_TEMP |
105 |
c - included some CPP flags ALLOW_BULKFORMULAE to make |
106 |
c sure ALLOW_ATM_TEMP, ALLOW_ATM_WIND are used only in |
107 |
c conjunction with defined ALLOW_BULKFORMULAE |
108 |
c - statement functions discarded |
109 |
c |
110 |
c Ralf.Giering@FastOpt.de 25-Mai-2000 |
111 |
c - total rewrite using new subroutines |
112 |
c |
113 |
c Detlef Stammer: include river run-off. Nov. 21, 2001 |
114 |
c |
115 |
c heimbach@mit.edu, 10-Jan-2002 |
116 |
c - changes to enable field swapping |
117 |
c |
118 |
c menemenlis@jpl.nasa.gov, 20-Dec-2002 |
119 |
c - Added EXF_READ_EVAP and EXF_NO_BULK_COMPUTATIONS. |
120 |
c |
121 |
c ================================================================== |
122 |
c SUBROUTINE exf_GetFFields |
123 |
c ================================================================== |
124 |
|
125 |
implicit none |
126 |
|
127 |
c == global variables == |
128 |
|
129 |
#include "EEPARAMS.h" |
130 |
#include "SIZE.h" |
131 |
#include "PARAMS.h" |
132 |
#include "DYNVARS.h" |
133 |
#include "GRID.h" |
134 |
|
135 |
#include "exf_fields.h" |
136 |
#include "exf_constants.h" |
137 |
|
138 |
#ifdef ALLOW_AUTODIFF_TAMC |
139 |
#include "tamc.h" |
140 |
#endif |
141 |
|
142 |
c == routine arguments == |
143 |
|
144 |
integer mythid |
145 |
integer mycurrentiter |
146 |
_RL mycurrenttime |
147 |
|
148 |
c == local variables == |
149 |
|
150 |
integer bi,bj |
151 |
integer i,j,k |
152 |
|
153 |
#ifdef ALLOW_BULKFORMULAE |
154 |
|
155 |
#ifdef ALLOW_ATM_TEMP |
156 |
integer iter |
157 |
_RL aln |
158 |
_RL delq |
159 |
_RL deltap |
160 |
_RL hqol |
161 |
_RL htol |
162 |
_RL huol |
163 |
_RL psimh |
164 |
_RL psixh |
165 |
_RL qstar |
166 |
_RL rd |
167 |
_RL re |
168 |
_RL rdn |
169 |
_RL rh |
170 |
_RL ssttmp |
171 |
_RL ssq |
172 |
_RL stable |
173 |
_RL tstar |
174 |
_RL t0 |
175 |
_RL ustar |
176 |
_RL uzn |
177 |
_RL shn |
178 |
_RL xsq |
179 |
_RL x |
180 |
_RL tau |
181 |
#ifdef ALLOW_AUTODIFF_TAMC |
182 |
integer ikey_1 |
183 |
integer ikey_2 |
184 |
#endif |
185 |
#endif /* ALLOW_ATM_TEMP */ |
186 |
|
187 |
_RL ustmp |
188 |
_RL us |
189 |
_RL cw |
190 |
_RL sw |
191 |
_RL sh |
192 |
_RL hs(1-olx:snx+olx,1-oly:sny+oly,nsx,nsy) |
193 |
_RL hl(1-olx:snx+olx,1-oly:sny+oly,nsx,nsy) |
194 |
_RL hfl |
195 |
|
196 |
#endif /* ALLOW_BULKFORMULAE */ |
197 |
|
198 |
c == external functions == |
199 |
|
200 |
integer ilnblnk |
201 |
external ilnblnk |
202 |
|
203 |
#ifdef ALLOW_BULKFORMULAE |
204 |
_RL exf_BulkqSat |
205 |
external exf_BulkqSat |
206 |
_RL exf_BulkCdn |
207 |
external exf_BulkCdn |
208 |
_RL exf_BulkRhn |
209 |
external exf_BulkRhn |
210 |
#endif /* ALLOW_BULKFORMULAE */ |
211 |
|
212 |
#ifndef ALLOW_ATM_WIND |
213 |
_RL TMP1 |
214 |
_RL TMP2 |
215 |
_RL TMP3 |
216 |
_RL TMP4 |
217 |
_RL TMP5 |
218 |
#endif |
219 |
|
220 |
c == end of interface == |
221 |
|
222 |
c determine forcing field records |
223 |
|
224 |
#ifdef EXF_READ_EVAP |
225 |
c Evaporation |
226 |
call exf_set_evap( mycurrenttime, mycurrentiter, mythid ) |
227 |
#endif EXF_READ_EVAP |
228 |
|
229 |
#ifdef ALLOW_BULKFORMULAE |
230 |
|
231 |
#if (defined (ALLOW_ATM_TEMP) || defined (ALLOW_ATM_WIND)) |
232 |
cph This statement cannot be a PARAMETER statement in the header, |
233 |
cph but must come here; it's not fortran77 standard |
234 |
aln = log(ht/zref) |
235 |
#endif |
236 |
|
237 |
c Determine where we are in time and set counters, flags and |
238 |
c the linear interpolation factors accordingly. |
239 |
#ifdef ALLOW_ATM_TEMP |
240 |
c Atmospheric temperature. |
241 |
call exf_set_atemp( mycurrenttime, mycurrentiter, mythid ) |
242 |
|
243 |
c Atmospheric humidity. |
244 |
call exf_set_aqh( mycurrenttime, mycurrentiter, mythid ) |
245 |
|
246 |
c Net long wave radiative flux. |
247 |
call exf_set_lwflux( mycurrenttime, mycurrentiter, mythid ) |
248 |
|
249 |
c Net short wave radiative flux. |
250 |
call exf_set_swflux( mycurrenttime, mycurrentiter, mythid ) |
251 |
|
252 |
c Precipitation. |
253 |
call exf_set_precip( mycurrenttime, mycurrentiter, mythid ) |
254 |
|
255 |
#ifdef ALLOW_ATEMP_CONTROL |
256 |
call ctrl_getatemp ( mycurrenttime, mycurrentiter, mythid ) |
257 |
#endif |
258 |
|
259 |
#ifdef ALLOW_AQH_CONTROL |
260 |
call ctrl_getaqh ( mycurrenttime, mycurrentiter, mythid ) |
261 |
#endif |
262 |
|
263 |
#else |
264 |
|
265 |
c Atmospheric heat flux. |
266 |
call exf_set_hflux( mycurrenttime, mycurrentiter, mythid ) |
267 |
|
268 |
c Salt flux. |
269 |
call exf_set_sflux( mycurrenttime, mycurrentiter, mythid ) |
270 |
|
271 |
#ifdef ALLOW_KPP |
272 |
c Net short wave radiative flux. |
273 |
call exf_set_swflux( mycurrenttime, mycurrentiter, mythid ) |
274 |
|
275 |
#endif /* ALLOW_KPP */ |
276 |
|
277 |
#endif /* ALLOW_ATM_TEMP */ |
278 |
|
279 |
#ifdef ALLOW_ATM_WIND |
280 |
c Zonal wind. |
281 |
call exf_set_uwind( mycurrenttime, mycurrentiter, mythid ) |
282 |
|
283 |
c Meridional wind. |
284 |
call exf_set_vwind( mycurrenttime, mycurrentiter, mythid ) |
285 |
|
286 |
#ifdef ALLOW_UWIND_CONTROL |
287 |
call ctrl_getuwind ( mycurrenttime, mycurrentiter, mythid ) |
288 |
#endif |
289 |
|
290 |
#ifdef ALLOW_VWIND_CONTROL |
291 |
call ctrl_getvwind ( mycurrenttime, mycurrentiter, mythid ) |
292 |
#endif |
293 |
|
294 |
#else |
295 |
c Zonal wind stress. |
296 |
call exf_set_ustress( mycurrenttime, mycurrentiter, mythid ) |
297 |
|
298 |
c Meridional wind stress. |
299 |
call exf_set_vstress( mycurrenttime, mycurrentiter, mythid ) |
300 |
|
301 |
#endif /* ALLOW_ATM_WIND */ |
302 |
|
303 |
#else /* ALLOW_BULKFORMULAE undefined */ |
304 |
c Atmospheric heat flux. |
305 |
call exf_set_hflux( mycurrenttime, mycurrentiter, mythid ) |
306 |
|
307 |
c Salt flux. |
308 |
call exf_set_sflux( mycurrenttime, mycurrentiter, mythid ) |
309 |
|
310 |
c Zonal wind stress. |
311 |
call exf_set_ustress( mycurrenttime, mycurrentiter, mythid ) |
312 |
|
313 |
c Meridional wind stress. |
314 |
call exf_set_vstress( mycurrenttime, mycurrentiter, mythid ) |
315 |
|
316 |
#ifdef ALLOW_KPP |
317 |
c Net short wave radiative flux. |
318 |
call exf_set_swflux( mycurrenttime, mycurrentiter, mythid ) |
319 |
#endif |
320 |
|
321 |
#endif /* ALLOW_BULKFORMULAE */ |
322 |
|
323 |
#if ~defined(EXF_NO_BULK_COMPUTATIONS) || ~defined(EXF_READ_EVAP) |
324 |
C-- Use atmospheric state to compute surace fluxes. |
325 |
|
326 |
c Loop over tiles. |
327 |
#ifdef ALLOW_AUTODIFF_TAMC |
328 |
C-- HPF directive to help TAMC |
329 |
CHPF$ INDEPENDENT |
330 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
331 |
do bj = mybylo(mythid),mybyhi(mythid) |
332 |
#ifdef ALLOW_AUTODIFF_TAMC |
333 |
C-- HPF directive to help TAMC |
334 |
CHPF$ INDEPENDENT |
335 |
#endif |
336 |
do bi = mybxlo(mythid),mybxhi(mythid) |
337 |
|
338 |
k = 1 |
339 |
|
340 |
do j = 1-oly,sny+oly |
341 |
do i = 1-olx,snx+olx |
342 |
|
343 |
#ifdef ALLOW_BULKFORMULAE |
344 |
|
345 |
#ifdef ALLOW_AUTODIFF_TAMC |
346 |
act1 = bi - myBxLo(myThid) |
347 |
max1 = myBxHi(myThid) - myBxLo(myThid) + 1 |
348 |
act2 = bj - myByLo(myThid) |
349 |
max2 = myByHi(myThid) - myByLo(myThid) + 1 |
350 |
act3 = myThid - 1 |
351 |
max3 = nTx*nTy |
352 |
act4 = ikey_dynamics - 1 |
353 |
|
354 |
ikey_1 = i |
355 |
& + sNx*(j-1) |
356 |
& + sNx*sNy*act1 |
357 |
& + sNx*sNy*max1*act2 |
358 |
& + sNx*sNy*max1*max2*act3 |
359 |
& + sNx*sNy*max1*max2*max3*act4 |
360 |
#endif |
361 |
|
362 |
c Compute the turbulent surface fluxes. |
363 |
c (Bulk formulae estimates) |
364 |
|
365 |
#ifdef ALLOW_ATM_WIND |
366 |
c Wind speed and direction. |
367 |
ustmp = uwind(i,j,bi,bj)*uwind(i,j,bi,bj) + |
368 |
& vwind(i,j,bi,bj)*vwind(i,j,bi,bj) |
369 |
if ( ustmp .ne. 0. _d 0 ) then |
370 |
us = sqrt(ustmp) |
371 |
cw = uwind(i,j,bi,bj)/us |
372 |
sw = vwind(i,j,bi,bj)/us |
373 |
else |
374 |
us = 0. _d 0 |
375 |
cw = 0. _d 0 |
376 |
sw = 0. _d 0 |
377 |
endif |
378 |
sh = max(us,umin) |
379 |
#else |
380 |
#ifdef ALLOW_ATM_TEMP |
381 |
|
382 |
c The variables us, sh and rdn have to be computed from |
383 |
c given wind stresses inverting relationship for neutral |
384 |
c drag coeff. cdn. |
385 |
c The inversion is based on linear and quadratic form of |
386 |
c cdn(umps); ustar can be directly computed from stress; |
387 |
|
388 |
ustmp = ustress(i,j,bi,bj)*ustress(i,j,bi,bj) + |
389 |
& vstress(i,j,bi,bj)*vstress(i,j,bi,bj) |
390 |
if ( ustmp .ne. 0. _d 0 ) then |
391 |
ustar = sqrt(ustmp/atmrho) |
392 |
cw = ustress(i,j,bi,bj)/sqrt(ustmp) |
393 |
sw = vstress(i,j,bi,bj)/sqrt(ustmp) |
394 |
else |
395 |
ustar = 0. _d 0 |
396 |
cw = 0. _d 0 |
397 |
sw = 0. _d 0 |
398 |
endif |
399 |
|
400 |
if ( ustar .eq. 0. _d 0 ) then |
401 |
us = 0. _d 0 |
402 |
else if ( ustar .lt. ustofu11 ) then |
403 |
tmp1 = -cquadrag_2/cquadrag_1/2 |
404 |
tmp2 = sqrt(tmp1*tmp1 + ustar*ustar/cquadrag_1) |
405 |
us = sqrt(tmp1 + tmp2) |
406 |
else |
407 |
tmp3 = clindrag_2/clindrag_1/3 |
408 |
tmp4 = ustar*ustar/clindrag_1/2 - tmp3**3 |
409 |
tmp5 = sqrt(ustar*ustar/clindrag_1* |
410 |
& (ustar*ustar/clindrag_1/4 - tmp3**3)) |
411 |
us = (tmp4 + tmp5)**(1/3) + |
412 |
& tmp3**2 * (tmp4 + tmp5)**(-1/3) - tmp3 |
413 |
endif |
414 |
|
415 |
if ( us .ne. 0 ) then |
416 |
rdn = ustar/us |
417 |
else |
418 |
rdn = 0. _d 0 |
419 |
end if |
420 |
|
421 |
sh = max(us,umin) |
422 |
#endif /* ALLOW_ATM_TEMP */ |
423 |
#endif /* ALLOW_ATM_WIND */ |
424 |
|
425 |
#ifdef ALLOW_ATM_TEMP |
426 |
|
427 |
c Initial guess: z/l=0.0; hu=ht=hq=z |
428 |
c Iterations: converge on z/l and hence the fluxes. |
429 |
c t0 : virtual temperature (K) |
430 |
c ssq : sea surface humidity (kg/kg) |
431 |
c deltap : potential temperature diff (K) |
432 |
|
433 |
if ( atemp(i,j,bi,bj) .ne. 0. _d 0 ) then |
434 |
t0 = atemp(i,j,bi,bj)* |
435 |
& (exf_one + humid_fac*aqh(i,j,bi,bj)) |
436 |
ssttmp = theta(i,j,k,bi,bj) |
437 |
ssq = saltsat* |
438 |
& exf_BulkqSat(ssttmp + cen2kel)/ |
439 |
& atmrho |
440 |
deltap = atemp(i,j,bi,bj) + gamma_blk*ht - |
441 |
& ssttmp - cen2kel |
442 |
delq = aqh(i,j,bi,bj) - ssq |
443 |
stable = exf_half + sign(exf_half, deltap) |
444 |
#ifdef ALLOW_AUTODIFF_TAMC |
445 |
CADJ STORE sh = comlev1_exf_1, key = ikey_1 |
446 |
#endif |
447 |
rdn = sqrt(exf_BulkCdn(sh)) |
448 |
ustar = rdn*sh |
449 |
tstar = exf_BulkRhn(stable)*deltap |
450 |
qstar = cdalton*delq |
451 |
|
452 |
do iter = 1,niter_bulk |
453 |
|
454 |
#ifdef ALLOW_AUTODIFF_TAMC |
455 |
ikey_2 = iter |
456 |
& + niter_bulk*(i-1) |
457 |
& + sNx*niter_bulk*(j-1) |
458 |
& + sNx*niter_bulk*sNy*act1 |
459 |
& + sNx*niter_bulk*sNy*max1*act2 |
460 |
& + sNx*niter_bulk*sNy*max1*max2*act3 |
461 |
& + sNx*niter_bulk*sNy*max1*max2*max3*act4 |
462 |
#endif |
463 |
|
464 |
#ifdef ALLOW_AUTODIFF_TAMC |
465 |
CADJ STORE rdn = comlev1_exf_2, key = ikey_2 |
466 |
CADJ STORE ustar = comlev1_exf_2, key = ikey_2 |
467 |
CADJ STORE qstar = comlev1_exf_2, key = ikey_2 |
468 |
CADJ STORE tstar = comlev1_exf_2, key = ikey_2 |
469 |
CADJ STORE sh = comlev1_exf_2, key = ikey_2 |
470 |
CADJ STORE us = comlev1_exf_2, key = ikey_2 |
471 |
#endif |
472 |
|
473 |
huol = czol*(tstar/t0 + |
474 |
& qstar/(exf_one/humid_fac+aqh(i,j,bi,bj)))/ |
475 |
& ustar**2 |
476 |
huol = max(huol,zolmin) |
477 |
stable = exf_half + sign(exf_half, huol) |
478 |
htol = huol*ht/hu |
479 |
hqol = huol*hq/hu |
480 |
|
481 |
c Evaluate all stability functions assuming hq = ht. |
482 |
xsq = max(sqrt(abs(exf_one - 16.*huol)),exf_one) |
483 |
x = sqrt(xsq) |
484 |
psimh = -psim_fac*huol*stable + |
485 |
& (exf_one - stable)* |
486 |
& log((exf_one + x*(exf_two + x))* |
487 |
& (exf_one + xsq)/8.) - exf_two*atan(x) + |
488 |
& pi*exf_half |
489 |
xsq = max(sqrt(abs(exf_one - 16.*htol)),exf_one) |
490 |
psixh = -psim_fac*htol*stable + (exf_one - stable)* |
491 |
& exf_two*log((exf_one + xsq)/exf_two) |
492 |
|
493 |
c Shift wind speed using old coefficient |
494 |
ccc rd = rdn/(exf_one + rdn/karman* |
495 |
ccc & (log(hu/zref) - psimh) ) |
496 |
rd = rdn/(exf_one - rdn/karman*psimh ) |
497 |
shn = sh*rd/rdn |
498 |
uzn = max(shn, umin) |
499 |
|
500 |
c Update the transfer coefficients at 10 meters |
501 |
c and neutral stability. |
502 |
|
503 |
rdn = sqrt(exf_BulkCdn(uzn)) |
504 |
|
505 |
c Shift all coefficients to the measurement height |
506 |
c and stability. |
507 |
c rd = rdn/(exf_one + rdn/karman*(log(hu/zref) - psimh)) |
508 |
rd = rdn/(exf_one - rdn/karman*psimh) |
509 |
rh = exf_BulkRhn(stable)/(exf_one + |
510 |
& exf_BulkRhn(stable)/ |
511 |
& karman*(aln - psixh)) |
512 |
re = cdalton/(exf_one + cdalton/karman*(aln - psixh)) |
513 |
|
514 |
c Update ustar, tstar, qstar using updated, shifted |
515 |
c coefficients. |
516 |
ustar = rd*sh |
517 |
qstar = re*delq |
518 |
tstar = rh*deltap |
519 |
tau = atmrho*ustar**2 |
520 |
tau = tau*us/sh |
521 |
|
522 |
enddo |
523 |
|
524 |
#ifdef ALLOW_AUTODIFF_TAMC |
525 |
CADJ STORE ustar = comlev1_exf_1, key = ikey_1 |
526 |
CADJ STORE qstar = comlev1_exf_1, key = ikey_1 |
527 |
CADJ STORE tstar = comlev1_exf_1, key = ikey_1 |
528 |
CADJ STORE tau = comlev1_exf_1, key = ikey_1 |
529 |
CADJ STORE cw = comlev1_exf_1, key = ikey_1 |
530 |
CADJ STORE sw = comlev1_exf_1, key = ikey_1 |
531 |
#endif |
532 |
|
533 |
hs(i,j,bi,bj) = atmcp*tau*tstar/ustar |
534 |
hl(i,j,bi,bj) = flamb*tau*qstar/ustar |
535 |
#ifndef EXF_READ_EVAP |
536 |
evap(i,j,bi,bj) = tau*qstar/ustar |
537 |
#endif EXF_READ_EVAP |
538 |
ustress(i,j,bi,bj) = tau*cw |
539 |
vstress(i,j,bi,bj) = tau*sw |
540 |
else |
541 |
ustress(i,j,bi,bj) = 0. _d 0 |
542 |
vstress(i,j,bi,bj) = 0. _d 0 |
543 |
hflux (i,j,bi,bj) = 0. _d 0 |
544 |
hs(i,j,bi,bj) = 0. _d 0 |
545 |
hl(i,j,bi,bj) = 0. _d 0 |
546 |
endif |
547 |
|
548 |
#else |
549 |
#ifdef ALLOW_ATM_WIND |
550 |
ustress(i,j,bi,bj) = atmrho*exf_BulkCdn(sh)*us* |
551 |
& uwind(i,j,bi,bj) |
552 |
vstress(i,j,bi,bj) = atmrho*exf_BulkCdn(sh)*us* |
553 |
& vwind(i,j,bi,bj) |
554 |
#endif /* ALLOW_ATM_WIND */ |
555 |
#endif /* ALLOW_ATM_TEMP */ |
556 |
enddo |
557 |
enddo |
558 |
enddo |
559 |
enddo |
560 |
|
561 |
c Add all contributions. |
562 |
do bj = mybylo(mythid),mybyhi(mythid) |
563 |
do bi = mybxlo(mythid),mybxhi(mythid) |
564 |
do j = 1,sny |
565 |
do i = 1,snx |
566 |
c Net surface heat flux. |
567 |
#ifdef ALLOW_ATM_TEMP |
568 |
hfl = 0. _d 0 |
569 |
hfl = hfl - hs(i,j,bi,bj) |
570 |
hfl = hfl - hl(i,j,bi,bj) |
571 |
hfl = hfl + lwflux(i,j,bi,bj) |
572 |
#ifndef ALLOW_KPP |
573 |
hfl = hfl + swflux(i,j,bi,bj) |
574 |
#endif /* ALLOW_KPP undef */ |
575 |
c Heat flux: |
576 |
hflux(i,j,bi,bj) = hfl*maskc(i,j,1,bi,bj) |
577 |
c Salt flux from Precipitation and Evaporation. |
578 |
sflux(i,j,bi,bj) = precip(i,j,bi,bj) - evap(i,j,bi,bj) |
579 |
#endif /* ALLOW_ATM_TEMP */ |
580 |
|
581 |
#else |
582 |
hflux(i,j,bi,bj) = hflux(i,j,bi,bj)*maskc(i,j,1,bi,bj) |
583 |
sflux(i,j,bi,bj) = sflux(i,j,bi,bj)*maskc(i,j,1,bi,bj) |
584 |
#endif /* ALLOW_BULKFORMULAE */ |
585 |
|
586 |
#ifdef ALLOW_RUNOFF |
587 |
sflux(i,j,bi,bj) = sflux(i,j,bi,bj) + |
588 |
& runoff(i,j,bi,bj)*maskc(i,j,1,bi,bj) |
589 |
#endif /* ALLOW_RUNOFF */ |
590 |
|
591 |
enddo |
592 |
enddo |
593 |
enddo |
594 |
enddo |
595 |
|
596 |
#endif EXF_NO_BULK_COMPUTATIONS |
597 |
|
598 |
c Update the tile edges. |
599 |
_EXCH_XY_R8(hflux, mythid) |
600 |
_EXCH_XY_R8(sflux, mythid) |
601 |
c _EXCH_XY_R8(ustress, mythid) |
602 |
c _EXCH_XY_R8(vstress, mythid) |
603 |
CALL EXCH_UV_XY_RS(ustress, vstress, .TRUE., myThid) |
604 |
|
605 |
#ifdef ALLOW_KPP |
606 |
_EXCH_XY_R8(swflux, mythid) |
607 |
#endif /* ALLOW_KPP */ |
608 |
|
609 |
end |