-
Notifications
You must be signed in to change notification settings - Fork 83
/
Copy pathshr_flux_mod.F90
2314 lines (1964 loc) · 99.4 KB
/
shr_flux_mod.F90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
module shr_flux_mod
! atm/ocn/flux calculations
! !USES:
use shr_kind_mod, only : R8=>SHR_KIND_R8, IN=>SHR_KIND_IN ! shared kinds
use shr_const_mod ! shared constants
use shr_sys_mod ! shared system routines
implicit none
private ! default private
! !PUBLIC TYPES:
! none
! !PUBLIC MEMBER FUNCTIONS:
public :: flux_atmOcn ! computes atm/ocn fluxes
public :: flux_atmOcn_diurnal ! computes atm/ocn fluxes with diurnal cycle
public :: flux_atmOcn_UA ! computes atm/ocn fluxes using University of Ariz algorithm (Zeng et al., 1998)
public :: flux_MOstability ! boundary layer stability scales/functions
public :: shr_flux_adjust_constants ! adjust constant values used in flux calculations. (used by CAM as well)
! !PRIVATE MEMBER FUNCTIONS:
private :: psi_ua
private :: qsat_ua
private :: rough_ua
private :: cuberoot
private :: cor30a
private :: psiuo
private :: psit_30
! !PUBLIC DATA MEMBERS:
integer(IN),parameter,public :: shr_flux_MOwScales = 1 ! w scales option
integer(IN),parameter,public :: shr_flux_MOfunctions = 2 ! functions option
real (R8),parameter,public :: shr_flux_MOgammaM = 3.59_R8
real (R8),parameter,public :: shr_flux_MOgammaS = 7.86_R8
!--- rename kinds for local readability only ---
integer,parameter :: debug = 0 ! internal debug level
! The follow variables are not declared as parameters so that they can be
! adjusted to support aquaplanet and potentially other simple model modes.
! The flux_adjust_constants subroutine is called to set the desired
! values. The default values are from shr_const_mod. Currently they are
! only used by the flux_atmocn routine.
real(R8) :: loc_zvir = shr_const_zvir
real(R8) :: loc_cpdair = shr_const_cpdair
real(R8) :: loc_cpvir = shr_const_cpvir
real(R8) :: loc_karman = shr_const_karman
real(R8) :: loc_g = shr_const_g
real(R8) :: loc_latvap = shr_const_latvap
real(R8) :: loc_latice = shr_const_latice
real(R8) :: loc_stebol = shr_const_stebol
real(R8) :: loc_tkfrz = shr_const_tkfrz
! These control convergence of the iterative flux calculation
! (For Large and Pond scheme only; not UA or COARE).
real(r8) :: flux_con_tol = 0.0_R8
integer(IN) :: flux_con_max_iter = 2
!--- cold air outbreak parameters (Mahrt & Sun 1995,MWR) -------------
logical :: use_coldair_outbreak_mod = .false.
real(R8),parameter :: alpha = 1.4_R8
real(R8),parameter :: maxscl =2._R8 ! maximum wind scaling for flux
real(R8),parameter :: td0 = -10._R8 ! start t-ts for scaling
character(len=*), parameter :: sourcefile = &
__FILE__
!===============================================================================
contains
!===============================================================================
subroutine shr_flux_adjust_constants( &
zvir, cpair, cpvir, karman, gravit, &
latvap, latice, stebol, flux_convergence_tolerance, &
flux_convergence_max_iteration, &
coldair_outbreak_mod)
! Adjust local constants. Used to support simple models.
real(R8), optional, intent(in) :: zvir
real(R8), optional, intent(in) :: cpair
real(R8), optional, intent(in) :: cpvir
real(R8), optional, intent(in) :: karman
real(R8), optional, intent(in) :: gravit
real(R8), optional, intent(in) :: latvap
real(R8), optional, intent(in) :: latice
real(R8), optional, intent(in) :: stebol
real(r8), optional, intent(in) :: flux_convergence_tolerance
integer(in), optional, intent(in) :: flux_convergence_max_iteration
logical, optional, intent(in) :: coldair_outbreak_mod
!----------------------------------------------------------------------------
if (present(zvir)) loc_zvir = zvir
if (present(cpair)) loc_cpdair = cpair
if (present(cpvir)) loc_cpvir = cpvir
if (present(karman)) loc_karman = karman
if (present(gravit)) loc_g = gravit
if (present(latvap)) loc_latvap = latvap
if (present(latice)) loc_latice = latice
if (present(stebol)) loc_stebol = stebol
if (present(flux_convergence_tolerance)) flux_con_tol = flux_convergence_tolerance
if (present(flux_convergence_max_iteration)) flux_con_max_iter = flux_convergence_max_iteration
if (present(coldair_outbreak_mod)) use_coldair_outbreak_mod = coldair_outbreak_mod
end subroutine shr_flux_adjust_constants
!===============================================================================
! !IROUTINE: flux_atmOcn -- internal atm/ocn flux calculation
!
! !DESCRIPTION:
!
! Internal atm/ocn flux calculation
!
! !REVISION HISTORY:
! 2002-Jun-10 - B. Kauffman - code migrated from cpl5 to cpl6
! 2003-Apr-02 - B. Kauffman - taux & tauy now utilize ocn velocity
! 2003-Apr-02 - B. Kauffman - tref,qref,duu10n mods as per Bill Large
! 2006-Nov-07 - B. Kauffman - code migrated from cpl6 to share
!
! 2011-Mar-13 - J. Nusbaumer - Water Isotope ocean flux added.
! 2019-May-16 - Jack Reeves Eyre (UA) and Kai Zhang (PNNL) -
! Added COARE/Fairall surface flux scheme option
! (ocn_surface_flux_scheme .eq. 1) based on code from
! Thomas Toniazzo (Bjerknes Centre, Bergen) ”
!===============================================================================
SUBROUTINE flux_atmOcn(logunit, nMax ,zbot ,ubot ,vbot ,thbot , &
& qbot ,s16O ,sHDO ,s18O ,rbot, &
& tbot ,us ,vs, pslv, &
& ts ,mask , seq_flux_atmocn_minwind, &
& sen ,lat ,lwup , &
& r16O, rhdo, r18O, &
& evap ,evap_16O, evap_HDO, evap_18O, &
& taux ,tauy ,tref ,qref , &
& ocn_surface_flux_scheme, &
& duu10n, ustar_sv ,re_sv ,ssq_sv, &
& missval)
! !USES:
use water_isotopes, only: wiso_flxoce !subroutine used to calculate water isotope fluxes.
implicit none
! !INPUT/OUTPUT PARAMETERS:
!--- input arguments --------------------------------
integer ,intent(in) :: logunit
integer(IN),intent(in) :: nMax ! data vector length
integer(IN),intent(in) :: mask (nMax) ! ocn domain mask 0 <=> out of domain
integer(IN),intent(in) :: ocn_surface_flux_scheme
real(R8) ,intent(in) :: zbot (nMax) ! atm level height (m)
real(R8) ,intent(in) :: ubot (nMax) ! atm u wind (m/s)
real(R8) ,intent(in) :: vbot (nMax) ! atm v wind (m/s)
real(R8) ,intent(in) :: thbot(nMax) ! atm potential T (K)
real(R8) ,intent(in) :: qbot (nMax) ! atm specific humidity (kg/kg)
real(R8) ,intent(in) :: s16O (nMax) ! atm H216O tracer conc. (kg/kg)
real(R8) ,intent(in) :: sHDO (nMax) ! atm HDO tracer conc. (kg/kg)
real(R8) ,intent(in) :: s18O (nMax) ! atm H218O tracer conc. (kg/kg)
real(R8) ,intent(in) :: r16O (nMax) ! ocn H216O tracer ratio/Rstd
real(R8) ,intent(in) :: rHDO (nMax) ! ocn HDO tracer ratio/Rstd
real(R8) ,intent(in) :: r18O (nMax) ! ocn H218O tracer ratio/Rstd
real(R8) ,intent(in) :: rbot (nMax) ! atm air density (kg/m^3)
real(R8) ,intent(in) :: tbot (nMax) ! atm T (K)
real(R8) ,intent(in) :: pslv (nMax) ! atm sea level pressure(Pa)
real(R8) ,intent(in) :: us (nMax) ! ocn u-velocity (m/s)
real(R8) ,intent(in) :: vs (nMax) ! ocn v-velocity (m/s)
real(R8) ,intent(in) :: ts (nMax) ! ocn temperature (K)
real(R8) ,intent(in) :: seq_flux_atmocn_minwind ! minimum wind speed for atmocn (m/s)
!--- output arguments -------------------------------
real(R8),intent(out) :: sen (nMax) ! heat flux: sensible (W/m^2)
real(R8),intent(out) :: lat (nMax) ! heat flux: latent (W/m^2)
real(R8),intent(out) :: lwup (nMax) ! heat flux: lw upward (W/m^2)
real(R8),intent(out) :: evap (nMax) ! water flux: evap ((kg/s)/m^2)
real(R8),intent(out) :: evap_16O (nMax) ! water flux: evap ((kg/s/m^2)
real(R8),intent(out) :: evap_HDO (nMax) ! water flux: evap ((kg/s)/m^2)
real(R8),intent(out) :: evap_18O (nMax) ! water flux: evap ((kg/s/m^2)
real(R8),intent(out) :: taux (nMax) ! surface stress, zonal (N)
real(R8),intent(out) :: tauy (nMax) ! surface stress, maridional (N)
real(R8),intent(out) :: tref (nMax) ! diag: 2m ref height T (K)
real(R8),intent(out) :: qref (nMax) ! diag: 2m ref humidity (kg/kg)
real(R8),intent(out) :: duu10n(nMax) ! diag: 10m wind speed squared (m/s)^2
real(R8),intent(out),optional :: ustar_sv(nMax) ! diag: ustar
real(R8),intent(out),optional :: re_sv (nMax) ! diag: sqrt of exchange coefficient (water)
real(R8),intent(out),optional :: ssq_sv (nMax) ! diag: sea surface humidity (kg/kg)
real(R8),intent(in) ,optional :: missval ! masked value
!--- local constants --------------------------------
real(R8),parameter :: zref = 10.0_R8 ! reference height (m)
real(R8),parameter :: ztref = 2.0_R8 ! reference height for air T (m)
!!++ Large only
!real(R8),parameter :: cexcd = 0.0346_R8 ! ratio Ch(water)/CD
!real(R8),parameter :: chxcds = 0.018_R8 ! ratio Ch(heat)/CD for stable case
!real(R8),parameter :: chxcdu = 0.0327_R8 ! ratio Ch(heat)/CD for unstable case
!!++ COARE only
real(R8),parameter :: zpbl =700.0_R8 ! PBL depth [m] for gustiness parametriz.
!--- local variables --------------------------------
integer(IN) :: n ! vector loop index
integer(IN) :: iter
real(R8) :: vmag ! surface wind magnitude (m/s)
real(R8) :: ssq ! sea surface humidity (kg/kg)
real(R8) :: delt ! potential T difference (K)
real(R8) :: delq ! humidity difference (kg/kg)
real(R8) :: stable ! stability factor
real(R8) :: rdn ! sqrt of neutral exchange coeff (momentum)
real(R8) :: rhn ! sqrt of neutral exchange coeff (heat)
real(R8) :: ren ! sqrt of neutral exchange coeff (water)
real(R8) :: rd ! sqrt of exchange coefficient (momentum)
real(R8) :: rh ! sqrt of exchange coefficient (heat)
real(R8) :: re ! sqrt of exchange coefficient (water)
real(R8) :: ustar ! ustar
real(r8) :: ustar_prev
real(R8) :: qstar ! qstar
real(R8) :: tstar ! tstar
real(R8) :: hol ! H (at zbot) over L
real(R8) :: xsq ! ?
real(R8) :: xqq ! ?
!!++ Large only
real(R8) :: psimh ! stability function at zbot (momentum)
real(R8) :: psixh ! stability function at zbot (heat and water)
real(R8) :: psix2 ! stability function at ztref reference height
real(R8) :: alz ! ln(zbot/zref)
real(R8) :: al2 ! ln(zref/ztref)
real(R8) :: u10n ! 10m neutral wind
real(R8) :: tau ! stress at zbot
real(R8) :: cp ! specific heat of moist air
real(R8) :: fac ! vertical interpolation factor
real(R8) :: spval ! local missing value
!!++ COARE only
real(R8) :: zo,zot,zoq ! roughness lengths
real(R8) :: hsb,hlb ! sens & lat heat flxs at zbot
real(R8) :: trf,qrf,urf,vrf ! reference-height quantities
!--- local functions --------------------------------
real(R8) :: qsat ! function: the saturation humididty of air (kg/m^3)
!!++ Large only (formula v*=[c4/U10+c5+c6*U10]*U10 in Large et al. 1994)
real(R8) :: cdn ! function: neutral drag coeff at 10m
!!++ Large only (stability functions)
real(R8) :: psimhu ! function: unstable part of psimh
real(R8) :: psixhu ! function: unstable part of psimx
real(R8) :: Umps ! dummy arg ~ wind velocity (m/s)
real(R8) :: Tk ! dummy arg ~ temperature (K)
real(R8) :: xd ! dummy arg ~ ?
!--- for cold air outbreak calc --------------------------------
real(R8) :: tdiff(nMax) ! tbot - ts
real(R8) :: vscl
qsat(Tk) = 640380.0_R8 / exp(5107.4_R8/Tk)
cdn(Umps) = 0.0027_R8 / Umps + 0.000142_R8 + 0.0000764_R8 * Umps
psimhu(xd) = log((1.0_R8+xd*(2.0_R8+xd))*(1.0_R8+xd*xd)/8.0_R8) - 2.0_R8*atan(xd) + 1.571_R8
psixhu(xd) = 2.0_R8 * log((1.0_R8 + xd*xd)/2.0_R8)
!--- formats ----------------------------------------
character(*),parameter :: subName = '(flux_atmOcn) '
character(*),parameter :: F00 = "('(flux_atmOcn) ',4a)"
!-------------------------------------------------------------------------------
! PURPOSE:
! computes atm/ocn surface fluxes
!
! NOTES:
! o all fluxes are positive downward
! o net heat flux = net sw + lw up + lw down + sen + lat
! o here, tstar = <WT>/U*, and qstar = <WQ>/U*.
! o wind speeds should all be above a minimum speed (eg. 1.0 m/s)
!
! ASSUMPTIONS:
! Large:
! o Neutral 10m drag coeff: cdn = .0027/U10 + .000142 + .0000764 U10
! o Neutral 10m stanton number: ctn = .0327 sqrt(cdn), unstable
! ctn = .0180 sqrt(cdn), stable
! o Neutral 10m dalton number: cen = .0346 sqrt(cdn)
! o The saturation humidity of air at T(K): qsat(T) (kg/m^3)
! COARE:
! o use COAREv3.0 function (tht 22/11/2013)
!-------------------------------------------------------------------------------
if (debug > 0) write(logunit,F00) "enter"
if (present(missval)) then
spval = missval
else
spval = shr_const_spval
endif
u10n = spval
rh = spval
psixh = spval
hol=spval
!--- for cold air outbreak calc --------------------------------
tdiff= tbot - ts
!!.................................................................
!! ocn_surface_flux_scheme = 0 : Default CESM1.2
!! = 1 : COARE algorithm
!! = 2 : UA algorithm (separate subroutine)
!!.................................................................
! Default flux scheme.
if (ocn_surface_flux_scheme .eq. 0) then
al2 = log(zref/ztref)
DO n=1,nMax
if (mask(n) /= 0) then
!--- compute some needed quantities ---
vmag = max(seq_flux_atmocn_minwind, sqrt( (ubot(n)-us(n))**2 + (vbot(n)-vs(n))**2) )
if (use_coldair_outbreak_mod) then
! Cold Air Outbreak Modification:
! Increase windspeed for negative tbot-ts
! based on Mahrt & Sun 1995,MWR
if (tdiff(n).lt.td0) then
vscl=min((1._R8+alpha*(abs(tdiff(n)-td0)**0.5_R8/abs(vmag))),maxscl)
vmag=vmag*vscl
endif
endif
ssq = 0.98_R8 * qsat(ts(n)) / rbot(n) ! sea surf hum (kg/kg)
delt = thbot(n) - ts(n) ! pot temp diff (K)
delq = qbot(n) - ssq ! spec hum dif (kg/kg)
alz = log(zbot(n)/zref)
cp = loc_cpdair*(1.0_R8 + loc_cpvir*ssq)
!------------------------------------------------------------
! first estimate of Z/L and ustar, tstar and qstar
!------------------------------------------------------------
!--- neutral coefficients, z/L = 0.0 ---
stable = 0.5_R8 + sign(0.5_R8 , delt)
rdn = sqrt(cdn(vmag))
rhn = (1.0_R8-stable) * 0.0327_R8 + stable * 0.018_R8
!(1.0_R8-stable) * chxcdu + stable * chxcds
ren = 0.0346_R8 !cexcd
!--- ustar, tstar, qstar ---
ustar = rdn * vmag
tstar = rhn * delt
qstar = ren * delq
ustar_prev = ustar*2.0_R8
iter = 0
do while( abs((ustar - ustar_prev)/ustar) > flux_con_tol .and. iter < flux_con_max_iter)
iter = iter + 1
ustar_prev = ustar
!--- compute stability & evaluate all stability functions ---
hol = loc_karman*loc_g*zbot(n)* &
(tstar/thbot(n)+qstar/(1.0_R8/loc_zvir+qbot(n)))/ustar**2
hol = sign( min(abs(hol),10.0_R8), hol )
stable = 0.5_R8 + sign(0.5_R8 , hol)
xsq = max(sqrt(abs(1.0_R8 - 16.0_R8*hol)) , 1.0_R8)
xqq = sqrt(xsq)
psimh = -5.0_R8*hol*stable + (1.0_R8-stable)*psimhu(xqq)
psixh = -5.0_R8*hol*stable + (1.0_R8-stable)*psixhu(xqq)
!--- shift wind speed using old coefficient ---
rd = rdn / (1.0_R8 + max(rdn/loc_karman*(alz-psimh), -0.5_r8))
u10n = vmag * rd / rdn
!--- update transfer coeffs at 10m and neutral stability ---
rdn = sqrt(cdn(u10n))
ren = 0.0346_R8 !cexcd
rhn = (1.0_R8-stable)*0.0327_R8 + stable * 0.018_R8
!(1.0_R8-stable) * chxcdu + stable * chxcds
!--- shift all coeffs to measurement height and stability ---
rd = rdn / (1.0_R8 + rdn/loc_karman*(alz-psimh))
rh = rhn / (1.0_R8 + rhn/loc_karman*(alz-psixh))
re = ren / (1.0_R8 + ren/loc_karman*(alz-psixh))
!--- update ustar, tstar, qstar using updated, shifted coeffs --
ustar = rd * vmag
tstar = rh * delt
qstar = re * delq
enddo
if (iter < 1) then
write(logunit,*) ustar,ustar_prev,flux_con_tol,flux_con_max_iter
call shr_sys_abort('No iterations performed in flux_atmocn_mod')
end if
!------------------------------------------------------------
! compute the fluxes
!------------------------------------------------------------
tau = rbot(n) * ustar * ustar
!--- momentum flux ---
taux(n) = tau * (ubot(n)-us(n)) / vmag
tauy(n) = tau * (vbot(n)-vs(n)) / vmag
!--- heat flux ---
sen (n) = cp * tau * tstar / ustar
lat (n) = loc_latvap * tau * qstar / ustar
lwup(n) = -loc_stebol * ts(n)**4
!--- water flux ---
evap(n) = lat(n)/loc_latvap
!---water isotope flux ---
call wiso_flxoce(2,rbot(n),zbot(n),s16O(n),ts(n),r16O(n),ustar,re,ssq,evap_16O(n), &
qbot(n),evap(n))
call wiso_flxoce(3,rbot(n),zbot(n),sHDO(n),ts(n),rHDO(n),ustar,re,ssq, evap_HDO(n),&
qbot(n),evap(n))
call wiso_flxoce(4,rbot(n),zbot(n),s18O(n),ts(n),r18O(n),ustar,re,ssq, evap_18O(n), &
qbot(n),evap(n))
!------------------------------------------------------------
! compute diagnositcs: 2m ref T & Q, 10m wind speed squared
!------------------------------------------------------------
hol = hol*ztref/zbot(n)
xsq = max( 1.0_R8, sqrt(abs(1.0_R8-16.0_R8*hol)) )
xqq = sqrt(xsq)
psix2 = -5.0_R8*hol*stable + (1.0_R8-stable)*psixhu(xqq)
fac = (rh/loc_karman) * (alz + al2 - psixh + psix2 )
tref(n) = thbot(n) - delt*fac
tref(n) = tref(n) - 0.01_R8*ztref ! pot temp to temp correction
fac = (re/loc_karman) * (alz + al2 - psixh + psix2 )
qref(n) = qbot(n) - delq*fac
duu10n(n) = u10n*u10n ! 10m wind speed squared
!------------------------------------------------------------
! optional diagnostics, needed for water tracer fluxes (dcn)
!------------------------------------------------------------
if (present(ustar_sv)) ustar_sv(n) = ustar
if (present(re_sv )) re_sv(n) = re
if (present(ssq_sv )) ssq_sv(n) = ssq
else
!------------------------------------------------------------
! no valid data here -- out of domain
!------------------------------------------------------------
sen (n) = spval ! sensible heat flux (W/m^2)
lat (n) = spval ! latent heat flux (W/m^2)
lwup (n) = spval ! long-wave upward heat flux (W/m^2)
evap (n) = spval ! evaporative water flux ((kg/s)/m^2)
evap_16O (n) = spval !water tracer flux (kg/s)/m^2)
evap_HDO (n) = spval !HDO tracer flux (kg/s)/m^2)
evap_18O (n) = spval !H218O tracer flux (kg/s)/m^2)
taux (n) = spval ! x surface stress (N)
tauy (n) = spval ! y surface stress (N)
tref (n) = spval ! 2m reference height temperature (K)
qref (n) = spval ! 2m reference height humidity (kg/kg)
duu10n(n) = spval ! 10m wind speed squared (m/s)^2
if (present(ustar_sv)) ustar_sv(n) = spval
if (present(re_sv )) re_sv (n) = spval
if (present(ssq_sv )) ssq_sv (n) = spval
endif
ENDDO
else if (ocn_surface_flux_scheme .eq. 1) then
!!.................................
!! use COARE algorithm
!!.................................
DO n=1,nMax
if (mask(n) /= 0) then
!--- compute some needed quantities ---
vmag = max(seq_flux_atmocn_minwind, sqrt( (ubot(n)-us(n))**2 + (vbot(n)-vs(n))**2) )
if (use_coldair_outbreak_mod) then
! Cold Air Outbreak Modification:
! Increase windspeed for negative tbot-ts
! based on Mahrt & Sun 1995,MWR
if (tdiff(n).lt.td0) then
vscl=min((1._R8+alpha*(abs(tdiff(n)-td0)**0.5_R8/abs(vmag))),maxscl)
vmag=vmag*vscl
endif
endif
ssq = 0.98_R8 * qsat(ts(n)) / rbot(n) ! sea surf hum (kg/kg)
call cor30a(ubot(n),vbot(n),tbot(n),qbot(n),rbot(n) & ! in atm params
& ,us(n),vs(n),ts(n),ssq & ! in surf params
& ,zpbl,zbot(n),zbot(n),zref,ztref,ztref & ! in heights
& ,tau,hsb,hlb & ! out: fluxes
& ,zo,zot,zoq,hol,ustar,tstar,qstar & ! out: ss scales
& ,rd,rh,re & ! out: exch. coeffs
& ,trf,qrf,urf,vrf) ! out: reference-height params
! for the sake of maintaining same defs
hol=zbot(n)/hol
rd=sqrt(rd)
rh=sqrt(rh)
re=sqrt(re)
!--- momentum flux ---
taux(n) = tau * (ubot(n)-us(n)) / vmag
tauy(n) = tau * (vbot(n)-vs(n)) / vmag
!--- heat flux ---
sen (n) = hsb
lat (n) = hlb
lwup(n) = -shr_const_stebol * ts(n)**4
!--- water flux ---
evap(n) = lat(n)/shr_const_latvap
!---water isotope flux ---
call wiso_flxoce(2,rbot(n),zbot(n),s16O(n),ts(n),r16O(n),ustar,re,ssq, evap_16O(n), &
qbot(n),evap(n))
call wiso_flxoce(3,rbot(n),zbot(n),sHDO(n),ts(n),rHDO(n),ustar,re,ssq, evap_HDO(n),&
qbot(n),evap(n))
call wiso_flxoce(4,rbot(n),zbot(n),s18O(n),ts(n),r18O(n),ustar,re,ssq, evap_18O(n), &
qbot(n),evap(n))
!------------------------------------------------------------
! compute diagnositcs: 2m ref T & Q, 10m wind speed squared
!------------------------------------------------------------
tref(n) = trf
qref(n) = qrf
duu10n(n) = urf**2+vrf**2
!------------------------------------------------------------
! optional diagnostics, needed for water tracer fluxes (dcn)
!------------------------------------------------------------
if (present(ustar_sv)) ustar_sv(n) = ustar
if (present(re_sv )) re_sv(n) = re
if (present(ssq_sv )) ssq_sv(n) = ssq
else
!------------------------------------------------------------
! no valid data here -- out of domain
!------------------------------------------------------------
sen (n) = spval ! sensible heat flux (W/m^2)
lat (n) = spval ! latent heat flux (W/m^2)
lwup (n) = spval ! long-wave upward heat flux (W/m^2)
evap (n) = spval ! evaporative water flux ((kg/s)/m^2)
evap_16O (n) = spval ! water tracer flux (kg/s)/m^2)
evap_HDO (n) = spval ! HDO tracer flux (kg/s)/m^2)
evap_18O (n) = spval ! H218O tracer flux (kg/s)/m^2)
taux (n) = spval ! x surface stress (N)
tauy (n) = spval ! y surface stress (N)
tref (n) = spval ! 2m reference height temperature (K)
qref (n) = spval ! 2m reference height humidity (kg/kg)
duu10n (n) = spval ! 10m wind speed squared (m/s)^2
if (present(ustar_sv)) ustar_sv(n) = spval
if (present(re_sv )) re_sv (n) = spval
if (present(ssq_sv )) ssq_sv (n) = spval
endif
ENDDO
else if (ocn_surface_flux_scheme .eq. 2) then
call flux_atmOcn_UA(logunit,&
nMax, zbot, ubot, vbot, thbot, &
qbot, s16O, sHDO, s18O, rbot, &
tbot, pslv, us, vs, &
ts, mask, sen, lat, lwup, &
r16O, rhdo, r18O, &
evap, evap_16O, evap_HDO, evap_18O, &
taux, tauy, tref, qref, &
duu10n, ustar_sv, re_sv, ssq_sv, &
missval)
else
call shr_sys_abort(subName//" subroutine flux_atmOcn requires ocn_surface_flux_scheme = 0, 1 or 2")
endif !! ocn_surface_flux_scheme
END subroutine flux_atmOcn
!===============================================================================
! !IROUTINE: flux_atmOcn_UA -- internal atm/ocn flux calculation
!
! !DESCRIPTION:
!
! Internal atm/ocn flux calculation
! using University of Arizona method.
!
! Reference:
! Zeng, X., M. Zhao, and R.E. Dickinson, 1998: Intercomparison of Bulk
! Aerodynamic Algorithms for the Computation of Sea Surface Fluxes
! Using TOGA COARE and TAO Data. J. Climate, 11, 2628–2644,
! https://doi.org/10.1175/1520-0442(1998)011<2628%3AIOBAAF>2.0.CO%3B2
!
! Equation numbers are from this paper.
!
! !REVISION HISTORY:
! 2017-Aug-28 - J. Reeves Eyre - code re-written for E3SM
! 2018-Oct-30 - J. Reeves Eyre - bug fix and add
! convective gustiness.
! 2019-May-08 - J. Reeves Eyre - remove convective gustiness
! and add cold air outbreak modification.
!===============================================================================
SUBROUTINE flux_atmOcn_UA(logunit, &
& nMax ,zbot ,ubot ,vbot ,thbot , &
& qbot ,s16O ,sHDO ,s18O ,rbot , &
& tbot , pslv ,us , vs , &
& ts ,mask ,sen ,lat ,lwup , &
& r16O, rhdo, r18O, &
& evap ,evap_16O, evap_HDO, evap_18O, &
& taux ,tauy ,tref ,qref , &
& duu10n, ustar_sv ,re_sv ,ssq_sv, &
& missval)
! !USES:
use water_isotopes, only: wiso_flxoce !subroutine used to calculate water isotope fluxes.
implicit none
! !INPUT/OUTPUT PARAMETERS:
!--- input arguments --------------------------------
integer ,intent(in) :: logunit
integer ,intent(in) :: nMax ! data vector length
integer ,intent(in) :: mask (nMax) ! ocn domain mask 0 <=> out of domain
real(R8) ,intent(in) :: zbot (nMax) ! atm level height (m)
real(R8) ,intent(in) :: ubot (nMax) ! atm u wind (m/s)
real(R8) ,intent(in) :: vbot (nMax) ! atm v wind (m/s)
real(R8) ,intent(in) :: thbot(nMax) ! atm potential T (K)
real(R8) ,intent(in) :: qbot (nMax) ! atm specific humidity (kg/kg)
real(R8) ,intent(in) :: s16O (nMax) ! atm H216O tracer conc. (kg/kg)
real(R8) ,intent(in) :: sHDO (nMax) ! atm HDO tracer conc. (kg/kg)
real(R8) ,intent(in) :: s18O (nMax) ! atm H218O tracer conc. (kg/kg)
real(R8) ,intent(in) :: r16O (nMax) ! ocn H216O tracer ratio/Rstd
real(R8) ,intent(in) :: rHDO (nMax) ! ocn HDO tracer ratio/Rstd
real(R8) ,intent(in) :: r18O (nMax) ! ocn H218O tracer ratio/Rstd
real(R8) ,intent(in) :: rbot (nMax) ! atm air density (kg/m^3)
real(R8) ,intent(in) :: tbot (nMax) ! atm T (K)
real(R8) ,intent(in) :: pslv (nMax) ! sea level pressure (Pa)
real(R8) ,intent(in) :: us (nMax) ! ocn u-velocity (m/s)
real(R8) ,intent(in) :: vs (nMax) ! ocn v-velocity (m/s)
real(R8) ,intent(in) :: ts (nMax) ! ocn temperature (K)
!--- output arguments -------------------------------
real(R8),intent(out) :: sen (nMax) ! heat flux: sensible (W/m^2)
real(R8),intent(out) :: lat (nMax) ! heat flux: latent (W/m^2)
real(R8),intent(out) :: lwup (nMax) ! heat flux: lw upward (W/m^2)
real(R8),intent(out) :: evap (nMax) ! water flux: evap ((kg/s)/m^2)
real(R8),intent(out) :: evap_16O (nMax) ! water flux: evap ((kg/s/m^2)
real(R8),intent(out) :: evap_HDO (nMax) ! water flux: evap ((kg/s)/m^2)
real(R8),intent(out) :: evap_18O (nMax) ! water flux: evap ((kg/s/m^2)
real(R8),intent(out) :: taux (nMax) ! surface stress, zonal (N)
real(R8),intent(out) :: tauy (nMax) ! surface stress, maridional (N)
real(R8),intent(out) :: tref (nMax) ! diag: 2m ref height T (K)
real(R8),intent(out) :: qref (nMax) ! diag: 2m ref humidity (kg/kg)
real(R8),intent(out) :: duu10n(nMax) ! diag: 10m wind speed squared (m/s)^2
real(R8),intent(out),optional :: ustar_sv(nMax) ! diag: ustar
real(R8),intent(out),optional :: re_sv (nMax) ! diag: sqrt of exchange coefficient (water)
real(R8),intent(out),optional :: ssq_sv (nMax) ! diag: sea surface humidity (kg/kg)
real(R8),intent(in) ,optional :: missval ! masked value
!--- local constants --------------------------------
real(R8),parameter :: zetam = -1.574_R8 ! Very unstable zeta cutoff for momentum (-)
real(R8),parameter :: zetat = -0.465_R8 ! Very unstable zeta cutoff for T/q (-)
real(R8),parameter :: umin = 0.1_R8 ! minimum wind speed (m/s)
real(R8),parameter :: zref = 10.0_R8 ! reference height (m)
real(R8),parameter :: ztref = 2.0_R8 ! reference height for air T (m)
real(R8),parameter :: beta = 1.0_R8 ! constant used in W* calculation (-)
real(R8),parameter :: zpbl = 1000.0_R8 ! PBL height used in W* calculation (m)
real(R8),parameter :: gamma = 0.0098_R8 ! Dry adiabatic lapse rate (K/m)
real(R8),parameter :: onethird = 1.0_R8/3.0_R8 ! Used repeatedly.
!--- local variables --------------------------------
integer(IN) :: n ! vector loop index
integer(IN) :: i ! iteration loop index
real(R8) :: vmag_abs ! surface wind magnitude (m s-1)
real(R8) :: vmag_rel ! surface wind magnitude relative to
! surface current (m s-1)
real(R8) :: vmag ! surface wind magnitude with large
! eddy correction and minimum value (m s-1)
! (This can change on each iteration.)
real(R8) :: thv ! virtual temperature (K)
real(R8) :: ssq ! sea surface humidity (kg/kg)
real(R8) :: delth ! potential T difference (K)
real(R8) :: delthv ! virtual potential T difference (K)
real(R8) :: delq ! humidity difference (kg/kg)
real(R8) :: ustar ! friction velocity (m s-1)
real(R8) :: qstar ! humidity scaling parameter (kg/kg)
real(R8) :: tstar ! temperature scaling parameter (K)
real(R8) :: thvstar ! virtual temperature scaling parameter (K)
real(R8) :: wstar ! convective velocity scale (m s-1)
real(R8) :: zeta ! dimensionless height (z / Obukhov length)
real(R8) :: obu ! Obukhov length (m)
real(R8) :: tau ! magnitude of wind stress (N m-2)
real(R8) :: cp ! specific heat of moist air (J kg-1 K-1)
real(R8) :: xlv ! Latent heat of vaporization (J kg-1)
real(R8) :: visa ! Kinematic viscosity of dry air (m2 s-1)
real(R8) :: tbot_oC ! Temperature used in visa (deg C)
real(R8) :: rb ! Bulk Richardson number (-)
real(R8) :: zo ! Roughness length for momentum (m)
real(R8) :: zoq ! Roughness length for moisture (m)
real(R8) :: zot ! Roughness length for heat (m)
real(R8) :: u10 ! 10-metre wind speed (m s-1)
real(R8) :: re ! Moisture exchange coefficient for compatibility
! with default algorithm.
real(R8) :: spval ! local missing value
real(R8) :: loc_epsilon ! Ratio of gas constants (-)
!--- for cold air outbreak calc --------------------------------
real(R8) :: tdiff(nMax) ! tbot - ts
real(R8) :: vscl
!--- formats ----------------------------------------
character(*),parameter :: subName = '(flux_atmOcn) '
character(*),parameter :: F00 = "('(flux_atmOcn) ',4a)"
!-----
! Straight from original subroutine.
if (debug > 0) write(logunit,F00) "enter"
if (present(missval)) then
spval = missval
else
spval = shr_const_spval
endif
!-----
! Evaluate loc_epsilon.
loc_epsilon = 1.0_R8 / (1.0_R8 + loc_zvir)
!--- for cold air outbreak calc --------------------------------
tdiff = tbot - ts
! Loop over grid points.
DO n=1,nMax
if (mask(n) /= 0) then
!-----Calculate some required near surface variables.---------
vmag_abs = sqrt( ubot(n)**2 + vbot(n)**2 )
vmag_rel = sqrt( (ubot(n)-us(n))**2 + (vbot(n)-vs(n))**2 )
! For Cold Air Outbreak Modification (based on Mahrt & Sun 1995,MWR):
if (use_coldair_outbreak_mod) then
! Increase windspeed for negative tbot-ts
if (tdiff(n).lt.td0) then
vscl=min((1._R8+alpha*(abs(tdiff(n)-td0)**0.5_R8/abs(vmag_rel))),maxscl)
vmag_rel=vmag_rel*vscl
endif
endif
delth = thbot(n) - ts(n) ! Pot. temp. difference with surface (K)
! Note this is equivalent to Zeng et al
! (1998) version = delt + 0.0098*zbot
thv = thbot(n)*(1.0_R8+0.61_R8*qbot(n)) ! Virtual potential temperature (K)
! EQN (17):
!ssq = 0.98_R8 * qsat_ua(ts(n),ps, & ! Surface specific humidity (kg kg-1)
! loc_epsilon)
ssq = 0.98_R8 * qsat_ua(ts(n),pslv(n), & ! Surface specific humidity (kg kg-1)
loc_epsilon)
delq = qbot(n) - ssq ! Difference to surface (kg kg-1)
delthv = delth*(1.0_R8+0.61_R8*qbot(n)) + & ! Difference of virtual potential
& 0.61_R8*thbot(n)*delq ! temperature with surface (K)
xlv = 1.0e+6_R8 * & ! Latent heat of vaporization (J kg-1)
& (2.501_R8 - 0.00237_R8 * (ts(n) - loc_tkfrz))
tbot_oC = tbot(n) - loc_tkfrz
visa = 1.326e-5_R8 * (1.0_R8 + & ! Kinematic viscosity of dry
& 6.542e-3_R8*tbot_oC + & ! air (m2 s-1) from Andreas (1989)
& 8.301e-6_R8*tbot_oC*tbot_oC - & ! CRREL Rep. 89-11
& 4.84e-9_R8*tbot_oC*tbot_oC*tbot_oC)
cp = loc_cpdair*(1.0_R8 + loc_cpvir*ssq) ! specific heat of moist air (J kg-1 K-1)
!-----Initial values of u* and convective velocity.-----------
ustar = 0.06_R8
wstar = 0.5_R8
! Update wind speed if unstable regime.
if (delthv.lt.0.0_R8) then
! EQN (19)
vmag = sqrt( vmag_rel**2 + beta*beta*wstar*wstar )
else
! EQN (18)
vmag = max(umin,vmag_rel)
endif
!-----Iterate to compute new u* and z0.-----------------------
do i = 1,5
! EQN (24)
zo = 0.013_R8*ustar*ustar/loc_g + 0.11_R8*visa/ustar
! EQN (9) assuming neutral
ustar = loc_karman*vmag/log(zbot(n)/zo)
enddo
!-----Assess stability.---------------------------------------
rb = loc_g*zbot(n)*delthv / (thv*vmag*vmag) ! bulk Richardson number
if(rb.ge.0.0_R8) then
! Neutral or stable: EQNs (4), (9), (13) and definition of rb.
zeta = rb*log(zbot(n)/zo) / &
& (1.0_R8 - 5.0_R8*min(rb,0.19_R8))
else
! Unstable: EQNs (4), (8), (12) and definition of rb.
zeta = rb*log(zbot(n)/zo)
endif
obu = zbot(n)/zeta ! Obukhov length
obu = sign(max(zbot(n)/10.0_R8, abs(obu)), obu)
!-----Main iterations (2-10 iterations would be fine).-------
do i=1,10
! Update roughness lengths.
call rough_ua(zo,zot,zoq,ustar,visa)
! Wind variables.
zeta = zbot(n) / obu
if (zeta.lt.zetam) then
! Very unstable regime
! EQN (7) with extra z0 term.
ustar = loc_karman * vmag / (log(zetam*obu/zo) - &
& psi_ua(1_IN, zetam) + &
& psi_ua(1_IN, zo/obu) + &
& 1.14_R8 * ((-zeta)**onethird - (-zetam)**onethird) )
else if (zeta.lt.0.0_R8) then
! Unstable regime
! EQN (8) with extra z0 term.
ustar = loc_karman * vmag / (log(zbot(n)/zo) - &
& psi_ua(1_IN,zeta) + psi_ua(1_IN,zo/obu) )
else if (zeta.le.1.0_R8) then
! Stable regime
! EQN (9) with extra z0 term.
ustar = loc_karman * vmag / (log(zbot(n)/zo) + &
& 5.0_R8*zeta - 5.0_R8*zo/obu)
else
! Very stable regime
! EQN (10) with extra z0 term.
ustar = loc_karman * vmag / (log(obu/zo) + 5.0_R8 - &
& 5.0_R8*zo/obu + &
& (5.0_R8*log(zeta) + zeta - 1.0_R8) )
endif
! Temperature variables.
if(zeta.lt.zetat) then
! Very unstable regime
! EQN (11) with extra z0 term.
tstar = loc_karman * delth / (log(zetat*obu/zot) - &
& psi_ua(2_IN, zetat) + &
& psi_ua(2_IN, zot/obu) + &
& 0.8_R8*((-zetat)**(-onethird) - (-zeta)**(-onethird)) )
else if (zeta.lt.0.0_R8) then
! Unstable regime
! EQN (12) with extra z0 term.
tstar = loc_karman * delth / &
& (log(zbot(n)/zot) - psi_ua(2_IN,zeta) + psi_ua(2_IN,zot/obu))
else if (zeta.le.1.0_R8) then
! Stable regime
! EQN (13) with extra z0 term.
tstar = loc_karman * delth / (log(zbot(n)/zot) + &
& 5.0_R8*zeta - 5.0_R8*zot/obu)
else
! Very stable regime
! EQN (14) with extra z0 term.
tstar = loc_karman * delth / (log(obu/zot) + &
& 5.0_R8 - 5.0_R8*zot/obu + &
& (5.0_R8*log(zeta) + zeta - 1.0_R8) )
endif
! Humidity variables.
! This is done with re to give variable to save out like
! in old algorithm.
if (zeta.lt.zetat) then
! Very unstable regime
! EQN (11) with extra z0 term.
re = loc_karman / (log(zetat*obu/zoq) - psi_ua(2_IN,zetat) + &
& psi_ua(2_IN,zoq/obu) + &
& 0.8_R8*((-zetat)**(-onethird) - (-zeta)**(-onethird)) )
else if (zeta.lt.0.0_R8) then
! Unstable regime
! EQN (12) with extra z0 term.
re = loc_karman / &
& (log(zbot(n)/zoq) - psi_ua(2_IN,zeta) + psi_ua(2_IN,zoq/obu))
else if (zeta.le.1.0_R8) then
! Stable regime
! EQN (13) with extra z0 term.
re = loc_karman / &
& (log(zbot(n)/zoq) + 5.0_R8*zeta - 5.0_R8*zoq/obu)
else
! Very stable regime
! EQN (14) with extra z0 term.
re = loc_karman / &
& (log(obu/zoq) + 5.0_R8 - 5.0_R8*zoq/obu + &
& (5.0_R8*log(zeta) + zeta - 1.0_R8) )
endif
qstar = re * delq
! Update Obukhov length.
thvstar = tstar*(1.0_R8 + 0.61_R8*qbot(n)) + 0.61_R8*thbot(n)*qstar
! EQN (4)
obu = ustar*ustar * thv / (loc_karman*loc_g*thvstar)
obu = sign( max(zbot(n)/10.0_R8, abs(obu)) ,obu)
! Update wind speed if in unstable regime.
if (delthv.lt.0.0_R8) then
! EQN (20)
wstar = beta * (-loc_g*ustar*thvstar*zpbl/thv)**onethird
! EQN (19)
vmag = sqrt(vmag_rel**2 + wstar*wstar)
else
! EQN (18)
vmag = max(umin,vmag_rel)
endif
enddo ! End of iterations for ustar, tstar, qstar etc.
!-----Calculate fluxes and wind stress.---------------------
!--- momentum flux ---
! This should ensure zero wind stress when (relative) wind speed is zero,
! components are consistent with total, and we don't ever divide by zero.
! EQN (21)
tau = rbot(n) * ustar * ustar
taux(n) = tau * (ubot(n)-us(n)) / max(umin, vmag_rel)
tauy(n) = tau * (vbot(n)-vs(n)) / max(umin, vmag_rel)
!--- heat flux ---
! EQNs (22) and (23)
sen (n) = cp * rbot(n) * tstar * ustar
lat (n) = xlv * rbot(n) * qstar * ustar
lwup(n) = -loc_stebol * ts(n)**4
!--- water flux ---
evap(n) = lat(n)/xlv
!---water isotope flux ---
call wiso_flxoce(2,rbot(n),zbot(n),s16O(n),ts(n),r16O(n),ustar,re,ssq,evap_16O(n), &
qbot(n),evap(n))
call wiso_flxoce(3,rbot(n),zbot(n),sHDO(n),ts(n),rHDO(n),ustar,re,ssq, evap_HDO(n),&
qbot(n),evap(n))
call wiso_flxoce(4,rbot(n),zbot(n),s18O(n),ts(n),r18O(n),ustar,re,ssq, evap_18O(n), &
qbot(n),evap(n))
!------------------------------------------------------------
! compute diagnositcs: 2m ref T & Q, 10m wind speed squared
!------------------------------------------------------------
zeta = zbot(n) / obu
if (zeta.lt.zetat) then
if (zeta.lt.zetam) then
! Very unstable regime for U.
! EQN (7)
u10 = vmag_abs + (ustar/loc_karman) * &
& 1.14_R8 * ((-zref/obu)**onethird - (-zeta)**onethird)
else
! Unstable regime for U.
! EQN (8)
u10 = vmag_abs + (ustar/loc_karman) * &
& (log(zref/zbot(n)) - (psi_ua(1_IN,zref/obu) - psi_ua(1_IN,zeta)) )
endif
! Very unstable regime for T and q.
! EQN (11)
tref(n) = thbot(n) + (tstar/loc_karman) * &
& 0.8_R8 * ((-zeta)**(-onethird) - (-ztref/obu)**(-onethird))
qref(n) = qbot(n) + (qstar/loc_karman) * &
& 0.8_R8 * ((-zeta)**(-onethird) - (-ztref/obu)**(-onethird))
else if (zeta.lt.0.0_R8) then
! Unstable regime.
! EQN (8)
u10 = vmag_abs + (ustar/loc_karman) * &
& (log(zref/zbot(n)) - (psi_ua(1_IN,zref/obu) - psi_ua(1_IN,zeta)) )
! EQN (12)
tref(n) = thbot(n) + (tstar/loc_karman) * &
& (log(ztref/zbot(n)) - (psi_ua(2_IN,ztref/obu) - psi_ua(2_IN,zeta)) )
qref(n) = qbot(n) + (qstar/loc_karman) * &
& (log(ztref/zbot(n)) - (psi_ua(2_IN,ztref/obu) - psi_ua(2_IN,zeta)) )
else if (zeta.le.1.0_R8) then
! Stable regime.
! EQN (9)
u10 = vmag_abs + (ustar/loc_karman) * &
& (log(zref/zbot(n)) + 5.0_R8*zref/obu - 5.0_R8*zeta)
! EQN (13)
tref(n) = thbot(n) + (tstar/loc_karman) * &
& (log(ztref/zbot(n)) + 5.0_R8*ztref/obu - 5.0_R8*zeta)
qref(n) = qbot(n) + (qstar/loc_karman) * &
& (log(ztref/zbot(n)) + 5.0_R8*ztref/obu - 5.0_R8*zeta)
else
! Very stable regime.
! EQN (10)
u10 = vmag_abs + (ustar/loc_karman) * &
& (5.0_R8*log(zref/zbot(n)) + zref/obu - zeta)
! EQN (14)
tref(n) = thbot(n) + (tstar/loc_karman) * &
& (5.0_R8*log(ztref/zbot(n)) + ztref/obu - zeta)
qref(n) = qbot(n) + (qstar/loc_karman) * &
& (5.0_R8*log(ztref/zbot(n)) + ztref/obu - zeta)
endif
tref(n) = tref(n) - gamma*ztref ! pot. temp to temp correction
duu10n(n) = u10*u10 ! 10m wind speed squared