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

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

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


Revision 1.3 - (hide annotations) (download)
Tue Jun 24 16:07:32 2003 UTC (21 years ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint51e_post, checkpoint51j_post, checkpoint51a_post, checkpoint51c_post, checkpoint51f_pre, checkpoint51, checkpoint51f_post, checkpoint51b_post, checkpoint51b_pre, branchpoint-genmake2, checkpoint51h_pre, checkpoint51g_post, checkpoint51d_post, checkpoint51i_pre
Branch point for: branch-genmake2
Changes since 1.2: +5 -13 lines
Merging for c51 vs. e34

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

  ViewVC Help
Powered by ViewVC 1.1.22