1 |
c $Header: /u/gcmpack/MITgcm/pkg/exf/exf_bulkformulae.F,v 1.9 2005/06/27 16:34:29 heimbach Exp $ |
2 |
|
3 |
#include "EXF_OPTIONS.h" |
4 |
|
5 |
subroutine exf_bulkformulae(mycurrenttime, mycurrentiter, mythid) |
6 |
|
7 |
c ================================================================== |
8 |
c SUBROUTINE exf_bulkformulae |
9 |
c ================================================================== |
10 |
c |
11 |
c o Read-in atmospheric state and/or surface fluxes from files. |
12 |
c |
13 |
c o Use bulk formulae to estimate turbulent and/or radiative |
14 |
c fluxes at the surface. |
15 |
c |
16 |
c NOTES: |
17 |
c ====== |
18 |
c |
19 |
c See EXF_OPTIONS.h for a description of the various possible |
20 |
c ocean-model forcing configurations. |
21 |
c |
22 |
c The bulk formulae of pkg/exf are not valid for sea-ice covered |
23 |
c oceans but they can be used in combination with a sea-ice model, |
24 |
c for example, pkg/seaice, to specify open water flux contributions. |
25 |
c |
26 |
c ================================================================== |
27 |
c |
28 |
c The calculation of the bulk surface fluxes has been adapted from |
29 |
c the NCOM model which uses the formulae given in Large and Pond |
30 |
c (1981 & 1982 ) |
31 |
c |
32 |
c |
33 |
c Header taken from NCOM version: ncom1.4.1 |
34 |
c ----------------------------------------- |
35 |
c |
36 |
c Following procedures and coefficients in Large and Pond |
37 |
c (1981 ; 1982) |
38 |
c |
39 |
c Output: Bulk estimates of the turbulent surface fluxes. |
40 |
c ------- |
41 |
c |
42 |
c hs - sensible heat flux (W/m^2), into ocean |
43 |
c hl - latent heat flux (W/m^2), into ocean |
44 |
c |
45 |
c Input: |
46 |
c ------ |
47 |
c |
48 |
c us - mean wind speed (m/s) at height hu (m) |
49 |
c th - mean air temperature (K) at height ht (m) |
50 |
c qh - mean air humidity (kg/kg) at height hq (m) |
51 |
c sst - sea surface temperature (K) |
52 |
c tk0 - Kelvin temperature at 0 Celsius (K) |
53 |
c |
54 |
c Assume 1) a neutral 10m drag coefficient = |
55 |
c |
56 |
c cdn = .0027/u10 + .000142 + .0000764 u10 |
57 |
c |
58 |
c 2) a neutral 10m stanton number = |
59 |
c |
60 |
c ctn = .0327 sqrt(cdn), unstable |
61 |
c ctn = .0180 sqrt(cdn), stable |
62 |
c |
63 |
c 3) a neutral 10m dalton number = |
64 |
c |
65 |
c cen = .0346 sqrt(cdn) |
66 |
c |
67 |
c 4) the saturation humidity of air at |
68 |
c |
69 |
c t(k) = exf_BulkqSat(t) (kg/m^3) |
70 |
c |
71 |
c Note: 1) here, tstar = <wt>/u*, and qstar = <wq>/u*. |
72 |
c 2) wind speeds should all be above a minimum speed, |
73 |
c say 0.5 m/s. |
74 |
c 3) with optional iteration loop, niter=3, should suffice. |
75 |
c 4) this version is for analyses inputs with hu = 10m and |
76 |
c ht = hq. |
77 |
c 5) sst enters in Celsius. |
78 |
c |
79 |
c ================================================================== |
80 |
c |
81 |
c started: Christian Eckert eckert@mit.edu 27-Aug-1999 |
82 |
c |
83 |
c changed: Christian Eckert eckert@mit.edu 14-Jan-2000 |
84 |
c - restructured the original version in order to have a |
85 |
c better interface to the MITgcmUV. |
86 |
c |
87 |
c Christian Eckert eckert@mit.edu 12-Feb-2000 |
88 |
c - Changed Routine names (package prefix: exf_) |
89 |
c |
90 |
c Patrick Heimbach, heimbach@mit.edu 04-May-2000 |
91 |
c - changed the handling of precip and sflux with respect |
92 |
c to CPP options ALLOW_BULKFORMULAE and ALLOW_ATM_TEMP |
93 |
c - included some CPP flags ALLOW_BULKFORMULAE to make |
94 |
c sure ALLOW_ATM_TEMP, ALLOW_ATM_WIND are used only in |
95 |
c conjunction with defined ALLOW_BULKFORMULAE |
96 |
c - statement functions discarded |
97 |
c |
98 |
c Ralf.Giering@FastOpt.de 25-Mai-2000 |
99 |
c - total rewrite using new subroutines |
100 |
c |
101 |
c Detlef Stammer: include river run-off. Nov. 21, 2001 |
102 |
c |
103 |
c heimbach@mit.edu, 10-Jan-2002 |
104 |
c - changes to enable field swapping |
105 |
c |
106 |
c mods for pkg/seaice: menemenlis@jpl.nasa.gov 20-Dec-2002 |
107 |
c |
108 |
c ================================================================== |
109 |
c SUBROUTINE exf_bulkformulae |
110 |
c ================================================================== |
111 |
|
112 |
implicit none |
113 |
|
114 |
c == global variables == |
115 |
|
116 |
#include "EEPARAMS.h" |
117 |
#include "SIZE.h" |
118 |
#include "PARAMS.h" |
119 |
#include "DYNVARS.h" |
120 |
#include "GRID.h" |
121 |
|
122 |
#include "exf_param.h" |
123 |
#include "exf_fields.h" |
124 |
#include "exf_constants.h" |
125 |
|
126 |
#ifdef ALLOW_AUTODIFF_TAMC |
127 |
#include "tamc.h" |
128 |
#endif |
129 |
|
130 |
c == routine arguments == |
131 |
|
132 |
integer mythid |
133 |
integer mycurrentiter |
134 |
_RL mycurrenttime |
135 |
|
136 |
#ifdef ALLOW_BULKFORMULAE |
137 |
|
138 |
c == local variables == |
139 |
|
140 |
integer bi,bj |
141 |
integer i,j,k |
142 |
|
143 |
_RL aln |
144 |
|
145 |
#ifdef ALLOW_ATM_TEMP |
146 |
integer iter |
147 |
_RL delq |
148 |
_RL deltap |
149 |
_RL hqol |
150 |
_RL htol |
151 |
_RL huol |
152 |
_RL psimh |
153 |
_RL psixh |
154 |
_RL qstar |
155 |
_RL rd |
156 |
_RL re |
157 |
_RL rdn |
158 |
_RL rh |
159 |
_RL ssttmp |
160 |
_RL ssq |
161 |
_RL stable |
162 |
_RL tstar |
163 |
_RL t0 |
164 |
_RL ustar |
165 |
_RL uzn |
166 |
_RL shn |
167 |
_RL xsq |
168 |
_RL x |
169 |
_RL tau |
170 |
#ifdef ALLOW_AUTODIFF_TAMC |
171 |
integer ikey_1 |
172 |
integer ikey_2 |
173 |
#endif |
174 |
#endif /* ALLOW_ATM_TEMP */ |
175 |
|
176 |
_RL ustmp |
177 |
_RL us |
178 |
_RL cw |
179 |
_RL sw |
180 |
_RL sh |
181 |
_RL hfl |
182 |
_RL tmpbulk |
183 |
|
184 |
c == external functions == |
185 |
|
186 |
integer ilnblnk |
187 |
external ilnblnk |
188 |
|
189 |
_RL exf_BulkqSat |
190 |
external exf_BulkqSat |
191 |
_RL exf_BulkCdn |
192 |
external exf_BulkCdn |
193 |
_RL exf_BulkRhn |
194 |
external exf_BulkRhn |
195 |
|
196 |
#ifndef ALLOW_ATM_WIND |
197 |
_RL TMP1 |
198 |
_RL TMP2 |
199 |
_RL TMP3 |
200 |
_RL TMP4 |
201 |
_RL TMP5 |
202 |
#endif |
203 |
|
204 |
c == end of interface == |
205 |
|
206 |
cph This statement cannot be a PARAMETER statement in the header, |
207 |
cph but must come here; it's not fortran77 standard |
208 |
aln = log(ht/zref) |
209 |
|
210 |
c-- Use atmospheric state to compute surface fluxes. |
211 |
|
212 |
c Loop over tiles. |
213 |
#ifdef ALLOW_AUTODIFF_TAMC |
214 |
C-- HPF directive to help TAMC |
215 |
CHPF$ INDEPENDENT |
216 |
#endif |
217 |
do bj = mybylo(mythid),mybyhi(mythid) |
218 |
#ifdef ALLOW_AUTODIFF_TAMC |
219 |
C-- HPF directive to help TAMC |
220 |
CHPF$ INDEPENDENT |
221 |
#endif |
222 |
do bi = mybxlo(mythid),mybxhi(mythid) |
223 |
|
224 |
k = 1 |
225 |
|
226 |
do j = 1,sny |
227 |
do i = 1,snx |
228 |
|
229 |
#ifdef ALLOW_AUTODIFF_TAMC |
230 |
act1 = bi - myBxLo(myThid) |
231 |
max1 = myBxHi(myThid) - myBxLo(myThid) + 1 |
232 |
act2 = bj - myByLo(myThid) |
233 |
max2 = myByHi(myThid) - myByLo(myThid) + 1 |
234 |
act3 = myThid - 1 |
235 |
max3 = nTx*nTy |
236 |
act4 = ikey_dynamics - 1 |
237 |
|
238 |
ikey_1 = i |
239 |
& + sNx*(j-1) |
240 |
& + sNx*sNy*act1 |
241 |
& + sNx*sNy*max1*act2 |
242 |
& + sNx*sNy*max1*max2*act3 |
243 |
& + sNx*sNy*max1*max2*max3*act4 |
244 |
#endif |
245 |
|
246 |
#ifdef ALLOW_DOWNWARD_RADIATION |
247 |
c-- Compute net from downward and downward from net longwave and |
248 |
c shortwave radiation, if needed. |
249 |
c lwflux = Stefan-Boltzman constant * emissivity * SST - lwdown |
250 |
c swflux = - ( 1 - albedo ) * swdown |
251 |
|
252 |
#ifdef ALLOW_ATM_TEMP |
253 |
if ( lwfluxfile .EQ. ' ' .AND. lwdownfile .NE. ' ' ) |
254 |
& lwflux(i,j,bi,bj) = 5.5 _d -08 * |
255 |
& ((theta(i,j,k,bi,bj)+cen2kel)**4) |
256 |
& - lwdown(i,j,bi,bj) |
257 |
if ( lwfluxfile .NE. ' ' .AND. lwdownfile .EQ. ' ' ) |
258 |
& lwdown(i,j,bi,bj) = 5.5 _d -08 * |
259 |
& ((theta(i,j,k,bi,bj)+cen2kel)**4) |
260 |
& - lwflux(i,j,bi,bj) |
261 |
#endif |
262 |
|
263 |
#if defined(ALLOW_ATM_TEMP) || defined(SHORTWAVE_HEATING) |
264 |
if ( swfluxfile .EQ. ' ' .AND. swdownfile .NE. ' ' ) |
265 |
& swflux(i,j,bi,bj) = -(1.0-exf_albedo) * swdown(i,j,bi,bj) |
266 |
if ( swfluxfile .NE. ' ' .AND. swdownfile .EQ. ' ' ) |
267 |
& swdown(i,j,bi,bj) = -swflux(i,j,bi,bj) / (1.0-exf_albedo) |
268 |
#endif |
269 |
|
270 |
#endif /* ALLOW_DOWNWARD_RADIATION */ |
271 |
|
272 |
c-- Compute the turbulent surface fluxes. |
273 |
|
274 |
#ifdef ALLOW_ATM_WIND |
275 |
c Wind speed and direction. |
276 |
ustmp = uwind(i,j,bi,bj)*uwind(i,j,bi,bj) + |
277 |
& vwind(i,j,bi,bj)*vwind(i,j,bi,bj) |
278 |
if ( ustmp .ne. 0. _d 0 ) then |
279 |
us = sqrt(ustmp) |
280 |
cw = uwind(i,j,bi,bj)/us |
281 |
sw = vwind(i,j,bi,bj)/us |
282 |
else |
283 |
us = 0. _d 0 |
284 |
cw = 0. _d 0 |
285 |
sw = 0. _d 0 |
286 |
endif |
287 |
sh = max(us,umin) |
288 |
#else /* ifndef ALLOW_ATM_WIND */ |
289 |
#ifdef ALLOW_ATM_TEMP |
290 |
|
291 |
c The variables us, sh and rdn have to be computed from |
292 |
c given wind stresses inverting relationship for neutral |
293 |
c drag coeff. cdn. |
294 |
c The inversion is based on linear and quadratic form of |
295 |
c cdn(umps); ustar can be directly computed from stress; |
296 |
|
297 |
ustmp = ustress(i,j,bi,bj)*ustress(i,j,bi,bj) + |
298 |
& vstress(i,j,bi,bj)*vstress(i,j,bi,bj) |
299 |
if ( ustmp .ne. 0. _d 0 ) then |
300 |
ustar = sqrt(ustmp/atmrho) |
301 |
cw = ustress(i,j,bi,bj)/sqrt(ustmp) |
302 |
sw = vstress(i,j,bi,bj)/sqrt(ustmp) |
303 |
else |
304 |
ustar = 0. _d 0 |
305 |
cw = 0. _d 0 |
306 |
sw = 0. _d 0 |
307 |
endif |
308 |
|
309 |
if ( ustar .eq. 0. _d 0 ) then |
310 |
us = 0. _d 0 |
311 |
else if ( ustar .lt. ustofu11 ) then |
312 |
tmp1 = -cquadrag_2/cquadrag_1/2 |
313 |
tmp2 = sqrt(tmp1*tmp1 + ustar*ustar/cquadrag_1) |
314 |
us = sqrt(tmp1 + tmp2) |
315 |
else |
316 |
tmp3 = clindrag_2/clindrag_1/3 |
317 |
tmp4 = ustar*ustar/clindrag_1/2 - tmp3**3 |
318 |
tmp5 = sqrt(ustar*ustar/clindrag_1* |
319 |
& (ustar*ustar/clindrag_1/4 - tmp3**3)) |
320 |
us = (tmp4 + tmp5)**(1/3) + |
321 |
& tmp3**2 * (tmp4 + tmp5)**(-1/3) - tmp3 |
322 |
endif |
323 |
|
324 |
if ( us .ne. 0 ) then |
325 |
rdn = ustar/us |
326 |
else |
327 |
rdn = 0. _d 0 |
328 |
end if |
329 |
|
330 |
sh = max(us,umin) |
331 |
#endif /* ALLOW_ATM_TEMP */ |
332 |
#endif /* ifndef ALLOW_ATM_WIND */ |
333 |
|
334 |
#ifdef ALLOW_ATM_TEMP |
335 |
|
336 |
c Initial guess: z/l=0.0; hu=ht=hq=z |
337 |
c Iterations: converge on z/l and hence the fluxes. |
338 |
c t0 : virtual temperature (K) |
339 |
c ssq : sea surface humidity (kg/kg) |
340 |
c deltap : potential temperature diff (K) |
341 |
|
342 |
if ( atemp(i,j,bi,bj) .ne. 0. _d 0 ) then |
343 |
t0 = atemp(i,j,bi,bj)* |
344 |
& (exf_one + humid_fac*aqh(i,j,bi,bj)) |
345 |
ssttmp = theta(i,j,k,bi,bj) |
346 |
tmpbulk= exf_BulkqSat(ssttmp + cen2kel) |
347 |
ssq = saltsat*tmpbulk/atmrho |
348 |
deltap = atemp(i,j,bi,bj) + gamma_blk*ht - |
349 |
& ssttmp - cen2kel |
350 |
delq = aqh(i,j,bi,bj) - ssq |
351 |
stable = exf_half + sign(exf_half, deltap) |
352 |
#ifdef ALLOW_AUTODIFF_TAMC |
353 |
CADJ STORE sh = comlev1_exf_1, key = ikey_1 |
354 |
#endif |
355 |
tmpbulk= exf_BulkCdn(sh) |
356 |
rdn = sqrt(tmpbulk) |
357 |
ustar = rdn*sh |
358 |
tmpbulk= exf_BulkRhn(stable) |
359 |
tstar = tmpbulk*deltap |
360 |
qstar = cdalton*delq |
361 |
|
362 |
do iter = 1,niter_bulk |
363 |
|
364 |
#ifdef ALLOW_AUTODIFF_TAMC |
365 |
ikey_2 = iter |
366 |
& + niter_bulk*(i-1) |
367 |
& + niter_bulk*sNx*(j-1) |
368 |
& + niter_bulk*sNx*sNy*act1 |
369 |
& + niter_bulk*sNx*sNy*max1*act2 |
370 |
& + niter_bulk*sNx*sNy*max1*max2*act3 |
371 |
& + niter_bulk*sNx*sNy*max1*max2*max3*act4 |
372 |
|
373 |
CADJ STORE rdn = comlev1_exf_2, key = ikey_2 |
374 |
CADJ STORE ustar = comlev1_exf_2, key = ikey_2 |
375 |
CADJ STORE qstar = comlev1_exf_2, key = ikey_2 |
376 |
CADJ STORE tstar = comlev1_exf_2, key = ikey_2 |
377 |
CADJ STORE sh = comlev1_exf_2, key = ikey_2 |
378 |
CADJ STORE us = comlev1_exf_2, key = ikey_2 |
379 |
#endif |
380 |
|
381 |
huol = czol*(tstar/t0 + |
382 |
& qstar/(exf_one/humid_fac+aqh(i,j,bi,bj)))/ |
383 |
& ustar**2 |
384 |
huol = max(huol,zolmin) |
385 |
stable = exf_half + sign(exf_half, huol) |
386 |
htol = huol*ht/hu |
387 |
hqol = huol*hq/hu |
388 |
|
389 |
c Evaluate all stability functions assuming hq = ht. |
390 |
xsq = max(sqrt(abs(exf_one - 16.*huol)),exf_one) |
391 |
x = sqrt(xsq) |
392 |
psimh = -psim_fac*huol*stable + |
393 |
& (exf_one - stable)* |
394 |
& (log((exf_one + x*(exf_two + x))* |
395 |
& (exf_one + xsq)/8.) - exf_two*atan(x) + |
396 |
& pi*exf_half) |
397 |
xsq = max(sqrt(abs(exf_one - 16.*htol)),exf_one) |
398 |
psixh = -psim_fac*htol*stable + (exf_one - stable)* |
399 |
& exf_two*log((exf_one + xsq)/exf_two) |
400 |
|
401 |
c Shift wind speed using old coefficient |
402 |
ccc rd = rdn/(exf_one + rdn/karman* |
403 |
ccc & (log(hu/zref) - psimh) ) |
404 |
rd = rdn/(exf_one - rdn/karman*psimh ) |
405 |
shn = sh*rd/rdn |
406 |
uzn = max(shn, umin) |
407 |
|
408 |
c Update the transfer coefficients at 10 meters |
409 |
c and neutral stability. |
410 |
|
411 |
tmpbulk= exf_BulkCdn(uzn) |
412 |
rdn = sqrt(tmpbulk) |
413 |
|
414 |
c Shift all coefficients to the measurement height |
415 |
c and stability. |
416 |
c rd = rdn/(exf_one + rdn/karman*(log(hu/zref) - psimh)) |
417 |
rd = rdn/(exf_one - rdn/karman*psimh) |
418 |
tmpbulk= exf_BulkRhn(stable) |
419 |
rh = tmpbulk/( exf_one + |
420 |
& tmpbulk/karman*(aln - psixh) ) |
421 |
re = cdalton/( exf_one + |
422 |
& cdalton/karman*(aln - psixh) ) |
423 |
|
424 |
c Update ustar, tstar, qstar using updated, shifted |
425 |
c coefficients. |
426 |
ustar = rd*sh |
427 |
qstar = re*delq |
428 |
tstar = rh*deltap |
429 |
tau = atmrho*ustar**2 |
430 |
tau = tau*us/sh |
431 |
|
432 |
enddo |
433 |
|
434 |
#ifdef ALLOW_AUTODIFF_TAMC |
435 |
CADJ STORE ustar = comlev1_exf_1, key = ikey_1 |
436 |
CADJ STORE qstar = comlev1_exf_1, key = ikey_1 |
437 |
CADJ STORE tstar = comlev1_exf_1, key = ikey_1 |
438 |
CADJ STORE tau = comlev1_exf_1, key = ikey_1 |
439 |
CADJ STORE cw = comlev1_exf_1, key = ikey_1 |
440 |
CADJ STORE sw = comlev1_exf_1, key = ikey_1 |
441 |
#endif |
442 |
|
443 |
hs(i,j,bi,bj) = atmcp*tau*tstar/ustar |
444 |
hl(i,j,bi,bj) = flamb*tau*qstar/ustar |
445 |
#ifndef EXF_READ_EVAP |
446 |
cdm evap(i,j,bi,bj) = tau*qstar/ustar |
447 |
cdm !!! need to change sign and to convert from kg/m^2/s to m/s !!! |
448 |
evap(i,j,bi,bj) = -recip_rhonil*tau*qstar/ustar |
449 |
#endif |
450 |
ustress(i,j,bi,bj) = tau*cw |
451 |
vstress(i,j,bi,bj) = tau*sw |
452 |
else |
453 |
ustress(i,j,bi,bj) = 0. _d 0 |
454 |
vstress(i,j,bi,bj) = 0. _d 0 |
455 |
hflux (i,j,bi,bj) = 0. _d 0 |
456 |
hs(i,j,bi,bj) = 0. _d 0 |
457 |
hl(i,j,bi,bj) = 0. _d 0 |
458 |
endif |
459 |
|
460 |
#else /* ifndef ALLOW_ATM_TEMP */ |
461 |
#ifdef ALLOW_ATM_WIND |
462 |
tmpbulk= exf_BulkCdn(sh) |
463 |
ustress(i,j,bi,bj) = atmrho*tmpbulk*us* |
464 |
& uwind(i,j,bi,bj) |
465 |
vstress(i,j,bi,bj) = atmrho*tmpbulk*us* |
466 |
& vwind(i,j,bi,bj) |
467 |
#endif |
468 |
#endif /* ifndef ALLOW_ATM_TEMP */ |
469 |
enddo |
470 |
enddo |
471 |
enddo |
472 |
enddo |
473 |
|
474 |
c Add all contributions. |
475 |
do bj = mybylo(mythid),mybyhi(mythid) |
476 |
do bi = mybxlo(mythid),mybxhi(mythid) |
477 |
do j = 1,sny |
478 |
do i = 1,snx |
479 |
c Net surface heat flux. |
480 |
#ifdef ALLOW_ATM_TEMP |
481 |
hfl = 0. _d 0 |
482 |
hfl = hfl - hs(i,j,bi,bj) |
483 |
hfl = hfl - hl(i,j,bi,bj) |
484 |
hfl = hfl + lwflux(i,j,bi,bj) |
485 |
#ifndef SHORTWAVE_HEATING |
486 |
hfl = hfl + swflux(i,j,bi,bj) |
487 |
#endif |
488 |
c Heat flux: |
489 |
hflux(i,j,bi,bj) = hfl |
490 |
c Salt flux from Precipitation and Evaporation. |
491 |
sflux(i,j,bi,bj) = evap(i,j,bi,bj) - precip(i,j,bi,bj) |
492 |
#endif /* ALLOW_ATM_TEMP */ |
493 |
|
494 |
enddo |
495 |
enddo |
496 |
enddo |
497 |
enddo |
498 |
|
499 |
#endif /* ALLOW_BULKFORMULAE */ |
500 |
|
501 |
end |