/[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.8 - (hide annotations) (download)
Tue Feb 18 05:33:54 2003 UTC (21 years, 4 months ago) by dimitri
Branch: MAIN
CVS Tags: checkpoint48f_post, checkpoint48g_post
Changes since 1.7: +103 -125 lines
Merging from release1_p12:
o Modifications for using pkg/exf with pkg/seaice
  - improved description of the various forcing configurations
  - added basic radiation bulk formulae to pkg/exf
  - units/sign fix for evap computation in exf_getffields.F
  - updated verification/global_with_exf/results/output.txt
o Added pkg/sbo for computing IERS Special Bureau for the Oceans
  (SBO) core products, including oceanic mass, center-of-mass,
  angular, and bottom pressure (see pkg/sbo/README.sbo).
o Lower bound for viscosity/diffusivity in pkg/kpp/kpp_routines.F
  to avoid negative values in shallow regions.
  - updated verification/natl_box/results/output.txt
  - updated verification/lab_sea/results/output.txt
o MPI gather, scatter: eesupp/src/gather_2d.F and scatter_2d.F
o Added useSingleCpuIO option (see PARAMS.h).
o Updated useSingleCpuIO option in mdsio_writefield.F to
  work with multi-field files, e.g., for single-file pickup.
o pkg/seaice:
  - bug fix in growth.F: QNET for no shortwave case
  - added HeffFile for specifying initial sea-ice thickness
  - changed SEAICE_EXTERNAL_FLUXES wind stress implementation
o Added missing /* */ to CPP comments in pkg/seaice, pkg/exf,
  kpp_transport_t.F, forward_step.F, and the_main_loop.F
o pkg/seaice:
  - adjoint-friendly modifications
  - added a SEAICE_WRITE_PICKUP at end of the_model_main.F

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

  ViewVC Help
Powered by ViewVC 1.1.22