forked from NCAR/ccpp-physics
-
Notifications
You must be signed in to change notification settings - Fork 37
/
Copy pathsatmedmfvdifq.F
2415 lines (2392 loc) · 80.3 KB
/
satmedmfvdifq.F
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
!> \file satmedmfvdifq.F
!! This file contains the CCPP-compliant SATMEDMF scheme (updated version) which
!! computes subgrid vertical turbulence mixing using scale-aware TKE-based moist
!! eddy-diffusion mass-flux (TKE-EDMF) parameterization (by Jongil Han).
module satmedmfvdifq
use mfpbltq_mod
use tridi_mod
use mfscuq_mod
contains
!> \defgroup module_satmedmfvdifq GFS TKE-EDMF PBL Module
!! This file contains the CCPP-compliant SATMEDMF scheme (updated version) which
!! computes subgrid vertical turbulence mixing using scale-aware TKE-based moist
!! eddy-diffusion mass-flux (TKE-EDMF) parameterization (by Jongil Han).
!> @{
!! \brief This subroutine contains all of the logic for the
!! scale-aware TKE-based moist eddy-diffusion mass-flux (TKE-EDMF, updated version) scheme.
!! For local turbulence mixing, a TKE closure model is used.
!! Updated version of satmedmfvdif.f (May 2019) to have better low level
!! inversion, to reduce the cold bias in lower troposphere,
!! and to reduce the negative wind speed bias in upper troposphere
!!
!! Incorporate the LES-based changes for TC simulation
!! (Chen et al.,2022, https://doi.org/10.1175/WAF-D-21-0168.1)
!! with additional improvements on MF working with Cu schemes
!! Xiaomin Chen, 5/2/2022
!!
!> \section arg_table_satmedmfvdifq_init Argument Table
!! \htmlinclude satmedmfvdifq_init.html
!!
subroutine satmedmfvdifq_init (satmedmf, &
& isatmedmf,isatmedmf_vdifq, &
& errmsg,errflg)
logical, intent(in ) :: satmedmf
integer, intent(in) :: isatmedmf,isatmedmf_vdifq
character(len=*), intent(out) :: errmsg
integer, intent(out) :: errflg
! Initialize CCPP error handling variables
errmsg = ''
errflg = 0
! Consistency checks
if (.not. satmedmf) then
write(errmsg,fmt='(*(a))') 'Logic error: satmedmf = .false.'
errflg = 1
return
end if
if (.not. isatmedmf==isatmedmf_vdifq) then
write(errmsg,fmt='(*(a))') 'Logic error: satmedmfvdif is ',
& 'called, but isatmedmf/=isatmedmf_vdifq.'
errflg = 1
return
end if
end subroutine satmedmfvdifq_init
!> \section arg_table_satmedmfvdifq_run Argument Table
!! \htmlinclude satmedmfvdifq_run.html
!!
!!\section gen_satmedmfvdifq GFS satmedmfvdifq General Algorithm
!! satmedmfvdifq_run() computes subgrid vertical turbulence mixing
!! using the scale-aware TKE-based moist eddy-diffusion mass-flux (EDMF) parameterization of
!! Han and Bretherton (2019) \cite Han_2019 .
!! -# The local turbulent mixing is represented by an eddy-diffusivity scheme which
!! is a function of a prognostic TKE.
!! -# For the convective boundary layer, nonlocal transport by large eddies
!! (mfpbltq.f), is represented using a mass flux approach (Siebesma et al.(2007) \cite Siebesma_2007 ).
!! -# A mass-flux approach is also used to represent the stratocumulus-top-induced turbulence
!! (mfscuq.f).
!! \section detail_satmedmfvidfq GFS satmedmfvdifq Detailed Algorithm
subroutine satmedmfvdifq_run(im,km,ntrac,ntcw,ntrw, &
& ntiw,ntke,grav,rd,cp,rv,hvap,hfus,fv,eps,epsm1, &
& dv,du,tdt,rtg,u1,v1,t1,q1,swh,hlw,xmu, &
& garea,zvfun,sigmaf, &
& psk,rbsoil,zorl,u10m,v10m,fm,fh, &
& tsea,heat,evap,stress,spd1,kpbl, &
& prsi,del,prsl,prslk,phii,phil,delt, &
& dspheat,dusfc,dvsfc,dtsfc,dqsfc,hpbl,dkt,dku, &
& kinver,xkzm_m,xkzm_h,xkzm_s,dspfac,bl_upfr,bl_dnfr, &
& rlmx,elmx,sfc_rlm,tc_pbl, &
& ntqv,dtend,dtidx,index_of_temperature,index_of_x_wind, &
& index_of_y_wind,index_of_process_pbl,gen_tend,ldiag3d, &
& errmsg,errflg)
!
use machine , only : kind_phys
use funcphys , only : fpvs
!
implicit none
!
!----------------------------------------------------------------------
integer, intent(in) :: im, km, ntrac, ntcw, ntrw, ntiw, &
& ntke, ntqv
integer, intent(in) :: sfc_rlm
integer, intent(in) :: tc_pbl
integer, intent(in) :: kinver(:)
integer, intent(out) :: kpbl(:)
logical, intent(in) :: gen_tend,ldiag3d
!
real(kind=kind_phys), intent(in) :: grav,rd,cp,rv,hvap,hfus,fv, &
& eps,epsm1
real(kind=kind_phys), intent(in) :: delt, xkzm_m, xkzm_h, xkzm_s
real(kind=kind_phys), intent(in) :: dspfac, bl_upfr, bl_dnfr
real(kind=kind_phys), intent(in) :: rlmx, elmx
real(kind=kind_phys), intent(inout) :: dv(:,:), du(:,:), &
& tdt(:,:), rtg(:,:,:)
real(kind=kind_phys), intent(in) :: &
& u1(:,:), v1(:,:), &
& t1(:,:), q1(:,:,:), &
& swh(:,:), hlw(:,:), &
& xmu(:), garea(:), &
& zvfun(:), sigmaf(:), &
& psk(:), rbsoil(:), &
& zorl(:), tsea(:), &
& u10m(:), v10m(:), &
& fm(:), fh(:), &
& evap(:), heat(:), &
& stress(:), spd1(:), &
& prsi(:,:), del(:,:), &
& prsl(:,:), prslk(:,:), &
& phii(:,:), phil(:,:)
real(kind=kind_phys), intent(inout), dimension(:,:,:) :: dtend
integer, intent(in) :: dtidx(:,:), index_of_temperature, &
& index_of_x_wind, index_of_y_wind, index_of_process_pbl
real(kind=kind_phys), intent(out) :: &
& dusfc(:), dvsfc(:), &
& dtsfc(:), dqsfc(:), &
& hpbl(:)
real(kind=kind_phys), intent(out) :: &
& dkt(:,:), dku(:,:)
!
logical, intent(in) :: dspheat
character(len=*), intent(out) :: errmsg
integer, intent(out) :: errflg
! flag for tke dissipative heating
!
!----------------------------------------------------------------------
!***
!*** local variables
!***
integer i,is,k,n,ndt,km1,kmpbl,kmscu,ntrac1,idtend
integer kps,kbx,kmx
integer lcld(im),kcld(im),krad(im),mrad(im)
integer kx1(im), kb1(im), kpblx(im)
!
real(kind=kind_phys) tke(im,km), tkeh(im,km-1), e2(im,0:km)
!
real(kind=kind_phys) theta(im,km),thvx(im,km), thlvx(im,km),
& qlx(im,km), thetae(im,km),thlx(im,km),
& slx(im,km), svx(im,km), qtx(im,km),
& tvx(im,km), pix(im,km), radx(im,km-1),
& dkq(im,km-1),cku(im,km-1), ckt(im,km-1)
!
real(kind=kind_phys) plyr(im,km), rhly(im,km), cfly(im,km),
& qstl(im,km)
!
real(kind=kind_phys) dtdz1(im), gdx(im),
& phih(im), phim(im),
& phims(im), prn(im,km-1),
& rbdn(im), rbup(im), thermal(im),
& ustar(im), wstar(im), hpblx(im),
& ust3(im), wst3(im), rho_a(im),
& z0(im), crb(im), tkemean(im),
& hgamt(im), hgamq(im),
& wscale(im),vpert(im),
& zol(im), sflux(im),
& sumx(im), tx1(im), tx2(im)
!
real(kind=kind_phys) radmin(im)
!
real(kind=kind_phys) zi(im,km+1), zl(im,km), zm(im,km),
& xkzo(im,km-1),xkzmo(im,km-1),
& xkzm_hx(im), xkzm_mx(im), tkmnz(im,km-1),
& rdzt(im,km-1),rlmnz(im,km),
& al(im,km-1), ad(im,km), au(im,km-1),
& f1(im,km), f2(im,km*(ntrac-1))
!
real(kind=kind_phys) elm(im,km), ele(im,km),
& ckz(im,km), chz(im,km),
& diss(im,km-1),prod(im,km-1),
& bf(im,km-1), shr2(im,km-1), wush(im,km),
& xlamue(im,km-1), xlamde(im,km-1),
& gotvx(im,km), rlam(im,km-1)
!
! variables for updrafts (thermals)
!
real(kind=kind_phys) tcko(im,km), qcko(im,km,ntrac),
& ucko(im,km), vcko(im,km),
& buou(im,km), xmf(im,km)
!
! variables for stratocumulus-top induced downdrafts
!
real(kind=kind_phys) tcdo(im,km), qcdo(im,km,ntrac),
& ucdo(im,km), vcdo(im,km),
& buod(im,km), xmfd(im,km)
!
! variables for Total Variation Diminishing (TVD) flux-limiter scheme
!
real(kind=kind_phys) e_half(im,km-1), e_diff(im,0:km-1),
& q_half(im,km-1,ntrac-1),
& qh(im,km-1,ntrac-1),
& q_diff(im,0:km-1,ntrac-1)
real(kind=kind_phys) rrkp, phkp
real(kind=kind_phys) tsumn(im), tsump(im), rtnp(im)
real(kind=kind_phys) sfcpbl(im), vez0fun(im)
!
logical pblflg(im), sfcflg(im), flg(im)
logical scuflg(im), pcnvflg(im)
logical mlenflg
!
! pcnvflg: true for unstable pbl
!
real(kind=kind_phys) aphi16, aphi5,
& wfac, cfac,
& gamcrt, gamcrq, sfcfrac,
! & conq, cont, conw,
& dsdz2, dsdzt, dkmax,
& dsig, dt2, dtodsd,
& dtodsu, g, factor, dz,
& gocp, gravi, zol1, zolcru,
& buop, shrp, dtn,
& prnum, prmax, prmin, prtke,
& prscu, pr0, ri,
& dw2, dw2min, zk,
& elmfac, elefac, dspmax,
& alp, clwt, cql,
& f0, robn, crbmin, crbmax,
& es, qs, value, onemrh,
& cfh, gamma, elocp, el2orc,
& epsi, beta, chx, cqx,
& rdt, rdz, qmin, qlmin,
& rimin, rbcr, rbint, tdzmin,
& rlmn, rlmn0, rlmn1, rlmn2,
& ttend, utend, vtend, qtend,
& zfac, zfmin, vk, spdk2,
& tkmin, tkbmx, xkgdx,
& xkinv1, xkinv2,
& zlup, zldn, cs0, csmf,
& tem, tem1, tem2, tem3,
& ptem, ptem0, ptem1, ptem2
!
real(kind=kind_phys) slfac
!
real(kind=kind_phys) vegflo, vegfup, z0lo, z0up, vc0, zc0
!
real(kind=kind_phys) ck0, ck1, ch0, ch1, ce0, rchck
!
real(kind=kind_phys) qlcr, zstblmax, hcrinv
!
real(kind=kind_phys) h1
real(kind=kind_phys) bfac, mffac
!!
parameter(bfac=100.)
parameter(wfac=7.0,cfac=4.5)
parameter(gamcrt=3.,gamcrq=0.,sfcfrac=0.1)
parameter(vk=0.4,rimin=-100.,slfac=0.1)
parameter(rbcr=0.25,zolcru=-0.02,tdzmin=1.e-3)
parameter(rlmn=30.,rlmn0=5.,rlmn1=5.,rlmn2=10.)
parameter(prmin=0.25,prmax=4.0)
parameter(pr0=1.0,prtke=1.0,prscu=0.67)
parameter(f0=1.e-4,crbmin=0.15,crbmax=0.35)
parameter(tkmin=1.e-9,tkbmx=0.2,dspmax=10.0)
parameter(qmin=1.e-8,qlmin=1.e-12,zfmin=1.e-8)
parameter(aphi5=5.,aphi16=16.)
parameter(elmfac=1.0,elefac=1.0,cql=100.)
parameter(dw2min=1.e-4,dkmax=1000.,xkgdx=1000.)
parameter(qlcr=3.5e-5,zstblmax=2500.)
parameter(xkinv1=0.15,xkinv2=0.3)
parameter(h1=0.33333333,hcrinv=250.)
parameter(vegflo=0.1,vegfup=1.0,z0lo=0.1,z0up=1.0)
parameter(vc0=1.0,zc0=1.0)
parameter(ck1=0.15,ch1=0.15)
parameter(cs0=0.4,csmf=0.5)
parameter(rchck=1.5,ndt=20)
if (tc_pbl == 0) then
ck0 = 0.4
ch0 = 0.4
ce0 = 0.4
else if (tc_pbl == 1) then
ck0 = 0.55
ch0 = 0.55
ce0 = 0.12
endif
gravi = 1.0 / grav
g = grav
gocp = g / cp
! cont=cp/g
! conq=hvap/g
! conw=1.0/g ! for del in pa
!! parameter(cont=1000.*cp/g,conq=1000.*hvap/g,conw=1000./g) !kpa
elocp = hvap / cp
el2orc = hvap * hvap / (rv * cp)
!
!************************************************************************
! Initialize CCPP error handling variables
errmsg = ''
errflg = 0
!> ## Compute preliminary variables from input arguments
dt2 = delt
rdt = 1. / dt2
!
! the code is written assuming ntke=ntrac
! if ntrac > ntke, the code needs to be modified
!
ntrac1 = ntrac - 1
km1 = km - 1
kmpbl = km / 2
kmscu = km / 2
!> - Compute physical height of the layer centers and interfaces from
!! the geopotential height (\p zi and \p zl)
do k=1,km
do i=1,im
zi(i,k) = phii(i,k) * gravi
zl(i,k) = phil(i,k) * gravi
xmf(i,k) = 0.
xmfd(i,k) = 0.
buou(i,k) = 0.
buod(i,k) = 0.
wush(i,k) = 0.
ckz(i,k) = ck1
chz(i,k) = ch1
rlmnz(i,k) = rlmn0
enddo
enddo
do i=1,im
zi(i,km+1) = phii(i,km+1) * gravi
enddo
do k=1,km
do i=1,im
zm(i,k) = zi(i,k+1)
enddo
enddo
!> - Compute horizontal grid size (\p gdx)
do i=1,im
gdx(i) = sqrt(garea(i))
enddo
!> - Initialize tke value at vertical layer centers and interfaces
!! from tracer (\p tke and \p tkeh)
do k=1,km
do i=1,im
tke(i,k) = max(q1(i,k,ntke), tkmin)
enddo
enddo
do k=1,km1
do i=1,im
tkeh(i,k) = 0.5 * (tke(i,k) + tke(i,k+1))
enddo
enddo
!> - Compute reciprocal of \f$ \Delta z \f$ (rdzt)
do k = 1,km1
do i=1,im
rdzt(i,k) = 1.0 / (zl(i,k+1) - zl(i,k))
prn(i,k) = pr0
enddo
enddo
!
!> - Compute reciprocal of pressure (tx1, tx2)
!> - Compute minimum turbulent mixing length (rlmnz)
!> - Compute background vertical diffusivities for scalars and momentum (xkzo and xkzmo)
!> - set background diffusivities with xkzm_h & xkzm_m for gdx >= xkgdx and
!! as a function of horizontal grid size for gdx < xkgdx
!! \n xkzm_hx = xkzm_h * (gdx / xkgdx)
!! \n xkzm_mx = xkzm_m * (gdx / xkgdx)
!
do i=1,im
kx1(i) = 1
tx1(i) = 1.0 / prsi(i,1)
tx2(i) = tx1(i)
if(gdx(i) >= xkgdx) then
xkzm_hx(i) = xkzm_h
xkzm_mx(i) = xkzm_m
else
tem = gdx(i) / xkgdx
xkzm_hx(i) = xkzm_h * tem
xkzm_mx(i) = xkzm_m * tem
endif
enddo
do k = 1,km1
do i=1,im
xkzo(i,k) = 0.0
xkzmo(i,k) = 0.0
if (k < kinver(i)) then
! minimum turbulent mixing length
ptem = prsi(i,k+1) * tx1(i)
tem1 = 1.0 - ptem
tem2 = tem1 * tem1 * 2.5
tem2 = min(1.0, exp(-tem2))
rlmnz(i,k)= rlmn * tem2
rlmnz(i,k)= max(rlmnz(i,k), rlmn0)
! vertical background diffusivity
tem2 = tem1 * tem1 * 10.0
tem2 = min(1.0, exp(-tem2))
xkzo(i,k) = xkzm_hx(i) * tem2
! vertical background diffusivity for
! momentum
if (ptem >= xkzm_s) then
xkzmo(i,k) = xkzm_mx(i)
kx1(i) = k + 1
else
if (k == kx1(i) .and. k > 1) tx2(i) = 1.0 / prsi(i,k)
tem1 = 1.0 - prsi(i,k+1) * tx2(i)
tem1 = tem1 * tem1 * 5.0
xkzmo(i,k) = xkzm_mx(i) * min(1.0, exp(-tem1))
endif
endif
enddo
enddo
!
!> - Some output variables and logical flags are initialized
do i = 1,im
z0(i) = 0.01 * zorl(i)
rho_a(i) = prsl(i,1)/(rd*t1(i,1)*(1.+fv*max(q1(i,1,1),qmin)))
dusfc(i) = 0.
dvsfc(i) = 0.
dtsfc(i) = 0.
dqsfc(i) = 0.
kpbl(i) = 1
hpbl(i) = 0.
kpblx(i) = 1
hpblx(i) = 0.
pblflg(i)= .true.
sfcflg(i)= .true.
if(rbsoil(i) > 0.) sfcflg(i) = .false.
pcnvflg(i)= .false.
scuflg(i)= .true.
if(scuflg(i)) then
radmin(i)= 0.
mrad(i) = km1
krad(i) = 1
lcld(i) = km1
kcld(i) = km1
endif
enddo
!
! compute a function for green vegetation fraction and surface roughness
!
do i = 1,im
tem = (sigmaf(i) - vegflo) / (vegfup - vegflo)
tem = min(max(tem, 0.), 1.)
tem1 = sqrt(tem)
ptem = (z0(i) - z0lo) / (z0up - z0lo)
ptem = min(max(ptem, 0.), 1.)
vez0fun(i) = (1. + vc0 * tem1) * (1. + zc0 * ptem)
enddo
!
!> - Compute \f$\theta\f$(theta), and \f$q_l\f$(qlx), \f$\theta_e\f$(thetae),
!! \f$\theta_v\f$(thvx),\f$\theta_{l,v}\f$ (thlvx) including ice water
do k=1,km
do i=1,im
pix(i,k) = psk(i) / prslk(i,k)
theta(i,k) = t1(i,k) * pix(i,k)
if(ntiw > 0) then
tem = max(q1(i,k,ntcw),qlmin)
tem1 = max(q1(i,k,ntiw),qlmin)
qlx(i,k) = tem + tem1
ptem = hvap*tem + (hvap+hfus)*tem1
slx(i,k) = cp * t1(i,k) + phil(i,k) - ptem
else
qlx(i,k) = max(q1(i,k,ntcw),qlmin)
slx(i,k) = cp * t1(i,k) + phil(i,k) - hvap*qlx(i,k)
endif
tem2 = 1.+fv*max(q1(i,k,1),qmin)-qlx(i,k)
thvx(i,k) = theta(i,k) * tem2
tvx(i,k) = t1(i,k) * tem2
qtx(i,k) = max(q1(i,k,1),qmin)+qlx(i,k)
thlx(i,k) = theta(i,k) - pix(i,k)*elocp*qlx(i,k)
thlvx(i,k) = thlx(i,k) * (1. + fv * qtx(i,k))
svx(i,k) = cp * tvx(i,k)
ptem1 = elocp * pix(i,k) * max(q1(i,k,1),qmin)
thetae(i,k)= theta(i,k) + ptem1
! gotvx(i,k) = g / tvx(i,k)
gotvx(i,k) = g / thvx(i,k)
enddo
enddo
!> - Compute an empirical cloud fraction based on
!! Xu and Randall (1996) \cite xu_and_randall_1996
do k = 1, km
do i = 1, im
plyr(i,k) = 0.01 * prsl(i,k) ! pa to mb (hpa)
! --- ... compute relative humidity
es = 0.01 * fpvs(t1(i,k)) ! fpvs in pa
qs = max(qmin, eps * es / (plyr(i,k) + epsm1*es))
rhly(i,k) = max(0.0, min(1.0, max(qmin, q1(i,k,1))/qs))
qstl(i,k) = qs
enddo
enddo
!
do k = 1, km
do i = 1, im
cfly(i,k) = 0.
clwt = 1.0e-6 * (plyr(i,k)*0.001)
if (qlx(i,k) > clwt) then
onemrh = max(1.e-10, 1.0-rhly(i,k))
tem1 = min(max((onemrh*qstl(i,k))**0.49,0.0001),1.0)
tem1 = cql / tem1
value = max(min( tem1*qlx(i,k), 50.0), 0.0)
tem2 = sqrt(sqrt(rhly(i,k)))
cfly(i,k) = min(max(tem2*(1.0-exp(-value)), 0.0), 1.0)
endif
enddo
enddo
!
!> - Compute buoyancy modified by clouds
!
do k = 1, km1
do i = 1, im
tem = 0.5 * (svx(i,k) + svx(i,k+1))
tem1 = 0.5 * (t1(i,k) + t1(i,k+1))
tem2 = 0.5 * (qstl(i,k) + qstl(i,k+1))
cfh = min(cfly(i,k+1),0.5*(cfly(i,k)+cfly(i,k+1)))
alp = g / tem
gamma = el2orc * tem2 / (tem1**2)
epsi = tem1 / elocp
beta = (1. + gamma*epsi*(1.+fv)) / (1. + gamma)
chx = cfh * alp * beta + (1. - cfh) * alp
cqx = cfh * alp * hvap * (beta - epsi)
cqx = cqx + (1. - cfh) * fv * g
ptem1 = (slx(i,k+1)-slx(i,k))*rdzt(i,k)
ptem2 = (qtx(i,k+1)-qtx(i,k))*rdzt(i,k)
bf(i,k) = chx * ptem1 + cqx * ptem2
enddo
enddo
!
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!
!> - Initialize diffusion coefficients to 0 and calculate the total
!! radiative heating rate (dku, dkt, radx)
do k=1,km
do i=1,im
dku(i,k) = 0.
dkt(i,k) = 0.
enddo
enddo
do k=1,km1
do i=1,im
dkq(i,k) = 0.
cku(i,k) = 0.
ckt(i,k) = 0.
tem = zi(i,k+1)-zi(i,k)
radx(i,k) = tem*(swh(i,k)*xmu(i)+hlw(i,k))
enddo
enddo
!> - Compute stable/unstable PBL flag (pblflg) based on the total
!! surface energy flux (\e false if the total surface energy flux
!! is into the surface)
do i = 1,im
sflux(i) = heat(i) + evap(i)*fv*theta(i,1)
if(.not.sfcflg(i) .or. sflux(i) <= 0.) pblflg(i)=.false.
enddo
!
!> ## Calculate the PBL height
!! The calculation of the boundary layer height follows Troen and Mahrt (1986) \cite troen_and_mahrt_1986 section 3. The approach is to find the level in the column where a modified bulk Richardson number exceeds a critical value.
!! - Compute critical bulk Richardson number (\f$Rb_{cr}\f$) (crb)
!! - For the unstable PBL, crb is a constant (0.25)
!! - For the stable boundary layer (SBL), \f$Rb_{cr}\f$ varies
!! with the surface Rossby number, \f$R_{0}\f$, as given by
!! Vickers and Mahrt (2004) \cite Vickers_2004
!! \f[
!! Rb_{cr}=0.16(10^{-7}R_{0})^{-0.18}
!! \f]
!! \f[
!! R_{0}=\frac{U_{10}}{f_{0}z_{0}}
!! \f]
!! where \f$U_{10}\f$ is the wind speed at 10m above the ground surface,
!! \f$f_0\f$ the Coriolis parameter, and \f$z_{0}\f$ the surface roughness
!! length. To avoid too much variation, we restrict \f$Rb_{cr}\f$ to vary
!! within the range of 0.15~0.35
do i = 1,im
if(pblflg(i)) then
! thermal(i) = thvx(i,1)
thermal(i) = thlvx(i,1)
crb(i) = rbcr
else
thermal(i) = tsea(i)*(1.+fv*max(q1(i,1,1),qmin))
tem = sqrt(u10m(i)**2+v10m(i)**2)
tem = max(tem, 1.)
robn = tem / (f0 * z0(i))
tem1 = 1.e-7 * robn
crb(i) = 0.16 * (tem1 ** (-0.18))
crb(i) = max(min(crb(i), crbmax), crbmin)
endif
enddo
!> - Compute \f$\frac{\Delta t}{\Delta z}\f$ , \f$u_*\f$
do i=1,im
dtdz1(i) = dt2 / (zi(i,2)-zi(i,1))
enddo
!
do i=1,im
ustar(i) = sqrt(stress(i))
enddo
!
!> - Compute buoyancy \f$\frac{\partial \theta_v}{\partial z}\f$ (bf)
!! and the wind shear squared (shr2)
!
do k = 1, km1
do i = 1, im
rdz = rdzt(i,k)
! bf(i,k) = gotvx(i,k)*(thvx(i,k+1)-thvx(i,k))*rdz
dw2 = (u1(i,k)-u1(i,k+1))**2
& + (v1(i,k)-v1(i,k+1))**2
shr2(i,k) = max(dw2,dw2min)*rdz*rdz
enddo
enddo
!
! Find first quess pbl height based on bulk richardson number (mrf pbl scheme)
! and also for diagnostic purpose
!
do i=1,im
flg(i) = .false.
rbup(i) = rbsoil(i)
enddo
!> - Given the thermal's properties and the critical Richardson number,
!! a loop is executed to find the first level above the surface (kpblx) where
!! the modified Richardson number is greater than the critical Richardson
!! number, using equation 10a from Troen and Mahrt (1996) \cite troen_and_mahrt_1986
!! (also equation 8 from Hong and Pan (1996) \cite hong_and_pan_1996):
do k = 1, kmpbl
do i = 1, im
if(.not.flg(i)) then
rbdn(i) = rbup(i)
if (tc_pbl == 0) then
spdk2 = max((u1(i,k)**2+v1(i,k)**2),1.)
! rbup(i) = (thvx(i,k)-thermal(i))*
! & (g*zl(i,k)/thvx(i,1))/spdk2
rbup(i) = (thlvx(i,k)-thermal(i))*
& (g*zl(i,k)/thlvx(i,1))/spdk2
else if (tc_pbl == 1) then
spdk2 = max((u1(i,k)-u1(i,1))**2+(v1(i,k)-v1(i,1))**2,1.)
& + bfac*ustar(i)**2
rbup(i) = (thlvx(i,k)-thermal(i))*
& (g*(zl(i,k)-zl(i,1))/thlvx(i,1))/spdk2
endif
kpblx(i) = k
flg(i) = rbup(i) > crb(i)
endif
enddo
enddo
!> - Once the level is found, some linear interpolation is performed to find
!! the exact height of the boundary layer top (where \f$R_{i} > Rb_{cr}\f$)
!! and the PBL height (hpbl and kpbl) and the PBL top index are saved.
do i = 1,im
if(kpblx(i) > 1) then
k = kpblx(i)
if(rbdn(i) >= crb(i)) then
rbint = 0.
elseif(rbup(i) <= crb(i)) then
rbint = 1.
else
rbint = (crb(i)-rbdn(i))/(rbup(i)-rbdn(i))
endif
hpblx(i) = zl(i,k-1) + rbint*(zl(i,k)-zl(i,k-1))
if(hpblx(i) < zi(i,kpblx(i))) kpblx(i)=kpblx(i)-1
else
hpblx(i) = zl(i,1)
kpblx(i) = 1
endif
hpbl(i) = hpblx(i)
kpbl(i) = kpblx(i)
if(kpbl(i) <= 1) pblflg(i)=.false.
enddo
!
! update thermal at a level of slfac*hpbl for unstable pbl
!
do i=1,im
sfcpbl(i) = slfac * hpbl(i)
kb1(i) = 1
flg(i) = .false.
if(pblflg(i)) then
flg(i) = .true.
endif
enddo
do k = 2, kmpbl
do i=1,im
if (flg(i) .and. zl(i,k) <= sfcpbl(i)) then
kb1(i) = k
else
flg(i) = .false.
endif
enddo
enddo
do i=1,im
if(pblflg(i)) kb1(i)=min(kb1(i),kpbl(i))
enddo
!
! re-compute pbl height with the updated thermal
!
do i=1,im
flg(i) = .true.
if(pblflg(i) .and. kb1(i) > 1) then
flg(i) = .false.
rbup(i) = rbsoil(i)
! thermal(i) = thvx(i,kb1(i))
thermal(i) = thlvx(i,kb1(i))
kpblx(i) = kb1(i)
hpblx(i) = zl(i,kb1(i))
endif
enddo
do k = 2, kmpbl
do i = 1, im
if(.not.flg(i) .and. k > kb1(i)) then
rbdn(i) = rbup(i)
if (tc_pbl == 0) then
spdk2 = max((u1(i,k)**2+v1(i,k)**2),1.)
! rbup(i) = (thvx(i,k)-thermal(i))*
! & (g*zl(i,k)/thvx(i,1))/spdk2
rbup(i) = (thlvx(i,k)-thermal(i))*
& (g*zl(i,k)/thlvx(i,1))/spdk2
else if (tc_pbl == 1) then
spdk2 = max((u1(i,k)-u1(i,1))**2+(v1(i,k)-v1(i,1))**2,1.)
& + bfac*ustar(i)**2
rbup(i) = (thlvx(i,k)-thermal(i))*
& (g*(zl(i,k)-zl(i,1))/thlvx(i,1))/spdk2
endif
kpblx(i) = k
flg(i) = rbup(i) > crb(i)
endif
enddo
enddo
do i = 1,im
if(pblflg(i) .and. kb1(i) > 1) then
k = kpblx(i)
if(rbdn(i) >= crb(i)) then
rbint = 0.
elseif(rbup(i) <= crb(i)) then
rbint = 1.
else
rbint = (crb(i)-rbdn(i))/(rbup(i)-rbdn(i))
endif
hpblx(i) = zl(i,k-1) + rbint*(zl(i,k)-zl(i,k-1))
if(hpblx(i) < zi(i,kpblx(i))) kpblx(i)=kpblx(i)-1
hpbl(i) = hpblx(i)
kpbl(i) = kpblx(i)
endif
enddo
!
! compute mean tke within pbl
!
do i = 1, im
sumx(i) = 0.
tkemean(i) = 0.
enddo
do k = 1, kmpbl
do i = 1, im
if(k < kpbl(i)) then
dz = zi(i,k+1) - zi(i,k)
tkemean(i) = tkemean(i) + tke(i,k) * dz
sumx(i) = sumx(i) + dz
endif
enddo
enddo
do i = 1, im
if(tkemean(i) > 0. .and. sumx(i) > 0.) then
tkemean(i) = tkemean(i) / sumx(i)
endif
enddo
!
! compute wind shear term as a sink term for updraft and downdraft
! velocity
!
kps = max(kmpbl, kmscu)
do k = 2, kps
do i = 1, im
dz = zi(i,k+1) - zi(i,k)
tem = (0.5*(u1(i,k-1)-u1(i,k+1))/dz)**2
tem1 = tem+(0.5*(v1(i,k-1)-v1(i,k+1))/dz)**2
wush(i,k) = csmf * sqrt(tem1)
enddo
enddo
!
!> ## Compute Monin-Obukhov similarity parameters
!! - Calculate the Monin-Obukhov nondimensional stability paramter, commonly
!! referred to as \f$\zeta\f$ using the following equation from Businger et al.(1971) \cite businger_et_al_1971
!! (eqn 28):
!! \f[
!! \zeta = Ri_{sfc}\frac{F_m^2}{F_h} = \frac{z}{L}
!! \f]
!! where \f$F_m\f$ and \f$F_h\f$ are surface Monin-Obukhov stability functions calculated in sfc_diff.f and
!! \f$L\f$ is the Obukhov length.
do i=1,im
zol(i) = max(rbsoil(i)*fm(i)*fm(i)/fh(i),rimin)
if(sfcflg(i)) then
zol(i) = min(zol(i),-zfmin)
else
zol(i) = max(zol(i),zfmin)
endif
!> - Calculate the nondimensional gradients of momentum and temperature (\f$\phi_m\f$ (phim) and \f$\phi_h\f$(phih)) are calculated using
!! eqns 5 and 6 from Hong and Pan (1996) \cite hong_and_pan_1996 depending on the surface layer stability:
!! - For the unstable and neutral conditions:
!! \f[
!! \phi_m=(1-16\frac{0.1h}{L})^{-1/4}
!! \phi_h=(1-16\frac{0.1h}{L})^{-1/2}
!! \f]
!! - For the stable regime
!! \f[
!! \phi_m=\phi_t=(1+5\frac{0.1h}{L})
!! \f]
zol1 = zol(i)*sfcfrac*hpbl(i)/zl(i,1)
if(sfcflg(i)) then
tem = 1.0 / (1. - aphi16*zol1)
phih(i) = sqrt(tem)
phim(i) = sqrt(phih(i))
tem1 = 1.0 / (1. - aphi16*zol(i))
phims(i) = sqrt(sqrt(tem1))
else
phim(i) = 1. + aphi5*zol1
phih(i) = phim(i)
phims(i) = 1. + aphi5*zol(i)
endif
enddo
!
!> - The \f$z/L\f$ (zol) is used as the stability criterion for the PBL.Currently,
!! strong unstable (convective) PBL for \f$z/L < -0.02\f$ and weakly and moderately
!! unstable PBL for \f$0>z/L>-0.02\f$
!> - Compute the velocity scale \f$w_s\f$ (wscale) (eqn 22 of Han et al. 2019). It
!! is represented by the value scaled at the top of the surface layer:
!! \f[
!! w_s=(u_*^3+7\alpha\kappa w_*^3)^{1/3}
!! \f]
!! where \f$u_*\f$ (ustar) is the surface friction velocity,\f$\alpha\f$ is the ratio
!! of the surface layer height to the PBL height (specified as sfcfrac =0.1),
!! \f$\kappa =0.4\f$ is the von Karman constant, and \f$w_*\f$ is the convective velocity
!! scale defined as eqn23 of Han et al.(2019):
!! \f[
!! w_{*}=[(g/T)\overline{(w'\theta_v^{'})}_0h]^{1/3}
!! \f]
do i=1,im
if(pblflg(i)) then
if(zol(i) < zolcru) then
pcnvflg(i) = .true.
endif
wst3(i) = gotvx(i,1)*sflux(i)*hpbl(i)
wstar(i)= wst3(i)**h1
ust3(i) = ustar(i)**3.
wscale(i)=(ust3(i)+wfac*vk*wst3(i)*sfcfrac)**h1
ptem = ustar(i)/aphi5
wscale(i) = max(wscale(i),ptem)
endif
enddo
!
!> ## The counter-gradient terms for temperature and humidity are calculated.
!! - Equation 4 of Hong and Pan (1996) \cite hong_and_pan_1996 and are used to calculate the "scaled virtual temperature excess near the surface" (equation 9 in Hong and Pan (1996) \cite hong_and_pan_1996) for use in the mass-flux algorithm.
!
do i = 1,im
if(pcnvflg(i)) then
hgamt(i) = heat(i)/wscale(i)
hgamq(i) = evap(i)/wscale(i)
vpert(i) = hgamt(i) + hgamq(i)*fv*theta(i,1)
vpert(i) = max(vpert(i),0.)
tem = min(cfac*vpert(i),gamcrt)
thermal(i)= thermal(i) + tem
endif
enddo
!
! enhance the pbl height by considering the thermal excess
! (overshoot pbl top)
!
do i=1,im
flg(i) = .true.
if(pcnvflg(i)) then
flg(i) = .false.
rbup(i) = rbsoil(i)
endif
enddo
do k = 2, kmpbl
do i = 1, im
if(.not.flg(i) .and. k > kb1(i)) then
rbdn(i) = rbup(i)
if (tc_pbl == 0) then
spdk2 = max((u1(i,k)**2+v1(i,k)**2),1.)
rbup(i) = (thlvx(i,k)-thermal(i))*
& (g*zl(i,k)/thlvx(i,1))/spdk2
else if (tc_pbl == 1) then
spdk2 = max((u1(i,k)-u1(i,1))**2+(v1(i,k)-v1(i,1))**2,1.)
& + bfac*ustar(i)**2
rbup(i) = (thlvx(i,k)-thermal(i))*
& (g*(zl(i,k)-zl(i,1))/thlvx(i,1))/spdk2
endif
kpbl(i) = k
flg(i) = rbup(i) > crb(i)
endif
enddo
enddo
do i = 1,im
if(pcnvflg(i)) then
k = kpbl(i)
if(rbdn(i) >= crb(i)) then
rbint = 0.
elseif(rbup(i) <= crb(i)) then
rbint = 1.
else
rbint = (crb(i)-rbdn(i))/(rbup(i)-rbdn(i))
endif
hpbl(i) = zl(i,k-1) + rbint*(zl(i,k)-zl(i,k-1))
if(hpbl(i) < zi(i,kpbl(i))) then
kpbl(i) = kpbl(i) - 1
endif
if(kpbl(i) <= 1) then
pcnvflg(i) = .false.
pblflg(i) = .false.
endif
endif
enddo
!
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! look for stratocumulus
!> ## Determine whether stratocumulus layers exist and compute quantities
!! - Starting at the PBL top and going downward, if the level is less than 2.5 km
!! and \f$q_l\geq q_{lcr}\f$ then set kcld = k (find the cloud top index in the PBL.
!! If no cloud water above the threshold is hound, \e scuflg is set to F.
do i=1,im
flg(i) = scuflg(i)
enddo
do k = 1, km1
do i=1,im
if(flg(i).and.zl(i,k) >= zstblmax) then
lcld(i)=k
flg(i)=.false.
endif
enddo
enddo
do i = 1, im
flg(i)=scuflg(i)
enddo
do k = kmscu,1,-1
do i = 1, im
if(flg(i) .and. k <= lcld(i)) then
if(qlx(i,k) >= qlcr) then
kcld(i)=k
flg(i)=.false.
endif
endif
enddo
enddo
do i = 1, im
if(scuflg(i) .and. kcld(i)==km1) scuflg(i)=.false.
enddo
!> - Starting at the PBL top and going downward, if the level is less
!! than the cloud top, find the level of the minimum radiative heating
!! rate wihin the cloud. If the level of the minimum is the lowest model
!! level or the minimum radiative heating rate is positive, then set
!! scuflg to F.
do i = 1, im
flg(i)=scuflg(i)
enddo
do k = kmscu,1,-1
do i = 1, im
if(flg(i) .and. k <= kcld(i)) then
if(qlx(i,k) >= qlcr) then
if(radx(i,k) < radmin(i)) then
radmin(i)=radx(i,k)
krad(i)=k
endif
else
flg(i)=.false.
endif
endif
enddo
enddo
do i = 1, im
if(scuflg(i) .and. krad(i) <= 1) scuflg(i)=.false.
if(scuflg(i) .and. radmin(i)>=0.) scuflg(i)=.false.
enddo
!
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!> ## Compute components for mass flux mixing by large thermals
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!> - If the PBL is convective, the updraft properties are initialized
!! to be the same as the state variables.
do k = 1, km
do i = 1, im
if(pcnvflg(i)) then
tcko(i,k) = t1(i,k)
ucko(i,k) = u1(i,k)
vcko(i,k) = v1(i,k)
endif
if(scuflg(i)) then
tcdo(i,k) = t1(i,k)
ucdo(i,k) = u1(i,k)
vcdo(i,k) = v1(i,k)
endif
enddo
enddo
do n = 1, ntrac1
do k = 1, km
do i = 1, im
if(pcnvflg(i)) then
qcko(i,k,n) = q1(i,k,n)
endif