/[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.7 - (hide annotations) (download)
Thu May 5 16:41:54 2005 UTC (19 years, 1 month ago) by dimitri
Branch: MAIN
Changes since 1.6: +10 -1 lines
added pkg/exf, SALTanom, and SALTSQan diagnostics

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

  ViewVC Help
Powered by ViewVC 1.1.22