1 |
|
c $Header$ |
2 |
|
|
3 |
|
#include "EXF_CPPOPTIONS.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_CPPOPTIONS.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 |
|
c == local variables == |
137 |
|
|
138 |
|
integer bi,bj |
139 |
|
integer i,j,k |
140 |
|
|
141 |
|
#ifdef ALLOW_BULKFORMULAE |
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 hs(1-olx:snx+olx,1-oly:sny+oly,nsx,nsy) |
182 |
|
_RL hl(1-olx:snx+olx,1-oly:sny+oly,nsx,nsy) |
183 |
|
_RL hfl |
184 |
|
|
185 |
|
#endif /* ALLOW_BULKFORMULAE */ |
186 |
|
|
187 |
|
c == external functions == |
188 |
|
|
189 |
|
integer ilnblnk |
190 |
|
external ilnblnk |
191 |
|
|
192 |
|
#ifdef ALLOW_BULKFORMULAE |
193 |
|
_RL exf_BulkqSat |
194 |
|
external exf_BulkqSat |
195 |
|
_RL exf_BulkCdn |
196 |
|
external exf_BulkCdn |
197 |
|
_RL exf_BulkRhn |
198 |
|
external exf_BulkRhn |
199 |
|
#endif /* ALLOW_BULKFORMULAE */ |
200 |
|
|
201 |
|
#ifndef ALLOW_ATM_WIND |
202 |
|
_RL TMP1 |
203 |
|
_RL TMP2 |
204 |
|
_RL TMP3 |
205 |
|
_RL TMP4 |
206 |
|
_RL TMP5 |
207 |
|
#endif |
208 |
|
|
209 |
|
c == end of interface == |
210 |
|
|
211 |
|
#ifdef ALLOW_BULKFORMULAE |
212 |
|
cph This statement cannot be a PARAMETER statement in the header, |
213 |
|
cph but must come here; it's not fortran77 standard |
214 |
|
aln = log(ht/zref) |
215 |
|
#endif |
216 |
|
|
217 |
|
c-- Use atmospheric state to compute surface fluxes. |
218 |
|
|
219 |
|
c Loop over tiles. |
220 |
|
#ifdef ALLOW_AUTODIFF_TAMC |
221 |
|
C-- HPF directive to help TAMC |
222 |
|
CHPF$ INDEPENDENT |
223 |
|
#endif |
224 |
|
do bj = mybylo(mythid),mybyhi(mythid) |
225 |
|
#ifdef ALLOW_AUTODIFF_TAMC |
226 |
|
C-- HPF directive to help TAMC |
227 |
|
CHPF$ INDEPENDENT |
228 |
|
#endif |
229 |
|
do bi = mybxlo(mythid),mybxhi(mythid) |
230 |
|
|
231 |
|
k = 1 |
232 |
|
|
233 |
|
do j = 1,sny |
234 |
|
do i = 1,snx |
235 |
|
|
236 |
|
#ifdef ALLOW_BULKFORMULAE |
237 |
|
|
238 |
|
#ifdef ALLOW_AUTODIFF_TAMC |
239 |
|
act1 = bi - myBxLo(myThid) |
240 |
|
max1 = myBxHi(myThid) - myBxLo(myThid) + 1 |
241 |
|
act2 = bj - myByLo(myThid) |
242 |
|
max2 = myByHi(myThid) - myByLo(myThid) + 1 |
243 |
|
act3 = myThid - 1 |
244 |
|
max3 = nTx*nTy |
245 |
|
act4 = ikey_dynamics - 1 |
246 |
|
|
247 |
|
ikey_1 = i |
248 |
|
& + sNx*(j-1) |
249 |
|
& + sNx*sNy*act1 |
250 |
|
& + sNx*sNy*max1*act2 |
251 |
|
& + sNx*sNy*max1*max2*act3 |
252 |
|
& + sNx*sNy*max1*max2*max3*act4 |
253 |
|
#endif |
254 |
|
|
255 |
|
#ifdef ALLOW_DOWNWARD_RADIATION |
256 |
|
c-- Compute net from downward and downward from net longwave and |
257 |
|
c shortwave radiation, if needed. |
258 |
|
c lwflux = Stefan-Boltzman constant * emissivity * SST - lwdown |
259 |
|
c swflux = - ( 1 - albedo ) * swdown |
260 |
|
|
261 |
|
#ifdef ALLOW_ATM_TEMP |
262 |
|
if ( lwfluxfile .EQ. ' ' .AND. lwdownfile .NE. ' ' ) |
263 |
|
& lwflux(i,j,bi,bj) = 5.5 _d -08 * |
264 |
|
& ((theta(i,j,k,bi,bj)+cen2kel)**4) |
265 |
|
& - lwdown(i,j,bi,bj) |
266 |
|
if ( lwfluxfile .NE. ' ' .AND. lwdownfile .EQ. ' ' ) |
267 |
|
& lwdown(i,j,bi,bj) = 5.5 _d -08 * |
268 |
|
& ((theta(i,j,k,bi,bj)+cen2kel)**4) |
269 |
|
& - lwflux(i,j,bi,bj) |
270 |
|
#endif |
271 |
|
|
272 |
|
#if defined(ALLOW_ATM_TEMP) || defined(SHORTWAVE_HEATING) |
273 |
|
if ( swfluxfile .EQ. ' ' .AND. swdownfile .NE. ' ' ) |
274 |
|
& swflux(i,j,bi,bj) = -0.9 _d 0 * swdown(i,j,bi,bj) |
275 |
|
if ( swfluxfile .NE. ' ' .AND. swdownfile .EQ. ' ' ) |
276 |
|
& swdown(i,j,bi,bj) = -1.111111 _d 0 * swflux(i,j,bi,bj) |
277 |
|
#endif |
278 |
|
|
279 |
|
#endif /* ALLOW_DOWNWARD_RADIATION */ |
280 |
|
|
281 |
|
c-- Compute the turbulent surface fluxes. |
282 |
|
|
283 |
|
#ifdef ALLOW_ATM_WIND |
284 |
|
c Wind speed and direction. |
285 |
|
ustmp = uwind(i,j,bi,bj)*uwind(i,j,bi,bj) + |
286 |
|
& vwind(i,j,bi,bj)*vwind(i,j,bi,bj) |
287 |
|
if ( ustmp .ne. 0. _d 0 ) then |
288 |
|
us = sqrt(ustmp) |
289 |
|
cw = uwind(i,j,bi,bj)/us |
290 |
|
sw = vwind(i,j,bi,bj)/us |
291 |
|
else |
292 |
|
us = 0. _d 0 |
293 |
|
cw = 0. _d 0 |
294 |
|
sw = 0. _d 0 |
295 |
|
endif |
296 |
|
sh = max(us,umin) |
297 |
|
#else /* ifndef ALLOW_ATM_WIND */ |
298 |
|
#ifdef ALLOW_ATM_TEMP |
299 |
|
|
300 |
|
c The variables us, sh and rdn have to be computed from |
301 |
|
c given wind stresses inverting relationship for neutral |
302 |
|
c drag coeff. cdn. |
303 |
|
c The inversion is based on linear and quadratic form of |
304 |
|
c cdn(umps); ustar can be directly computed from stress; |
305 |
|
|
306 |
|
ustmp = ustress(i,j,bi,bj)*ustress(i,j,bi,bj) + |
307 |
|
& vstress(i,j,bi,bj)*vstress(i,j,bi,bj) |
308 |
|
if ( ustmp .ne. 0. _d 0 ) then |
309 |
|
ustar = sqrt(ustmp/atmrho) |
310 |
|
cw = ustress(i,j,bi,bj)/sqrt(ustmp) |
311 |
|
sw = vstress(i,j,bi,bj)/sqrt(ustmp) |
312 |
|
else |
313 |
|
ustar = 0. _d 0 |
314 |
|
cw = 0. _d 0 |
315 |
|
sw = 0. _d 0 |
316 |
|
endif |
317 |
|
|
318 |
|
if ( ustar .eq. 0. _d 0 ) then |
319 |
|
us = 0. _d 0 |
320 |
|
else if ( ustar .lt. ustofu11 ) then |
321 |
|
tmp1 = -cquadrag_2/cquadrag_1/2 |
322 |
|
tmp2 = sqrt(tmp1*tmp1 + ustar*ustar/cquadrag_1) |
323 |
|
us = sqrt(tmp1 + tmp2) |
324 |
|
else |
325 |
|
tmp3 = clindrag_2/clindrag_1/3 |
326 |
|
tmp4 = ustar*ustar/clindrag_1/2 - tmp3**3 |
327 |
|
tmp5 = sqrt(ustar*ustar/clindrag_1* |
328 |
|
& (ustar*ustar/clindrag_1/4 - tmp3**3)) |
329 |
|
us = (tmp4 + tmp5)**(1/3) + |
330 |
|
& tmp3**2 * (tmp4 + tmp5)**(-1/3) - tmp3 |
331 |
|
endif |
332 |
|
|
333 |
|
if ( us .ne. 0 ) then |
334 |
|
rdn = ustar/us |
335 |
|
else |
336 |
|
rdn = 0. _d 0 |
337 |
|
end if |
338 |
|
|
339 |
|
sh = max(us,umin) |
340 |
|
#endif /* ALLOW_ATM_TEMP */ |
341 |
|
#endif /* ifndef ALLOW_ATM_WIND */ |
342 |
|
|
343 |
|
#ifdef ALLOW_ATM_TEMP |
344 |
|
|
345 |
|
c Initial guess: z/l=0.0; hu=ht=hq=z |
346 |
|
c Iterations: converge on z/l and hence the fluxes. |
347 |
|
c t0 : virtual temperature (K) |
348 |
|
c ssq : sea surface humidity (kg/kg) |
349 |
|
c deltap : potential temperature diff (K) |
350 |
|
|
351 |
|
if ( atemp(i,j,bi,bj) .ne. 0. _d 0 ) then |
352 |
|
t0 = atemp(i,j,bi,bj)* |
353 |
|
& (exf_one + humid_fac*aqh(i,j,bi,bj)) |
354 |
|
ssttmp = theta(i,j,k,bi,bj) |
355 |
|
ssq = saltsat* |
356 |
|
& exf_BulkqSat(ssttmp + cen2kel)/ |
357 |
|
& atmrho |
358 |
|
deltap = atemp(i,j,bi,bj) + gamma_blk*ht - |
359 |
|
& ssttmp - cen2kel |
360 |
|
delq = aqh(i,j,bi,bj) - ssq |
361 |
|
stable = exf_half + sign(exf_half, deltap) |
362 |
|
#ifdef ALLOW_AUTODIFF_TAMC |
363 |
|
CADJ STORE sh = comlev1_exf_1, key = ikey_1 |
364 |
|
#endif |
365 |
|
rdn = sqrt(exf_BulkCdn(sh)) |
366 |
|
ustar = rdn*sh |
367 |
|
tstar = exf_BulkRhn(stable)*deltap |
368 |
|
qstar = cdalton*delq |
369 |
|
|
370 |
|
do iter = 1,niter_bulk |
371 |
|
|
372 |
|
#ifdef ALLOW_AUTODIFF_TAMC |
373 |
|
ikey_2 = iter |
374 |
|
& + niter_bulk*(i-1) |
375 |
|
& + sNx*niter_bulk*(j-1) |
376 |
|
& + sNx*niter_bulk*sNy*act1 |
377 |
|
& + sNx*niter_bulk*sNy*max1*act2 |
378 |
|
& + sNx*niter_bulk*sNy*max1*max2*act3 |
379 |
|
& + sNx*niter_bulk*sNy*max1*max2*max3*act4 |
380 |
|
|
381 |
|
CADJ STORE rdn = comlev1_exf_2, key = ikey_2 |
382 |
|
CADJ STORE ustar = comlev1_exf_2, key = ikey_2 |
383 |
|
CADJ STORE qstar = comlev1_exf_2, key = ikey_2 |
384 |
|
CADJ STORE tstar = comlev1_exf_2, key = ikey_2 |
385 |
|
CADJ STORE sh = comlev1_exf_2, key = ikey_2 |
386 |
|
CADJ STORE us = comlev1_exf_2, key = ikey_2 |
387 |
|
#endif |
388 |
|
|
389 |
|
huol = czol*(tstar/t0 + |
390 |
|
& qstar/(exf_one/humid_fac+aqh(i,j,bi,bj)))/ |
391 |
|
& ustar**2 |
392 |
|
huol = max(huol,zolmin) |
393 |
|
stable = exf_half + sign(exf_half, huol) |
394 |
|
htol = huol*ht/hu |
395 |
|
hqol = huol*hq/hu |
396 |
|
|
397 |
|
c Evaluate all stability functions assuming hq = ht. |
398 |
|
xsq = max(sqrt(abs(exf_one - 16.*huol)),exf_one) |
399 |
|
x = sqrt(xsq) |
400 |
|
psimh = -psim_fac*huol*stable + |
401 |
|
& (exf_one - stable)* |
402 |
|
& log((exf_one + x*(exf_two + x))* |
403 |
|
& (exf_one + xsq)/8.) - exf_two*atan(x) + |
404 |
|
& pi*exf_half |
405 |
|
xsq = max(sqrt(abs(exf_one - 16.*htol)),exf_one) |
406 |
|
psixh = -psim_fac*htol*stable + (exf_one - stable)* |
407 |
|
& exf_two*log((exf_one + xsq)/exf_two) |
408 |
|
|
409 |
|
c Shift wind speed using old coefficient |
410 |
|
ccc rd = rdn/(exf_one + rdn/karman* |
411 |
|
ccc & (log(hu/zref) - psimh) ) |
412 |
|
rd = rdn/(exf_one - rdn/karman*psimh ) |
413 |
|
shn = sh*rd/rdn |
414 |
|
uzn = max(shn, umin) |
415 |
|
|
416 |
|
c Update the transfer coefficients at 10 meters |
417 |
|
c and neutral stability. |
418 |
|
|
419 |
|
rdn = sqrt(exf_BulkCdn(uzn)) |
420 |
|
|
421 |
|
c Shift all coefficients to the measurement height |
422 |
|
c and stability. |
423 |
|
c rd = rdn/(exf_one + rdn/karman*(log(hu/zref) - psimh)) |
424 |
|
rd = rdn/(exf_one - rdn/karman*psimh) |
425 |
|
rh = exf_BulkRhn(stable)/(exf_one + |
426 |
|
& exf_BulkRhn(stable)/ |
427 |
|
& karman*(aln - psixh)) |
428 |
|
re = cdalton/(exf_one + cdalton/karman*(aln - psixh)) |
429 |
|
|
430 |
|
c Update ustar, tstar, qstar using updated, shifted |
431 |
|
c coefficients. |
432 |
|
ustar = rd*sh |
433 |
|
qstar = re*delq |
434 |
|
tstar = rh*deltap |
435 |
|
tau = atmrho*ustar**2 |
436 |
|
tau = tau*us/sh |
437 |
|
|
438 |
|
enddo |
439 |
|
|
440 |
|
#ifdef ALLOW_AUTODIFF_TAMC |
441 |
|
CADJ STORE ustar = comlev1_exf_1, key = ikey_1 |
442 |
|
CADJ STORE qstar = comlev1_exf_1, key = ikey_1 |
443 |
|
CADJ STORE tstar = comlev1_exf_1, key = ikey_1 |
444 |
|
CADJ STORE tau = comlev1_exf_1, key = ikey_1 |
445 |
|
CADJ STORE cw = comlev1_exf_1, key = ikey_1 |
446 |
|
CADJ STORE sw = comlev1_exf_1, key = ikey_1 |
447 |
|
#endif |
448 |
|
|
449 |
|
hs(i,j,bi,bj) = atmcp*tau*tstar/ustar |
450 |
|
hl(i,j,bi,bj) = flamb*tau*qstar/ustar |
451 |
|
#ifndef EXF_READ_EVAP |
452 |
|
cdm evap(i,j,bi,bj) = tau*qstar/ustar |
453 |
|
cdm !!! need to change sign and to convert from kg/m^2/s to m/s !!! |
454 |
|
evap(i,j,bi,bj) = -recip_rhonil*tau*qstar/ustar |
455 |
|
#endif |
456 |
|
ustress(i,j,bi,bj) = tau*cw |
457 |
|
vstress(i,j,bi,bj) = tau*sw |
458 |
|
else |
459 |
|
ustress(i,j,bi,bj) = 0. _d 0 |
460 |
|
vstress(i,j,bi,bj) = 0. _d 0 |
461 |
|
hflux (i,j,bi,bj) = 0. _d 0 |
462 |
|
hs(i,j,bi,bj) = 0. _d 0 |
463 |
|
hl(i,j,bi,bj) = 0. _d 0 |
464 |
|
endif |
465 |
|
|
466 |
|
#else /* ifndef ALLOW_ATM_TEMP */ |
467 |
|
#ifdef ALLOW_ATM_WIND |
468 |
|
ustress(i,j,bi,bj) = atmrho*exf_BulkCdn(sh)*us* |
469 |
|
& uwind(i,j,bi,bj) |
470 |
|
vstress(i,j,bi,bj) = atmrho*exf_BulkCdn(sh)*us* |
471 |
|
& vwind(i,j,bi,bj) |
472 |
|
#endif |
473 |
|
#endif /* ifndef ALLOW_ATM_TEMP */ |
474 |
|
enddo |
475 |
|
enddo |
476 |
|
enddo |
477 |
|
enddo |
478 |
|
|
479 |
|
c Add all contributions. |
480 |
|
do bj = mybylo(mythid),mybyhi(mythid) |
481 |
|
do bi = mybxlo(mythid),mybxhi(mythid) |
482 |
|
do j = 1,sny |
483 |
|
do i = 1,snx |
484 |
|
c Net surface heat flux. |
485 |
|
#ifdef ALLOW_ATM_TEMP |
486 |
|
hfl = 0. _d 0 |
487 |
|
hfl = hfl - hs(i,j,bi,bj) |
488 |
|
hfl = hfl - hl(i,j,bi,bj) |
489 |
|
hfl = hfl + lwflux(i,j,bi,bj) |
490 |
|
#ifndef SHORTWAVE_HEATING |
491 |
|
hfl = hfl + swflux(i,j,bi,bj) |
492 |
|
#endif |
493 |
|
c Heat flux: |
494 |
|
hflux(i,j,bi,bj) = hfl |
495 |
|
c Salt flux from Precipitation and Evaporation. |
496 |
|
sflux(i,j,bi,bj) = evap(i,j,bi,bj) - precip(i,j,bi,bj) |
497 |
|
#endif /* ALLOW_ATM_TEMP */ |
498 |
|
|
499 |
|
#endif /* ALLOW_BULKFORMULAE */ |
500 |
|
|
501 |
|
enddo |
502 |
|
enddo |
503 |
|
enddo |
504 |
|
enddo |
505 |
|
|
506 |
|
end |