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

Contents of /MITgcm/pkg/exf/exf_getffields.F

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


Revision 1.9 - (show annotations) (download)
Tue Feb 25 06:35:46 2003 UTC (21 years, 3 months ago) by dimitri
Branch: MAIN
CVS Tags: checkpoint48i_post, checkpoint48h_post
Changes since 1.8: +23 -7 lines
o Added net flux to downward flux conversion to pkg/exf/exf_getffields.F
o Added SEAICE_initialHEFF to pkg/seaice

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

  ViewVC Help
Powered by ViewVC 1.1.22