forked from wrf-model/WRF
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathmodule_bl_temf.F
1459 lines (1371 loc) · 64 KB
/
module_bl_temf.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
!wrf:model_layer:physics
!
!
!
!
module module_bl_temf
contains
!
!-------------------------------------------------------------------
!
subroutine temfpbl(u3d,v3d,th3d,t3d,qv3d,qc3d,qi3d,p3d,p3di,pi3d,rho, &
rublten,rvblten,rthblten, &
rqvblten,rqcblten,rqiblten,flag_qi, &
g,cp,rcp,r_d,r_v,cpv, &
z,xlv,psfc, &
p_top, &
znt,ht,ust,zol,hol,hpbl,psim,psih, &
xland,hfx,qfx,tsk,qsfc,gz1oz0,wspd,br, &
dt,dtmin,kpbl2d, &
svp1,svp2,svp3,svpt0,ep1,ep2,karman,eomeg,stbolt, &
kh_temf,km_temf, &
u10,v10,t2, &
te_temf,shf_temf,qf_temf,uw_temf,vw_temf, &
wupd_temf,mf_temf,thup_temf,qtup_temf,qlup_temf, &
cf3d_temf,cfm_temf, &
hd_temf,lcl_temf,hct_temf, &
flhc,flqc,exch_temf, &
fCor, &
ids,ide, jds,jde, kds,kde, &
ims,ime, jms,jme, kms,kme, &
its,ite, jts,jte, kts,kte &
)
!-------------------------------------------------------------------
implicit none
!-------------------------------------------------------------------
! New variables for TEMF
!-- te_temf Total energy from this scheme
!-- shf_temf Sensible heat flux profile from this scheme (kinematic)
!-- qf_temf Moisture flux profile from this scheme (kinematic)
!-- uw_temf U momentum flux component from this scheme
!-- vw_temf V momentum flux component from this scheme
!-- kh_temf Exchange coefficient for heat (3D)
!-- km_temf Exchange coefficient for momentum (3D)
!-- wupd_temf Updraft velocity from TEMF BL scheme
!-- mf_temf Mass flux from TEMF BL scheme
!-- thup_temf Updraft thetal from TEMF BL scheme
!-- qtup_temf Updraft qt from TEMF BL scheme
!-- qlup_temf Updraft ql from TEMF BL scheme
!-- cf3d_temf 3D cloud fraction from TEMF BL scheme
!-- cfm_temf Column cloud fraction from TEMF BL scheme
!-- exch_temf Surface exchange coefficient (as for moisture) from TEMF surface layer scheme
!-- flhc Surface exchange coefficient for heat (needed by surface scheme)
!-- flqc Surface exchange coefficient for moisture (including moisture availablity)
!-- fCor Coriolis parameter (from grid%f)
!
!-- u3d 3d u-velocity interpolated to theta points (m/s)
!-- v3d 3d v-velocity interpolated to theta points (m/s)
!-- th3d 3d potential temperature (k)
!-- t3d temperature (k)
!-- qv3d 3d water vapor mixing ratio (kg/kg)
!-- qc3d 3d cloud mixing ratio (kg/kg)
!-- qi3d 3d ice mixing ratio (kg/kg)
! (note: if P_QI<PARAM_FIRST_SCALAR this should be zero filled)
!-- p3d 3d pressure (pa)
!-- p3di 3d pressure (pa) at interface level
!-- pi3d 3d exner function (dimensionless)
!-- rho 3d dry air density (kg/m^3)
!-- rublten u tendency due to
! pbl parameterization (m/s/s)
!-- rvblten v tendency due to
! pbl parameterization (m/s/s)
!-- rthblten theta tendency due to
! pbl parameterization (K/s)
!-- rqvblten qv tendency due to
! pbl parameterization (kg/kg/s)
!-- rqcblten qc tendency due to
! pbl parameterization (kg/kg/s)
!-- rqiblten qi tendency due to
! pbl parameterization (kg/kg/s)
!-- cp heat capacity at constant pressure for dry air (j/kg/k)
!-- g acceleration due to gravity (m/s^2)
!-- rovcp r/cp
!-- r_d gas constant for dry air (j/kg/k)
!-- rovg r/g
!-- z height above sea level (m)
!-- xlv latent heat of vaporization (j/kg)
!-- r_v gas constant for water vapor (j/kg/k)
!-- psfc pressure at the surface (pa)
!-- znt roughness length (m)
!-- ht terrain height ASL (m)
!-- ust u* in similarity theory (m/s)
!-- zol z/l height over monin-obukhov length
!-- hol pbl height over monin-obukhov length
!-- hpbl pbl height (m)
!-- psim similarity stability function for momentum
!-- psih similarity stability function for heat
!-- xland land mask (1 for land, 2 for water)
!-- hfx upward heat flux at the surface (w/m^2)
!-- qfx upward moisture flux at the surface (kg/m^2/s)
!-- tsk surface temperature (k)
!-- qsfc surface specific humidity (kg/kg)
!-- gz1oz0 log(z/z0) where z0 is roughness length
!-- wspd wind speed at lowest model level (m/s)
!-- u10 u-wind speed at 10 m (m/s)
!-- v10 v-wind speed at 10 m (m/s)
!-- br bulk richardson number in surface layer
!-- dt time step (s)
!-- dtmin time step (minute)
!-- rvovrd r_v divided by r_d (dimensionless)
!-- svp1 constant for saturation vapor pressure (kpa)
!-- svp2 constant for saturation vapor pressure (dimensionless)
!-- svp3 constant for saturation vapor pressure (k)
!-- svpt0 constant for saturation vapor pressure (k)
!-- ep1 constant for virtual temperature (r_v/r_d - 1) (dimensionless)
!-- ep2 constant for specific humidity calculation
!-- karman von karman constant
!-- eomeg angular velocity of earths rotation (rad/s)
!-- stbolt stefan-boltzmann constant (w/m^2/k^4)
!-- ids start index for i in domain
!-- ide end index for i in domain
!-- jds start index for j in domain
!-- jde end index for j in domain
!-- kds start index for k in domain
!-- kde end index for k in domain
!-- ims start index for i in memory
!-- ime end index for i in memory
!-- jms start index for j in memory
!-- jme end index for j in memory
!-- kms start index for k in memory
!-- kme end index for k in memory
!-- its start index for i in tile
!-- ite end index for i in tile
!-- jts start index for j in tile
!-- jte end index for j in tile
!-- kts start index for k in tile
!-- kte end index for k in tile
!-------------------------------------------------------------------
! Arguments
!
integer, intent(in ) :: ids,ide, jds,jde, kds,kde, &
ims,ime, jms,jme, kms,kme, &
its,ite, jts,jte, kts,kte
!
real, intent(in ) :: dt,dtmin,g,cp,rcp,r_d,r_v,xlv,cpv
!
real, intent(in ) :: svp1,svp2,svp3,svpt0
real, intent(in ) :: ep1,ep2,karman,eomeg,stbolt
!
real, dimension( ims:ime, kms:kme, jms:jme ) , &
intent(in ) :: qv3d, qc3d, qi3d, &
p3d, pi3d, th3d, t3d, &
z, rho
!
real, dimension( ims:ime, kms:kme, jms:jme ) , &
intent(inout) :: te_temf
real, dimension( ims:ime, kms:kme, jms:jme ) , &
intent( out) :: shf_temf, qf_temf, uw_temf, vw_temf , &
wupd_temf, mf_temf, thup_temf, qtup_temf, &
qlup_temf,cf3d_temf
real, dimension( ims:ime, jms:jme ) , &
intent(inout) :: flhc, flqc, exch_temf
real, dimension( ims:ime, jms:jme ) , &
intent(in ) :: fCor
real, dimension( ims:ime, jms:jme ) , &
intent( out) :: hd_temf, lcl_temf, hct_temf, cfm_temf
!
real, dimension( ims:ime, kms:kme, jms:jme ) , &
intent(in ) :: p3di
!
real, dimension( ims:ime, kms:kme, jms:jme ) , &
intent(inout) :: rublten, rvblten, &
rthblten, &
rqvblten, rqcblten, rqiblten
!
real, dimension( ims:ime, kms:kme, jms:jme ) , &
intent(inout) :: kh_temf, km_temf
real, dimension( ims:ime, jms:jme ) , &
intent(inout) :: u10, v10, t2
!
real, dimension( ims:ime, jms:jme ) , &
intent(in ) :: xland, &
psim, psih, gz1oz0, br, &
psfc, tsk, qsfc
!
real, dimension( ims:ime, jms:jme ) , &
intent(inout) :: hfx, qfx
real, dimension( ims:ime, jms:jme ) , &
intent(inout) :: hol, ust, hpbl, znt, wspd, zol
real, dimension( ims:ime, jms:jme ) , &
intent(in ) :: ht
!
real, dimension( ims:ime, kms:kme, jms:jme ) , &
intent(in ) :: u3d, v3d
!
integer, dimension( ims:ime, jms:jme ) , &
intent(out ) :: kpbl2d
!
logical, intent(in) :: flag_qi
!
! real, dimension( ims:ime, kms:kme, jms:jme ), &
! optional , &
! intent(inout) :: rqiblten
!
real, optional, intent(in ) :: p_top
!
!-------------------------------------------------------
! Local variables
integer :: j
do j = jts,jte
call temf2d(J=j,ux=u3d(ims,kms,j),vx=v3d(ims,kms,j) &
,tx=t3d(ims,kms,j),thx=th3d(ims,kms,j) &
,qvx=qv3d(ims,kms,j),qcx=qc3d(ims,kms,j) &
,qix=qi3d(ims,kms,j) &
,p2d=p3d(ims,kms,j),p2di=p3di(ims,kms,j) &
,pi2d=pi3d(ims,kms,j),rho=rho(ims,kms,j) &
,rubltenx=rublten(ims,kms,j),rvbltenx=rvblten(ims,kms,j) &
,rthbltenx=rthblten(ims,kms,j),rqvbltenx=rqvblten(ims,kms,j) &
,rqcbltenx=rqcblten(ims,kms,j),rqibltenx=rqiblten(ims,kms,j) &
,g=g,cp=cp,rcp=rcp,r_d=r_d,r_v=r_v,cpv=cpv &
,z2d=z(ims,kms,j) &
,xlv=xlv &
,psfcpa=psfc(ims,j),znt=znt(ims,j),zsrf=ht(ims,j),ust=ust(ims,j) &
,zol=zol(ims,j),hol=hol(ims,j),hpbl=hpbl(ims,j) &
,psim=psim(ims,j) &
,psih=psih(ims,j),xland=xland(ims,j) &
,hfx=hfx(ims,j),qfx=qfx(ims,j) &
,tsk=tsk(ims,j),qsfc=qsfc(ims,j),gz1oz0=gz1oz0(ims,j) &
,wspd=wspd(ims,j),br=br(ims,j) &
,dt=dt,dtmin=dtmin,kpbl1d=kpbl2d(ims,j) &
,svp1=svp1,svp2=svp2,svp3=svp3,svpt0=svpt0 &
,ep1=ep1,ep2=ep2,karman=karman,eomeg=eomeg &
,stbolt=stbolt &
,kh_temfx=kh_temf(ims,kms,j),km_temfx=km_temf(ims,kms,j) &
,u10=u10(ims,j),v10=v10(ims,j),t2=t2(ims,j) &
,te_temfx=te_temf(ims,kms,j) &
,shf_temfx=shf_temf(ims,kms,j),qf_temfx=qf_temf(ims,kms,j) &
,uw_temfx=uw_temf(ims,kms,j),vw_temfx=vw_temf(ims,kms,j) &
,wupd_temfx=wupd_temf(ims,kms,j),mf_temfx=mf_temf(ims,kms,j) &
,thup_temfx=thup_temf(ims,kms,j),qtup_temfx=qtup_temf(ims,kms,j) &
,qlup_temfx=qlup_temf(ims,kms,j) &
,cf3d_temfx=cf3d_temf(ims,kms,j),cfm_temfx=cfm_temf(ims,j) &
,hd_temfx=hd_temf(ims,j),lcl_temfx=lcl_temf(ims,j) &
,hct_temfx=hct_temf(ims,j),exch_temfx=exch_temf(ims,j) &
,flhc=flhc(ims,j),flqc=flqc(ims,j) &
,fCor=fCor(ims,j) &
,ids=ids,ide=ide, jds=jds,jde=jde, kds=kds,kde=kde &
,ims=ims,ime=ime, jms=jms,jme=jme, kms=kms,kme=kme &
,its=its,ite=ite, jts=jts,jte=jte, kts=kts,kte=kte )
enddo
!
end subroutine temfpbl
!
!-------------------------------------------------------------------
!
subroutine temf2d(j,ux,vx,tx,thx,qvx,qcx,qix,p2d,p2di,pi2d,rho, &
rubltenx,rvbltenx,rthbltenx, &
rqvbltenx,rqcbltenx,rqibltenx, &
g,cp,rcp,r_d,r_v,cpv, &
z2d, &
xlv,psfcpa, &
znt,zsrf,ust,zol,hol,hpbl,psim,psih, &
xland,hfx,qfx,tsk,qsfc,gz1oz0,wspd,br, &
dt,dtmin,kpbl1d, &
svp1,svp2,svp3,svpt0,ep1,ep2,karman,eomeg,stbolt, &
kh_temfx,km_temfx, &
u10,v10,t2, &
te_temfx,shf_temfx,qf_temfx,uw_temfx,vw_temfx, &
wupd_temfx,mf_temfx,thup_temfx,qtup_temfx,qlup_temfx, &
cf3d_temfx,cfm_temfx, &
hd_temfx,lcl_temfx,hct_temfx,exch_temfx, &
flhc,flqc, &
fCor, &
ids,ide, jds,jde, kds,kde, &
ims,ime, jms,jme, kms,kme, &
its,ite, jts,jte, kts,kte &
)
!-------------------------------------------------------------------
implicit none
!-------------------------------------------------------------------
!
! This is the Total Energy - Mass Flux (TEMF) PBL scheme.
! Initial implementation 2010 by Wayne Angevine, CIRES/NOAA ESRL.
! References:
! Angevine et al., 2010, MWR
! Angevine, 2005, JAM
! Mauritsen et al., 2007, JAS
!
!-------------------------------------------------------------------
!
integer, intent(in ) :: ids,ide, jds,jde, kds,kde, &
ims,ime, jms,jme, kms,kme, &
its,ite, jts,jte, kts,kte, j
!
real, intent(in ) :: dt,dtmin,g,cp,rcp,r_d,r_v,cpv,xlv
!
real, intent(in ) :: svp1,svp2,svp3,svpt0
real, intent(in ) :: ep1,ep2,karman,eomeg,stbolt
!
real, dimension( ims:ime, kms:kme ), &
intent(in) :: z2d
!
real, dimension( ims:ime, kms:kme ) , &
intent(in ) :: ux, vx
real, dimension( ims:ime, kms:kme ) , &
intent(inout) :: te_temfx
real, dimension( ims:ime, kms:kme ) , &
intent( out) :: shf_temfx, qf_temfx, uw_temfx, vw_temfx , &
wupd_temfx, mf_temfx,thup_temfx, &
qtup_temfx, qlup_temfx, cf3d_temfx
real, dimension( ims:ime ) , &
intent( out) :: hd_temfx, lcl_temfx, hct_temfx, cfm_temfx
real, dimension( ims:ime ) , &
intent(in ) :: fCor
real, dimension( ims:ime ) , &
intent(inout) :: flhc, flqc, exch_temfx
real, dimension( ims:ime, kms:kme ) , &
intent(in ) :: tx, thx, qvx, qcx, qix, pi2d, rho
real, dimension( ims:ime, kms:kme ) , &
intent(in ) :: p2di, p2d
!
real, dimension( ims:ime, kms:kme ) , &
intent(inout) :: rubltenx, rvbltenx, rthbltenx, &
rqvbltenx, rqcbltenx, rqibltenx
!
real, dimension( ims:ime ) , &
intent(inout) :: hol, ust, hpbl, znt
real, dimension( ims:ime ) , &
intent(in ) :: xland, zsrf
real, dimension( ims:ime ) , &
intent(inout) :: hfx, qfx
!
real, dimension( ims:ime ), intent(inout) :: wspd
real, dimension( ims:ime ), intent(in ) :: br
!
real, dimension( ims:ime ), intent(in ) :: psim, psih
real, dimension( ims:ime ), intent(in ) :: gz1oz0
!
real, dimension( ims:ime ), intent(in ) :: psfcpa
real, dimension( ims:ime ), intent(in ) :: tsk, qsfc
real, dimension( ims:ime ), intent(inout) :: zol
integer, dimension( ims:ime ), intent(out ) :: kpbl1d
real, dimension( ims:ime, kms:kme ) , &
intent(inout) :: kh_temfx, km_temfx
!
real, dimension( ims:ime ) , &
intent(inout) :: u10, v10, t2
!
!
!-----------------------------------------------------------
! Local variables
!
! TE model constants
logical, parameter :: MFopt = .true. ! Use mass flux or not
real, parameter :: visc_temf = 1.57e-4 ! WA TEST bigger minimum K
real, parameter :: conduc_temf = 1.57e-4 / 0.733
real, parameter :: Pr_temf = 0.733
real, parameter :: TEmin = 1e-3
real, parameter :: ftau0 = 0.17
real, parameter :: fth0 = 0.145
real, parameter :: critRi = 0.25
real, parameter :: Cf = 0.185
real, parameter :: CN = 2.0
real, parameter :: Ceps = 0.070
real, parameter :: Cgamma = Ceps
real, parameter :: Cphi = Ceps
real, parameter :: PrT0 = Cphi/Ceps * ftau0**2 / 2. / fth0**2
! EDMF constants
real, parameter :: CM = 0.03 ! Proportionality constant for subcloud MF
real, parameter :: Cdelt = 0.006 ! Prefactor for detrainment rate
real, parameter :: Cw = 0.5 ! Prefactor for surface wUPD
real, parameter :: Cc = 3.0 ! Prefactor for convective length scale
real, parameter :: lasymp = 200.0 ! Asymptotic length scale WA 11/20/09
real, parameter :: hmax = 4000.0 ! Max hd,hct WA 11/20/09
!
integer :: i, k, kt ! Loop variable
integer, dimension( its:ite) :: h0idx
real, dimension( its:ite) :: h0
real, dimension( its:ite) :: wstr, ang, wm
real, dimension( its:ite) :: hd,lcl,hct,ht
real, dimension( its:ite) :: convection_TKE_surface_src, sfcFTE
real, dimension( its:ite) :: sfcTHVF
real, dimension( its:ite) :: z0t
integer, dimension( its:ite) :: hdidx,lclidx,hctidx,htidx
integer, dimension( its:ite) :: hmax_idx
integer, dimension( its:ite) :: tval
real, dimension( its:ite, kts:kte) :: thetal, qt
real, dimension( its:ite, kts:kte) :: u_temf, v_temf
real, dimension( its:ite, kts:kte) :: rv, rl, rt
real, dimension( its:ite, kts:kte) :: chi_poisson, gam
real, dimension( its:ite, kts:kte) :: dthdz, dqtdz, dudz, dvdz
real, dimension( its:ite, kts:kte) :: lepsmin
real, dimension( its:ite, kts:kte) :: thetav
real, dimension( its:ite, kts:kte) :: MFCth, MFCq, MFCu, MFCv
real, dimension( its:ite, kts:kte) :: MFCql, MFCthv, MFCTE
real, dimension( its:ite, kts:kte) :: epsmf, deltmf, dMdz
real, dimension( its:ite, kts:kte) :: UUPD, VUPD
real, dimension( its:ite, kts:kte) :: thetavUPD, qlUPD, TEUPD
real, dimension( its:ite, kts:kte) :: thetavUPDmoist, wupd_dry
real, dimension( its:ite, kts:kte) :: B, Bmoist
real, dimension( its:ite, kts:kte) :: zm, zt, dzm, dzt
real, dimension( its:ite, kts:kte) :: dthUPDdz, dqtup_temfxdz, dwUPDdz
real, dimension( its:ite, kts:kte) :: dwUPDmoistdz
real, dimension( its:ite, kts:kte) :: dUUPDdz, dVUPDdz, dTEUPDdz
real, dimension( its:ite, kts:kte) :: TUPD, rstUPD, rUPD, rlUPD, qstUPD
real, dimension( its:ite, kts:kte) :: N2, S, Ri, beta, ftau, fth, ratio
real, dimension( its:ite, kts:kte) :: TKE, TE2
real, dimension( its:ite, kts:kte) :: ustrtilde, linv, leps
real, dimension( its:ite, kts:kte) :: km, kh
real, dimension( its:ite, kts:kte) :: Fz, QFK, uwk, vwk
real, dimension( its:ite, kts:kte) :: km_conv, kh_conv, lconv
real, dimension( its:ite, kts:kte) :: alpha2, beta2 ! For thetav flux calculation
real, dimension( its:ite, kts:kte) :: THVF, buoy_src, srcs
real, dimension( its:ite, kts:kte) :: u_new, v_new
real, dimension( its:ite, kts:kte) :: thx_new, qvx_new, qcx_new
real, dimension( its:ite, kts:kte) :: thup_new, qvup_new
real, dimension( its:ite, kts:kte) :: beta1 ! For saturation humidity calculations
real Cepsmf ! Prefactor for entrainment rate
real red_fact ! for reducing MF components
logical is_convective
! Vars for cloud fraction calculation
real, dimension( its:ite, kts:kte) :: au, sigq, qst, satdef
real sigq2, rst
!----------------------------------------------------------------------
! Grid staggering: Matlab version has mass and turbulence levels.
! WRF has full levels (with w) and half levels (u,v,theta,q*). Both
! sets of levels use the same indices (kts:kte). See pbl_driver or
! WRF Physics doc for (a few) details.
! So *mass levels correspond to half levels.*
! WRF full levels are ignored, we define our own turbulence levels
! in order to put the first one below the first half level.
! Another difference is that
! the Matlab version (and the Mauritsen et al. paper) consider the
! first mass level to be at z0 (effectively the surface). WRF considers
! the first half level to be above the effective surface. The first half
! level, at k=1, has nonzero values of u,v for example. Here we convert
! all incoming variables to internal ones with the correct indexing
! in order to make the code consistent with the Matlab version. We
! already had to do this for thetal and qt anyway, so the only additional
! overhead is for u and v.
! I use suffixes m for mass and t for turbulence as in Matlab for things
! like indices.
! Note that zsrf is the terrain height ASL, from Registry variable ht.
! Translations (Matlab to WRF):
! dzt -> calculated below
! dzm -> not supplied, calculated below
! k -> karman
! z0 -> znt
! z0t -> not in WRF, calculated below
! zt -> calculated below
! zm -> (z2d - zsrf) but NOTE zm(1) is now z0 (znt) and zm(2) is
! z2d(1) - zsrf
!
! WA I take the temperature at z0 to be
! TSK. This isn't exactly robust.
! WA 2/16/11 removed calculation of flhc, flqc which are not needed here.
! These should be removed from the calling sequence someday.
!
! Other notes:
! - I have often used 1 instead of kts below, because the scheme demands
! to know where the surface is. It won't work if kts .NE. 1.
do i = its,ite ! Main loop
! Get incoming surface theta from TSK (WA for now)
thetal(i,1) = tsk(i) / pi2d(i,1) ! WA really should use Exner func. at z0
qt(i,1) = qvx(i,1)
rv(i,1) = qt(i,1) / (1.-qt(i,1)) ! Water vapor
rl(i,1) = 0.
rt(i,1) = rv(i,1) + rl(i,1) ! Total water (without ice)
chi_poisson(i,1) = rcp * (1.+rv(i,1)/ep2) / (1.+rv(i,1)*cpv/cp)
gam(i,1) = rv(i,1) * r_v / (cp + rv(i,1)*cpv)
thetav(i,1) = thetal(i,1) * (1. + 0.608*qt(i,1) - qcx(i,1)) ! WA 4/6/10 allow environment liquid
z0t(i) = znt(i)
! Convert incoming theta to thetal and qv,qc to qt
! NOTE this is where the indexing gets changed from WRF to TEMF basis
do k = kts+1,kte
! Convert specific humidities to mixing ratios
rv(i,k) = qvx(i,k-1) / (1.-qvx(i,k-1)) ! Water vapor
rl(i,k) = qcx(i,k-1) / (1.-qcx(i,k-1)) ! Liquid water
rt(i,k) = rv(i,k) + rl(i,k) ! Total water (without ice)
chi_poisson(i,k) = rcp * (1.+rv(i,k)/ep2) / (1.+rv(i,k)*cpv/cp)
gam(i,k) = rt(i,k) * r_v / (cp + rt(i,k)*cpv)
thetal(i,k) = thx(i,k-1) * &
((ep2+rv(i,k))/(ep2+rt(i,k)))**chi_poisson(i,k) * &
(rv(i,k)/rt(i,k))**(-gam(i,k)) * exp( -xlv*rl(i,k) / &
((cp + rt(i,k)*cpv) * tx(i,k)))
qt(i,k) = qvx(i,k-1) + qcx(i,k-1)
thetav(i,k) = thetal(i,k) * (1. + 0.608*qt(i,k) - qcx(i,k-1)) ! WA 4/6/10 allow environment liquid
end do
! Convert incoming u,v to internal u_temf, v_temf
! NOTE this is where the indexing gets changed from WRF to TEMF basis
u_temf(i,1) = 0. ! zero winds at z0
v_temf(i,1) = 0.
do k = kts+1,kte
u_temf(i,k) = ux(i,k-1)
v_temf(i,k) = vx(i,k-1)
end do
! Get delta height at half (mass) levels
zm(i,1) = znt(i)
dzt(i,1) = z2d(i,1) - zsrf(i) - zm(i,1)
! Get height and delta at turbulence levels
zt(i,1) = (z2d(i,1) - zsrf(i) - znt(i)) / 2.
do kt = kts+1,kte
zm(i,kt) = z2d(i,kt-1) - zsrf(i) ! Convert indexing from WRF to TEMF
zt(i,kt) = (zm(i,kt) + z2d(i,kt) - zsrf(i)) / 2.
dzm(i,kt) = zt(i,kt) - zt(i,kt-1)
dzt(i,kt) = z2d(i,kt+1) - z2d(i,kt)
end do
dzm(i,1) = dzm(i,2) ! WA why?
dzt(i,kte) = dzt(i,kte-1) ! WA 12/23/09
! Gradients at first level
dthdz(i,1) = (thetal(i,2)-thetal(i,1)) / (zt(i,1) * log10(zm(i,2)/z0t(i)))
dqtdz(i,1) = (qt(i,2)-qt(i,1)) / (zt(i,1) * log10(zm(i,2)/z0t(i)))
dudz(i,1) = (u_temf(i,2)-u_temf(i,1)) / (zt(i,1) * log10(zm(i,2)/znt(i)))
dvdz(i,1) = (v_temf(i,2)-v_temf(i,1)) / (zt(i,1) * log10(zm(i,2)/znt(i)))
! Surface thetaV flux from Stull p.147
sfcTHVF(i) = hfx(i)/(rho(i,1)*cp) * (1.+0.608*(qvx(i,1)+qcx(i,1))) + 0.608*thetav(i,1)*qf_temfx(i,1)
! WA use hd_temf to calculate w* instead of finding h0 here????
! Watch initialization!
h0idx(i) = 1
h0(i) = zm(i,1)
lepsmin(i,kts) = 0.
! WA 2/11/13 find index just above hmax for use below
hmax_idx(i) = kte-1
do k = kts+1,kte-1
lepsmin(i,k) = 0.
! Mean gradients
dthdz(i,k) = (thetal(i,k+1) - thetal(i,k)) / dzt(i,k)
dqtdz(i,k) = (qt(i,k+1) - qt(i,k)) / dzt(i,k)
dudz(i,k) = (u_temf(i,k+1) - u_temf(i,k)) / dzt(i,k)
dvdz(i,k) = (v_temf(i,k+1) - v_temf(i,k)) / dzt(i,k)
! Find h0 (should eventually be interpolated for smoothness)
if (thetav(i,k) > thetav(i,1) .AND. h0idx(i) .EQ. 1) then
! WA 9/28/11 limit h0 as for hd and hct
if (zm(i,k) < hmax) then
h0idx(i) = k
h0(i) = zm(i,k)
else
h0idx(i) = k
h0(i) = hmax
end if
end if
! WA 2/11/13 find index just above hmax for use below
if (zm(i,k) > hmax) then
hmax_idx(i) = min(hmax_idx(i),k)
end if
end do
! Gradients at top level
dthdz(i,kte) = dthdz(i,kte-1)
dqtdz(i,kte) = dqtdz(i,kte-1)
dudz(i,kte) = dudz(i,kte-1)
dvdz(i,kte) = dvdz(i,kte-1)
if ( hfx(i) > 0.) then
! wstr(i) = (g * h0(i) / thetav(i,2) * shf_temfx(i,1) ) ** (1./3.)
wstr(i) = (g * h0(i) / thetav(i,2) * hfx(i)/(rho(i,1)*cp) ) ** (1./3.)
else
wstr(i) = 0.
end if
! Set flag convective or not for use below
is_convective = wstr(i) > 0. .AND. MFopt .AND. dthdz(i,1)<0. .AND. dthdz(i,2)<0. ! WA 12/16/09 require two levels of negative (unstable) gradient
! Find stability parameters and length scale (on turbulence levels)
do kt = 1,kte-1
N2(i,kt) = 2. * g / (thetav(i,kt) + thetav(i,kt+1))*dthdz(i,kt)
S(i,kt) = sqrt(dudz(i,kt)**2. + dvdz(i,kt)**2.)
Ri(i,kt) = N2(i,kt) / S(i,kt)**2.
if (S(i,kt) < 1e-15) then
if (N2(i,kt) >= 0) then
Ri(i,kt) = 10.
else
Ri(i,kt) = -1.
end if
end if
beta(i,kt) = 2. * g / (thetav(i,kt)+thetav(i,kt+1))
if (Ri(i,kt) > 0) then
ratio(i,kt) = Ri(i,kt)/(Cphi**2.*ftau0**2./(2.*Ceps**2.*fth0**2.)+3.*Ri(i,kt))
ftau(i,kt) = ftau0 * ((3./4.) / (1.+4.*Ri(i,kt)) + 1./4.)
fth(i,kt) = fth0 / (1.+4.*Ri(i,kt))
TE2(i,kt) = 2. * te_temfx(i,kt) * ratio(i,kt) * N2(i,kt) / beta(i,kt)**2.
else
ratio(i,kt) = Ri(i,kt)/(Cphi**2.*ftau0**2./(-2.*Ceps**2.*fth0**2.)+2.*Ri(i,kt))
ftau(i,kt) = ftau0
fth(i,kt) = fth0
TE2(i,kt) = 0.
end if
TKE(i,kt) = te_temfx(i,kt) * (1. - ratio(i,kt))
ustrtilde(i,kt) = sqrt(ftau(i,kt) * TKE(i,kt))
if (N2(i,kt) > 0.) then
linv(i,kt) = 1./karman / zt(i,kt) + abs(fCor(i)) / &
(Cf*ustrtilde(i,kt)) + &
sqrt(N2(i,kt))/(CN*ustrtilde(i,kt)) + 1./lasymp
else
linv(i,kt) = 1./karman / zt(i,kt) + abs(fCor(i)) / &
(Cf*ustrtilde(i,kt)) + 1./lasymp
end if
leps(i,kt) = 1./linv(i,kt)
leps(i,kt) = max(leps(i,kt),lepsmin(i,kt))
end do
S(i,kte) = 0.0
N2(i,kte) = 0.0
TKE(i,kte) = 0.0
linv(i,kte) = linv(i,kte-1)
leps(i,kte) = leps(i,kte-1)
! Find diffusion coefficients
! First use basic formulae for stable and neutral cases,
! then for convective conditions, and finally choose the larger
! WA 12/23/09 use convective form up to hd/2 always
! WA 12/28/09 after restructuring, this block is above MF block,
! so hd is not yet available for this timestep, must use h0,
! or use hd from previous timestep but be careful about initialization.
do kt = 1,kte-1 ! WA 12/22/09
! WA 4/8/10 remove beta term to avoid negative and huge values
! of km due to very small denominator. This is an interim fix
! until we find something better (more theoretically sound).
! km(i,kt) = TKE(i,kt)**1.5 * ftau(i,kt)**2. / (-beta(i,kt) * fth(i,kt) * sqrt(TE2(i,kt)) + Ceps * sqrt(TKE(i,kt)*te_temfx(i,kt)) / leps(i,kt))
km(i,kt) = TKE(i,kt)**1.5 * ftau(i,kt)**2. / (Ceps * sqrt(TKE(i,kt)*te_temfx(i,kt)) / leps(i,kt))
kh(i,kt) = 2. * leps(i,kt) * fth(i,kt)**2. * TKE(i,kt) / sqrt(te_temfx(i,kt)) / Cphi
if ( is_convective) then
! WA 2/20/14 trap rare "equality" of h0 and zt (only when h0 = hmax)
if (kt <= h0idx(i) .AND. h0(i)-zt(i,kt) > 1e-15) then
lconv(i,kt) = 1. / (1. / (karman*zt(i,kt)) + Cc / (karman * (h0(i) - zt(i,kt))))
else
lconv(i,kt) = 0.
end if
! WA 12/15/09 use appropriate coeffs to match kh_conv and kh at neutral
kh_conv(i,kt) = ftau0**2. / Ceps / PrT0 * sqrt(TKE(i,kt)) * lconv(i,kt)
if (kh_conv(i,kt) < 0.) then
kh_conv(i,kt) = 0.
end if
km_conv(i,kt) = PrT0 * kh_conv(i,kt)
if (zt(i,kt) <= h0(i)/2.) then
km(i,kt) = km_conv(i,kt)
kh(i,kt) = kh_conv(i,kt)
end if
if (zt(i,kt) > h0(i)/2. .AND. kt <= h0idx(i)) then
km(i,kt) = max(km(i,kt),km_conv(i,kt),visc_temf)
kh(i,kt) = max(kh(i,kt),kh_conv(i,kt),conduc_temf)
end if
end if ! is_convective
km(i,kt) = max(km(i,kt),visc_temf)
kh(i,kt) = max(kh(i,kt),conduc_temf)
Fz(i,kt) = -kh(i,kt) * dthdz(i,kt) ! Diffusive heat flux
end do
km(i,kte) = km(i,kte-1)
kh(i,kte) = kh(i,kte-1)
Fz(i,kte) = 0.0
!*** Mass flux block starts here ***
if ( is_convective) then
Cepsmf = 2. / max(200.,h0(i))
Cepsmf = max(Cepsmf,0.002)
do k = kts,kte
! Calculate lateral entrainment fraction for subcloud layer
! epsilon and delta are defined on mass grid (half levels)
epsmf(i,k) = Cepsmf
end do
! Initialize updraft
thup_temfx(i,1) = thetal(i,1) ! No excess
qtup_temfx(i,1) = qt(i,1) ! No excess
rUPD(i,1) = qtup_temfx(i,1) / (1. - qtup_temfx(i,1))
wupd_temfx(i,1) = Cw * wstr(i)
wupd_dry(i,1) = Cw * wstr(i)
UUPD(i,1) = u_temf(i,1)
VUPD(i,1) = v_temf(i,1)
thetavUPD(i,1) = thup_temfx(i,1) * (1. + 0.608*qtup_temfx(i,1)) ! WA Assumes no liquid
thetavUPDmoist(i,1) = thup_temfx(i,1) * (1. + 0.608*qtup_temfx(i,1)) ! WA Assumes no liquid
TEUPD(i,1) = te_temfx(i,1) + g / thetav(i,1) * sfcTHVF(i)
qlUPD(i,1) = qcx(i,1) ! WA allow environment liquid
TUPD(i,1) = thup_temfx(i,1) * pi2d(i,1)
rstUPD(i,1) = rsat(p2d(i,1),TUPD(i,1),ep2)
rlUPD(i,1) = 0.
! Calculate updraft parameters counting up
do k = 2,kte
! WA 2/11/13 use hmax index to prevent oddness high up
if ( k < hmax_idx(i)) then
dthUPDdz(i,k-1) = -epsmf(i,k) * (thup_temfx(i,k-1) - thetal(i,k-1))
thup_temfx(i,k) = thup_temfx(i,k-1) + dthUPDdz(i,k-1) * dzm(i,k-1)
dqtup_temfxdz(i,k-1) = -epsmf(i,k) * (qtup_temfx(i,k-1) - qt(i,k-1))
qtup_temfx(i,k) = qtup_temfx(i,k-1) + dqtup_temfxdz(i,k-1) * dzm(i,k-1)
thetavUPD(i,k) = thup_temfx(i,k) * (1. + 0.608*qtup_temfx(i,k)) ! WA Assumes no liquid
B(i,k-1) = g * (thetavUPD(i,k) - thetav(i,k)) / thetav(i,k)
if ( wupd_dry(i,k-1) < 1e-15 ) then
wupd_dry(i,k) = 0.
else
dwUPDdz(i,k-1) = -2. *epsmf(i,k)*wupd_dry(i,k-1) + 0.33*B(i,k-1)/wupd_dry(i,k-1)
wupd_dry(i,k) = wupd_dry(i,k-1) + dwUPDdz(i,k-1) * dzm(i,k-1)
end if
dUUPDdz(i,k-1) = -epsmf(i,k) * (UUPD(i,k-1) - u_temf(i,k-1))
UUPD(i,k) = UUPD(i,k-1) + dUUPDdz(i,k-1) * dzm(i,k-1)
dVUPDdz(i,k-1) = -epsmf(i,k) * (VUPD(i,k-1) - v_temf(i,k-1))
VUPD(i,k) = VUPD(i,k-1) + dVUPDdz(i,k-1) * dzm(i,k-1)
dTEUPDdz(i,k-1) = -epsmf(i,k) * (TEUPD(i,k-1) - te_temfx(i,k-1))
TEUPD(i,k) = TEUPD(i,k-1) + dTEUPDdz(i,k-1) * dzm(i,k-1)
! Alternative updraft velocity based on moist thetav
! Need thetavUPDmoist, qlUPD
rUPD(i,k) = qtup_temfx(i,k) / (1. - qtup_temfx(i,k))
! WA Updraft temperature assuming no liquid
TUPD(i,k) = thup_temfx(i,k) * pi2d(i,k)
! Updraft saturation mixing ratio
rstUPD(i,k) = rsat(p2d(i,k-1),TUPD(i,k),ep2)
! Correct to actual temperature (Sommeria & Deardorff 1977)
beta1(i,k) = 0.622 * (xlv/(r_d*TUPD(i,k))) * (xlv/(cp*TUPD(i,k)))
rstUPD(i,k) = rstUPD(i,k) * (1.0+beta1(i,k)*rUPD(i,k)) / (1.0+beta1(i,k)*rstUPD(i,k))
qstUPD(i,k) = rstUPD(i,k) / (1. + rstUPD(i,k))
if (rUPD(i,k) > rstUPD(i,k)) then
rlUPD(i,k) = rUPD(i,k) - rstUPD(i,k)
qlUPD(i,k) = rlUPD(i,k) / (1. + rlUPD(i,k))
thetavUPDmoist(i,k) = (thup_temfx(i,k) + ((xlv/cp)*qlUPD(i,k)/pi2d(i,k))) * &
(1. + 0.608*qstUPD(i,k) - qlUPD(i,k))
else
rlUPD(i,k) = 0.
qlUPD(i,k) = qcx(i,k-1) ! WA 4/6/10 allow environment liquid
thetavUPDmoist(i,k) = thup_temfx(i,k) * (1. + 0.608*qtup_temfx(i,k))
end if
Bmoist(i,k-1) = g * (thetavUPDmoist(i,k) - thetav(i,k)) / thetav(i,k)
if ( wupd_temfx(i,k-1) < 1e-15 ) then
wupd_temfx(i,k) = 0.
else
dwUPDmoistdz(i,k-1) = -2. *epsmf(i,k)*wupd_temfx(i,k-1) + 0.33*Bmoist(i,k-1)/wupd_temfx(i,k-1)
wupd_temfx(i,k) = wupd_temfx(i,k-1) + dwUPDmoistdz(i,k-1) * dzm(i,k-1)
end if
else
thup_temfx(i,k) = thetal(i,k)
qtup_temfx(i,k) = qt(i,k)
wupd_dry(i,k) = 0.
UUPD(i,k) = u_temf(i,k)
VUPD(i,k) = v_temf(i,k)
TEUPD(i,k) = te_temfx(i,k)
qlUPD(i,k) = qcx(i,k-1)
wupd_temfx(i,k) = 0.
end if
end do
! Find hd based on wUPD
if (wupd_dry(i,1) == 0.) then
hdidx(i) = 1
else
hdidx(i) = kte ! In case wUPD <= 0 not found
do k = 2,kte
! if (wupd_dry(i,k) <= 0.) then
if (wupd_dry(i,k) <= 0. .OR. zm(i,k) > hmax) then
hdidx(i) = k
goto 100 ! FORTRAN made me do it!
end if
end do
end if
100 hd(i) = zm(i,hdidx(i))
kpbl1d(i) = hdidx(i)
hpbl(i) = hd(i) ! WA 5/11/10 hpbl is height. Should still be replaced by something that works whether convective or not.
! Find LCL, hct, and ht
lclidx(i) = kte ! In case LCL not found
do k = kts,kte
if ( k < hmax_idx(i) .AND. rUPD(i,k) > rstUPD(i,k)) then
lclidx(i) = k
goto 200
end if
end do
200 lcl(i) = zm(i,lclidx(i))
if (hd(i) > lcl(i)) then ! Forced cloud (at least) occurs
! Find hct based on wUPDmoist
if (wupd_temfx(i,1) == 0.) then
hctidx(i) = 1
else
hctidx(i) = kte ! In case wUPD <= 0 not found
do k = 2,kte
if (wupd_temfx(i,k) <= 0. .OR. zm(i,k) > hmax) then
hctidx(i) = k
goto 300 ! FORTRAN made me do it!
end if
end do
end if
300 hct(i) = zm(i,hctidx(i))
if (hctidx(i) <= hdidx(i)+1) then ! No active cloud
hct(i) = hd(i)
hctidx(i) = hdidx(i)
else
end if
else ! No cloud
hct(i) = hd(i)
hctidx(i) = hdidx(i)
end if
ht(i) = max(hd(i),hct(i))
htidx(i) = max(hdidx(i),hctidx(i))
! Now truncate updraft at ht with taper
do k = 1,kte
if (zm(i,k) < 0.9*ht(i)) then ! Below taper region
tval(i) = 1
else if (zm(i,k) >= 0.9*ht(i) .AND. zm(i,k) <= 1.0*ht(i)) then
! Within taper region
tval(i) = 1. - ((zm(i,k) - 0.9*ht(i)) / (1.0*ht(i) - 0.9*ht(i)))
else ! Above taper region
tval(i) = 0.
end if
thup_temfx(i,k) = tval(i) * thup_temfx(i,k) + (1-tval(i))*thetal(i,k)
thetavUPD(i,k) = tval(i) * thetavUPD(i,k) + (1-tval(i))*thetav(i,k)
qtup_temfx(i,k) = tval(i) * qtup_temfx(i,k) + (1-tval(i)) * qt(i,k)
! WA 6/21/13 was a subscript error when k=1
if (k > 1) then
qlUPD(i,k) = tval(i) * qlUPD(i,k) + (1-tval(i)) * qcx(i,k-1)
end if
UUPD(i,k) = tval(i) * UUPD(i,k) + (1-tval(i)) * u_temf(i,k)
VUPD(i,k) = tval(i) * VUPD(i,k) + (1-tval(i)) * v_temf(i,k)
TEUPD(i,k) = tval(i) * TEUPD(i,k) + (1-tval(i)) * te_temfx(i,k)
if (zm(i,k) > ht(i)) then ! WA this is just for cleanliness
wupd_temfx(i,k) = 0.
dwUPDmoistdz(i,k) = 0.
wupd_dry(i,k) = 0.
dwUPDdz(i,k) = 0.
end if
end do
! Calculate lateral detrainment rate for cloud layer
deltmf(i,1) = Cepsmf
do k = 2,kte-1
if (hctidx(i) > hdidx(i)+1) then ! Some cloud
deltmf(i,k) = 0.9 * Cepsmf + Cdelt * (atan((zm(i,k)-(lcl(i)+(hct(i)-lcl(i))/1.5))/ &
((hct(i)-lcl(i))/8))+(3.1415926/2))/3.1415926
else if (k < hdidx(i)) then ! No cloud, below hd
deltmf(i,k) = Cepsmf + 0.05 * 1. / (hd(i) - zm(i,k))
else if (k >= hdidx(i)) then ! No cloud, above hd
deltmf(i,k) = deltmf(i,k-1)
end if
end do
! Calculate mass flux (defined on turbulence levels)
mf_temfx(i,1) = CM * wstr(i)
do kt = 2,kte-1
dMdz(i,kt) = (epsmf(i,kt) - deltmf(i,kt)) * mf_temfx(i,kt-1) * dzt(i,kt)
mf_temfx(i,kt) = mf_temfx(i,kt-1) + dMdz(i,kt)
end do
! WA 12/28/09 If mass flux component > diffusive
! component at second level,
! reduce M to prevent a stable layer
MFCth(i,2) = mf_temfx(i,2) * (thup_temfx(i,2)-thetal(i,2) + thup_temfx(i,3)-thetal(i,3)) / 2.
if (MFCth(i,2) > Fz(i,2)) then
red_fact = Fz(i,2) / MFCth(i,2)
do kt = 1,kte
mf_temfx(i,kt) = mf_temfx(i,kt) * red_fact
end do
end if ! Reduce M to prevent stable layer at second level
! Calculate mass flux contributions to fluxes (defined on turb levels)
! Use log interpolation at first level
MFCth(i,1) = mf_temfx(i,1) * (thup_temfx(i,1)-thetal(i,1) &
+ (thup_temfx(i,2)-thetal(i,2) - &
(thup_temfx(i,1)-thetal(i,1))) * log(zt(i,1)/znt(i))/log(zm(i,2)/znt(i)))
MFCq(i,1) = mf_temfx(i,1) * (qtup_temfx(i,1)-qt(i,1) &
+ (qtup_temfx(i,2)-qt(i,2) - &
(qtup_temfx(i,1)-qt(i,1))) * log(zt(i,1)/znt(i))/log(zm(i,2)/znt(i)))
MFCu(i,1) = mf_temfx(i,1) * (UUPD(i,1)-u_temf(i,1) &
+ (UUPD(i,2)-u_temf(i,2) - &
(UUPD(i,1)-u_temf(i,1))) * log(zt(i,1)/znt(i))/log(zm(i,2)/znt(i)))
MFCv(i,1) = mf_temfx(i,1) * (VUPD(i,1)-v_temf(i,1) &
+ (VUPD(i,2)-v_temf(i,2) - &
(VUPD(i,1)-v_temf(i,1))) * log(zt(i,1)/znt(i))/log(zm(i,2)/znt(i)))
MFCql(i,1) = mf_temfx(i,1) * (qlUPD(i,1)-qcx(i,1) &
+ (qlUPD(i,2)-qcx(i,2) - &
(qlUPD(i,1)-qcx(i,1))) * log(zt(i,1)/znt(i))/log(zm(i,2)/znt(i)))
MFCTE(i,1) = mf_temfx(i,1) * (TEUPD(i,1)-te_temfx(i,1) &
+ (TEUPD(i,2)-te_temfx(i,2) - &
(TEUPD(i,1)-te_temfx(i,1))) * log(zt(i,1)/znt(i))/log(zm(i,2)/znt(i))) ! WA Check this
do kt = 2,kte-1
MFCth(i,kt) = mf_temfx(i,kt) * (thup_temfx(i,kt)-thetal(i,kt) + thup_temfx(i,kt+1)-thetal(i,kt+1)) / 2.
MFCq(i,kt) = mf_temfx(i,kt) * (qtup_temfx(i,kt)-qt(i,kt) + qtup_temfx(i,kt+1)-qt(i,kt+1)) / 2.
MFCu(i,kt) = mf_temfx(i,kt) * (UUPD(i,kt)-u_temf(i,kt) + UUPD(i,kt+1)-u_temf(i,kt+1)) / 2.
MFCv(i,kt) = mf_temfx(i,kt) * (VUPD(i,kt)-v_temf(i,kt) + VUPD(i,kt+1)-v_temf(i,kt+1)) / 2.
MFCql(i,kt) = mf_temfx(i,kt) * (qlUPD(i,kt)-qcx(i,kt-1) + qlUPD(i,kt+1)-qcx(i,kt)) / 2.
MFCTE(i,kt) = mf_temfx(i,kt) * (TEUPD(i,kt)-te_temfx(i,kt)) ! TE is on turb levels
end do
MFCth(i,kte) = 0
MFCq(i,kte) = 0
MFCu(i,kte) = 0
MFCv(i,kte) = 0
MFCql(i,kte) = 0
MFCTE(i,kte) = 0
! Calculate cloud fraction (on mass levels)
cf3d_temfx(i,1) = 0.0
cfm_temfx(i) = 0.0
do k = 2,kte
! if (wupd_temfx(i,k-1) >= 1.0e-15 .AND. wupd_temfx(i,k) >= 1.0e-15 .AND. .NOT. isnan(wupd_temfx(i,k-1)) .AND. .NOT. isnan(wupd_temfx(i,k))) then
if (wupd_temfx(i,k-1) >= 1.0e-15 .AND. wupd_temfx(i,k) >= 1.0e-15) then
au(i,k) = ((mf_temfx(i,k-1)+mf_temfx(i,k))/2.0) / ((wupd_temfx(i,k-1)+wupd_temfx(i,k))/2.0) ! WA average before divide, is that best?
else
au(i,k) = 0.0
end if
sigq2 = au(i,k) * (qtup_temfx(i,k)-qt(i,k))
if (sigq2 > 0.0) then
sigq(i,k) = sqrt(sigq2)
else
sigq(i,k) = 0.0
end if
! rst = rsat(p2d(i,k),thx(i,k)*pi2d(i,k),ep2)
rst = rsat(p2d(i,k-1),thx(i,k-1)*pi2d(i,k-1),ep2)
qst(i,k) = rst / (1. + rst)
satdef(i,k) = qt(i,k) - qst(i,k)
if (satdef(i,k) <= 0.0) then
if (sigq(i,k) > 1.0e-15) then
cf3d_temfx(i,k) = max(0.5 + 0.36 * atan(1.55*(satdef(i,k)/sigq(i,k))),0.0)
else
cf3d_temfx(i,k) = 0.0
end if
else
cf3d_temfx(i,k) = 1.0
end if
if (zm(i,k) < lcl(i)) then
cf3d_temfx(i,k) = 0.0
end if
! Put max value so far into cfm
if (zt(i,k) <= hmax) then
cfm_temfx(i) = max(cf3d_temfx(i,k),cfm_temfx(i))
end if
end do
else ! not is_convective, no MF components
do kt = 1,kte
MFCth(i,kt) = 0
MFCq(i,kt) = 0
MFCu(i,kt) = 0
MFCv(i,kt) = 0
MFCql(i,kt) = 0
MFCTE(i,kt) = 0
end do
lcl(i) = zm(i,kte-1)
hct(i) = zm(i,1)
hctidx(i) = 1
hd(i) = zm(i,1)
hdidx(i) = 1
ht(i) = hd(i)
! Cloud fraction calculations
cf3d_temfx(i,1) = 0.0
cfm_temfx(i) = 0.0
do k = 2,kte
if (qcx(i,k-1) > 1.0e-15) then
cf3d_temfx(i,k) = 1.0
else
cf3d_temfx(i,k) = 0.0
end if
! Put max value so far into cfm
if (zt(i,k) <= hmax) then
cfm_temfx(i) = max(cf3d_temfx(i,k),cfm_temfx(i))
end if
end do
end if ! MF components or not
cf3d_temfx(i,kte) = 0.0
! Mass flux block ends here
! Flux profiles
do kt = 2,kte
! Fz(i,kt) = -kh(i,kt) * dthdz(i,kt)
shf_temfx(i,kt) = Fz(i,kt) + MFCth(i,kt)
QFK(i,kt) = -kh(i,kt) * dqtdz(i,kt)
qf_temfx(i,kt) = QFK(i,kt) + MFCq(i,kt)
uwk(i,kt) = -km(i,kt) * dudz(i,kt)
uw_temfx(i,kt) = uwk(i,kt) + MFCu(i,kt)
vwk(i,kt) = -km(i,kt) * dvdz(i,kt)
vw_temfx(i,kt) = vwk(i,kt) + MFCv(i,kt)
end do
! Surface momentum fluxes
! WA TEST 11/7/13 use w* as a component of the mean wind inside the
! u* calculation instead of in the velocity scale below (Felix)
! ust(i) = sqrt(ftau(i,1)/ftau0) * sqrt(u_temf(i,2)**2. + v_temf(i,2)**2.) * leps(i,1) / log(zm(i,2)/znt(i)) / zt(i,1)
ust(i) = sqrt(ftau(i,1)/ftau0) * sqrt(u_temf(i,2)**2. + v_temf(i,2)**2. + (0.5*wstr(i))**2.) * leps(i,1) / log(zm(i,2)/znt(i)) / zt(i,1)
ang(i) = atan2(v_temf(i,2),u_temf(i,2))
uw_temfx(i,1) = -cos(ang(i)) * ust(i)**2.
vw_temfx(i,1) = -sin(ang(i)) * ust(i)**2.