-
Notifications
You must be signed in to change notification settings - Fork 5
/
cu_ntiedtke.F90
3594 lines (3276 loc) · 136 KB
/
cu_ntiedtke.F90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
!=================================================================================================================
module cu_ntiedtke_common
use ccpp_kind_types,only: kind_phys
implicit none
save
real(kind=kind_phys):: alf
real(kind=kind_phys):: als
real(kind=kind_phys):: alv
real(kind=kind_phys):: cpd
real(kind=kind_phys):: g
real(kind=kind_phys):: rd
real(kind=kind_phys):: rv
real(kind=kind_phys),parameter:: t13 = 1.0/3.0
real(kind=kind_phys),parameter:: tmelt = 273.16
real(kind=kind_phys),parameter:: c1es = 610.78
real(kind=kind_phys),parameter:: c3les = 17.2693882
real(kind=kind_phys),parameter:: c3ies = 21.875
real(kind=kind_phys),parameter:: c4les = 35.86
real(kind=kind_phys),parameter:: c4ies = 7.66
real(kind=kind_phys),parameter:: rtwat = tmelt
real(kind=kind_phys),parameter:: rtber = tmelt-5.
real(kind=kind_phys),parameter:: rtice = tmelt-23.
integer,parameter:: momtrans = 2
real(kind=kind_phys),parameter:: entrdd = 2.0e-4
real(kind=kind_phys),parameter:: cmfcmax = 1.0
real(kind=kind_phys),parameter:: cmfcmin = 1.e-10
real(kind=kind_phys),parameter:: cmfdeps = 0.30
real(kind=kind_phys),parameter:: zdnoprc = 2.0e4
real(kind=kind_phys),parameter:: cprcon = 1.4e-3
real(kind=kind_phys),parameter:: pgcoef = 0.7
real(kind=kind_phys):: rcpd
real(kind=kind_phys):: c2es
real(kind=kind_phys):: c5les
real(kind=kind_phys):: c5ies
real(kind=kind_phys):: r5alvcp
real(kind=kind_phys):: r5alscp
real(kind=kind_phys):: ralvdcp
real(kind=kind_phys):: ralsdcp
real(kind=kind_phys):: ralfdcp
real(kind=kind_phys):: vtmpc1
real(kind=kind_phys):: zrg
logical,parameter:: nonequil = .true.
logical,parameter:: lmfpen = .true.
logical,parameter:: lmfmid = .true.
logical,parameter:: lmfscv = .true.
logical,parameter:: lmfdd = .true.
logical,parameter:: lmfdudv = .true.
!=================================================================================================================
end module cu_ntiedtke_common
!=================================================================================================================
module cu_ntiedtke
use ccpp_kind_types,only: kind_phys
use cu_ntiedtke_common
implicit none
private
public:: cu_ntiedtke_run, &
cu_ntiedtke_init, &
cu_ntiedtke_finalize
contains
!=================================================================================================================
!>\section arg_table_cu_ntiedtke_init
!!\html\include cu_ntiedtke_init.html
!!
subroutine cu_ntiedtke_init(con_cp,con_rd,con_rv,con_xlv,con_xls,con_xlf,con_grav,errmsg,errflg)
!=================================================================================================================
!input arguments:
real(kind=kind_phys),intent(in):: &
con_cp, &
con_rd, &
con_rv, &
con_xlv, &
con_xls, &
con_xlf, &
con_grav
!--- output arguments:
character(len=*),intent(out):: errmsg
integer,intent(out):: errflg
!-----------------------------------------------------------------------------------------------------------------
alf = con_xlf
als = con_xls
alv = con_xlv
cpd = con_cp
g = con_grav
rd = con_rd
rv = con_rv
rcpd = 1.0/con_cp
c2es = c1es*con_rd/con_rv
c5les = c3les*(tmelt-c4les)
c5ies = c3ies*(tmelt-c4ies)
r5alvcp = c5les*con_xlv*rcpd
r5alscp = c5ies*con_xls*rcpd
ralvdcp = con_xlv*rcpd
ralsdcp = con_xls*rcpd
ralfdcp = con_xlf*rcpd
vtmpc1 = con_rv/con_rd-1.0
zrg = 1.0/con_grav
errmsg = 'cu_ntiedtke_init OK'
errflg = 0
end subroutine cu_ntiedtke_init
!=================================================================================================================
!>\section arg_table_cu_ntiedtke_finalize
!!\html\include cu_ntiedtke_finalize.html
!!
subroutine cu_ntiedtke_finalize(errmsg,errflg)
!=================================================================================================================
!--- output arguments:
character(len=*),intent(out):: errmsg
integer,intent(out):: errflg
!-----------------------------------------------------------------------------------------------------------------
errmsg = 'cu_ntiedtke_finalize OK'
errflg = 0
end subroutine cu_ntiedtke_finalize
!=================================================================================================================
!>\section arg_table_cu_ntiedtke_run
!!\html\include cu_ntiedtke_run.html
!!
! level 1 subroutine 'cu_ntiedkte_run'
subroutine cu_ntiedtke_run(pu,pv,pt,pqv,pqc,pqi,pqvf,ptf,poz,pzz,pomg, &
& pap,paph,evap,hfx,zprecc,lndj,lq,km,km1,dt,dx,errmsg,errflg)
!=================================================================================================================
! this is the interface between the model and the mass flux convection module
! m.tiedtke e.c.m.w.f. 1989
! j.morcrette 1992
!--------------------------------------------
! modifications
! C. zhang & Yuqing Wang 2011-2017
!
! modified from IPRC IRAM - yuqing wang, university of hawaii (ICTP REGCM4.4).
!
! The current version is stable. There are many updates to the old Tiedtke scheme (cu_physics=6)
! update notes:
! the new Tiedtke scheme is similar to the Tiedtke scheme used in REGCM4 and ECMWF cy40r1.
! the major differences to the old Tiedtke (cu_physics=6) scheme are,
! (a) New trigger functions for deep and shallow convections (Jakob and Siebesma 2003;
! Bechtold et al. 2004, 2008, 2014).
! (b) Non-equilibrium situations are considered in the closure for deep convection
! (Bechtold et al. 2014).
! (c) New convection time scale for the deep convection closure (Bechtold et al. 2008).
! (d) New entrainment and detrainment rates for all convection types (Bechtold et al. 2008).
! (e) New formula for the conversion from cloud water/ice to rain/snow (Sundqvist 1978).
! (f) Different way to include cloud scale pressure gradients (Gregory et al. 1997;
! Wu and Yanai 1994)
!
! other reference: tiedtke (1989, mwr, 117, 1779-1800)
! IFS documentation - cy33r1, cy37r2, cy38r1, cy40r1
!
! Note for climate simulation of Tropical Cyclones
! This version of Tiedtke scheme was tested with YSU PBL scheme, RRTMG radation
! schemes, and WSM6 microphysics schemes, at horizontal resolution around 20 km
! Set: momtrans = 2.
! pgcoef = 0.7 to 1.0 is good depends on the basin
! nonequil = .false.
! Note for the diurnal simulation of precipitaton
! When nonequil = .true., the CAPE is relaxed toward to a value from PBL
! It can improve the diurnal precipitation over land.
!--- input arguments:
integer,intent(in):: lq,km,km1
integer,intent(in),dimension(:):: lndj
real(kind=kind_phys),intent(in):: dt
real(kind=kind_phys),intent(in),dimension(:):: dx
real(kind=kind_phys),intent(in),dimension(:):: evap,hfx
real(kind=kind_phys),intent(in),dimension(:,:):: pqvf,ptf
real(kind=kind_phys),intent(in),dimension(:,:):: poz,pomg,pap
real(kind=kind_phys),intent(in),dimension(:,:):: pzz,paph
!--- inout arguments:
real(kind=kind_phys),intent(inout),dimension(:):: zprecc
real(kind=kind_phys),intent(inout),dimension(:,:):: pu,pv,pt,pqv,pqc,pqi
!--- output arguments:
character(len=*),intent(out):: errmsg
integer,intent(out):: errflg
!--- local variables and arrays:
logical,dimension(lq):: locum
integer:: i,j,k
integer,dimension(lq):: icbot,ictop,ktype
real(kind=kind_phys):: ztmst,fliq,fice,ztc,zalf,tt
real(kind=kind_phys):: ztpp1,zew,zqs,zcor
real(kind=kind_phys):: dxref
real(kind=kind_phys),dimension(lq):: pqhfl,prsfc,pssfc,phhfl,zrain
real(kind=kind_phys),dimension(lq):: scale_fac,scale_fac2
real(kind=kind_phys),dimension(lq,km):: pum1,pvm1,ztt,ptte,pqte,pvom,pvol,pverv,pgeo
real(kind=kind_phys),dimension(lq,km):: zqq,pcte
real(kind=kind_phys),dimension(lq,km):: ztp1,zqp1,ztu,zqu,zlu,zlude,zmfu,zmfd,zqsat
real(kind=kind_phys),dimension(lq,km1):: pgeoh
!-----------------------------------------------------------------------------------------------------------------
!
ztmst=dt
!
! set scale-dependency factor when dx is < 15 km
!
dxref = 15000.
do j=1,lq
if (dx(j).lt.dxref) then
scale_fac(j) = (1.06133+log(dxref/dx(j)))**3
scale_fac2(j) = scale_fac(j)**0.5
else
scale_fac(j) = 1.+1.33e-5*dx(j)
scale_fac2(j) = 1.
end if
end do
!
! 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, &
& ktype, icbot, ictop, ztu, zqu, &
& zlu, zlude, zmfu, zmfd, zrain, &
& pcte, phhfl, lndj, pgeoh, dx, &
& scale_fac, scale_fac2)
!
! 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
!
errmsg = 'cu_ntiedtke_run OK'
errflg = 0
!
return
end subroutine cu_ntiedtke_run
!#############################################################
!
! 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, &
& ktype, kcbot, kctop, ptu, pqu, &
& plu, plude, pmfu, pmfd, prain, &
& pcte, phhfl, lndj, zgeoh, dx, &
& scale_fac, scale_fac2)
implicit none
!
!***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)
!-----------------------------------------------------------------
!--- input arguments:
integer,intent(in):: klev,klon,klevp1,klevm1
integer,intent(in),dimension(klon):: lndj
real(kind=kind_phys),intent(in):: ztmst
real(kind=kind_phys),intent(in),dimension(klon):: dx
real(kind=kind_phys),intent(in),dimension(klon):: pqhfl,phhfl
real(kind=kind_phys),intent(in),dimension(klon):: scale_fac,scale_fac2
real(kind=kind_phys),intent(in),dimension(klon,klev):: pten,pqen,puen,pven,pverv
real(kind=kind_phys),intent(in),dimension(klon,klev):: pap,pgeo
real(kind=kind_phys),intent(in),dimension(klon,klevp1):: paph,zgeoh
!--- inout arguments:
integer,intent(inout),dimension(klon):: ktype,kcbot,kctop
logical,intent(inout),dimension(klon):: ldcum
real(kind=kind_phys),intent(inout),dimension(klon):: pqsen
real(kind=kind_phys),intent(inout),dimension(klon):: prsfc,pssfc,prain
real(kind=kind_phys),intent(inout),dimension(klon,klev):: pcte,ptte,pqte,pvom,pvol
real(kind=kind_phys),intent(inout),dimension(klon,klev):: ptu,pqu,plu,plude,pmfu,pmfd
!--- local variables and arrays:
logical:: llo1
logical,dimension(klon):: loddraf,llo2
integer:: jl,jk,ik
integer:: ikb,ikt,icum,itopm2
integer,dimension(klon):: kdpl,idtop,ictop0,ilwmin
integer,dimension(klon,klev):: ilab
real(kind=kind_phys):: zcons,zcons2,zqumqe,zdqmin,zdh,zmfmax
real(kind=kind_phys):: zalfaw,zalv,zqalv,zc5ldcp,zc4les,zhsat,zgam,zzz,zhhat
real(kind=kind_phys):: zpbmpt,zro,zdz,zdp,zeps,zfac,wspeed
real(kind=kind_phys):: zduten,zdvten,ztdis,pgf_u,pgf_v
real(kind=kind_phys):: zlon
real(kind=kind_phys):: ztau,zerate,zderate,zmfa
real(kind=kind_phys),dimension(klon):: zmfs
real(kind=kind_phys),dimension(klon):: zsfl,zcape,zcape1,zcape2,ztauc,ztaubl,zheat
real(kind=kind_phys),dimension(klon):: wup,zdqcv
real(kind=kind_phys),dimension(klon):: wbase,zmfuub
real(kind=kind_phys),dimension(klon):: upbl
real(kind=kind_phys),dimension(klon):: zhcbase,zmfub,zmfub1,zdhpbl
real(kind=kind_phys),dimension(klon):: zmfuvb,zsum12,zsum22
real(kind=kind_phys),dimension(klon):: zrfl
real(kind=kind_phys),dimension(klev):: pmean
real(kind=kind_phys),dimension(klon,klev):: pmfude_rate,pmfdde_rate
real(kind=kind_phys),dimension(klon,klev):: zdpmel
real(kind=kind_phys),dimension(klon,klev):: zmfuus,zmfdus,zuv2,ztenu,ztenv
real(kind=kind_phys),dimension(klon,klev):: ztenh,zqenh,zqsenh,ztd,zqd
real(kind=kind_phys),dimension(klon,klev):: zmfus,zmfds,zmfuq,zmfdq,zdmfup,zdmfdp,zmful
real(kind=kind_phys),dimension(klon,klev):: zuu,zvu,zud,zvd,zlglac
real(kind=kind_phys),dimension(klon,klevp1):: pmflxr,pmflxs
!-------------------------------------------
! 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)
itopm2 = kctop(jl)
zpbmpt = paph(jl,ikb) - paph(jl,itopm2)
if ( ktype(jl) == 1 .and. zpbmpt < zdnoprc ) ktype(jl) = 2
if ( ktype(jl) == 2 .and. zpbmpt >= zdnoprc ) ktype(jl) = 1
ictop0(jl) = kctop(jl)
end if
zrfl(jl)=zdmfup(jl,1)
end do
do jk=2,klev
do jl=1,klon
zrfl(jl)=zrfl(jl)+zdmfup(jl,jk)
end do
end do
do jk = 1,klev
do jl = 1,klon
pmfd(jl,jk) = 0.
zmfds(jl,jk) = 0.
zmfdq(jl,jk) = 0.
zdmfdp(jl,jk) = 0.
zdpmel(jl,jk) = 0.
end do
end do
!-----------------------------------------
!* 6.0 cumulus downdraft calculations
!-----------------------------------------
if(lmfdd) then
!* (a) determine lfs in 'cudlfsn'
!--------------------------------------
call cudlfsn &
& (klon, klev,&
& kcbot, kctop, lndj, ldcum, &
& ztenh, zqenh, puen, pven, &
& pten, pqsen, pgeo, &
& zgeoh, paph, ptu, pqu, plu, &
& zuu, zvu, zmfub, zrfl, &
& ztd, zqd, zud, zvd, &
& pmfd, zmfds, zmfdq, zdmfdp, &
& idtop, loddraf)
!* (b) determine downdraft t,q and fluxes in 'cuddrafn'
!------------------------------------------------------------
call cuddrafn &
& (klon, klev, loddraf, &
& ztenh, zqenh, puen, pven, &
& pgeo, zgeoh, paph, zrfl, &
& ztd, zqd, zud, zvd, pmfu, &
& pmfd, zmfds, zmfdq, zdmfdp, pmfdde_rate)
!-----------------------------------------------------------
end if
!
!-----------------------------------------------------------------------
!* 6.0 closure and clean work
! ------
!-- 6.1 recalculate cloud base massflux from a cape closure
! for deep convection (ktype=1)
!
do jl=1,klon
if(ldcum(jl) .and. ktype(jl) .eq. 1) then
ikb = kcbot(jl)
ikt = kctop(jl)
zheat(jl)=0.0
zcape(jl)=0.0
zcape1(jl)=0.0
zcape2(jl)=0.0
zmfub1(jl)=zmfub(jl)
ztauc(jl) = (zgeoh(jl,ikt)-zgeoh(jl,ikb)) / &
((2.+ min(15.0,wup(jl)))*g)
if(lndj(jl) .eq. 0) then
upbl(jl) = 2.+ upbl(jl)/(paph(jl,klev+1)-paph(jl,ikb))
ztaubl(jl) = (zgeoh(jl,ikb)-zgeoh(jl,klev+1))/(g*upbl(jl))
ztaubl(jl) = min(300., ztaubl(jl))
else
ztaubl(jl) = ztauc(jl)
end if
end if
end do
!
do jk = 1 , klev
do jl = 1 , klon
llo1 = ldcum(jl) .and. ktype(jl) .eq. 1
if ( llo1 .and. jk <= kcbot(jl) .and. jk > kctop(jl) ) then
ikb = kcbot(jl)
zdz = pgeo(jl,jk-1)-pgeo(jl,jk)
zdp = pap(jl,jk)-pap(jl,jk-1)
zheat(jl) = zheat(jl) + ((pten(jl,jk-1)-pten(jl,jk)+zdz*rcpd) / &
ztenh(jl,jk)+vtmpc1*(pqen(jl,jk-1)-pqen(jl,jk))) * &
(g*(pmfu(jl,jk)+pmfd(jl,jk)))
zcape1(jl) = zcape1(jl) + ((ptu(jl,jk)-ztenh(jl,jk))/ztenh(jl,jk) + &
vtmpc1*(pqu(jl,jk)-zqenh(jl,jk))-plu(jl,jk))*zdp
end if
if ( llo1 .and. jk >= kcbot(jl) ) then
if((paph(jl,klev+1)-paph(jl,kdpl(jl)))<50.e2) then
zdp = paph(jl,jk+1)-paph(jl,jk)
zcape2(jl) = zcape2(jl) + ztaubl(jl)* &
((1.+vtmpc1*pqen(jl,jk))*ptte(jl,jk)+vtmpc1*pten(jl,jk)*pqte(jl,jk))*zdp
end if
end if
end do
end do
do jl=1,klon
if(ldcum(jl).and.ktype(jl).eq.1) then
ikb = kcbot(jl)
ikt = kctop(jl)
ztauc(jl) = max(ztmst,ztauc(jl))
ztauc(jl) = max(360.,ztauc(jl))
ztauc(jl) = min(10800.,ztauc(jl))
ztau = ztauc(jl) * scale_fac(jl)
if(nonequil) then
zcape2(jl)= max(0.,zcape2(jl))
zcape(jl) = max(0.,min(zcape1(jl)-zcape2(jl),5000.))
else
zcape(jl) = max(0.,min(zcape1(jl),5000.))
end if
zheat(jl) = max(1.e-4,zheat(jl))
zmfub1(jl) = (zcape(jl)*zmfub(jl))/(zheat(jl)*ztau)
zmfub1(jl) = max(zmfub1(jl),0.001)
zmfmax=(paph(jl,ikb)-paph(jl,ikb-1))*zcons2
zmfub1(jl)=min(zmfub1(jl),zmfmax)
end if
end do
!
!* 6.2 recalculate convective fluxes due to effect of
! downdrafts on boundary layer moist static energy budget (ktype=2)
!--------------------------------------------------------
do jl=1,klon
if(ldcum(jl) .and. ktype(jl) .eq. 2) then
ikb=kcbot(jl)
if(pmfd(jl,ikb).lt.0.0 .and. loddraf(jl)) then
zeps=-pmfd(jl,ikb)/max(zmfub(jl),cmfcmin)
else
zeps=0.
endif
zqumqe=pqu(jl,ikb)+plu(jl,ikb)- &
& zeps*zqd(jl,ikb)-(1.-zeps)*zqenh(jl,ikb)
zdqmin=max(0.01*zqenh(jl,ikb),cmfcmin)
zmfmax=(paph(jl,ikb)-paph(jl,ikb-1))*zcons2
! using moist static engergy closure instead of moisture closure
zdh=cpd*(ptu(jl,ikb)-zeps*ztd(jl,ikb)- &
& (1.-zeps)*ztenh(jl,ikb))+alv*zqumqe
zdh=g*max(zdh,1.e5*zdqmin)
if(zdhpbl(jl).gt.0.)then
zmfub1(jl)=zdhpbl(jl)/zdh
else
zmfub1(jl) = zmfub(jl)
end if
zmfub1(jl) = zmfub1(jl)/scale_fac2(jl)
zmfub1(jl) = min(zmfub1(jl),zmfmax)
end if
!* 6.3 mid-level convection - nothing special
!---------------------------------------------------------
if(ldcum(jl) .and. ktype(jl) .eq. 3 ) then
zmfub1(jl) = zmfub(jl)
end if
end do
!* 6.4 scaling the downdraft mass flux
!---------------------------------------------------------
do jk=1,klev
do jl=1,klon
if( ldcum(jl) ) then
zfac=zmfub1(jl)/max(zmfub(jl),cmfcmin)
pmfd(jl,jk)=pmfd(jl,jk)*zfac
zmfds(jl,jk)=zmfds(jl,jk)*zfac
zmfdq(jl,jk)=zmfdq(jl,jk)*zfac
zdmfdp(jl,jk)=zdmfdp(jl,jk)*zfac
pmfdde_rate(jl,jk) = pmfdde_rate(jl,jk)*zfac
end if
end do
end do
!* 6.5 scaling the updraft mass flux
! --------------------------------------------------------
do jl = 1,klon
if ( ldcum(jl) ) zmfs(jl) = zmfub1(jl)/max(cmfcmin,zmfub(jl))
end do
do jk = 2 , klev
do jl = 1,klon
if ( ldcum(jl) .and. jk >= kctop(jl)-1 ) then
ikb = kcbot(jl)
if ( jk>ikb ) then
zdz = ((paph(jl,klev+1)-paph(jl,jk))/(paph(jl,klev+1)-paph(jl,ikb)))
pmfu(jl,jk) = pmfu(jl,ikb)*zdz
end if
zmfmax = (paph(jl,jk)-paph(jl,jk-1))*zcons2
if ( pmfu(jl,jk)*zmfs(jl) > zmfmax ) then
zmfs(jl) = min(zmfs(jl),zmfmax/pmfu(jl,jk))
end if
end if
end do
end do
do jk = 2 , klev
do jl = 1,klon
if ( ldcum(jl) .and. jk <= kcbot(jl) .and. jk >= kctop(jl)-1 ) then
pmfu(jl,jk) = pmfu(jl,jk)*zmfs(jl)
zmfus(jl,jk) = zmfus(jl,jk)*zmfs(jl)
zmfuq(jl,jk) = zmfuq(jl,jk)*zmfs(jl)
zmful(jl,jk) = zmful(jl,jk)*zmfs(jl)
zdmfup(jl,jk) = zdmfup(jl,jk)*zmfs(jl)
plude(jl,jk) = plude(jl,jk)*zmfs(jl)
pmfude_rate(jl,jk) = pmfude_rate(jl,jk)*zmfs(jl)
end if
end do
end do
!* 6.6 if ktype = 2, kcbot=kctop is not allowed
! ---------------------------------------------------
do jl = 1,klon
if ( ktype(jl) == 2 .and. &
kcbot(jl) == kctop(jl) .and. kcbot(jl) >= klev-1 ) then
ldcum(jl) = .false.
ktype(jl) = 0
end if
end do
if ( .not. lmfscv .or. .not. lmfpen ) then
do jl = 1,klon
llo2(jl) = .false.
if ( (.not. lmfscv .and. ktype(jl) == 2) .or. &
(.not. lmfpen .and. ktype(jl) == 1) ) then
llo2(jl) = .true.
ldcum(jl) = .false.
end if
end do
end if
!* 6.7 set downdraft mass fluxes to zero above cloud top
!----------------------------------------------------
do jl = 1,klon
if ( loddraf(jl) .and. idtop(jl) <= kctop(jl) ) then
idtop(jl) = kctop(jl) + 1
end if
end do
do jk = 2 , klev
do jl = 1,klon
if ( loddraf(jl) ) then
if ( jk < idtop(jl) ) then
pmfd(jl,jk) = 0.
zmfds(jl,jk) = 0.
zmfdq(jl,jk) = 0.
pmfdde_rate(jl,jk) = 0.
zdmfdp(jl,jk) = 0.
else if ( jk == idtop(jl) ) then
pmfdde_rate(jl,jk) = 0.
end if
end if
end do
end do
!----------------------------------------------------------
!* 7.0 determine final convective fluxes in 'cuflx'
!----------------------------------------------------------
call cuflxn &
& ( klon, klev, ztmst &
& , pten, pqen, pqsen, ztenh, zqenh &
& , paph, pap, zgeoh, lndj, ldcum &
& , kcbot, kctop, idtop, itopm2 &
& , ktype, loddraf &
& , pmfu, pmfd, zmfus, zmfds &
& , zmfuq, zmfdq, zmful, plude &
& , zdmfup, zdmfdp, zdpmel, zlglac &
& , prain, pmfdde_rate, pmflxr, pmflxs )
! some adjustments needed
do jl=1,klon
zmfs(jl) = 1.
zmfuub(jl)=0.
end do
do jk = 2 , klev
do jl = 1,klon
if ( loddraf(jl) .and. jk >= idtop(jl)-1 ) then
zmfmax = pmfu(jl,jk)*0.98
if ( pmfd(jl,jk)+zmfmax+1.e-15 < 0. ) then
zmfs(jl) = min(zmfs(jl),-zmfmax/pmfd(jl,jk))
end if
end if
end do
end do
do jk = 2 , klev
do jl = 1 , klon
if ( zmfs(jl) < 1. .and. jk >= idtop(jl)-1 ) then
pmfd(jl,jk) = pmfd(jl,jk)*zmfs(jl)
zmfds(jl,jk) = zmfds(jl,jk)*zmfs(jl)
zmfdq(jl,jk) = zmfdq(jl,jk)*zmfs(jl)
pmfdde_rate(jl,jk) = pmfdde_rate(jl,jk)*zmfs(jl)
zmfuub(jl) = zmfuub(jl) - (1.-zmfs(jl))*zdmfdp(jl,jk)
pmflxr(jl,jk+1) = pmflxr(jl,jk+1) + zmfuub(jl)
zdmfdp(jl,jk) = zdmfdp(jl,jk)*zmfs(jl)
end if
end do
end do
do jk = 2 , klev - 1
do jl = 1, klon
if ( loddraf(jl) .and. jk >= idtop(jl)-1 ) then
zerate = -pmfd(jl,jk) + pmfd(jl,jk-1) + pmfdde_rate(jl,jk)
if ( zerate < 0. ) then
pmfdde_rate(jl,jk) = pmfdde_rate(jl,jk) - zerate
end if
end if
if ( ldcum(jl) .and. jk >= kctop(jl)-1 ) then
zerate = pmfu(jl,jk) - pmfu(jl,jk+1) + pmfude_rate(jl,jk)
if ( zerate < 0. ) then
pmfude_rate(jl,jk) = pmfude_rate(jl,jk) - zerate
end if
zdmfup(jl,jk) = pmflxr(jl,jk+1) + pmflxs(jl,jk+1) - &
pmflxr(jl,jk) - pmflxs(jl,jk)
zdmfdp(jl,jk) = 0.
end if
end do
end do
! avoid negative humidities at ddraught top
do jl = 1,klon
if ( loddraf(jl) ) then
jk = idtop(jl)
ik = min(jk+1,klev)
if ( zmfdq(jl,jk) < 0.3*zmfdq(jl,ik) ) then
zmfdq(jl,jk) = 0.3*zmfdq(jl,ik)
end if
end if
end do
! avoid negative humidities near cloud top because gradient of precip flux
! and detrainment / liquid water flux are too large
do jk = 2 , klev
do jl = 1, klon
if ( ldcum(jl) .and. jk >= kctop(jl)-1 .and. jk < kcbot(jl) ) then
zdz = ztmst*g/(paph(jl,jk+1)-paph(jl,jk))
zmfa = zmfuq(jl,jk+1) + zmfdq(jl,jk+1) - &
zmfuq(jl,jk) - zmfdq(jl,jk) + &
zmful(jl,jk+1) - zmful(jl,jk) + zdmfup(jl,jk)
zmfa = (zmfa-plude(jl,jk))*zdz
if ( pqen(jl,jk)+zmfa < 0. ) then
plude(jl,jk) = plude(jl,jk) + 2.*(pqen(jl,jk)+zmfa)/zdz
end if
if ( plude(jl,jk) < 0. ) plude(jl,jk) = 0.
end if
if ( .not. ldcum(jl) ) pmfude_rate(jl,jk) = 0.
if ( abs(pmfd(jl,jk-1)) < 1.0e-20 ) pmfdde_rate(jl,jk) = 0.
end do
end do
do jl=1,klon
prsfc(jl) = pmflxr(jl,klev+1)
pssfc(jl) = pmflxs(jl,klev+1)
end do
!----------------------------------------------------------------
!* 8.0 update tendencies for t and q in subroutine cudtdq
!----------------------------------------------------------------
call cudtdqn(klon,klev,itopm2,kctop,idtop,ldcum,loddraf, &
ztmst,paph,zgeoh,pgeo,pten,ztenh,pqen,zqenh,pqsen, &
zlglac,plude,pmfu,pmfd,zmfus,zmfds,zmfuq,zmfdq,zmful, &
zdmfup,zdmfdp,zdpmel,ptte,pqte,pcte)
!----------------------------------------------------------------
!* 9.0 update tendencies for u and u in subroutine cududv
!----------------------------------------------------------------
if(lmfdudv) then
do jk = klev-1 , 2 , -1
ik = jk + 1
do jl = 1,klon
if ( ldcum(jl) ) then
if ( jk == kcbot(jl) .and. ktype(jl) < 3 ) then
ikb = kdpl(jl)
zuu(jl,jk) = puen(jl,ikb-1)
zvu(jl,jk) = pven(jl,ikb-1)
else if ( jk == kcbot(jl) .and. ktype(jl) == 3 ) then
zuu(jl,jk) = puen(jl,jk-1)
zvu(jl,jk) = pven(jl,jk-1)
end if
if ( jk < kcbot(jl) .and. jk >= kctop(jl) ) then
if(momtrans .eq. 1)then
zfac = 0.
if ( ktype(jl) == 1 .or. ktype(jl) == 3 ) zfac = 2.
if ( ktype(jl) == 1 .and. jk <= kctop(jl)+2 ) zfac = 3.
zerate = pmfu(jl,jk) - pmfu(jl,ik) + &
(1.+zfac)*pmfude_rate(jl,jk)
zderate = (1.+zfac)*pmfude_rate(jl,jk)
zmfa = 1./max(cmfcmin,pmfu(jl,jk))
zuu(jl,jk) = (zuu(jl,ik)*pmfu(jl,ik) + &
zerate*puen(jl,jk)-zderate*zuu(jl,ik))*zmfa
zvu(jl,jk) = (zvu(jl,ik)*pmfu(jl,ik) + &
zerate*pven(jl,jk)-zderate*zvu(jl,ik))*zmfa
else
pgf_u = -pgcoef*0.5*(pmfu(jl,ik)*(puen(jl,ik)-puen(jl,jk))+&
pmfu(jl,jk)*(puen(jl,jk)-puen(jl,jk-1)))
pgf_v = -pgcoef*0.5*(pmfu(jl,ik)*(pven(jl,ik)-pven(jl,jk))+&
pmfu(jl,jk)*(pven(jl,jk)-pven(jl,jk-1)))
zerate = pmfu(jl,jk) - pmfu(jl,ik) + pmfude_rate(jl,jk)
zderate = pmfude_rate(jl,jk)
zmfa = 1./max(cmfcmin,pmfu(jl,jk))
zuu(jl,jk) = (zuu(jl,ik)*pmfu(jl,ik) + &
zerate*puen(jl,jk)-zderate*zuu(jl,ik)+pgf_u)*zmfa
zvu(jl,jk) = (zvu(jl,ik)*pmfu(jl,ik) + &
zerate*pven(jl,jk)-zderate*zvu(jl,ik)+pgf_v)*zmfa
end if
end if
end if
end do
end do
if(lmfdd) then
do jk = 3 , klev
ik = jk - 1
do jl = 1,klon
if ( ldcum(jl) ) then
if ( jk == idtop(jl) ) then
zud(jl,jk) = 0.5*(zuu(jl,jk)+puen(jl,ik))
zvd(jl,jk) = 0.5*(zvu(jl,jk)+pven(jl,ik))
else if ( jk > idtop(jl) ) then
zerate = -pmfd(jl,jk) + pmfd(jl,ik) + pmfdde_rate(jl,jk)
zmfa = 1./min(-cmfcmin,pmfd(jl,jk))
zud(jl,jk) = (zud(jl,ik)*pmfd(jl,ik) - &
zerate*puen(jl,ik)+pmfdde_rate(jl,jk)*zud(jl,ik))*zmfa
zvd(jl,jk) = (zvd(jl,ik)*pmfd(jl,ik) - &
zerate*pven(jl,ik)+pmfdde_rate(jl,jk)*zvd(jl,ik))*zmfa
end if
end if
end do
end do
end if
! --------------------------------------------------
! rescale massfluxes for stability in Momentum
!------------------------------------------------------------------------
zmfs(:) = 1.
do jk = 2 , klev
do jl = 1, klon
if ( ldcum(jl) .and. jk >= kctop(jl)-1 ) then
zmfmax = (paph(jl,jk)-paph(jl,jk-1))*zcons