/[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.6 - (hide annotations) (download)
Wed Dec 18 19:23:29 2002 UTC (21 years, 6 months ago) by cheisey
Branch: MAIN
CVS Tags: checkpoint47e_post
Branch point for: branch-exfmods-curt
Changes since 1.5: +9 -1 lines
Some TMP variables should be declared.

1 cheisey 1.6 c $Header: /u/u0/gcmpack/MITgcm/pkg/exf/exf_getffields.F,v 1.5 2002/12/17 19:47:41 cheisey 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     c o Get the surface fluxes either from file or as derived from bulk
12     c formulae that use the atmospheric state.
13     c
14     c The calculation of the bulk surface fluxes has been adapted from
15     c the NCOM model which uses the formulae given in Large and Pond
16     c (1981 & 1982 )
17     c
18     c
19     c Header taken from NCOM version: ncom1.4.1
20     c -----------------------------------------
21     c
22     c Following procedures and coefficients in Large and Pond
23     c (1981 ; 1982)
24     c
25     c Output: Bulk estimates of the turbulent surface fluxes.
26     c -------
27     c
28     c hs - sensible heat flux (W/m^2), into ocean
29     c hl - latent heat flux (W/m^2), into ocean
30     c
31     c Input:
32     c ------
33     c
34     c us - mean wind speed (m/s) at height hu (m)
35     c th - mean air temperature (K) at height ht (m)
36     c qh - mean air humidity (kg/kg) at height hq (m)
37     c sst - sea surface temperature (K)
38     c tk0 - Kelvin temperature at 0 Celsius (K)
39     c
40     c Assume 1) a neutral 10m drag coefficient =
41     c
42     c cdn = .0027/u10 + .000142 + .0000764 u10
43     c
44     c 2) a neutral 10m stanton number =
45     c
46     c ctn = .0327 sqrt(cdn), unstable
47     c ctn = .0180 sqrt(cdn), stable
48     c
49     c 3) a neutral 10m dalton number =
50     c
51     c cen = .0346 sqrt(cdn)
52     c
53     c 4) the saturation humidity of air at
54     c
55     c t(k) = exf_BulkqSat(t) (kg/m^3)
56     c
57     c Note: 1) here, tstar = <wt>/u*, and qstar = <wq>/u*.
58     c 2) wind speeds should all be above a minimum speed,
59     c say 0.5 m/s.
60     c 3) with optional iteration loop, niter=3, should suffice.
61     c 4) this version is for analyses inputs with hu = 10m and
62     c ht = hq.
63     c 5) sst enters in Celsius.
64     c
65     c ------------------------------------
66     c
67     c By setting CPP options in the header file *EXF_CPPOPTIONS.h* it
68     c is possible to combine data sets in four different ways:
69     c
70     c The following options are available:
71     c
72     c ALLOW_ATM_TEMP (UAT)
73     c ALLOW_ATM_WIND (UAW)
74     c
75     c
76     c UAT | UAW | action
77     c ----------------------------------------------------
78     c undefined | undefined | Use surface fluxes.
79     c undefined | defined | Assume cdn(u) given to
80     c | | infer the wind stress.
81     c defined | undefined | Compute wind field from
82     c | | given stress assuming a
83     c | | linear relation.
84     c defined | defined | Use the bulk formulae.
85     c ----------------------------------------------------
86     c
87     c Implementations of the bulk formulae exist for the follwing
88     c versions of the MITgcm:
89     c
90     c MITgcm : Patrick Heimbach
91     c MITgcmUV: Christian Eckert
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 heimbach 1.4 c - total rewrite using new subroutines
112 heimbach 1.1 c
113 heimbach 1.4 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 heimbach 1.1 c
118     c ==================================================================
119     c SUBROUTINE exf_GetFFields
120     c ==================================================================
121    
122     implicit none
123    
124     c == global variables ==
125    
126     #include "EEPARAMS.h"
127     #include "SIZE.h"
128     #include "PARAMS.h"
129     #include "DYNVARS.h"
130     #include "GRID.h"
131    
132     #include "exf_fields.h"
133     #include "exf_constants.h"
134    
135 heimbach 1.3 #ifdef ALLOW_AUTODIFF_TAMC
136     #include "tamc.h"
137     #endif
138    
139 heimbach 1.1 c == routine arguments ==
140    
141     integer mythid
142     integer mycurrentiter
143     _RL mycurrenttime
144    
145     c == local variables ==
146    
147     integer bi,bj
148     integer i,j,k
149    
150     #ifdef ALLOW_BULKFORMULAE
151    
152     #ifdef ALLOW_ATM_TEMP
153     integer iter
154     _RL aln
155     _RL delq
156     _RL deltap
157     _RL hqol
158     _RL htol
159     _RL huol
160     _RL psimh
161     _RL psixh
162     _RL qstar
163     _RL rd
164     _RL re
165     _RL rdn
166     _RL rh
167 heimbach 1.3 _RL ssttmp
168 heimbach 1.1 _RL ssq
169     _RL stable
170     _RL tstar
171     _RL t0
172     _RL ustar
173     _RL uzn
174 heimbach 1.3 _RL shn
175 heimbach 1.1 _RL xsq
176     _RL x
177 heimbach 1.3 _RL tau
178 heimbach 1.1 _RL evap(1-olx:snx+olx,1-oly:sny+oly,nsx,nsy)
179 heimbach 1.3 #ifdef ALLOW_AUTODIFF_TAMC
180     integer ikey_1
181     integer ikey_2
182     #endif
183 heimbach 1.1 #endif /* ALLOW_ATM_TEMP */
184    
185 heimbach 1.3 _RL ustmp
186 heimbach 1.1 _RL us
187     _RL cw
188     _RL sw
189     _RL sh
190     _RL hs(1-olx:snx+olx,1-oly:sny+oly,nsx,nsy)
191     _RL hl(1-olx:snx+olx,1-oly:sny+oly,nsx,nsy)
192     _RL hfl
193    
194     #endif /* ALLOW_BULKFORMULAE */
195    
196     c == external functions ==
197    
198     integer ilnblnk
199     external ilnblnk
200    
201     #ifdef ALLOW_BULKFORMULAE
202     _RL exf_BulkqSat
203     external exf_BulkqSat
204     _RL exf_BulkCdn
205     external exf_BulkCdn
206     _RL exf_BulkRhn
207     external exf_BulkRhn
208     #endif /* ALLOW_BULKFORMULAE */
209 cheisey 1.6
210     #ifndef ALLOW_ATM_WIND
211     _RL TMP1
212     _RL TMP2
213     _RL TMP3
214     _RL TMP4
215     _RL TMP5
216     #endif
217 heimbach 1.1
218     c == end of interface ==
219    
220     c determine forcing field records
221    
222     #ifdef ALLOW_BULKFORMULAE
223    
224 heimbach 1.3 #if (defined (ALLOW_ATM_TEMP) || defined (ALLOW_ATM_WIND))
225 heimbach 1.4 cph This statement cannot be a PARAMETER statement in the header,
226     cph but must come here; it's not fortran77 standard
227 heimbach 1.3 aln = log(ht/zref)
228     #endif
229    
230 heimbach 1.1 c Determine where we are in time and set counters, flags and
231     c the linear interpolation factors accordingly.
232     #ifdef ALLOW_ATM_TEMP
233     c Atmospheric temperature.
234 heimbach 1.3 call exf_set_atemp( mycurrenttime, mycurrentiter, mythid )
235 heimbach 1.1
236     c Atmospheric humidity.
237 heimbach 1.3 call exf_set_aqh( mycurrenttime, mycurrentiter, mythid )
238 heimbach 1.1
239     c Net long wave radiative flux.
240 heimbach 1.3 call exf_set_lwflux( mycurrenttime, mycurrentiter, mythid )
241 heimbach 1.1
242     c Net short wave radiative flux.
243 heimbach 1.3 call exf_set_swflux( mycurrenttime, mycurrentiter, mythid )
244 heimbach 1.1
245     c Precipitation.
246 heimbach 1.3 call exf_set_precip( mycurrenttime, mycurrentiter, mythid )
247 heimbach 1.1
248 heimbach 1.3 #ifdef ALLOW_ATEMP_CONTROL
249     call ctrl_getatemp ( mycurrenttime, mycurrentiter, mythid )
250     #endif
251    
252     #ifdef ALLOW_AQH_CONTROL
253     call ctrl_getaqh ( mycurrenttime, mycurrentiter, mythid )
254     #endif
255 heimbach 1.1
256     #else
257    
258     c Atmospheric heat flux.
259 heimbach 1.3 call exf_set_hflux( mycurrenttime, mycurrentiter, mythid )
260 heimbach 1.1
261     c Salt flux.
262 heimbach 1.3 call exf_set_sflux( mycurrenttime, mycurrentiter, mythid )
263 heimbach 1.1
264     #ifdef ALLOW_KPP
265     c Net short wave radiative flux.
266 heimbach 1.3 call exf_set_swflux( mycurrenttime, mycurrentiter, mythid )
267 heimbach 1.1
268     #endif /* ALLOW_KPP */
269    
270     #endif /* ALLOW_ATM_TEMP */
271    
272     #ifdef ALLOW_ATM_WIND
273     c Zonal wind.
274 heimbach 1.3 call exf_set_uwind( mycurrenttime, mycurrentiter, mythid )
275 heimbach 1.1
276     c Meridional wind.
277 heimbach 1.3 call exf_set_vwind( mycurrenttime, mycurrentiter, mythid )
278    
279     #ifdef ALLOW_UWIND_CONTROL
280     call ctrl_getuwind ( mycurrenttime, mycurrentiter, mythid )
281     #endif
282    
283     #ifdef ALLOW_VWIND_CONTROL
284     call ctrl_getvwind ( mycurrenttime, mycurrentiter, mythid )
285     #endif
286 heimbach 1.1
287     #else
288     c Zonal wind stress.
289 heimbach 1.3 call exf_set_ustress( mycurrenttime, mycurrentiter, mythid )
290 heimbach 1.1
291     c Meridional wind stress.
292 heimbach 1.3 call exf_set_vstress( mycurrenttime, mycurrentiter, mythid )
293 heimbach 1.1
294     #endif /* ALLOW_ATM_WIND */
295    
296     #else /* ALLOW_BULKFORMULAE undefined */
297     c Atmospheric heat flux.
298 heimbach 1.3 call exf_set_hflux( mycurrenttime, mycurrentiter, mythid )
299 heimbach 1.1
300     c Salt flux.
301 heimbach 1.3 call exf_set_sflux( mycurrenttime, mycurrentiter, mythid )
302 heimbach 1.1
303     c Zonal wind stress.
304 heimbach 1.3 call exf_set_ustress( mycurrenttime, mycurrentiter, mythid )
305 heimbach 1.1
306     c Meridional wind stress.
307 heimbach 1.3 call exf_set_vstress( mycurrenttime, mycurrentiter, mythid )
308 heimbach 1.1
309     #ifdef ALLOW_KPP
310     c Net short wave radiative flux.
311 heimbach 1.3 call exf_set_swflux( mycurrenttime, mycurrentiter, mythid )
312     #endif
313 heimbach 1.1
314     #endif /* ALLOW_BULKFORMULAE */
315    
316     c Loop over tiles.
317 heimbach 1.3 #ifdef ALLOW_AUTODIFF_TAMC
318     C-- HPF directive to help TAMC
319     CHPF$ INDEPENDENT
320     #endif /* ALLOW_AUTODIFF_TAMC */
321 heimbach 1.1 do bj = mybylo(mythid),mybyhi(mythid)
322 heimbach 1.3 #ifdef ALLOW_AUTODIFF_TAMC
323     C-- HPF directive to help TAMC
324     CHPF$ INDEPENDENT
325     #endif
326 heimbach 1.1 do bi = mybxlo(mythid),mybxhi(mythid)
327 heimbach 1.3
328 heimbach 1.1 k = 1
329    
330     do j = 1-oly,sny+oly
331     do i = 1-olx,snx+olx
332    
333     #ifdef ALLOW_BULKFORMULAE
334    
335 heimbach 1.3 #ifdef ALLOW_AUTODIFF_TAMC
336     act1 = bi - myBxLo(myThid)
337     max1 = myBxHi(myThid) - myBxLo(myThid) + 1
338     act2 = bj - myByLo(myThid)
339     max2 = myByHi(myThid) - myByLo(myThid) + 1
340     act3 = myThid - 1
341     max3 = nTx*nTy
342     act4 = ikey_dynamics - 1
343    
344     ikey_1 = i
345     & + sNx*(j-1)
346     & + sNx*sNy*act1
347     & + sNx*sNy*max1*act2
348     & + sNx*sNy*max1*max2*act3
349     & + sNx*sNy*max1*max2*max3*act4
350     #endif
351    
352 heimbach 1.1 c Compute the turbulent surface fluxes.
353     c (Bulk formulae estimates)
354    
355     #ifdef ALLOW_ATM_WIND
356     c Wind speed and direction.
357 heimbach 1.3 ustmp = uwind(i,j,bi,bj)*uwind(i,j,bi,bj) +
358     & vwind(i,j,bi,bj)*vwind(i,j,bi,bj)
359     if ( ustmp .ne. 0. _d 0 ) then
360     us = sqrt(ustmp)
361 heimbach 1.1 cw = uwind(i,j,bi,bj)/us
362     sw = vwind(i,j,bi,bj)/us
363     else
364 heimbach 1.3 us = 0. _d 0
365 heimbach 1.1 cw = 0. _d 0
366     sw = 0. _d 0
367     endif
368     sh = max(us,umin)
369     #else
370     #ifdef ALLOW_ATM_TEMP
371    
372 heimbach 1.3 c The variables us, sh and rdn have to be computed from
373     c given wind stresses inverting relationship for neutral
374     c drag coeff. cdn.
375 heimbach 1.1 c The inversion is based on linear and quadratic form of
376     c cdn(umps); ustar can be directly computed from stress;
377    
378 heimbach 1.3 ustmp = ustress(i,j,bi,bj)*ustress(i,j,bi,bj) +
379     & vstress(i,j,bi,bj)*vstress(i,j,bi,bj)
380     if ( ustmp .ne. 0. _d 0 ) then
381     ustar = sqrt(ustmp/atmrho)
382     cw = ustress(i,j,bi,bj)/sqrt(ustmp)
383     sw = vstress(i,j,bi,bj)/sqrt(ustmp)
384     else
385     ustar = 0. _d 0
386     cw = 0. _d 0
387     sw = 0. _d 0
388     endif
389 heimbach 1.1
390     if ( ustar .eq. 0. _d 0 ) then
391     us = 0. _d 0
392     else if ( ustar .lt. ustofu11 ) then
393     tmp1 = -cquadrag_2/cquadrag_1/2
394     tmp2 = sqrt(tmp1*tmp1 + ustar*ustar/cquadrag_1)
395     us = sqrt(tmp1 + tmp2)
396     else
397     tmp3 = clindrag_2/clindrag_1/3
398     tmp4 = ustar*ustar/clindrag_1/2 - tmp3**3
399     tmp5 = sqrt(ustar*ustar/clindrag_1*
400     & (ustar*ustar/clindrag_1/4 - tmp3**3))
401     us = (tmp4 + tmp5)**(1/3) +
402     & tmp3**2 * (tmp4 + tmp5)**(-1/3) - tmp3
403     endif
404    
405     if ( us .ne. 0 ) then
406     rdn = ustar/us
407     else
408     rdn = 0. _d 0
409     end if
410    
411     sh = max(us,umin)
412     #endif /* ALLOW_ATM_TEMP */
413     #endif /* ALLOW_ATM_WIND */
414    
415     #ifdef ALLOW_ATM_TEMP
416 heimbach 1.3
417 heimbach 1.1 c Initial guess: z/l=0.0; hu=ht=hq=z
418     c Iterations: converge on z/l and hence the fluxes.
419     c t0 : virtual temperature (K)
420     c ssq : sea surface humidity (kg/kg)
421     c deltap : potential temperature diff (K)
422    
423     if ( atemp(i,j,bi,bj) .ne. 0. _d 0 ) then
424     t0 = atemp(i,j,bi,bj)*
425     & (exf_one + humid_fac*aqh(i,j,bi,bj))
426 heimbach 1.3 ssttmp = theta(i,j,k,bi,bj)
427 heimbach 1.1 ssq = saltsat*
428 heimbach 1.3 & exf_BulkqSat(ssttmp + cen2kel)/
429 heimbach 1.1 & atmrho
430     deltap = atemp(i,j,bi,bj) + gamma_blk*ht -
431 heimbach 1.3 & ssttmp - cen2kel
432 heimbach 1.1 delq = aqh(i,j,bi,bj) - ssq
433     stable = exf_half + sign(exf_half, deltap)
434 heimbach 1.3 #ifdef ALLOW_AUTODIFF_TAMC
435     CADJ STORE sh = comlev1_exf_1, key = ikey_1
436     #endif
437 heimbach 1.1 rdn = sqrt(exf_BulkCdn(sh))
438     ustar = rdn*sh
439     tstar = exf_BulkRhn(stable)*deltap
440     qstar = cdalton*delq
441    
442     do iter = 1,niter_bulk
443    
444 heimbach 1.3 #ifdef ALLOW_AUTODIFF_TAMC
445     ikey_2 = iter
446     & + niter_bulk*(i-1)
447     & + sNx*niter_bulk*(j-1)
448     & + sNx*niter_bulk*sNy*act1
449     & + sNx*niter_bulk*sNy*max1*act2
450     & + sNx*niter_bulk*sNy*max1*max2*act3
451     & + sNx*niter_bulk*sNy*max1*max2*max3*act4
452     #endif
453    
454     #ifdef ALLOW_AUTODIFF_TAMC
455     CADJ STORE rdn = comlev1_exf_2, key = ikey_2
456     CADJ STORE ustar = comlev1_exf_2, key = ikey_2
457     CADJ STORE qstar = comlev1_exf_2, key = ikey_2
458     CADJ STORE tstar = comlev1_exf_2, key = ikey_2
459     CADJ STORE sh = comlev1_exf_2, key = ikey_2
460     CADJ STORE us = comlev1_exf_2, key = ikey_2
461     #endif
462    
463 heimbach 1.1 huol = czol*(tstar/t0 +
464     & qstar/(exf_one/humid_fac+aqh(i,j,bi,bj)))/
465     & ustar**2
466     huol = max(huol,zolmin)
467     stable = exf_half + sign(exf_half, huol)
468     htol = huol*ht/hu
469     hqol = huol*hq/hu
470    
471     c Evaluate all stability functions assuming hq = ht.
472     xsq = max(sqrt(abs(exf_one - 16.*huol)),exf_one)
473     x = sqrt(xsq)
474     psimh = -psim_fac*huol*stable +
475     & (exf_one - stable)*
476     & log((exf_one + x*(exf_two + x))*
477     & (exf_one + xsq)/8.) - exf_two*atan(x) +
478     & pi*exf_half
479     xsq = max(sqrt(abs(exf_one - 16.*htol)),exf_one)
480     psixh = -psim_fac*htol*stable + (exf_one - stable)*
481     & exf_two*log((exf_one + xsq)/exf_two)
482    
483     c Shift wind speed using old coefficient
484 heimbach 1.3 ccc rd = rdn/(exf_one + rdn/karman*
485     ccc & (log(hu/zref) - psimh) )
486 heimbach 1.1 rd = rdn/(exf_one - rdn/karman*psimh )
487 heimbach 1.3 shn = sh*rd/rdn
488     uzn = max(shn, umin)
489 heimbach 1.1
490     c Update the transfer coefficients at 10 meters
491     c and neutral stability.
492 heimbach 1.3
493 heimbach 1.1 rdn = sqrt(exf_BulkCdn(uzn))
494    
495     c Shift all coefficients to the measurement height
496     c and stability.
497     c rd = rdn/(exf_one + rdn/karman*(log(hu/zref) - psimh))
498     rd = rdn/(exf_one - rdn/karman*psimh)
499     rh = exf_BulkRhn(stable)/(exf_one +
500     & exf_BulkRhn(stable)/
501     & karman*(aln - psixh))
502     re = cdalton/(exf_one + cdalton/karman*(aln - psixh))
503    
504     c Update ustar, tstar, qstar using updated, shifted
505     c coefficients.
506     ustar = rd*sh
507     qstar = re*delq
508     tstar = rh*deltap
509 heimbach 1.3 tau = atmrho*ustar**2
510     tau = tau*us/sh
511 heimbach 1.1
512     enddo
513    
514 heimbach 1.3 #ifdef ALLOW_AUTODIFF_TAMC
515     CADJ STORE ustar = comlev1_exf_1, key = ikey_1
516     CADJ STORE qstar = comlev1_exf_1, key = ikey_1
517     CADJ STORE tstar = comlev1_exf_1, key = ikey_1
518     CADJ STORE tau = comlev1_exf_1, key = ikey_1
519     CADJ STORE cw = comlev1_exf_1, key = ikey_1
520     CADJ STORE sw = comlev1_exf_1, key = ikey_1
521     #endif
522    
523 heimbach 1.1 hs(i,j,bi,bj) = atmcp*tau*tstar/ustar
524     hl(i,j,bi,bj) = flamb*tau*qstar/ustar
525     evap(i,j,bi,bj) = tau*qstar/ustar
526     ustress(i,j,bi,bj) = tau*cw
527     vstress(i,j,bi,bj) = tau*sw
528     else
529     ustress(i,j,bi,bj) = 0. _d 0
530     vstress(i,j,bi,bj) = 0. _d 0
531     hflux (i,j,bi,bj) = 0. _d 0
532     hs(i,j,bi,bj) = 0. _d 0
533     hl(i,j,bi,bj) = 0. _d 0
534     endif
535    
536     #else
537     #ifdef ALLOW_ATM_WIND
538     ustress(i,j,bi,bj) = atmrho*exf_BulkCdn(sh)*us*
539     & uwind(i,j,bi,bj)
540     vstress(i,j,bi,bj) = atmrho*exf_BulkCdn(sh)*us*
541     & vwind(i,j,bi,bj)
542     #endif /* ALLOW_ATM_WIND */
543     #endif /* ALLOW_ATM_TEMP */
544     enddo
545     enddo
546     enddo
547     enddo
548    
549     c Add all contributions.
550     do bj = mybylo(mythid),mybyhi(mythid)
551     do bi = mybxlo(mythid),mybxhi(mythid)
552     do j = 1,sny
553     do i = 1,snx
554     c Net surface heat flux.
555     #ifdef ALLOW_ATM_TEMP
556     hfl = 0. _d 0
557     hfl = hfl - hs(i,j,bi,bj)
558     hfl = hfl - hl(i,j,bi,bj)
559     hfl = hfl + lwflux(i,j,bi,bj)
560     #ifndef ALLOW_KPP
561     hfl = hfl + swflux(i,j,bi,bj)
562     #endif /* ALLOW_KPP undef */
563     c Heat flux:
564 heimbach 1.3 hflux(i,j,bi,bj) = hfl*maskc(i,j,1,bi,bj)
565 heimbach 1.1 c Salt flux from Precipitation and Evaporation.
566     sflux(i,j,bi,bj) = precip(i,j,bi,bj) - evap(i,j,bi,bj)
567     #endif /* ALLOW_ATM_TEMP */
568    
569     #else
570 heimbach 1.3 hflux(i,j,bi,bj) = hflux(i,j,bi,bj)*maskc(i,j,1,bi,bj)
571     sflux(i,j,bi,bj) = sflux(i,j,bi,bj)*maskc(i,j,1,bi,bj)
572 heimbach 1.1 #endif /* ALLOW_BULKFORMULAE */
573 heimbach 1.4
574     #ifdef ALLOW_RUNOFF
575     sflux(i,j,bi,bj) = sflux(i,j,bi,bj) +
576     & runoff(i,j,bi,bj)*maskc(i,j,1,bi,bj)
577     #endif /* ALLOW_RUNOFF */
578    
579 heimbach 1.1 enddo
580     enddo
581     enddo
582     enddo
583    
584     c Update the tile edges.
585     _EXCH_XY_R8(hflux, mythid)
586     _EXCH_XY_R8(sflux, mythid)
587 cheisey 1.5 c _EXCH_XY_R8(ustress, mythid)
588     c _EXCH_XY_R8(vstress, mythid)
589     CALL EXCH_UV_XY_RS(ustress, vstress, .TRUE., myThid)
590 heimbach 1.1
591     #ifdef ALLOW_KPP
592     _EXCH_XY_R8(swflux, mythid)
593     #endif /* ALLOW_KPP */
594    
595     end

  ViewVC Help
Powered by ViewVC 1.1.22