-
Notifications
You must be signed in to change notification settings - Fork 5
/
module_convection_prep.F
3960 lines (3685 loc) · 167 KB
/
module_convection_prep.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
module module_convection_prep
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!
! !MODULE: module_convection_prep.F
!
! !DESCRIPTION: This module performs the diagnostics of variables needed
! by GEOS-CHEM convection\_mod from WRF cumulus
! parameterization. (xfeng, 11/14/2018)
!
!----------------------------NTIEDEKE SCHEME----------------------------
!
use module_model_constants, only:rd=>r_d, rv=>r_v, &
& cpd=>cp, alv=>xlv, als=>xls, alf=>xlf, t13=>Prandtl, g
implicit none
real,private :: rcpd,vtmpc1,tmelt, &
c1es,c2es,c3les,c3ies,c4les,c4ies,c5les,c5ies,zrg
real,private :: r5alvcp,r5alscp,ralvdcp,ralsdcp,ralfdcp,rtwat,rtber,rtice
real,private :: entrdd,cmfcmax,cmfcmin,cmfdeps,zdnoprc,cprcon,pgcoef
integer,private :: momtrans
parameter( &
rcpd=1.0/cpd, &
tmelt=273.16, &
zrg=1.0/g, &
c1es=610.78, &
c2es=c1es*rd/rv, &
c3les=17.2693882, &
c3ies=21.875, &
c4les=35.86, &
c4ies=7.66, &
c5les=c3les*(tmelt-c4les), &
c5ies=c3ies*(tmelt-c4ies), &
r5alvcp=c5les*alv*rcpd, &
r5alscp=c5ies*als*rcpd, &
ralvdcp=alv*rcpd, &
ralsdcp=als*rcpd, &
ralfdcp=alf*rcpd, &
rtwat=tmelt, &
rtber=tmelt-5., &
rtice=tmelt-23., &
vtmpc1=rv/rd-1.0 )
!
! entrdd: average entrainment & detrainment rate for downdrafts
! ------
!
parameter(entrdd = 2.0e-4)
!
! cmfcmax: maximum massflux value allowed for updrafts etc
! -------
!
parameter(cmfcmax = 1.0)
!
! cmfcmin: minimum massflux value (for safety)
! -------
!
parameter(cmfcmin = 1.e-10)
!
! cmfdeps: fractional massflux for downdrafts at lfs
! -------
!
parameter(cmfdeps = 0.30)
! zdnoprc: deep cloud is thicker than this height (Unit:Pa)
!
parameter(zdnoprc = 2.0e4)
! -------
!
! cprcon: coefficient from cloud water to rain water
!
parameter(cprcon = 1.4e-3)
! -------
!
! momtrans: momentum transport method
! ( 1 = IFS40r1 method; 2 = new method )
!
parameter(momtrans = 2 )
! -------
!
! coefficient for pressure gradient intensity
! (0.7 - 1.0 is recommended in this vesion of Tiedtke scheme)
parameter(pgcoef=0.7)
! -------
!
logical :: nonequil
! nonequil: representing equilibrium and nonequilibrium convection
! ( .false. [equilibrium: removing all CAPE]; .true. [nonequilibrium: relaxing CAPE toward CAPE from PBL].
! Ref. Bechtold et al. 2014 JAS )
!
parameter(nonequil = .true. )
!
!--------------------
! switches for deep, mid, shallow convections, downdraft, and momentum transport
! ------------------
logical :: lmfpen,lmfmid,lmfscv,lmfdd,lmfdudv
parameter(lmfpen=.true.,lmfmid=.true.,lmfscv=.true.,lmfdd=.true.,lmfdudv=.true.)
!--------------------------------------------------------------------------------------------------------------
contains
!-----------------------------------------------------------------------
subroutine cu_ntiedtke( &
dt,itimestep,stepcu &
,raincv,pratec,qfx,hfx &
,u3d,v3d,w,t3d,qv3d,qc3d,qi3d,pi3d,rho3d &
,qvften,thften &
,dz8w,pcps,p8w,xland,cu_act_flag,dx &
,ids,ide, jds,jde, kds,kde &
,ims,ime, jms,jme, kms,kme &
,its,ite, jts,jte, kts,kte &
,rthcuten,rqvcuten,rqccuten,rqicuten &
,rucuten, rvcuten &
! ,f_qv ,f_qc ,f_qr ,f_qi ,f_qs &
,pmflxrain, pmflxsnow &
,cmfmc, dqrcu, cloud_bot, cloud_top &
,dtrainu, dtraind, dtrain, reevapcn &
)
!-------------------------------------------------------------------
implicit none
!-------------------------------------------------------------------
!-- 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)
!-- rho3d 3d air density (kg/m^3)
!-- p8w 3d hydrostatic pressure at full levels (pa)
!-- pcps 3d hydrostatic pressure at half levels (pa)
!-- pi3d 3d exner function (dimensionless)
!-- qvften 3d total advective + PBL moisture tendency (kg kg-1 s-1)
!-- thften 3d total advective + PBL + radiative temperature tendency (k s-1)
!-- rthcuten theta tendency due to
! cumulus scheme precipitation (k/s)
!-- rucuten u wind tendency due to
! cumulus scheme precipitation (k/s)
!-- rvcuten v wind tendency due to
! cumulus scheme precipitation (k/s)
!-- rqvcuten qv tendency due to
! cumulus scheme precipitation (kg/kg/s)
!-- rqccuten qc tendency due to
! cumulus scheme precipitation (kg/kg/s)
!-- rqicuten qi tendency due to
! cumulus scheme precipitation (kg/kg/s)
!-- rainc accumulated total cumulus scheme precipitation (mm)
!-- raincv cumulus scheme precipitation (mm)
!-- pratec precipitiation rate from cumulus scheme (mm/s)
!-- pmflxrain convective rain flux (kg/m^2/s) ***fox change***
!-- pmflxsnow convective snow flux (kg/m^2/s) ***fox change***
!-- cmfmc mass flux in updrafts (kg/m^2/s) ***fox change***
!-- dqrcu to determin the cloud base level ***fox change***
!-- dtrain detrainment flux (kg/m^2/s) ***fox change***
!-- reevapcn flux difference of precip in downdrafts (kg/m^2/s) ***fox change***
!-- dz8w dz between full levels (m)
!-- qfx upward moisture flux at the surface (kg/m^2/s)
!-- hfx upward heat flux at the surface (w/m^2)
!-- dt time step (s)
!-- 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
!-------------------------------------------------------------------
integer, intent(in) :: ids,ide, jds,jde, kds,kde, &
ims,ime, jms,jme, kms,kme, &
its,ite, jts,jte, kts,kte, &
itimestep, &
stepcu
real, intent(in) :: &
dt, &
dx
real, dimension(ims:ime, jms:jme), intent(in) :: &
xland
real, dimension(ims:ime, jms:jme), intent(inout) :: &
raincv, pratec
real, dimension(ims:ime, kms:kme, jms:jme), intent(inout) :: &
pmflxrain, pmflxsnow, cmfmc
real, dimension(ims:ime, kms:kme, jms:jme), intent(inout) :: &
dqrcu
integer, dimension(ims:ime, jms:jme), intent(inout) :: &
cloud_bot, cloud_top
real, dimension(ims:ime, kms:kme, jms:jme), intent(inout) :: &
dtrainu, dtraind, dtrain
real, dimension(ims:ime, kms:kme, jms:jme), intent(inout) :: &
reevapcn
!***fox change***
! add output variables pmflxrain, pmflxsnow [convective rain/snow flux kg/(m^2*s)]
! add output variables cmfmc [mass flux in updrafts kg/(m^2*s)]
! add output variables dqrcu cloud_bot cloud_top
! add output variables dtrainu dtraind dtrain [up/downdrafts detrainment flux kg/(m^2*s)]
! add output variables reevapcn [flux difference of precip in downdrafts kg/(m^2*s)]
logical, dimension(ims:ime,jms:jme), intent(inout) :: &
cu_act_flag
real, dimension(ims:ime, kms:kme, jms:jme), intent(in) :: &
dz8w, &
pcps, &
p8w, &
pi3d, &
qc3d, &
qvften, &
thften, &
qi3d, &
qv3d, &
rho3d, &
t3d, &
u3d, &
v3d, &
w
real, dimension(ims:ime, jms:jme) :: &
qfx, &
hfx
!--------------------------- optional vars ----------------------------
real, dimension(ims:ime, kms:kme, jms:jme), &
optional, intent(inout) :: &
rqccuten, &
rqicuten, &
rqvcuten, &
rthcuten, &
rucuten, &
rvcuten
!
! flags relating to the optional tendency arrays declared above
! models that carry the optional tendencies will provdide the
! optional arguments at compile time; these flags all the model
! to determine at run-time whether a particular tracer is in
! use or not.
!
! logical, optional :: &
! f_qv &
! ,f_qc &
! ,f_qr &
! ,f_qi &
! ,f_qs
!--------------------------- local vars ------------------------------
real :: &
delt, &
rdelt
real , dimension(its:ite) :: &
rcs, &
rn, &
evap, &
heatflux
integer , dimension(its:ite) :: slimsk
real , dimension(its:ite, kts:kte+1) :: &
pmflxr1, &
pmflxs1
real , dimension(its:ite, kts:kte) :: &
zmfu1, &
zmfd1
integer , dimension(its:ite) :: icbot1, &
ictop1
real , dimension(its:ite, kts:kte) :: &
pmfude_rate1, &
pmfdde_rate1
real , dimension(its:ite, kts:kte) :: &
zdmfdp1
!***fox change***
! add local variables pmflxr1, pmflxs1, zmfu1, zmfd1, icbot1, ictop1, pmfude_rate1, pmfdde_rate1, zdmfdp1
real , dimension(its:ite, kts:kte+1) :: &
prsi, &
ghti, &
zi
real , dimension(its:ite, kts:kte) :: &
dot, &
prsl, &
q1, &
q2, &
q3, &
q1b, &
t1b, &
q11, &
q12, &
t1, &
u1, &
v1, &
zl, &
omg, &
ghtl
integer, dimension(its:ite) :: &
kbot, &
ktop
integer :: &
i, &
im, &
j, &
k, &
km, &
kp, &
kx, &
kx1
!-------other local variables----
integer :: zz, pp
!-----------------------------------------------------------------------
!
!
!*** check to see if this is a convection timestep
!
!-----------------------------------------------------------------------
do j=jts,jte
do i=its,ite
cu_act_flag(i,j)=.true.
enddo
enddo
im=ite-its+1
kx=kte-kts+1
kx1=kx+1
delt=dt*stepcu
rdelt=1./delt
!------------- j loop (outer) --------------------------------------------------
do j=jts,jte
! --------------- compute zi and zl -----------------------------------------
do i=its,ite
zi(i,kts)=0.0
enddo
!
do k=kts,kte
do i=its,ite
zi(i,k+1)=zi(i,k)+dz8w(i,k,j)
enddo
enddo
!
do k=kts,kte
do i=its,ite
zl(i,k)=0.5*(zi(i,k)+zi(i,k+1))
enddo
enddo
! --------------- end compute zi and zl -------------------------------------
do i=its,ite
slimsk(i)=int(abs(xland(i,j)-2.))
enddo
do k=kts,kte
kp=k+1
do i=its,ite
dot(i,k)=-0.5*g*rho3d(i,k,j)*(w(i,k,j)+w(i,kp,j))
enddo
enddo
pp = 0
do k=kts,kte
zz = kte-pp
do i=its,ite
u1(i,zz)=u3d(i,k,j)
v1(i,zz)=v3d(i,k,j)
t1(i,zz)=t3d(i,k,j)
q1(i,zz)=qv3d(i,k,j)
if(itimestep == 1) then
q1b(i,zz)=0.
t1b(i,zz)=0.
else
q1b(i,zz)=qvften(i,k,j)
t1b(i,zz)=thften(i,k,j)
endif
q2(i,zz)=qc3d(i,k,j)
q3(i,zz)=qi3d(i,k,j)
omg(i,zz)=dot(i,k)
ghtl(i,zz)=zl(i,k)
prsl(i,zz) = pcps(i,k,j)
enddo
pp = pp + 1
enddo
pp = 0
do k=kts,kte+1
zz = kte+1-pp
do i=its,ite
ghti(i,zz) = zi(i,k)
prsi(i,zz) = p8w(i,k,j)
enddo
pp = pp + 1
enddo
!
do i=its,ite
evap(i) = qfx(i,j)
heatflux(i)= hfx(i,j)
enddo
!----------------------------------------------------------------------------------------
call tiecnvn(u1,v1,t1,q1,q2,q3,q1b,t1b,ghtl,ghti,omg,prsl,prsi,evap,heatflux, &
rn,pmflxr1,pmflxs1,zmfu1,zmfd1,icbot1,ictop1,pmfude_rate1,pmfdde_rate1, &
zdmfdp1,slimsk,im,kx,kx1,delt,dx)
do i=its,ite
raincv(i,j)=rn(i)/stepcu
pratec(i,j)=rn(i)/(stepcu * dt)
enddo
! add pmflxr1, pmflxs1
pp = 0
do k=kts,kte+1
zz = kte+1-pp
do i=its,ite
pmflxrain(i,k,j)=pmflxr1(i,zz)
pmflxsnow(i,k,j)=pmflxs1(i,zz)
enddo
pp = pp + 1
enddo
! add icbot1 and ictop1
do i=its,ite
cloud_bot(i,j)=kte - icbot1(i)
cloud_top(i,j)=kte - ictop1(i)
if (cloud_top(i,j) .gt. kte)then
cloud_top(i,j) = 0
endif
enddo
! add zmfu1, zmfd1
pp=0
do k=kts,kte
zz = kte-pp
do i=its,ite
cmfmc(i,k,j)=zmfu1(i,zz) + zmfd1(i,zz)
if (k <= cloud_bot(i,j)) then
cmfmc(i,k,j)=0.0
endif
cmfmc(i,k,j)=min(cmfmc(i,k,j),0.2)
cmfmc(i,k,j)=max(cmfmc(i,k,j),0.0)
enddo
pp = pp + 1
enddo
! add pmfude_rate1, pmfdde_rate1, reevapcn
pp=0
do k=kts,kte
zz = kte-pp
do i=its,ite
dtrainu(i,k,j)=pmfude_rate1(i,zz)
dtraind(i,k,j)=pmfdde_rate1(i,zz)
reevapcn(i,k,j)=zdmfdp1(i,zz)
enddo
pp = pp + 1
enddo
pp = 0
do k=kts,kte
zz = kte-pp
do i=its,ite
rthcuten(i,k,j)=(t1(i,zz)-t3d(i,k,j))/pi3d(i,k,j)*rdelt
rqvcuten(i,k,j)=(q1(i,zz)-qv3d(i,k,j))*rdelt
rucuten(i,k,j) =(u1(i,zz)-u3d(i,k,j))*rdelt
rvcuten(i,k,j) =(v1(i,zz)-v3d(i,k,j))*rdelt
enddo
pp = pp + 1
enddo
if(present(rqccuten))then
! if ( f_qc ) then
pp = 0
do k=kts,kte
zz = kte-pp
do i=its,ite
rqccuten(i,k,j)=(q2(i,zz)-qc3d(i,k,j))*rdelt
enddo
pp = pp + 1
enddo
endif
! endif
if(present(rqicuten))then
! if ( f_qi ) then
pp = 0
do k=kts,kte
zz = kte-pp
do i=its,ite
rqicuten(i,k,j)=(q3(i,zz)-qi3d(i,k,j))*rdelt
enddo
pp = pp + 1
enddo
endif
! endif
enddo
! add dqrcu
do j=jts,jte
do k=kts,kte
do i=its,ite
dqrcu(i, k, j) = 0.
enddo
enddo
enddo
! for GEOS-CHEM DO_CLOUD_CONVECTION MOD
do j=jts,jte
do k=kts,kte
do i=its,ite
if(k .eq. cloud_bot(i,j)) then
dqrcu(i, k, j) = 1.0
endif
enddo
enddo
enddo
! add dtrain
do j=jts,jte
do k=kts,kte
do i=its,ite
dtrain(i,k,j) = dtrainu(i,k,j) - dtraind(i,k,j)
if (k <= cloud_bot(i,j)) then
dtrain(i,k,j)=0.0
endif
dtrain(i,k,j) = min(dtrain(i,k,j),0.2)
dtrain(i,k,j) = max(dtrain(i,k,j),0.0)
enddo
enddo
enddo
end subroutine cu_ntiedtke
!-----------------------------------------------------------------------
subroutine ntiedtkeinit(rthcuten,rqvcuten,rqccuten,rqicuten, &
rucuten,rvcuten,rthften,rqvften, &
restart,p_qc,p_qi,p_first_scalar, &
allowed_to_read, &
ids, ide, jds, jde, kds, kde, &
ims, ime, jms, jme, kms, kme, &
its, ite, jts, jte, kts, kte)
!--------------------------------------------------------------------
implicit none
!--------------------------------------------------------------------
logical , intent(in) :: allowed_to_read,restart
integer , intent(in) :: ids, ide, jds, jde, kds, kde, &
ims, ime, jms, jme, kms, kme, &
its, ite, jts, jte, kts, kte
integer , intent(in) :: p_first_scalar, p_qi, p_qc
real, dimension( ims:ime , kms:kme , jms:jme ) , intent(out) :: &
rthcuten, &
rqvcuten, &
rqccuten, &
rqicuten, &
rucuten,rvcuten,&
rthften,rqvften
integer :: i, j, k, itf, jtf, ktf
jtf=min0(jte,jde-1)
ktf=min0(kte,kde-1)
itf=min0(ite,ide-1)
if(.not.restart)then
do j=jts,jtf
do k=kts,ktf
do i=its,itf
rthcuten(i,k,j)=0.
rqvcuten(i,k,j)=0.
rucuten(i,k,j)=0.
rvcuten(i,k,j)=0.
enddo
enddo
enddo
DO j=jts,jtf
DO k=kts,ktf
DO i=its,itf
rthften(i,k,j)=0.
rqvften(i,k,j)=0.
ENDDO
ENDDO
ENDDO
if (p_qc .ge. p_first_scalar) then
do j=jts,jtf
do k=kts,ktf
do i=its,itf
rqccuten(i,k,j)=0.
enddo
enddo
enddo
endif
if (p_qi .ge. p_first_scalar) then
do j=jts,jtf
do k=kts,ktf
do i=its,itf
rqicuten(i,k,j)=0.
enddo
enddo
enddo
endif
endif
end subroutine ntiedtkeinit
!-----------------------------------------------------------------
! level 1 subroutine 'tiecnvn'
!-----------------------------------------------------------------
subroutine tiecnvn(pu,pv,pt,pqv,pqc,pqi,pqvf,ptf,poz,pzz,pomg, &
& pap,paph,evap,hfx,zprecc,pmflxr,pmflxs,zmfu,zmfd,icbot, &
& ictop,pmfude_rate,pmfdde_rate,zdmfdp,lndj,lq,km,km1,dt,dx)
!-----------------------------------------------------------------
! this is the interface between the model and the mass
! flux convection module
!-----------------------------------------------------------------
implicit none
!
real pu(lq,km), pv(lq,km), pt(lq,km), pqv(lq,km)
real poz(lq,km), pomg(lq,km), evap(lq), zprecc(lq)
real pzz(lq,km1)
real pum1(lq,km), pvm1(lq,km), ztt(lq,km), &
& ptte(lq,km), pqte(lq,km), pvom(lq,km), pvol(lq,km), &
& pverv(lq,km), pgeo(lq,km), pap(lq,km), paph(lq,km1)
real pqhfl(lq), zqq(lq,km), &
& prsfc(lq), pssfc(lq), pcte(lq,km), &
& phhfl(lq), hfx(lq), pgeoh(lq,km1)
real ztp1(lq,km), zqp1(lq,km), ztu(lq,km), zqu(lq,km), &
& zlu(lq,km), zlude(lq,km), zmfu(lq,km), zmfd(lq,km), &
& zqsat(lq,km), pqc(lq,km), pqi(lq,km), zrain(lq)
real pqvf(lq,km), ptf(lq,km)
real pmflxr(lq,km1), pmflxs(lq,km1)
real pmfude_rate(lq,km), pmfdde_rate(lq,km)
real zdmfdp(lq,km)
! ***fox change***
! add variables pmflxr, pmflxs, zmfu, zmfd, icbot, ictop
! add variables pmfude_rate, pmfdde_rate, zdmfdp
integer icbot(lq), ictop(lq), ktype(lq), lndj(lq)
logical locum(lq)
!
real ztmst,fliq,fice,ztc,zalf,tt
integer i,j,k,lq,km,km1
real dt,dx,ztpp1
real zew,zqs,zcor
!
ztmst=dt
!
! masv flux diagnostics.
!
do j=1,lq
zrain(j)=0.0
locum(j)=.false.
prsfc(j)=0.0
pssfc(j)=0.0
pqhfl(j)=evap(j)
phhfl(j)=hfx(j)
pgeoh(j,km1)=g*pzz(j,km1)
end do
!
! convert model variables for mflux scheme
!
do k=1,km
do j=1,lq
pcte(j,k)=0.0
pvom(j,k)=0.0
pvol(j,k)=0.0
ztp1(j,k)=pt(j,k)
zqp1(j,k)=pqv(j,k)/(1.0+pqv(j,k))
pum1(j,k)=pu(j,k)
pvm1(j,k)=pv(j,k)
pverv(j,k)=pomg(j,k)
pgeo(j,k)=g*poz(j,k)
pgeoh(j,k)=g*pzz(j,k)
tt=ztp1(j,k)
zew = foeewm(tt)
zqs = zew/pap(j,k)
zqs = min(0.5,zqs)
zcor = 1./(1.-vtmpc1*zqs)
zqsat(j,k)=zqs*zcor
pqte(j,k)=pqvf(j,k)
zqq(j,k) =pqte(j,k)
ptte(j,k)=ptf(j,k)
ztt(j,k) =ptte(j,k)
end do
end do
!
!-----------------------------------------------------------------------
!* 2. call 'cumastrn'(master-routine for cumulus parameterization)
!
call cumastrn &
& (lq, km, km1, km-1, ztp1, &
& zqp1, pum1, pvm1, pverv, zqsat,&
& pqhfl, ztmst, pap, paph, pgeo, &
& ptte, pqte, pvom, pvol, prsfc,&
& pssfc, locum, pmflxr, pmflxs, &
& pmfude_rate, pmfdde_rate, zdmfdp, &
& ktype, icbot, ictop, ztu, zqu, &
& zlu, zlude, zmfu, zmfd, zrain,&
& pcte, phhfl, lndj, pgeoh, dx)
!
! to include the cloud water and cloud ice detrained from convection
!
do k=1,km
do j=1,lq
if(pcte(j,k).gt.0.) then
fliq=foealfa(ztp1(j,k))
fice=1.0-fliq
pqc(j,k)=pqc(j,k)+fliq*pcte(j,k)*ztmst
pqi(j,k)=pqi(j,k)+fice*pcte(j,k)*ztmst
endif
end do
end do
!
do k=1,km
do j=1,lq
pt(j,k)= ztp1(j,k)+(ptte(j,k)-ztt(j,k))*ztmst
zqp1(j,k)=zqp1(j,k)+(pqte(j,k)-zqq(j,k))*ztmst
pqv(j,k)=zqp1(j,k)/(1.0-zqp1(j,k))
end do
end do
do j=1,lq
zprecc(j)=amax1(0.0,(prsfc(j)+pssfc(j))*ztmst)
end do
if (lmfdudv) then
do k=1,km
do j=1,lq
pu(j,k)=pu(j,k)+pvom(j,k)*ztmst
pv(j,k)=pv(j,k)+pvol(j,k)*ztmst
end do
end do
endif
!
return
end subroutine tiecnvn
!--------------------------------------------------------------
!
! level 2 subroutines
!--------------------------------------------------------------
! subroutine cumastrn
subroutine cumastrn &
& (klon, klev, klevp1, klevm1, pten, &
& pqen, puen, pven, pverv, pqsen,&
& pqhfl, ztmst, pap, paph, pgeo, &
& ptte, pqte, pvom, pvol, prsfc,&
& pssfc, ldcum, pmflxr, pmflxs, &
& pmfude_rate, pmfdde_rate, zdmfdp, &
& ktype, kcbot, kctop, ptu, pqu, &
& plu, plude, pmfu, pmfd, prain,&
& pcte, phhfl, lndj, zgeoh, dx)
implicit none
!***fox change***
! add output variables pmflxr, pmflxs, pmfude_rate, pmfdde_rate, zdmfdp
!***cumastrn* master routine for cumulus massflux-scheme
! m.tiedtke e.c.m.w.f. 1986/1987/1989
! modifications
! y.wang i.p.r.c 2001
! c.zhang 2012
!***purpose
! -------
! this routine computes the physical tendencies of the
! prognostic variables t,q,u and v due to convective processes.
! processes considered are: convective fluxes, formation of
! precipitation, evaporation of falling rain below cloud base,
! saturated cumulus downdrafts.
!***method
! ------
! parameterization is done using a massflux-scheme.
! (1) define constants and parameters
! (2) specify values (t,q,qs...) at half levels and
! initialize updraft- and downdraft-values in 'cuinin'
! (3) calculate cloud base in 'cutypen', calculate cloud types in cutypen,
! and specify cloud base massflux
! (4) do cloud ascent in 'cuascn' in absence of downdrafts
! (5) do downdraft calculations:
! (a) determine values at lfs in 'cudlfsn'
! (b) determine moist descent in 'cuddrafn'
! (c) recalculate cloud base massflux considering the
! effect of cu-downdrafts
! (6) do final adjusments to convective fluxes in 'cuflxn',
! do evaporation in subcloud layer
! (7) calculate increments of t and q in 'cudtdqn'
! (8) calculate increments of u and v in 'cududvn'
!***externals.
! ----------
! cuinin: initializes values at vertical grid used in cu-parametr.
! cutypen: cloud bypes, 1: deep cumulus 2: shallow cumulus
! cuascn: cloud ascent for entraining plume
! cudlfsn: determines values at lfs for downdrafts
! cuddrafn:does moist descent for cumulus downdrafts
! cuflxn: final adjustments to convective fluxes (also in pbl)
! cudqdtn: updates tendencies for t and q
! cududvn: updates tendencies for u and v
!***switches.
! --------
! lmfmid=.t. midlevel convection is switched on
! lmfdd=.t. cumulus downdrafts switched on
! lmfdudv=.t. cumulus friction switched on
!***
! model parameters (defined in subroutine cuparam)
! ------------------------------------------------
! entrdd entrainment rate for cumulus downdrafts
! cmfcmax maximum massflux value allowed for
! cmfcmin minimum massflux value (for safety)
! cmfdeps fractional massflux for downdrafts at lfs
! cprcon coefficient for conversion from cloud water to rain
!***reference.
! ----------
! paper on massflux scheme (tiedtke,1989)
!-----------------------------------------------------------------
integer klev,klon,klevp1,klevm1
real pten(klon,klev), pqen(klon,klev),&
& puen(klon,klev), pven(klon,klev),&
& ptte(klon,klev), pqte(klon,klev),&
& pvom(klon,klev), pvol(klon,klev),&
& pqsen(klon,klev), pgeo(klon,klev),&
& pap(klon,klev), paph(klon,klevp1),&
& pverv(klon,klev), pqhfl(klon),&
& phhfl(klon)
real ptu(klon,klev), pqu(klon,klev),&
& plu(klon,klev), plude(klon,klev),&
& pmfu(klon,klev), pmfd(klon,klev),&
& prain(klon),&
& prsfc(klon), pssfc(klon)
real ztenh(klon,klev), zqenh(klon,klev),&
& zgeoh(klon,klevp1), zqsenh(klon,klev),&
& ztd(klon,klev), zqd(klon,klev),&
& zmfus(klon,klev), zmfds(klon,klev),&
& zmfuq(klon,klev), zmfdq(klon,klev),&
& zdmfup(klon,klev), zdmfdp(klon,klev),&
& zmful(klon,klev), zrfl(klon),&
& zuu(klon,klev), zvu(klon,klev),&
& zud(klon,klev), zvd(klon,klev),&
& zlglac(klon,klev)
real pmflxr(klon,klevp1), pmflxs(klon,klevp1)
real zhcbase(klon),&
& zmfub(klon), zmfub1(klon),&
& zdhpbl(klon)
real zsfl(klon), zdpmel(klon,klev),&
& pcte(klon,klev), zcape(klon),&
& zcape1(klon), zcape2(klon),&
& ztauc(klon), ztaubl(klon),&
& zheat(klon)
real wup(klon), zdqcv(klon)
real wbase(klon), zmfuub(klon)
real upbl(klon)
real dx
real pmfude_rate(klon,klev), pmfdde_rate(klon,klev)
real zmfuus(klon,klev), zmfdus(klon,klev)
real zuv2(klon,klev),ztenu(klon,klev),ztenv(klon,klev)
real zmfuvb(klon),zsum12(klon),zsum22(klon)
integer ilab(klon,klev), idtop(klon),&
& ictop0(klon), ilwmin(klon)
integer kdpl(klon)
integer kcbot(klon), kctop(klon),&
& ktype(klon), lndj(klon)
logical ldcum(klon)
logical loddraf(klon), llo1, llo2(klon)
! local varaiables
real zcons,zcons2,zqumqe,zdqmin,zdh,zmfmax
real zalfaw,zalv,zqalv,zc5ldcp,zc4les,zhsat,zgam,zzz,zhhat
real zpbmpt,zro,zdz,zdp,zeps,zfac,wspeed
integer jl,jk,ik
integer ikb,ikt,icum,itopm2
real ztmst,ztau,zerate,zderate,zmfa
real zmfs(klon),pmean(klev),zlon
real zduten,zdvten,ztdis,pgf_u,pgf_v
!-------------------------------------------
! 1. specify constants and parameters
!-------------------------------------------
zcons=1./(g*ztmst)
zcons2=3./(g*ztmst)
!--------------------------------------------------------------
!* 2. initialize values at vertical grid points in 'cuini'
!--------------------------------------------------------------
call cuinin &
& (klon, klev, klevp1, klevm1, pten, &
& pqen, pqsen, puen, pven, pverv,&
& pgeo, paph, zgeoh, ztenh, zqenh,&
& zqsenh, ilwmin, ptu, pqu, ztd, &
& zqd, zuu, zvu, zud, zvd, &
& pmfu, pmfd, zmfus, zmfds, zmfuq,&
& zmfdq, zdmfup, zdmfdp, zdpmel, plu, &
& plude, ilab)
!----------------------------------
!* 3.0 cloud base calculations
!----------------------------------
!* (a) determine cloud base values in 'cutypen',
! and the cumulus type 1 or 2
! -------------------------------------------
call cutypen &
& ( klon, klev, klevp1, klevm1, pqen,&
& ztenh, zqenh, zqsenh, zgeoh, paph,&
& phhfl, pqhfl, pgeo, pqsen, pap,&
& pten, lndj, ptu, pqu, ilab,&
& ldcum, kcbot, ictop0, ktype, wbase, plu, kdpl)
!* (b) assign the first guess mass flux at cloud base
! ------------------------------------------
do jl=1,klon
zdhpbl(jl)=0.0
upbl(jl) = 0.0
idtop(jl)=0
end do
do jk=2,klev
do jl=1,klon
if(jk.ge.kcbot(jl) .and. ldcum(jl)) then
zdhpbl(jl)=zdhpbl(jl)+(alv*pqte(jl,jk)+cpd*ptte(jl,jk))&
& *(paph(jl,jk+1)-paph(jl,jk))
if(lndj(jl) .eq. 0) then
wspeed = sqrt(puen(jl,jk)**2 + pven(jl,jk)**2)
upbl(jl) = upbl(jl) + wspeed*(paph(jl,jk+1)-paph(jl,jk))
end if
end if
end do
end do
do jl=1,klon
if(ldcum(jl)) then
ikb=kcbot(jl)
zmfmax = (paph(jl,ikb)-paph(jl,ikb-1))*zcons2
if(ktype(jl) == 1) then
zmfub(jl)= 0.1*zmfmax
else if ( ktype(jl) == 2 ) then
zqumqe = pqu(jl,ikb) + plu(jl,ikb) - zqenh(jl,ikb)
zdqmin = max(0.01*zqenh(jl,ikb),1.e-10)
zdh = cpd*(ptu(jl,ikb)-ztenh(jl,ikb)) + alv*zqumqe
zdh = g*max(zdh,1.e5*zdqmin)
if ( zdhpbl(jl) > 0. ) then
zmfub(jl) = zdhpbl(jl)/zdh
zmfub(jl) = min(zmfub(jl),zmfmax)
else
zmfub(jl) = 0.1*zmfmax
ldcum(jl) = .false.
end if
end if
else
zmfub(jl) = 0.
end if
end do
!------------------------------------------------------
!* 4.0 determine cloud ascent for entraining plume
!------------------------------------------------------
!* (a) do ascent in 'cuasc'in absence of downdrafts
!----------------------------------------------------------
call cuascn &
& (klon, klev, klevp1, klevm1, ztenh,&
& zqenh, puen, pven, pten, pqen,&
& pqsen, pgeo, zgeoh, pap, paph,&
& pqte, pverv, ilwmin, ldcum, zhcbase,&
& ktype, ilab, ptu, pqu, plu,&
& zuu, zvu, pmfu, zmfub,&
& zmfus, zmfuq, zmful, plude, zdmfup,&
& kcbot, kctop, ictop0, icum, ztmst,&
& zqsenh, zlglac, lndj, wup, wbase, kdpl, pmfude_rate )
!* (b) check cloud depth and change entrainment rate accordingly
! calculate precipitation rate (for downdraft calculation)
!------------------------------------------------------------------
do jl=1,klon
if ( ldcum(jl) ) then
ikb = kcbot(jl)