/[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.31 - (hide annotations) (download)
Mon Oct 20 03:13:32 2014 UTC (9 years, 7 months ago) by gforget
Branch: MAIN
CVS Tags: checkpoint66c, checkpoint66b, checkpoint66a, checkpoint65z, checkpoint65x, checkpoint65y, checkpoint65r, checkpoint65s, checkpoint65p, checkpoint65q, checkpoint65v, checkpoint65w, checkpoint65t, checkpoint65u, checkpoint65j, checkpoint65k, checkpoint65h, checkpoint65i, checkpoint65n, checkpoint65o, checkpoint65l, checkpoint65m, checkpoint65g
Changes since 1.30: +4 -1 lines
- ECCO_OPTIONS.h is needed when including ecco_cost.h, ecco.h
- AUTODIFF_OPTIONS.h is needed when including tamc.h, tamc_keys.h
- CTRL_OPTIONS.h is needed when including ctrl.h

1 gforget 1.31 C $Header: /u/gcmpack/MITgcm/pkg/exf/exf_bulkformulae.F,v 1.30 2014/05/02 22:14:52 jmc Exp $
2 jmc 1.13 C $Name: $
3 heimbach 1.2
4 edhill 1.4 #include "EXF_OPTIONS.h"
5 gforget 1.31 #ifdef ALLOW_AUTODIFF
6     # include "AUTODIFF_OPTIONS.h"
7     #endif
8 heimbach 1.2
9 mlosch 1.23 CBOP
10     C !ROUTINE: EXF_BULKFORMULAE
11     C !INTERFACE:
12     SUBROUTINE EXF_BULKFORMULAE( myTime, myIter, myThid )
13    
14     C !DESCRIPTION: \bv
15     C *==========================================================*
16     C | SUBROUTINE EXF_BULKFORMULAE
17     C | o Calculate bulk formula fluxes over open ocean
18     C | either following
19 jmc 1.25 C | Large and Pond, 1981 & 1982
20 mlosch 1.23 C | or (if defined ALLOW_BULK_LARGEYEAGER04)
21     C | Large and Yeager, 2004, NCAR/TN-460+STR.
22     C *==========================================================*
23     C \ev
24     C
25     C
26     C NOTES:
27     C ======
28     C
29     C See EXF_OPTIONS.h for a description of the various possible
30     C ocean-model forcing configurations.
31     C
32     C The bulk formulae of pkg/exf are not valid for sea-ice covered
33     C oceans but they can be used in combination with a sea-ice model,
34     C for example, pkg/seaice, to specify open water flux contributions.
35     C
36     C ==================================================================
37     C
38     C for Large and Pond, 1981 & 1982
39     C
40     C The calculation of the bulk surface fluxes has been adapted from
41     C the NCOM model which uses the formulae given in Large and Pond
42     C (1981 & 1982 )
43     C
44     C
45     C Header taken from NCOM version: ncom1.4.1
46     C -----------------------------------------
47     C
48     C Following procedures and coefficients in Large and Pond
49     C (1981 ; 1982)
50     C
51     C Output: Bulk estimates of the turbulent surface fluxes.
52     C -------
53     C
54     C hs - sensible heat flux (W/m^2), into ocean
55     C hl - latent heat flux (W/m^2), into ocean
56     C
57     C Input:
58     C ------
59     C
60     C us - mean wind speed (m/s) at height hu (m)
61     C th - mean air temperature (K) at height ht (m)
62     C qh - mean air humidity (kg/kg) at height hq (m)
63     C sst - sea surface temperature (K)
64     C tk0 - Kelvin temperature at 0 Celsius (K)
65     C
66     C Assume 1) a neutral 10m drag coefficient =
67     C
68     C cdn = .0027/u10 + .000142 + .0000764 u10
69     C
70     C 2) a neutral 10m stanton number =
71     C
72     C ctn = .0327 sqrt(cdn), unstable
73     C ctn = .0180 sqrt(cdn), stable
74     C
75     C 3) a neutral 10m dalton number =
76     C
77     C cen = .0346 sqrt(cdn)
78     C
79     C 4) the saturation humidity of air at
80     C
81     C t(k) = exf_BulkqSat(t) (kg/m^3)
82     C
83     C Note: 1) here, tstar = <wt>/u*, and qstar = <wq>/u*.
84     C 2) wind speeds should all be above a minimum speed,
85     C say 0.5 m/s.
86     C 3) with optional iteration loop, niter=3, should suffice.
87     C 4) this version is for analyses inputs with hu = 10m and
88     C ht = hq.
89     C 5) sst enters in Celsius.
90     C
91     C ==================================================================
92     C
93     C started: Christian Eckert eckert@mit.edu 27-Aug-1999
94     C
95     C changed: Christian Eckert eckert@mit.edu 14-Jan-2000
96     C - restructured the original version in order to have a
97     C better interface to the MITgcmUV.
98     C
99     C Christian Eckert eckert@mit.edu 12-Feb-2000
100     C - Changed Routine names (package prefix: exf_)
101     C
102     C Patrick Heimbach, heimbach@mit.edu 04-May-2000
103     C - changed the handling of precip and sflux with respect
104     C to CPP options ALLOW_BULKFORMULAE and ALLOW_ATM_TEMP
105     C - included some CPP flags ALLOW_BULKFORMULAE to make
106     C sure ALLOW_ATM_TEMP, ALLOW_ATM_WIND are used only in
107     C conjunction with defined ALLOW_BULKFORMULAE
108     C - statement functions discarded
109     C
110     C Ralf.Giering@FastOpt.de 25-Mai-2000
111     C - total rewrite using new subroutines
112     C
113     C Detlef Stammer: include river run-off. Nov. 21, 2001
114     C
115     C heimbach@mit.edu, 10-Jan-2002
116     C - changes to enable field swapping
117     C
118     C mods for pkg/seaice: menemenlis@jpl.nasa.gov 20-Dec-2002
119     C
120     C martin.losch@awi.de: merged with exf_bulk_largeyeager04, 21-May-2010
121     C
122     C ==================================================================
123     C
124     C for Large and Yeager, 2004
125     C
126     C === Turbulent Fluxes :
127     C * use the approach "B": shift coeff to height & stability of the
128     C atmosphere state (instead of "C": shift temp & humid to the height
129     C of wind, then shift the coeff to this height & stability of the atmos).
130     C * similar to EXF (except over sea-ice) ; default parameter values
131     C taken from Large & Yeager.
132     C * assume that Qair & Tair inputs are from the same height (zq=zt)
133     C * formulae in short:
134     C wind stress = (ust,vst) = rhoA * Cd * Ws * (del.u,del.v)
135     C Sensib Heat flux = fsha = rhoA * Ch * Ws * del.T * CpAir
136     C Latent Heat flux = flha = rhoA * Ce * Ws * del.Q * Lvap
137     C = -Evap * Lvap
138     C with Ws = wind speed = sqrt(del.u^2 +del.v^2) ;
139     C del.T = Tair - Tsurf ; del.Q = Qair - Qsurf ;
140     C Cd,Ch,Ce = drag coefficient, Stanton number and Dalton number
141     C respectively [no-units], function of height & stability
142    
143     C !USES:
144     IMPLICIT NONE
145     C === Global variables ===
146 heimbach 1.2 #include "EEPARAMS.h"
147     #include "SIZE.h"
148     #include "PARAMS.h"
149     #include "DYNVARS.h"
150 mlosch 1.23 #include "GRID.h"
151 heimbach 1.2
152 jmc 1.15 #include "EXF_PARAM.h"
153     #include "EXF_FIELDS.h"
154     #include "EXF_CONSTANTS.h"
155 heimbach 1.2
156     #ifdef ALLOW_AUTODIFF_TAMC
157     #include "tamc.h"
158     #endif
159    
160 mlosch 1.23 C !INPUT/OUTPUT PARAMETERS:
161     C input:
162     C myTime :: Current time in simulation
163     C myIter :: Current iteration number in simulation
164     C myThid :: My Thread Id number
165     _RL myTime
166     INTEGER myIter
167     INTEGER myThid
168     C output:
169     CEOP
170 heimbach 1.2
171 heimbach 1.3 #ifdef ALLOW_BULKFORMULAE
172 mlosch 1.23 C == external Functions
173 heimbach 1.3
174 mlosch 1.23 C == Local variables ==
175     C i,j :: grid point indices
176     C bi,bj :: tile indices
177     C ssq :: saturation specific humidity [kg/kg] in eq. with Sea-Surface water
178     INTEGER i,j,bi,bj
179 heimbach 1.2
180 heimbach 1.16 _RL czol
181 mlosch 1.23 _RL Tsf ! surface temperature [K]
182     _RL wsm ! limited wind speed [m/s] (> umin)
183     _RL t0 ! virtual temperature [K]
184     C these need to be 2D-arrays for vectorizing code
185     _RL tstar (1:sNx,1:sNy) ! turbulent temperature scale [K]
186     _RL qstar (1:sNx,1:sNy) ! turbulent humidity scale [kg/kg]
187     _RL ustar (1:sNx,1:sNy) ! friction velocity [m/s]
188     _RL tau (1:sNx,1:sNy) ! surface stress coef = rhoA * Ws * sqrt(Cd)
189     _RL rdn (1:sNx,1:sNy) ! neutral, zref (=10m) values of rd
190     _RL rd (1:sNx,1:sNy) ! = sqrt(Cd) [-]
191     _RL delq (1:sNx,1:sNy) ! specific humidity difference [kg/kg]
192 mlosch 1.19 _RL deltap(1:sNx,1:sNy)
193 jmc 1.30 #ifdef EXF_CALC_ATMRHO
194     _RL atmrho_loc(1:sNx,1:sNy) ! local atmospheric density [kg/m^3]
195     #endif
196    
197 jmc 1.25 #ifdef ALLOW_BULK_LARGEYEAGER04
198 mlosch 1.23 _RL dzTmp
199 jmc 1.25 #endif
200 mlosch 1.23 _RL SSTtmp
201     _RL ssq
202     _RL re ! = Ce / sqrt(Cd) [-]
203     _RL rh ! = Ch / sqrt(Cd) [-]
204     _RL ren, rhn ! neutral, zref (=10m) values of re, rh
205     _RL usn, usm
206     _RL stable ! = 1 if stable ; = 0 if unstable
207     _RL huol ! stability parameter at zwd [-] (=z/Monin-Obuklov length)
208     _RL htol ! stability parameter at zth [-]
209     _RL hqol
210     _RL x ! stability function [-]
211     _RL xsq ! = x^2 [-]
212     _RL psimh ! momentum stability function
213     _RL psixh ! latent & sensib. stability function
214     _RL zwln ! = log(zwd/zref)
215     _RL ztln ! = log(zth/zref)
216     _RL tmpbulk
217     _RL recip_rhoConstFresh
218     INTEGER ksrf, ksrfp1
219     INTEGER iter
220     C solve4Stress :: if F, by-pass momentum turb. flux (wind-stress) calculation
221     LOGICAL solve4Stress
222 jmc 1.25 _RL windStress ! surface wind-stress (@ grid cell center)
223 mlosch 1.23
224 heimbach 1.2 #ifdef ALLOW_AUTODIFF_TAMC
225     integer ikey_1
226     integer ikey_2
227     #endif
228    
229 jmc 1.25 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
230    
231 gforget 1.26 if ( useAtmWind ) then
232     solve4Stress = .TRUE.
233     else
234 mlosch 1.23 #ifdef ALLOW_BULK_LARGEYEAGER04
235 gforget 1.26 solve4Stress = wspeedfile .NE. ' '
236 mlosch 1.23 #else
237 gforget 1.26 solve4Stress = .FALSE.
238 mlosch 1.23 #endif
239 gforget 1.26 endif
240 heimbach 1.2
241 mlosch 1.23 C-- Set surface parameters :
242 mlosch 1.22 zwln = LOG(hu/zref)
243     ztln = LOG(ht/zref)
244 heimbach 1.16 czol = hu*karman*gravity_mks
245 mlosch 1.23 C-- abbreviation
246 mlosch 1.22 recip_rhoConstFresh = 1. _d 0/rhoConstFresh
247 jmc 1.25
248 jmc 1.30 C Loop over tiles.
249 heimbach 1.2 #ifdef ALLOW_AUTODIFF_TAMC
250     C-- HPF directive to help TAMC
251     CHPF$ INDEPENDENT
252     #endif
253 mlosch 1.23 DO bj = myByLo(myThid),myByHi(myThid)
254 heimbach 1.2 #ifdef ALLOW_AUTODIFF_TAMC
255     C-- HPF directive to help TAMC
256     CHPF$ INDEPENDENT
257     #endif
258 mlosch 1.23 DO bi = myBxLo(myThid),myBxHi(myThid)
259     ksrf = 1
260     ksrfp1 = 2
261     DO j = 1,sNy
262     DO i = 1,sNx
263 heimbach 1.2
264     #ifdef ALLOW_AUTODIFF_TAMC
265 heimbach 1.12 act1 = bi - myBxLo(myThid)
266     max1 = myBxHi(myThid) - myBxLo(myThid) + 1
267     act2 = bj - myByLo(myThid)
268     max2 = myByHi(myThid) - myByLo(myThid) + 1
269     act3 = myThid - 1
270     max3 = nTx*nTy
271     act4 = ikey_dynamics - 1
272    
273     ikey_1 = i
274     & + sNx*(j-1)
275     & + sNx*sNy*act1
276     & + sNx*sNy*max1*act2
277     & + sNx*sNy*max1*max2*act3
278     & + sNx*sNy*max1*max2*max3*act4
279 heimbach 1.2 #endif
280    
281     #ifdef ALLOW_ATM_TEMP
282    
283 heimbach 1.20 #ifdef ALLOW_AUTODIFF_TAMC
284     deltap(i,j) = 0. _d 0
285     delq(i,j) = 0. _d 0
286     #endif
287 mlosch 1.23 C--- Compute turbulent surface fluxes
288     C- Pot. Temp and saturated specific humidity
289     IF ( atemp(i,j,bi,bj) .NE. 0. _d 0 ) THEN
290     C- Surface Temp.
291     Tsf = theta(i,j,ksrf,bi,bj) + cen2kel
292 jmc 1.29 IF ( Nr.GE.2 .AND. sstExtrapol.GT.0. _d 0 ) THEN
293 mlosch 1.23 SSTtmp = sstExtrapol
294     & *( theta(i,j,ksrf,bi,bj)-theta(i,j,ksrfp1,bi,bj) )
295     & * maskC(i,j,ksrfp1,bi,bj)
296     Tsf = Tsf + MAX( SSTtmp, 0. _d 0 )
297     ENDIF
298 jmc 1.30
299 mlosch 1.23 tmpbulk = cvapor_fac*exp(-cvapor_exp/Tsf)
300 jmc 1.30 #ifdef EXF_CALC_ATMRHO
301     atmrho_loc(i,j) = apressure(i,j,bi,bj) /
302     & (287.04 _d 0*atemp(i,j,bi,bj)
303     & *(1. _d 0 + humid_fac*aqh(i,j,bi,bj)))
304     ssq = saltsat*tmpbulk/atmrho_loc(i,j)
305     #else
306 mlosch 1.23 ssq = saltsat*tmpbulk/atmrho
307 jmc 1.30 #endif
308 mlosch 1.23 deltap(i,j) = atemp(i,j,bi,bj) + gamma_blk*ht - Tsf
309 mlosch 1.19 delq(i,j) = aqh(i,j,bi,bj) - ssq
310 mlosch 1.23 C-- no negative evaporation over ocean:
311     IF ( noNegativeEvap ) delq(i,j) = MIN( 0. _d 0, delq(i,j) )
312    
313     C-- initial guess for exchange coefficients:
314     C take U_N = del.U ; stability from del.Theta ;
315     stable = exf_half + SIGN(exf_half, deltap(i,j))
316    
317     IF ( solve4Stress ) THEN
318 heimbach 1.2 #ifdef ALLOW_AUTODIFF_TAMC
319 heimbach 1.12 CADJ STORE sh(i,j,bi,bj) = comlev1_exf_1, key = ikey_1
320 mlosch 1.23 #endif /* ALLOW_AUTODIFF_TAMC */
321     C-- Wind speed
322     wsm = sh(i,j,bi,bj)
323     tmpbulk = exf_scal_BulkCdn *
324     & ( cdrag_1/wsm + cdrag_2 + cdrag_3*wsm )
325     rdn(i,j) = SQRT(tmpbulk)
326     ustar(i,j) = rdn(i,j)*wsm
327     ELSE
328 gforget 1.24 rdn(i,j) = 0. _d 0
329 mlosch 1.23 windStress = wStress(i,j,bi,bj)
330 jmc 1.30 #ifdef EXF_CALC_ATMRHO
331     ustar(i,j) = sqrt(windStress/atmrho_loc(i,j))
332     tau(i,j) = sqrt(windStress*atmrho_loc(i,j))
333     #else
334 mlosch 1.23 ustar(i,j) = sqrt(windStress/atmrho)
335     c tau(i,j) = windStress/ustar(i,j)
336     tau(i,j) = sqrt(windStress*atmrho)
337 jmc 1.30 #endif
338 mlosch 1.23 ENDIF
339 mlosch 1.21
340     C-- initial guess for exchange other coefficients:
341 mlosch 1.22 rhn = (exf_one-stable)*cstanton_1 + stable*cstanton_2
342 mlosch 1.23 ren = cDalton
343 mlosch 1.21 C-- calculate turbulent scales
344 mlosch 1.23 tstar(i,j)=rhn*deltap(i,j)
345     qstar(i,j)=ren*delq(i,j)
346 mlosch 1.19
347 mlosch 1.23 ELSE
348 jmc 1.25 C atemp(i,j,bi,bj) .EQ. 0.
349 mlosch 1.19 tstar (i,j) = 0. _d 0
350     qstar (i,j) = 0. _d 0
351     ustar (i,j) = 0. _d 0
352     tau (i,j) = 0. _d 0
353     rdn (i,j) = 0. _d 0
354 mlosch 1.23 ENDIF
355     ENDDO
356     ENDDO
357     DO iter=1,niter_bulk
358     DO j = 1,sNy
359     DO i = 1,sNx
360     IF ( atemp(i,j,bi,bj) .NE. 0. _d 0 ) THEN
361     C--- iterate with psi-functions to find transfer coefficients
362    
363 heimbach 1.2 #ifdef ALLOW_AUTODIFF_TAMC
364 heimbach 1.20 ikey_2 = i
365     & + sNx*(j-1)
366     & + sNx*sNy*(iter-1)
367     & + sNx*sNy*niter_bulk*act1
368     & + sNx*sNy*niter_bulk*max1*act2
369     & + sNx*sNy*niter_bulk*max1*max2*act3
370     & + sNx*sNy*niter_bulk*max1*max2*max3*act4
371 mlosch 1.23 CADJ STORE rdn (i,j) = comlev1_exf_2, key = ikey_2
372     CADJ STORE ustar (i,j) = comlev1_exf_2, key = ikey_2
373     CADJ STORE qstar (i,j) = comlev1_exf_2, key = ikey_2
374     CADJ STORE tstar (i,j) = comlev1_exf_2, key = ikey_2
375     CADJ STORE sh (i,j,bi,bj) = comlev1_exf_2, key = ikey_2
376     CADJ STORE wspeed(i,j,bi,bj) = comlev1_exf_2, key = ikey_2
377 heimbach 1.2 #endif
378    
379 mlosch 1.23 t0 = atemp(i,j,bi,bj)*
380 mlosch 1.19 & (exf_one + humid_fac*aqh(i,j,bi,bj))
381 mlosch 1.23 huol = ( tstar(i,j)/t0 +
382     & qstar(i,j)/(exf_one/humid_fac+aqh(i,j,bi,bj))
383     & )*czol/(ustar(i,j)*ustar(i,j))
384     #ifdef ALLOW_BULK_LARGEYEAGER04
385     tmpbulk = MIN(abs(huol),10. _d 0)
386     huol = SIGN(tmpbulk , huol)
387     #else
388     C-- Large&Pond1981:
389 mlosch 1.19 huol = max(huol,zolmin)
390 mlosch 1.23 #endif /* ALLOW_BULK_LARGEYEAGER04 */
391 mlosch 1.19 htol = huol*ht/hu
392     hqol = huol*hq/hu
393 mlosch 1.22 stable = exf_half + sign(exf_half, huol)
394 mlosch 1.23
395 jmc 1.30 C Evaluate all stability functions assuming hq = ht.
396 mlosch 1.23 IF ( solve4Stress ) THEN
397     #ifdef ALLOW_BULK_LARGEYEAGER04
398     C-- Large&Yeager04:
399     xsq = SQRT( ABS(exf_one - huol*16. _d 0) )
400     #else
401     C-- Large&Pond1981:
402     xsq = max(sqrt(abs(exf_one - 16.*huol)),exf_one)
403     #endif /* ALLOW_BULK_LARGEYEAGER04 */
404     x = SQRT(xsq)
405     psimh = -psim_fac*huol*stable
406     & +(exf_one-stable)*
407     & ( LOG( (exf_one + exf_two*x + xsq)*
408     & (exf_one+xsq)*.125 _d 0 )
409     & -exf_two*ATAN(x) + exf_half*pi )
410 gforget 1.24 ELSE
411     psimh = 0. _d 0
412 mlosch 1.23 ENDIF
413     #ifdef ALLOW_BULK_LARGEYEAGER04
414     C-- Large&Yeager04:
415     xsq = SQRT( ABS(exf_one - htol*16. _d 0) )
416     #else
417     C-- Large&Pond1981:
418     xsq = max(sqrt(abs(exf_one - 16.*htol)),exf_one)
419     #endif /* ALLOW_BULK_LARGEYEAGER04 */
420 mlosch 1.22 psixh = -psim_fac*htol*stable + (exf_one-stable)*
421     & ( exf_two*LOG(exf_half*(exf_one+xsq)) )
422 mlosch 1.23
423     IF ( solve4Stress ) THEN
424 heimbach 1.28 CADJ STORE rdn(i,j) = comlev1_exf_2, key = ikey_2
425 mlosch 1.23 C- Shift wind speed using old coefficient
426     #ifdef ALLOW_BULK_LARGEYEAGER04
427     C-- Large&Yeager04:
428     dzTmp = (zwln-psimh)/karman
429     usn = wspeed(i,j,bi,bj)/(exf_one + rdn(i,j)*dzTmp )
430     #else
431     C-- Large&Pond1981:
432 mlosch 1.22 c rd(i,j)= rdn(i,j)/(exf_one - rdn(i,j)/karman*psimh )
433 mlosch 1.23 c usn = sh(i,j,bi,bj)*rd(i,j)/rdn(i,j)
434 jmc 1.30 C ML: the original formulation above is replaced below to be
435     C similar to largeyeager04, but does not give the same results, strange
436 mlosch 1.23 usn = sh(i,j,bi,bj)/(exf_one - rdn(i,j)/karman*psimh)
437     #endif /* ALLOW_BULK_LARGEYEAGER04 */
438     usm = MAX(usn, umin)
439    
440     C- Update the 10m, neutral stability transfer coefficients (momentum)
441     tmpbulk = exf_scal_BulkCdn *
442     & ( cdrag_1/usm + cdrag_2 + cdrag_3*usm )
443     rdn(i,j) = SQRT(tmpbulk)
444     #ifdef ALLOW_BULK_LARGEYEAGER04
445     rd(i,j) = rdn(i,j)/(exf_one + rdn(i,j)*dzTmp)
446     #else
447     rd(i,j) = rdn(i,j)/(exf_one - rdn(i,j)/karman*psimh)
448     #endif /* ALLOW_BULK_LARGEYEAGER04 */
449     ustar(i,j) = rd(i,j)*sh(i,j,bi,bj)
450     C- Coeff:
451 jmc 1.30 #ifdef EXF_CALC_ATMRHO
452     tau(i,j) = atmrho_loc(i,j)*rd(i,j)*wspeed(i,j,bi,bj)
453     #else
454 mlosch 1.23 tau(i,j) = atmrho*rd(i,j)*wspeed(i,j,bi,bj)
455 jmc 1.30 #endif
456 mlosch 1.23 ENDIF
457    
458     C- Update the 10m, neutral stability transfer coefficients (sens&evap)
459 mlosch 1.22 rhn = (exf_one-stable)*cstanton_1 + stable*cstanton_2
460 mlosch 1.23 ren = cDalton
461    
462     C- Shift all coefficients to the measurement height and stability.
463     rh = rhn/(exf_one + rhn*(ztln-psixh)/karman)
464     re = ren/(exf_one + ren*(ztln-psixh)/karman)
465 heimbach 1.2
466 mlosch 1.23 C-- Update ustar, tstar, qstar using updated, shifted coefficients.
467 mlosch 1.19 qstar(i,j) = re*delq(i,j)
468     tstar(i,j) = rh*deltap(i,j)
469 heimbach 1.2
470 mlosch 1.23 ENDIF
471     ENDDO
472     ENDDO
473 jmc 1.30 C end of iteration loop
474 mlosch 1.23 ENDDO
475     DO j = 1,sNy
476     DO i = 1,sNx
477     IF ( atemp(i,j,bi,bj) .NE. 0. _d 0 ) THEN
478 heimbach 1.2
479     #ifdef ALLOW_AUTODIFF_TAMC
480 mlosch 1.23 ikey_1 = i
481     & + sNx*(j-1)
482     & + sNx*sNy*act1
483     & + sNx*sNy*max1*act2
484     & + sNx*sNy*max1*max2*act3
485     & + sNx*sNy*max1*max2*max3*act4
486 mlosch 1.21 CADJ STORE qstar(i,j) = comlev1_exf_1, key = ikey_1
487     CADJ STORE tstar(i,j) = comlev1_exf_1, key = ikey_1
488 mlosch 1.23 CADJ STORE tau (i,j) = comlev1_exf_1, key = ikey_1
489     CADJ STORE rd (i,j) = comlev1_exf_1, key = ikey_1
490     #endif /* ALLOW_AUTODIFF_TAMC */
491    
492     C- Turbulent Fluxes
493 mlosch 1.22 hs(i,j,bi,bj) = atmcp*tau(i,j)*tstar(i,j)
494     hl(i,j,bi,bj) = flamb*tau(i,j)*qstar(i,j)
495 heimbach 1.2 #ifndef EXF_READ_EVAP
496 mlosch 1.23 C change sign and convert from kg/m^2/s to m/s via rhoConstFresh
497     c evap(i,j,bi,bj) = -recip_rhonil*tau(i,j)*qstar(i,j)
498 mlosch 1.22 evap(i,j,bi,bj) = -recip_rhoConstFresh*tau(i,j)*qstar(i,j)
499 mlosch 1.23 C but older version was using rhonil instead:
500     c evap(i,j,bi,bj) = -recip_rhonil*tau(i,j)*qstar(i,j)
501 heimbach 1.2 #endif
502 gforget 1.26 if ( useAtmWind ) then
503 mlosch 1.23 c ustress(i,j,bi,bj) = tau(i,j)*rd(i,j)*ws*cw(i,j,bi,bj)
504     c vstress(i,j,bi,bj) = tau(i,j)*rd(i,j)*ws*sw(i,j,bi,bj)
505     C- jmc: below is how it should be written ; different from above when
506     C both wind-speed & 2 compon. of the wind are specified, and in
507     C- this case, formula below is better.
508 mlosch 1.22 tmpbulk = tau(i,j)*rd(i,j)
509     ustress(i,j,bi,bj) = tmpbulk*uwind(i,j,bi,bj)
510     vstress(i,j,bi,bj) = tmpbulk*vwind(i,j,bi,bj)
511 gforget 1.26 endif
512 mlosch 1.21
513     c IF ( ATEMP(i,j,bi,bj) .NE. 0. )
514 mlosch 1.19 else
515 gforget 1.26 if ( useAtmWind ) ustress(i,j,bi,bj) = 0. _d 0
516     if ( useAtmWind ) vstress(i,j,bi,bj) = 0. _d 0
517 mlosch 1.19 hflux (i,j,bi,bj) = 0. _d 0
518 mlosch 1.23 evap (i,j,bi,bj) = 0. _d 0
519     hs (i,j,bi,bj) = 0. _d 0
520     hl (i,j,bi,bj) = 0. _d 0
521 mlosch 1.21
522     c IF ( ATEMP(i,j,bi,bj) .NE. 0. )
523 mlosch 1.19 endif
524 heimbach 1.2
525     #else /* ifndef ALLOW_ATM_TEMP */
526 gforget 1.26 if ( useAtmWind ) then
527 mlosch 1.23 wsm = sh(i,j,bi,bj)
528     tmpbulk = exf_scal_BulkCdn *
529 mlosch 1.21 & ( cdrag_1/wsm + cdrag_2 + cdrag_3*wsm )
530 mlosch 1.19 ustress(i,j,bi,bj) = atmrho*tmpbulk*wspeed(i,j,bi,bj)*
531 mlosch 1.23 & uwind(i,j,bi,bj)
532 mlosch 1.19 vstress(i,j,bi,bj) = atmrho*tmpbulk*wspeed(i,j,bi,bj)*
533 mlosch 1.23 & vwind(i,j,bi,bj)
534 gforget 1.26 endif
535 heimbach 1.2 #endif /* ifndef ALLOW_ATM_TEMP */
536 mlosch 1.23 ENDDO
537     ENDDO
538     ENDDO
539     ENDDO
540 heimbach 1.2
541 heimbach 1.3 #endif /* ALLOW_BULKFORMULAE */
542 heimbach 1.2
543 mlosch 1.23 RETURN
544     END

  ViewVC Help
Powered by ViewVC 1.1.22