/[MITgcm]/MITgcm/pkg/exf/exf_getffields.F
ViewVC logotype

Annotation of /MITgcm/pkg/exf/exf_getffields.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.2.6.7 - (hide annotations) (download)
Mon Mar 3 17:27:31 2003 UTC (21 years, 3 months ago) by heimbach
Branch: ecco-branch
CVS Tags: icebear5, ecco_c44_e27
Branch point for: icebear
Changes since 1.2.6.6: +7 -8 lines
o in exf_getffields.F
  Reduce i-/j-loop to interior domain, discarding overlaps.
  That also fixes wrong TAF-key computations for key_1, key_2
  with bulf formulae.
o exf_init.F
  modify #ifdef for exf_init_evap
o exf_getffieldrec.F, ctrl_getrec.F
  The following INT-usages are not safe:
      fldsecs  = int(fldsecs/fldperiod)*fldperiod
      fldcount = int(fldsecs/fldperiod) + 1
  and were modified.

1 heimbach 1.2.6.7 c $Header: /u/gcmpack/MITgcm/pkg/exf/exf_getffields.F,v 1.2.6.6 2003/02/28 09:34:05 dimitri Exp $
2 heimbach 1.1
3     #include "EXF_CPPOPTIONS.h"
4    
5 heimbach 1.2.6.1 subroutine exf_GetFFields( mycurrenttime, mycurrentiter, mythid )
6 heimbach 1.1
7     c ==================================================================
8     c SUBROUTINE exf_GetFFields
9     c ==================================================================
10     c
11 dimitri 1.2.6.3 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 heimbach 1.1 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 dimitri 1.2.6.3 c ==================================================================
80 heimbach 1.1 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 heimbach 1.2.6.2 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 dimitri 1.2.6.3 c mods for pkg/seaice: menemenlis@jpl.nasa.gov 20-Dec-2002
107     c
108 heimbach 1.1 c ==================================================================
109     c SUBROUTINE exf_GetFFields
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 dimitri 1.2.6.5 #include "exf_param.h"
123 heimbach 1.1 #include "exf_fields.h"
124     #include "exf_constants.h"
125    
126 heimbach 1.2.6.1 #ifdef ALLOW_AUTODIFF_TAMC
127     #include "tamc.h"
128     #endif
129    
130 heimbach 1.1 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     #ifdef ALLOW_ATM_TEMP
144     integer iter
145     _RL aln
146     _RL delq
147     _RL deltap
148     _RL hqol
149     _RL htol
150     _RL huol
151     _RL psimh
152     _RL psixh
153     _RL qstar
154     _RL rd
155     _RL re
156     _RL rdn
157     _RL rh
158 heimbach 1.2.6.1 _RL ssttmp
159 heimbach 1.1 _RL ssq
160     _RL stable
161     _RL tstar
162     _RL t0
163     _RL ustar
164     _RL uzn
165 heimbach 1.2.6.1 _RL shn
166 heimbach 1.1 _RL xsq
167     _RL x
168 heimbach 1.2.6.1 _RL tau
169     #ifdef ALLOW_AUTODIFF_TAMC
170     integer ikey_1
171     integer ikey_2
172     #endif
173 heimbach 1.1 #endif /* ALLOW_ATM_TEMP */
174    
175 heimbach 1.2.6.1 _RL ustmp
176 heimbach 1.1 _RL us
177     _RL cw
178     _RL sw
179     _RL sh
180     _RL hs(1-olx:snx+olx,1-oly:sny+oly,nsx,nsy)
181     _RL hl(1-olx:snx+olx,1-oly:sny+oly,nsx,nsy)
182     _RL hfl
183    
184     #endif /* ALLOW_BULKFORMULAE */
185    
186     c == external functions ==
187    
188     integer ilnblnk
189     external ilnblnk
190    
191     #ifdef ALLOW_BULKFORMULAE
192     _RL exf_BulkqSat
193     external exf_BulkqSat
194     _RL exf_BulkCdn
195     external exf_BulkCdn
196     _RL exf_BulkRhn
197     external exf_BulkRhn
198     #endif /* ALLOW_BULKFORMULAE */
199    
200 dimitri 1.2.6.3 #ifndef ALLOW_ATM_WIND
201     _RL TMP1
202     _RL TMP2
203     _RL TMP3
204     _RL TMP4
205     _RL TMP5
206     #endif
207 heimbach 1.1
208 dimitri 1.2.6.3 c == end of interface ==
209 heimbach 1.1
210     #ifdef ALLOW_BULKFORMULAE
211 heimbach 1.2.6.2 cph This statement cannot be a PARAMETER statement in the header,
212     cph but must come here; it's not fortran77 standard
213 heimbach 1.2.6.1 aln = log(ht/zref)
214     #endif
215    
216 dimitri 1.2.6.3 c-- read forcing fields from files and temporal interpolation
217 heimbach 1.1
218 dimitri 1.2.6.3 #ifdef ALLOW_ATM_WIND
219 heimbach 1.1
220 dimitri 1.2.6.3 c Zonal wind.
221     call exf_set_uwind ( mycurrenttime, mycurrentiter, mythid )
222 heimbach 1.1
223 dimitri 1.2.6.3 c Meridional wind.
224     call exf_set_vwind ( mycurrenttime, mycurrentiter, mythid )
225 heimbach 1.1
226 dimitri 1.2.6.3 #ifdef ALLOW_UWIND_CONTROL
227     call ctrl_getuwind ( mycurrenttime, mycurrentiter, mythid )
228 heimbach 1.2.6.1 #endif
229    
230 dimitri 1.2.6.3 #ifdef ALLOW_VWIND_CONTROL
231     call ctrl_getvwind ( mycurrenttime, mycurrentiter, mythid )
232 heimbach 1.2.6.1 #endif
233 heimbach 1.1
234 dimitri 1.2.6.3 #else /* ifndef ALLOW_ATM_WIND */
235 heimbach 1.1
236 dimitri 1.2.6.3 c Zonal wind stress.
237     call exf_set_ustress( mycurrenttime, mycurrentiter, mythid )
238 heimbach 1.1
239 dimitri 1.2.6.3 c Meridional wind stress.
240     call exf_set_vstress( mycurrenttime, mycurrentiter, mythid )
241 heimbach 1.1
242 dimitri 1.2.6.3 #endif /* ifndef ALLOW_ATM_WIND */
243 heimbach 1.1
244 dimitri 1.2.6.3 #ifdef ALLOW_ATM_TEMP
245 heimbach 1.1
246 dimitri 1.2.6.3 c Atmospheric temperature.
247     call exf_set_atemp ( mycurrenttime, mycurrentiter, mythid )
248 heimbach 1.1
249 dimitri 1.2.6.3 c Atmospheric humidity.
250     call exf_set_aqh ( mycurrenttime, mycurrentiter, mythid )
251 heimbach 1.1
252 dimitri 1.2.6.3 c Net long wave radiative flux.
253     call exf_set_lwflux ( mycurrenttime, mycurrentiter, mythid )
254 heimbach 1.2.6.1
255 dimitri 1.2.6.3 c Precipitation.
256     call exf_set_precip ( mycurrenttime, mycurrentiter, mythid )
257 heimbach 1.2.6.1
258 dimitri 1.2.6.3 #ifdef ALLOW_ATEMP_CONTROL
259     call ctrl_getatemp ( mycurrenttime, mycurrentiter, mythid )
260 heimbach 1.2.6.1 #endif
261 heimbach 1.1
262 dimitri 1.2.6.3 #ifdef ALLOW_AQH_CONTROL
263     call ctrl_getaqh ( mycurrenttime, mycurrentiter, mythid )
264     #endif
265 heimbach 1.1
266 dimitri 1.2.6.3 #else /* ifndef ALLOW_ATM_TEMP */
267 heimbach 1.1
268     c Atmospheric heat flux.
269 dimitri 1.2.6.3 call exf_set_hflux ( mycurrenttime, mycurrentiter, mythid )
270 heimbach 1.1
271     c Salt flux.
272 dimitri 1.2.6.3 call exf_set_sflux ( mycurrenttime, mycurrentiter, mythid )
273 heimbach 1.1
274 dimitri 1.2.6.3 #endif /* ifndef ALLOW_ATM_TEMP */
275 heimbach 1.1
276 dimitri 1.2.6.3 #if defined(ALLOW_ATM_TEMP) || defined(SHORTWAVE_HEATING)
277 heimbach 1.1 c Net short wave radiative flux.
278 dimitri 1.2.6.3 call exf_set_swflux ( mycurrenttime, mycurrentiter, mythid )
279 heimbach 1.2.6.1 #endif
280 heimbach 1.1
281 dimitri 1.2.6.3 #ifdef EXF_READ_EVAP
282     c Evaporation
283     call exf_set_evap ( mycurrenttime, mycurrentiter, mythid )
284 dimitri 1.2.6.4 #endif
285 dimitri 1.2.6.3
286     #ifdef ALLOW_DOWNWARD_RADIATION
287    
288     c Downward shortwave radiation.
289     call exf_set_swdown ( mycurrenttime, mycurrentiter, mythid )
290    
291     c Downward longwave radiation.
292     call exf_set_lwdown ( mycurrenttime, mycurrentiter, mythid )
293    
294     #endif
295    
296     c-- Use atmospheric state to compute surface fluxes.
297 heimbach 1.1
298     c Loop over tiles.
299 heimbach 1.2.6.1 #ifdef ALLOW_AUTODIFF_TAMC
300     C-- HPF directive to help TAMC
301     CHPF$ INDEPENDENT
302 dimitri 1.2.6.3 #endif
303 heimbach 1.1 do bj = mybylo(mythid),mybyhi(mythid)
304 heimbach 1.2.6.1 #ifdef ALLOW_AUTODIFF_TAMC
305     C-- HPF directive to help TAMC
306     CHPF$ INDEPENDENT
307     #endif
308 heimbach 1.1 do bi = mybxlo(mythid),mybxhi(mythid)
309 heimbach 1.2.6.1
310 heimbach 1.1 k = 1
311    
312 heimbach 1.2.6.7 do j = 1,sny
313     do i = 1,snx
314 heimbach 1.1
315     #ifdef ALLOW_BULKFORMULAE
316    
317 heimbach 1.2.6.1 #ifdef ALLOW_AUTODIFF_TAMC
318     act1 = bi - myBxLo(myThid)
319     max1 = myBxHi(myThid) - myBxLo(myThid) + 1
320     act2 = bj - myByLo(myThid)
321     max2 = myByHi(myThid) - myByLo(myThid) + 1
322     act3 = myThid - 1
323     max3 = nTx*nTy
324     act4 = ikey_dynamics - 1
325    
326 heimbach 1.2.6.7 ikey_1 = i
327     & + sNx*(j-1)
328 heimbach 1.2.6.1 & + sNx*sNy*act1
329     & + sNx*sNy*max1*act2
330     & + sNx*sNy*max1*max2*act3
331     & + sNx*sNy*max1*max2*max3*act4
332     #endif
333    
334 dimitri 1.2.6.3 #ifdef ALLOW_DOWNWARD_RADIATION
335 dimitri 1.2.6.5 c-- Compute net from downward and downward from net longwave and
336     c shortwave radiation, if needed.
337 dimitri 1.2.6.3 c lwflux = Stefan-Boltzman constant * emissivity * SST - lwdown
338     c swflux = - ( 1 - albedo ) * swdown
339 dimitri 1.2.6.5
340     #ifdef ALLOW_ATM_TEMP
341     if ( lwfluxfile .EQ. ' ' .AND. lwdownfile .NE. ' ' )
342     & lwflux(i,j,bi,bj) = 5.5 _d -08 *
343 dimitri 1.2.6.3 & ((theta(i,j,k,bi,bj)+cen2kel)**4)
344     & - lwdown(i,j,bi,bj)
345 dimitri 1.2.6.5 if ( lwfluxfile .NE. ' ' .AND. lwdownfile .EQ. ' ' )
346     & lwdown(i,j,bi,bj) = 5.5 _d -08 *
347     & ((theta(i,j,k,bi,bj)+cen2kel)**4)
348     & - lwflux(i,j,bi,bj)
349 dimitri 1.2.6.3 #endif
350 dimitri 1.2.6.5
351     #if defined(ALLOW_ATM_TEMP) || defined(SHORTWAVE_HEATING)
352     if ( swfluxfile .EQ. ' ' .AND. swdownfile .NE. ' ' )
353     & swflux(i,j,bi,bj) = -0.9 _d 0 * swdown(i,j,bi,bj)
354     if ( swfluxfile .NE. ' ' .AND. swdownfile .EQ. ' ' )
355     & swdown(i,j,bi,bj) = -1.111111 _d 0 * swflux(i,j,bi,bj)
356     #endif
357    
358     #endif /* ALLOW_DOWNWARD_RADIATION */
359 dimitri 1.2.6.3
360     c-- Compute the turbulent surface fluxes.
361 heimbach 1.1
362     #ifdef ALLOW_ATM_WIND
363     c Wind speed and direction.
364 heimbach 1.2.6.1 ustmp = uwind(i,j,bi,bj)*uwind(i,j,bi,bj) +
365     & vwind(i,j,bi,bj)*vwind(i,j,bi,bj)
366     if ( ustmp .ne. 0. _d 0 ) then
367     us = sqrt(ustmp)
368 heimbach 1.1 cw = uwind(i,j,bi,bj)/us
369     sw = vwind(i,j,bi,bj)/us
370     else
371 heimbach 1.2.6.1 us = 0. _d 0
372 heimbach 1.1 cw = 0. _d 0
373     sw = 0. _d 0
374     endif
375     sh = max(us,umin)
376 dimitri 1.2.6.3 #else /* ifndef ALLOW_ATM_WIND */
377 heimbach 1.1 #ifdef ALLOW_ATM_TEMP
378    
379 heimbach 1.2.6.1 c The variables us, sh and rdn have to be computed from
380     c given wind stresses inverting relationship for neutral
381     c drag coeff. cdn.
382 heimbach 1.1 c The inversion is based on linear and quadratic form of
383     c cdn(umps); ustar can be directly computed from stress;
384    
385 heimbach 1.2.6.1 ustmp = ustress(i,j,bi,bj)*ustress(i,j,bi,bj) +
386     & vstress(i,j,bi,bj)*vstress(i,j,bi,bj)
387     if ( ustmp .ne. 0. _d 0 ) then
388     ustar = sqrt(ustmp/atmrho)
389     cw = ustress(i,j,bi,bj)/sqrt(ustmp)
390     sw = vstress(i,j,bi,bj)/sqrt(ustmp)
391     else
392     ustar = 0. _d 0
393     cw = 0. _d 0
394     sw = 0. _d 0
395     endif
396 heimbach 1.1
397     if ( ustar .eq. 0. _d 0 ) then
398     us = 0. _d 0
399     else if ( ustar .lt. ustofu11 ) then
400     tmp1 = -cquadrag_2/cquadrag_1/2
401     tmp2 = sqrt(tmp1*tmp1 + ustar*ustar/cquadrag_1)
402     us = sqrt(tmp1 + tmp2)
403     else
404     tmp3 = clindrag_2/clindrag_1/3
405     tmp4 = ustar*ustar/clindrag_1/2 - tmp3**3
406     tmp5 = sqrt(ustar*ustar/clindrag_1*
407     & (ustar*ustar/clindrag_1/4 - tmp3**3))
408     us = (tmp4 + tmp5)**(1/3) +
409     & tmp3**2 * (tmp4 + tmp5)**(-1/3) - tmp3
410     endif
411    
412     if ( us .ne. 0 ) then
413     rdn = ustar/us
414     else
415     rdn = 0. _d 0
416     end if
417    
418     sh = max(us,umin)
419     #endif /* ALLOW_ATM_TEMP */
420 dimitri 1.2.6.3 #endif /* ifndef ALLOW_ATM_WIND */
421 heimbach 1.1
422     #ifdef ALLOW_ATM_TEMP
423 heimbach 1.2.6.1
424 heimbach 1.1 c Initial guess: z/l=0.0; hu=ht=hq=z
425     c Iterations: converge on z/l and hence the fluxes.
426     c t0 : virtual temperature (K)
427     c ssq : sea surface humidity (kg/kg)
428     c deltap : potential temperature diff (K)
429    
430     if ( atemp(i,j,bi,bj) .ne. 0. _d 0 ) then
431     t0 = atemp(i,j,bi,bj)*
432     & (exf_one + humid_fac*aqh(i,j,bi,bj))
433 heimbach 1.2.6.1 ssttmp = theta(i,j,k,bi,bj)
434 heimbach 1.1 ssq = saltsat*
435 heimbach 1.2.6.1 & exf_BulkqSat(ssttmp + cen2kel)/
436 heimbach 1.1 & atmrho
437     deltap = atemp(i,j,bi,bj) + gamma_blk*ht -
438 heimbach 1.2.6.1 & ssttmp - cen2kel
439 heimbach 1.1 delq = aqh(i,j,bi,bj) - ssq
440     stable = exf_half + sign(exf_half, deltap)
441 heimbach 1.2.6.1 #ifdef ALLOW_AUTODIFF_TAMC
442     CADJ STORE sh = comlev1_exf_1, key = ikey_1
443     #endif
444 heimbach 1.1 rdn = sqrt(exf_BulkCdn(sh))
445     ustar = rdn*sh
446     tstar = exf_BulkRhn(stable)*deltap
447     qstar = cdalton*delq
448    
449     do iter = 1,niter_bulk
450    
451 heimbach 1.2.6.1 #ifdef ALLOW_AUTODIFF_TAMC
452     ikey_2 = iter
453 heimbach 1.2.6.7 & + niter_bulk*(i-1)
454     & + sNx*niter_bulk*(j-1)
455 heimbach 1.2.6.1 & + sNx*niter_bulk*sNy*act1
456     & + sNx*niter_bulk*sNy*max1*act2
457     & + sNx*niter_bulk*sNy*max1*max2*act3
458     & + sNx*niter_bulk*sNy*max1*max2*max3*act4
459    
460     CADJ STORE rdn = comlev1_exf_2, key = ikey_2
461     CADJ STORE ustar = comlev1_exf_2, key = ikey_2
462     CADJ STORE qstar = comlev1_exf_2, key = ikey_2
463     CADJ STORE tstar = comlev1_exf_2, key = ikey_2
464     CADJ STORE sh = comlev1_exf_2, key = ikey_2
465     CADJ STORE us = comlev1_exf_2, key = ikey_2
466     #endif
467    
468 heimbach 1.1 huol = czol*(tstar/t0 +
469     & qstar/(exf_one/humid_fac+aqh(i,j,bi,bj)))/
470     & ustar**2
471     huol = max(huol,zolmin)
472     stable = exf_half + sign(exf_half, huol)
473     htol = huol*ht/hu
474     hqol = huol*hq/hu
475    
476     c Evaluate all stability functions assuming hq = ht.
477     xsq = max(sqrt(abs(exf_one - 16.*huol)),exf_one)
478     x = sqrt(xsq)
479     psimh = -psim_fac*huol*stable +
480     & (exf_one - stable)*
481     & log((exf_one + x*(exf_two + x))*
482     & (exf_one + xsq)/8.) - exf_two*atan(x) +
483     & pi*exf_half
484     xsq = max(sqrt(abs(exf_one - 16.*htol)),exf_one)
485     psixh = -psim_fac*htol*stable + (exf_one - stable)*
486     & exf_two*log((exf_one + xsq)/exf_two)
487    
488     c Shift wind speed using old coefficient
489 heimbach 1.2.6.1 ccc rd = rdn/(exf_one + rdn/karman*
490     ccc & (log(hu/zref) - psimh) )
491 heimbach 1.1 rd = rdn/(exf_one - rdn/karman*psimh )
492 heimbach 1.2.6.1 shn = sh*rd/rdn
493     uzn = max(shn, umin)
494 heimbach 1.1
495     c Update the transfer coefficients at 10 meters
496     c and neutral stability.
497 heimbach 1.2.6.1
498 heimbach 1.1 rdn = sqrt(exf_BulkCdn(uzn))
499    
500     c Shift all coefficients to the measurement height
501     c and stability.
502     c rd = rdn/(exf_one + rdn/karman*(log(hu/zref) - psimh))
503     rd = rdn/(exf_one - rdn/karman*psimh)
504     rh = exf_BulkRhn(stable)/(exf_one +
505     & exf_BulkRhn(stable)/
506     & karman*(aln - psixh))
507     re = cdalton/(exf_one + cdalton/karman*(aln - psixh))
508    
509     c Update ustar, tstar, qstar using updated, shifted
510     c coefficients.
511     ustar = rd*sh
512     qstar = re*delq
513     tstar = rh*deltap
514 heimbach 1.2.6.1 tau = atmrho*ustar**2
515     tau = tau*us/sh
516 heimbach 1.1
517     enddo
518    
519 heimbach 1.2.6.1 #ifdef ALLOW_AUTODIFF_TAMC
520     CADJ STORE ustar = comlev1_exf_1, key = ikey_1
521     CADJ STORE qstar = comlev1_exf_1, key = ikey_1
522     CADJ STORE tstar = comlev1_exf_1, key = ikey_1
523     CADJ STORE tau = comlev1_exf_1, key = ikey_1
524     CADJ STORE cw = comlev1_exf_1, key = ikey_1
525     CADJ STORE sw = comlev1_exf_1, key = ikey_1
526     #endif
527    
528 heimbach 1.1 hs(i,j,bi,bj) = atmcp*tau*tstar/ustar
529     hl(i,j,bi,bj) = flamb*tau*qstar/ustar
530 dimitri 1.2.6.3 #ifndef EXF_READ_EVAP
531     cdm evap(i,j,bi,bj) = tau*qstar/ustar
532     cdm !!! need to change sign and to convert from kg/m^2/s to m/s !!!
533     evap(i,j,bi,bj) = -recip_rhonil*tau*qstar/ustar
534     #endif
535 heimbach 1.1 ustress(i,j,bi,bj) = tau*cw
536     vstress(i,j,bi,bj) = tau*sw
537     else
538     ustress(i,j,bi,bj) = 0. _d 0
539     vstress(i,j,bi,bj) = 0. _d 0
540     hflux (i,j,bi,bj) = 0. _d 0
541     hs(i,j,bi,bj) = 0. _d 0
542     hl(i,j,bi,bj) = 0. _d 0
543     endif
544    
545 dimitri 1.2.6.3 #else /* ifndef ALLOW_ATM_TEMP */
546 heimbach 1.1 #ifdef ALLOW_ATM_WIND
547     ustress(i,j,bi,bj) = atmrho*exf_BulkCdn(sh)*us*
548     & uwind(i,j,bi,bj)
549     vstress(i,j,bi,bj) = atmrho*exf_BulkCdn(sh)*us*
550     & vwind(i,j,bi,bj)
551 dimitri 1.2.6.3 #endif
552     #endif /* ifndef ALLOW_ATM_TEMP */
553 heimbach 1.1 enddo
554     enddo
555     enddo
556     enddo
557    
558     c Add all contributions.
559     do bj = mybylo(mythid),mybyhi(mythid)
560     do bi = mybxlo(mythid),mybxhi(mythid)
561     do j = 1,sny
562     do i = 1,snx
563     c Net surface heat flux.
564     #ifdef ALLOW_ATM_TEMP
565     hfl = 0. _d 0
566     hfl = hfl - hs(i,j,bi,bj)
567     hfl = hfl - hl(i,j,bi,bj)
568     hfl = hfl + lwflux(i,j,bi,bj)
569 dimitri 1.2.6.3 #ifndef SHORTWAVE_HEATING
570 heimbach 1.1 hfl = hfl + swflux(i,j,bi,bj)
571 dimitri 1.2.6.3 #endif
572 heimbach 1.1 c Heat flux:
573 dimitri 1.2.6.3 hflux(i,j,bi,bj) = hfl
574 heimbach 1.1 c Salt flux from Precipitation and Evaporation.
575 dimitri 1.2.6.3 sflux(i,j,bi,bj) = evap(i,j,bi,bj) - precip(i,j,bi,bj)
576 heimbach 1.1 #endif /* ALLOW_ATM_TEMP */
577    
578     #endif /* ALLOW_BULKFORMULAE */
579 heimbach 1.2.6.2
580     #ifdef ALLOW_RUNOFF
581 dimitri 1.2.6.3 sflux(i,j,bi,bj) = sflux(i,j,bi,bj) - runoff(i,j,bi,bj)
582     #endif
583    
584     hflux(i,j,bi,bj) = hflux(i,j,bi,bj)*maskc(i,j,1,bi,bj)
585     sflux(i,j,bi,bj) = sflux(i,j,bi,bj)*maskc(i,j,1,bi,bj)
586 heimbach 1.2.6.2
587 heimbach 1.1 enddo
588     enddo
589     enddo
590     enddo
591    
592     c Update the tile edges.
593     _EXCH_XY_R8(hflux, mythid)
594     _EXCH_XY_R8(sflux, mythid)
595     _EXCH_XY_R8(ustress, mythid)
596     _EXCH_XY_R8(vstress, mythid)
597    
598 dimitri 1.2.6.3 #ifdef SHORTWAVE_HEATING
599 heimbach 1.1 _EXCH_XY_R8(swflux, mythid)
600 dimitri 1.2.6.3 #endif
601 heimbach 1.1
602     end

  ViewVC Help
Powered by ViewVC 1.1.22