-
Notifications
You must be signed in to change notification settings - Fork 3
/
mp_wsm6.F90
2449 lines (2311 loc) · 90.6 KB
/
mp_wsm6.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 mp_wsm6
use ccpp_kind_types,only: kind_phys
use module_libmassv,only: vrec,vsqrt
use mp_radar
implicit none
private
public:: mp_wsm6_run, &
mp_wsm6_init, &
mp_wsm6_finalize, &
refl10cm_wsm6
real(kind=kind_phys),parameter,private:: dtcldcr = 120. ! maximum time step for minor loops
real(kind=kind_phys),parameter,private:: n0r = 8.e6 ! intercept parameter rain
!real(kind=kind_phys),parameter,private:: n0g = 4.e6 ! intercept parameter graupel
real(kind=kind_phys),parameter,private:: avtr = 841.9 ! a constant for terminal velocity of rain
real(kind=kind_phys),parameter,private:: bvtr = 0.8 ! a constant for terminal velocity of rain
real(kind=kind_phys),parameter,private:: r0 = .8e-5 ! 8 microm in contrast to 10 micro m
real(kind=kind_phys),parameter,private:: peaut = .55 ! collection efficiency
real(kind=kind_phys),parameter,private:: xncr = 3.e8 ! maritime cloud in contrast to 3.e8 in tc80
real(kind=kind_phys),parameter,private:: xmyu = 1.718e-5 ! the dynamic viscosity kgm-1s-1
real(kind=kind_phys),parameter,private:: avts = 11.72 ! a constant for terminal velocity of snow
real(kind=kind_phys),parameter,private:: bvts = .41 ! a constant for terminal velocity of snow
!real(kind=kind_phys),parameter,private:: avtg = 330. ! a constant for terminal velocity of graupel
!real(kind=kind_phys),parameter,private:: bvtg = 0.8 ! a constant for terminal velocity of graupel
!real(kind=kind_phys),parameter,private:: deng = 500. ! density of graupel ! set later with hail_opt
real(kind=kind_phys),parameter,private:: lamdarmax = 8.e4 ! limited maximum value for slope parameter of rain
real(kind=kind_phys),parameter,private:: lamdasmax = 1.e5 ! limited maximum value for slope parameter of snow
!real(kind=kind_phys),parameter,private:: lamdagmax = 6.e4 ! limited maximum value for slope parameter of graupel
real(kind=kind_phys),parameter,private:: dicon = 11.9 ! constant for the cloud-ice diamter
real(kind=kind_phys),parameter,private:: dimax = 500.e-6 ! limited maximum value for the cloud-ice diamter
real(kind=kind_phys),parameter,private:: pfrz1 = 100. ! constant in Biggs freezing
real(kind=kind_phys),parameter,private:: pfrz2 = 0.66 ! constant in Biggs freezing
real(kind=kind_phys),parameter,private:: qcrmin = 1.e-9 ! minimun values for qr, qs, and qg
real(kind=kind_phys),parameter,private:: eacrc = 1.0 ! Snow/cloud-water collection efficiency
real(kind=kind_phys),parameter,private:: dens = 100.0 ! Density of snow
real(kind=kind_phys),parameter,private:: qs0 = 6.e-4 ! threshold amount for aggretion to occur
real(kind=kind_phys),parameter,public :: n0smax = 1.e11 ! maximum n0s (t=-90C unlimited)
real(kind=kind_phys),parameter,public :: n0s = 2.e6 ! temperature dependent intercept parameter snow
real(kind=kind_phys),parameter,public :: alpha = .12 ! .122 exponen factor for n0s
real(kind=kind_phys),save:: &
qc0,qck1, &
bvtr1,bvtr2,bvtr3,bvtr4,g1pbr, &
g3pbr,g4pbr,g5pbro2,pvtr,eacrr,pacrr, &
bvtr6,g6pbr, &
precr1,precr2,roqimax,bvts1, &
bvts2,bvts3,bvts4,g1pbs,g3pbs,g4pbs, &
n0g,avtg,bvtg,deng,lamdagmax, & !RAS13.3 - set these in wsm6init
g5pbso2,pvts,pacrs,precs1,precs2,pidn0r, &
xlv1,pacrc,pi, &
bvtg1,bvtg2,bvtg3,bvtg4,g1pbg, &
g3pbg,g4pbg,g5pbgo2,pvtg,pacrg, &
precg1,precg2,pidn0g, &
rslopermax,rslopesmax,rslopegmax, &
rsloperbmax,rslopesbmax,rslopegbmax, &
rsloper2max,rslopes2max,rslopeg2max, &
rsloper3max,rslopes3max,rslopeg3max
real(kind=kind_phys),public,save:: pidn0s,pidnc
contains
!=================================================================================================================
!>\section arg_table_mp_wsm6_init
!!\html\include mp_wsm6_init.html
!!
subroutine mp_wsm6_init(den0,denr,dens,cl,cpv,hail_opt,errmsg,errflg)
!=================================================================================================================
!input arguments:
integer,intent(in):: hail_opt ! RAS
real(kind=kind_phys),intent(in):: den0,denr,dens,cl,cpv
!output arguments:
character(len=*),intent(out):: errmsg
integer,intent(out):: errflg
!-----------------------------------------------------------------------------------------------------------------
! RAS13.1 define graupel parameters as graupel-like or hail-like,
! depending on namelist option
if(hail_opt .eq. 1) then !Hail!
n0g = 4.e4
deng = 700.
avtg = 285.0
bvtg = 0.8
lamdagmax = 2.e4
else !Graupel!
n0g = 4.e6
deng = 500
avtg = 330.0
bvtg = 0.8
lamdagmax = 6.e4
endif
!
pi = 4.*atan(1.)
xlv1 = cl-cpv
!
qc0 = 4./3.*pi*denr*r0**3*xncr/den0 ! 0.419e-3 -- .61e-3
qck1 = .104*9.8*peaut/(xncr*denr)**(1./3.)/xmyu*den0**(4./3.) ! 7.03
pidnc = pi*denr/6. ! syb
!
bvtr1 = 1.+bvtr
bvtr2 = 2.5+.5*bvtr
bvtr3 = 3.+bvtr
bvtr4 = 4.+bvtr
bvtr6 = 6.+bvtr
g1pbr = rgmma(bvtr1)
g3pbr = rgmma(bvtr3)
g4pbr = rgmma(bvtr4) ! 17.837825
g6pbr = rgmma(bvtr6)
g5pbro2 = rgmma(bvtr2) ! 1.8273
pvtr = avtr*g4pbr/6.
eacrr = 1.0
pacrr = pi*n0r*avtr*g3pbr*.25*eacrr
precr1 = 2.*pi*n0r*.78
precr2 = 2.*pi*n0r*.31*avtr**.5*g5pbro2
roqimax = 2.08e22*dimax**8
!
bvts1 = 1.+bvts
bvts2 = 2.5+.5*bvts
bvts3 = 3.+bvts
bvts4 = 4.+bvts
g1pbs = rgmma(bvts1) !.8875
g3pbs = rgmma(bvts3)
g4pbs = rgmma(bvts4) ! 12.0786
g5pbso2 = rgmma(bvts2)
pvts = avts*g4pbs/6.
pacrs = pi*n0s*avts*g3pbs*.25
precs1 = 4.*n0s*.65
precs2 = 4.*n0s*.44*avts**.5*g5pbso2
pidn0r = pi*denr*n0r
pidn0s = pi*dens*n0s
!
pacrc = pi*n0s*avts*g3pbs*.25*eacrc
!
bvtg1 = 1.+bvtg
bvtg2 = 2.5+.5*bvtg
bvtg3 = 3.+bvtg
bvtg4 = 4.+bvtg
g1pbg = rgmma(bvtg1)
g3pbg = rgmma(bvtg3)
g4pbg = rgmma(bvtg4)
pacrg = pi*n0g*avtg*g3pbg*.25
g5pbgo2 = rgmma(bvtg2)
pvtg = avtg*g4pbg/6.
precg1 = 2.*pi*n0g*.78
precg2 = 2.*pi*n0g*.31*avtg**.5*g5pbgo2
pidn0g = pi*deng*n0g
!
rslopermax = 1./lamdarmax
rslopesmax = 1./lamdasmax
rslopegmax = 1./lamdagmax
rsloperbmax = rslopermax ** bvtr
rslopesbmax = rslopesmax ** bvts
rslopegbmax = rslopegmax ** bvtg
rsloper2max = rslopermax * rslopermax
rslopes2max = rslopesmax * rslopesmax
rslopeg2max = rslopegmax * rslopegmax
rsloper3max = rsloper2max * rslopermax
rslopes3max = rslopes2max * rslopesmax
rslopeg3max = rslopeg2max * rslopegmax
!+---+-----------------------------------------------------------------+
!.. Set these variables needed for computing radar reflectivity. These
!.. get used within radar_init to create other variables used in the
!.. radar module.
xam_r = PI*denr/6.
xbm_r = 3.
xmu_r = 0.
xam_s = PI*dens/6.
xbm_s = 3.
xmu_s = 0.
xam_g = PI*deng/6.
xbm_g = 3.
xmu_g = 0.
call radar_init
errmsg = 'mp_wsm6_init OK'
errflg = 0
end subroutine mp_wsm6_init
!=================================================================================================================
!>\section arg_table_mp_wsm6_finalize
!!\html\include mp_wsm6_finalize.html
!!
subroutine mp_wsm6_finalize(errmsg,errflg)
!=================================================================================================================
!--- output arguments:
character(len=*),intent(out):: errmsg
integer,intent(out):: errflg
!-----------------------------------------------------------------------------------------------------------------
errmsg = 'mp_wsm6_finalize OK'
errflg = 0
end subroutine mp_wsm6_finalize
!=================================================================================================================
!>\section arg_table_mp_wsm6_run
!!\html\include mp_wsm6_run.html
!!
subroutine mp_wsm6_run(t,q,qc,qi,qr,qs,qg,den,p,delz,delt, &
g,cpd,cpv,rd,rv,t0c,ep1,ep2,qmin,xls, &
xlv0,xlf0,den0,denr,cliq,cice,psat, &
rain,rainncv,sr,snow,snowncv,graupel, &
graupelncv,rainprod2d,evapprod2d, &
its,ite,kts,kte,errmsg,errflg &
)
!=================================================================================================================!
! This code is a 6-class GRAUPEL phase microphyiscs scheme (WSM6) of the
! Single-Moment MicroPhyiscs (WSMMP). The WSMMP assumes that ice nuclei
! number concentration is a function of temperature, and seperate assumption
! is developed, in which ice crystal number concentration is a function
! of ice amount. A theoretical background of the ice-microphysics and related
! processes in the WSMMPs are described in Hong et al. (2004).
! All production terms in the WSM6 scheme are described in Hong and Lim (2006).
! All units are in m.k.s. and source/sink terms in kgkg-1s-1.
!
! WSM6 cloud scheme
!
! Coded by Song-You Hong and Jeong-Ock Jade Lim (Yonsei Univ.)
! Summer 2003
!
! Implemented by Song-You Hong (Yonsei Univ.) and Jimy Dudhia (NCAR)
! Summer 2004
!
! further modifications :
! semi-lagrangian sedimentation (JH,2010),hong, aug 2009
! ==> higher accuracy and efficient at lower resolutions
! reflectivity computation from greg thompson, lim, jun 2011
! ==> only diagnostic, but with removal of too large drops
! add hail option from afwa, aug 2014
! ==> switch graupel or hail by changing no, den, fall vel.
! effective radius of hydrometeors, bae from kiaps, jan 2015
! ==> consistency in solar insolation of rrtmg radiation
! bug fix in melting terms, bae from kiaps, nov 2015
! ==> density of air is divided, which has not been
!
! Reference) Hong, Dudhia, Chen (HDC, 2004) Mon. Wea. Rev.
! Hong and Lim (HL, 2006) J. Korean Meteor. Soc.
! Dudhia, Hong and Lim (DHL, 2008) J. Meteor. Soc. Japan
! Lin, Farley, Orville (LFO, 1983) J. Appl. Meteor.
! Rutledge, Hobbs (RH83, 1983) J. Atmos. Sci.
! Rutledge, Hobbs (RH84, 1984) J. Atmos. Sci.
! Juang and Hong (JH, 2010) Mon. Wea. Rev.
!
!input arguments:
integer,intent(in):: its,ite,kts,kte
real(kind=kind_phys),intent(in),dimension(its:,:):: &
den, &
p, &
delz
real(kind=kind_phys),intent(in):: &
delt, &
g, &
cpd, &
cpv, &
t0c, &
den0, &
rd, &
rv, &
ep1, &
ep2, &
qmin, &
xls, &
xlv0, &
xlf0, &
cliq, &
cice, &
psat, &
denr
!inout arguments:
real(kind=kind_phys),intent(inout),dimension(its:,:):: &
t
real(kind=kind_phys),intent(inout),dimension(its:,:):: &
q, &
qc, &
qi, &
qr, &
qs, &
qg
real(kind=kind_phys),intent(inout),dimension(its:):: &
rain, &
rainncv, &
sr
real(kind=kind_phys),intent(inout),dimension(its:),optional:: &
snow, &
snowncv
real(kind=kind_phys),intent(inout),dimension(its:),optional:: &
graupel, &
graupelncv
real(kind=kind_phys),intent(inout),dimension(its:,:),optional:: &
rainprod2d, &
evapprod2d
!output arguments:
character(len=*),intent(out):: errmsg
integer,intent(out):: errflg
!local variables and arrays:
real(kind=kind_phys),dimension(its:ite,kts:kte,3):: &
rh, &
qsat, &
rslope, &
rslope2, &
rslope3, &
rslopeb, &
qrs_tmp, &
falk, &
fall, &
work1
real(kind=kind_phys),dimension(its:ite,kts:kte):: &
fallc, &
falkc, &
work1c, &
work2c, &
workr, &
worka
real(kind=kind_phys),dimension(its:ite,kts:kte):: &
den_tmp, &
delz_tmp
real(kind=kind_phys),dimension(its:ite,kts:kte):: &
pigen, &
pidep, &
pcond, &
prevp, &
psevp, &
pgevp, &
psdep, &
pgdep, &
praut, &
psaut, &
pgaut, &
piacr, &
pracw, &
praci, &
pracs, &
psacw, &
psaci, &
psacr, &
pgacw, &
pgaci, &
pgacr, &
pgacs, &
paacw, &
psmlt, &
pgmlt, &
pseml, &
pgeml
real(kind=kind_phys),dimension(its:ite,kts:kte):: &
qsum, &
xl, &
cpm, &
work2, &
denfac, &
xni, &
denqrs1, &
denqrs2, &
denqrs3, &
denqci, &
n0sfac
real(kind=kind_phys),dimension(its:ite):: &
delqrs1, &
delqrs2, &
delqrs3, &
delqi
real(kind=kind_phys),dimension(its:ite):: &
tstepsnow, &
tstepgraup
integer,dimension(its:ite):: &
mstep, &
numdt
logical,dimension(its:ite):: flgcld
real(kind=kind_phys):: &
cpmcal, xlcal, diffus, &
viscos, xka, venfac, conden, diffac, &
x, y, z, a, b, c, d, e, &
qdt, holdrr, holdrs, holdrg, supcol, supcolt, pvt, &
coeres, supsat, dtcld, xmi, eacrs, satdt, &
qimax, diameter, xni0, roqi0, &
fallsum, fallsum_qsi, fallsum_qg, &
vt2i,vt2r,vt2s,vt2g,acrfac,egs,egi, &
xlwork2, factor, source, value, &
xlf, pfrzdtc, pfrzdtr, supice, alpha2, delta2, delta3
real(kind=kind_phys):: vt2ave
real(kind=kind_phys):: holdc, holdci
integer:: i, j, k, mstepmax, &
iprt, latd, lond, loop, loops, ifsat, n, idim, kdim
!Temporaries used for inlining fpvs function
real(kind=kind_phys):: dldti, xb, xai, tr, xbi, xa, hvap, cvap, hsub, dldt, ttp
! variables for optimization
real(kind=kind_phys),dimension(its:ite):: dvec1,tvec1
real(kind=kind_phys):: temp
!-----------------------------------------------------------------------------------------------------------------
! compute internal functions
!
cpmcal(x) = cpd*(1.-max(x,qmin))+max(x,qmin)*cpv
xlcal(x) = xlv0-xlv1*(x-t0c)
!----------------------------------------------------------------
! diffus: diffusion coefficient of the water vapor
! viscos: kinematic viscosity(m2s-1)
! Optimizatin : A**B => exp(log(A)*(B))
!
diffus(x,y) = 8.794e-5 * exp(log(x)*(1.81)) / y ! 8.794e-5*x**1.81/y
viscos(x,y) = 1.496e-6 * (x*sqrt(x)) /(x+120.)/y ! 1.496e-6*x**1.5/(x+120.)/y
xka(x,y) = 1.414e3*viscos(x,y)*y
diffac(a,b,c,d,e) = d*a*a/(xka(c,d)*rv*c*c)+1./(e*diffus(c,b))
venfac(a,b,c) = exp(log((viscos(b,c)/diffus(b,a)))*((.3333333))) &
/sqrt(viscos(b,c))*sqrt(sqrt(den0/c))
conden(a,b,c,d,e) = (max(b,qmin)-c)/(1.+d*d/(rv*e)*c/(a*a))
!
!
idim = ite-its+1
kdim = kte-kts+1
!
!----------------------------------------------------------------
! paddint 0 for negative values generated by dynamics
!
do k = kts, kte
do i = its, ite
qc(i,k) = max(qc(i,k),0.0)
qr(i,k) = max(qr(i,k),0.0)
qi(i,k) = max(qi(i,k),0.0)
qs(i,k) = max(qs(i,k),0.0)
qg(i,k) = max(qg(i,k),0.0)
enddo
enddo
!
!----------------------------------------------------------------
! latent heat for phase changes and heat capacity. neglect the
! changes during microphysical process calculation emanuel(1994)
!
do k = kts, kte
do i = its, ite
cpm(i,k) = cpmcal(q(i,k))
xl(i,k) = xlcal(t(i,k))
enddo
enddo
do k = kts, kte
do i = its, ite
delz_tmp(i,k) = delz(i,k)
den_tmp(i,k) = den(i,k)
enddo
enddo
!
!----------------------------------------------------------------
! initialize the surface rain, snow, graupel
!
do i = its, ite
rainncv(i) = 0.
if(present(snowncv) .and. present(snow)) snowncv(i) = 0.
if(present(graupelncv) .and. present(graupel)) graupelncv(i) = 0.
sr(i) = 0.
! new local array to catch step snow and graupel
tstepsnow(i) = 0.
tstepgraup(i) = 0.
enddo
!
!----------------------------------------------------------------
! compute the minor time steps.
!
loops = max(nint(delt/dtcldcr),1)
dtcld = delt/loops
if(delt.le.dtcldcr) dtcld = delt
!
do loop = 1,loops
!
!----------------------------------------------------------------
! initialize the large scale variables
!
do i = its, ite
mstep(i) = 1
flgcld(i) = .true.
enddo
!
! do k = kts, kte
! do i = its, ite
! denfac(i,k) = sqrt(den0/den(i,k))
! enddo
! enddo
do k = kts, kte
do i = its,ite
dvec1(i) = den(i,k)
enddo
call vrec(tvec1,dvec1,ite-its+1)
do i = its, ite
tvec1(i) = tvec1(i)*den0
enddo
call vsqrt(dvec1,tvec1,ite-its+1)
do i = its,ite
denfac(i,k) = dvec1(i)
enddo
enddo
!
! Inline expansion for fpvs
! qsat(i,k,1) = fpvs(t(i,k),0,rd,rv,cpv,cliq,cice,xlv0,xls,psat,t0c)
! qsat(i,k,2) = fpvs(t(i,k),1,rd,rv,cpv,cliq,cice,xlv0,xls,psat,t0c)
hsub = xls
hvap = xlv0
cvap = cpv
ttp=t0c+0.01
dldt=cvap-cliq
xa=-dldt/rv
xb=xa+hvap/(rv*ttp)
dldti=cvap-cice
xai=-dldti/rv
xbi=xai+hsub/(rv*ttp)
do k = kts, kte
do i = its, ite
tr=ttp/t(i,k)
qsat(i,k,1)=psat*exp(log(tr)*(xa))*exp(xb*(1.-tr))
qsat(i,k,1) = min(qsat(i,k,1),0.99*p(i,k))
qsat(i,k,1) = ep2 * qsat(i,k,1) / (p(i,k) - qsat(i,k,1))
qsat(i,k,1) = max(qsat(i,k,1),qmin)
rh(i,k,1) = max(q(i,k) / qsat(i,k,1),qmin)
tr=ttp/t(i,k)
if(t(i,k).lt.ttp) then
qsat(i,k,2)=psat*exp(log(tr)*(xai))*exp(xbi*(1.-tr))
else
qsat(i,k,2)=psat*exp(log(tr)*(xa))*exp(xb*(1.-tr))
endif
qsat(i,k,2) = min(qsat(i,k,2),0.99*p(i,k))
qsat(i,k,2) = ep2 * qsat(i,k,2) / (p(i,k) - qsat(i,k,2))
qsat(i,k,2) = max(qsat(i,k,2),qmin)
rh(i,k,2) = max(q(i,k) / qsat(i,k,2),qmin)
enddo
enddo
!
!----------------------------------------------------------------
! initialize the variables for microphysical physics
!
!
do k = kts, kte
do i = its, ite
prevp(i,k) = 0.
psdep(i,k) = 0.
pgdep(i,k) = 0.
praut(i,k) = 0.
psaut(i,k) = 0.
pgaut(i,k) = 0.
pracw(i,k) = 0.
praci(i,k) = 0.
piacr(i,k) = 0.
psaci(i,k) = 0.
psacw(i,k) = 0.
pracs(i,k) = 0.
psacr(i,k) = 0.
pgacw(i,k) = 0.
paacw(i,k) = 0.
pgaci(i,k) = 0.
pgacr(i,k) = 0.
pgacs(i,k) = 0.
pigen(i,k) = 0.
pidep(i,k) = 0.
pcond(i,k) = 0.
psmlt(i,k) = 0.
pgmlt(i,k) = 0.
pseml(i,k) = 0.
pgeml(i,k) = 0.
psevp(i,k) = 0.
pgevp(i,k) = 0.
falk(i,k,1) = 0.
falk(i,k,2) = 0.
falk(i,k,3) = 0.
fall(i,k,1) = 0.
fall(i,k,2) = 0.
fall(i,k,3) = 0.
fallc(i,k) = 0.
falkc(i,k) = 0.
xni(i,k) = 1.e3
enddo
enddo
!-------------------------------------------------------------
! Ni: ice crystal number concentraiton [HDC 5c]
!-------------------------------------------------------------
do k = kts, kte
do i = its, ite
temp = (den(i,k)*max(qi(i,k),qmin))
temp = sqrt(sqrt(temp*temp*temp))
xni(i,k) = min(max(5.38e7*temp,1.e3),1.e6)
enddo
enddo
!
!----------------------------------------------------------------
! compute the fallout term:
! first, vertical terminal velosity for minor loops
!----------------------------------------------------------------
do k = kts, kte
do i = its, ite
qrs_tmp(i,k,1) = qr(i,k)
qrs_tmp(i,k,2) = qs(i,k)
qrs_tmp(i,k,3) = qg(i,k)
enddo
enddo
call slope_wsm6(qrs_tmp,den_tmp,denfac,t,rslope,rslopeb,rslope2,rslope3, &
work1,its,ite,kts,kte)
!
do k = kte, kts, -1
do i = its, ite
workr(i,k) = work1(i,k,1)
qsum(i,k) = max( (qs(i,k)+qg(i,k)), 1.E-15)
if( qsum(i,k) .gt. 1.e-15 ) then
worka(i,k) = (work1(i,k,2)*qs(i,k) + work1(i,k,3)*qg(i,k)) &
/ qsum(i,k)
else
worka(i,k) = 0.
endif
denqrs1(i,k) = den(i,k)*qr(i,k)
denqrs2(i,k) = den(i,k)*qs(i,k)
denqrs3(i,k) = den(i,k)*qg(i,k)
if(qr(i,k).le.0.0) workr(i,k) = 0.0
enddo
enddo
call nislfv_rain_plm(idim,kdim,den_tmp,denfac,t,delz_tmp,workr,denqrs1, &
delqrs1,dtcld,1,1)
call nislfv_rain_plm6(idim,kdim,den_tmp,denfac,t,delz_tmp,worka, &
denqrs2,denqrs3,delqrs2,delqrs3,dtcld,1,1)
do k = kts, kte
do i = its, ite
qr(i,k) = max(denqrs1(i,k)/den(i,k),0.)
qs(i,k) = max(denqrs2(i,k)/den(i,k),0.)
qg(i,k) = max(denqrs3(i,k)/den(i,k),0.)
fall(i,k,1) = denqrs1(i,k)*workr(i,k)/delz(i,k)
fall(i,k,2) = denqrs2(i,k)*worka(i,k)/delz(i,k)
fall(i,k,3) = denqrs3(i,k)*worka(i,k)/delz(i,k)
enddo
enddo
do i = its, ite
fall(i,1,1) = delqrs1(i)/delz(i,1)/dtcld
fall(i,1,2) = delqrs2(i)/delz(i,1)/dtcld
fall(i,1,3) = delqrs3(i)/delz(i,1)/dtcld
enddo
do k = kts, kte
do i = its, ite
qrs_tmp(i,k,1) = qr(i,k)
qrs_tmp(i,k,2) = qs(i,k)
qrs_tmp(i,k,3) = qg(i,k)
enddo
enddo
call slope_wsm6(qrs_tmp,den_tmp,denfac,t,rslope,rslopeb,rslope2,rslope3, &
work1,its,ite,kts,kte)
!
do k = kte, kts, -1
do i = its, ite
supcol = t0c-t(i,k)
n0sfac(i,k) = max(min(exp(alpha*supcol),n0smax/n0s),1.)
if(t(i,k).gt.t0c) then
!---------------------------------------------------------------
! psmlt: melting of snow [HL A33] [RH83 A25]
! (T>T0: S->R)
!---------------------------------------------------------------
xlf = xlf0
work2(i,k) = venfac(p(i,k),t(i,k),den(i,k))
if(qs(i,k).gt.0.) then
coeres = rslope2(i,k,2)*sqrt(rslope(i,k,2)*rslopeb(i,k,2))
psmlt(i,k) = xka(t(i,k),den(i,k))/xlf*(t0c-t(i,k))*pi/2. &
*n0sfac(i,k)*(precs1*rslope2(i,k,2) &
+precs2*work2(i,k)*coeres)/den(i,k)
psmlt(i,k) = min(max(psmlt(i,k)*dtcld/mstep(i), &
-qs(i,k)/mstep(i)),0.)
qs(i,k) = qs(i,k) + psmlt(i,k)
qr(i,k) = qr(i,k) - psmlt(i,k)
t(i,k) = t(i,k) + xlf/cpm(i,k)*psmlt(i,k)
endif
!---------------------------------------------------------------
! pgmlt: melting of graupel [HL A23] [LFO 47]
! (T>T0: G->R)
!---------------------------------------------------------------
if(qg(i,k).gt.0.) then
coeres = rslope2(i,k,3)*sqrt(rslope(i,k,3)*rslopeb(i,k,3))
pgmlt(i,k) = xka(t(i,k),den(i,k))/xlf &
*(t0c-t(i,k))*(precg1*rslope2(i,k,3) &
+precg2*work2(i,k)*coeres)/den(i,k)
pgmlt(i,k) = min(max(pgmlt(i,k)*dtcld/mstep(i), &
-qg(i,k)/mstep(i)),0.)
qg(i,k) = qg(i,k) + pgmlt(i,k)
qr(i,k) = qr(i,k) - pgmlt(i,k)
t(i,k) = t(i,k) + xlf/cpm(i,k)*pgmlt(i,k)
endif
endif
enddo
enddo
!---------------------------------------------------------------
! Vice [ms-1] : fallout of ice crystal [HDC 5a]
!---------------------------------------------------------------
do k = kte, kts, -1
do i = its, ite
if(qi(i,k).le.0.) then
work1c(i,k) = 0.
else
xmi = den(i,k)*qi(i,k)/xni(i,k)
diameter = max(min(dicon * sqrt(xmi),dimax), 1.e-25)
work1c(i,k) = 1.49e4*exp(log(diameter)*(1.31))
endif
enddo
enddo
!
! forward semi-laglangian scheme (JH), PCM (piecewise constant), (linear)
!
do k = kte, kts, -1
do i = its, ite
denqci(i,k) = den(i,k)*qi(i,k)
enddo
enddo
call nislfv_rain_plm(idim,kdim,den_tmp,denfac,t,delz_tmp,work1c,denqci, &
delqi,dtcld,1,0)
do k = kts, kte
do i = its, ite
qi(i,k) = max(denqci(i,k)/den(i,k),0.)
enddo
enddo
do i = its, ite
fallc(i,1) = delqi(i)/delz(i,1)/dtcld
enddo
!
!----------------------------------------------------------------
! rain (unit is mm/sec;kgm-2s-1: /1000*delt ===> m)==> mm for wrf
!
do i = its, ite
fallsum = fall(i,kts,1)+fall(i,kts,2)+fall(i,kts,3)+fallc(i,kts)
fallsum_qsi = fall(i,kts,2)+fallc(i,kts)
fallsum_qg = fall(i,kts,3)
if(fallsum.gt.0.) then
rainncv(i) = fallsum*delz(i,kts)/denr*dtcld*1000. + rainncv(i)
rain(i) = fallsum*delz(i,kts)/denr*dtcld*1000. + rain(i)
endif
if(fallsum_qsi.gt.0.) then
tstepsnow(i) = fallsum_qsi*delz(i,kts)/denr*dtcld*1000. &
+ tstepsnow(i)
if(present(snowncv) .and. present(snow)) then
snowncv(i) = fallsum_qsi*delz(i,kts)/denr*dtcld*1000. &
+ snowncv(i)
snow(i) = fallsum_qsi*delz(i,kts)/denr*dtcld*1000. + snow(i)
endif
endif
if(fallsum_qg.gt.0.) then
tstepgraup(i) = fallsum_qg*delz(i,kts)/denr*dtcld*1000. &
+ tstepgraup(i)
if(present (graupelncv) .and. present (graupel)) then
graupelncv(i) = fallsum_qg*delz(i,kts)/denr*dtcld*1000. &
+ graupelncv(i)
graupel(i) = fallsum_qg*delz(i,kts)/denr*dtcld*1000. + graupel(i)
endif
endif
if(present (snowncv)) then
if(fallsum.gt.0.)sr(i)=(snowncv(i) + graupelncv(i))/(rainncv(i)+1.e-12)
else
if(fallsum.gt.0.)sr(i)=(tstepsnow(i) + tstepgraup(i))/(rainncv(i)+1.e-12)
endif
enddo
!
!---------------------------------------------------------------
! pimlt: instantaneous melting of cloud ice [HL A47] [RH83 A28]
! (T>T0: I->C)
!---------------------------------------------------------------
do k = kts, kte
do i = its, ite
supcol = t0c-t(i,k)
xlf = xls-xl(i,k)
if(supcol.lt.0.) xlf = xlf0
if(supcol.lt.0.and.qi(i,k).gt.0.) then
qc(i,k) = qc(i,k) + qi(i,k)
t(i,k) = t(i,k) - xlf/cpm(i,k)*qi(i,k)
qi(i,k) = 0.
endif
!---------------------------------------------------------------
! pihmf: homogeneous freezing of cloud water below -40c [HL A45]
! (T<-40C: C->I)
!---------------------------------------------------------------
if(supcol.gt.40..and.qc(i,k).gt.0.) then
qi(i,k) = qi(i,k) + qc(i,k)
t(i,k) = t(i,k) + xlf/cpm(i,k)*qc(i,k)
qc(i,k) = 0.
endif
!---------------------------------------------------------------
! pihtf: heterogeneous freezing of cloud water [HL A44]
! (T0>T>-40C: C->I)
!---------------------------------------------------------------
if(supcol.gt.0..and.qc(i,k).gt.qmin) then
! pfrzdtc = min(pfrz1*(exp(pfrz2*supcol)-1.) &
! * den(i,k)/denr/xncr*qc(i,k)**2*dtcld,qc(i,k))
supcolt=min(supcol,50.)
pfrzdtc = min(pfrz1*(exp(pfrz2*supcolt)-1.) &
* den(i,k)/denr/xncr*qc(i,k)*qc(i,k)*dtcld,qc(i,k))
qi(i,k) = qi(i,k) + pfrzdtc
t(i,k) = t(i,k) + xlf/cpm(i,k)*pfrzdtc
qc(i,k) = qc(i,k)-pfrzdtc
endif
!---------------------------------------------------------------
! pgfrz: freezing of rain water [HL A20] [LFO 45]
! (T<T0, R->G)
!---------------------------------------------------------------
if(supcol.gt.0..and.qr(i,k).gt.0.) then
! pfrzdtr = min(20.*pi**2*pfrz1*n0r*denr/den(i,k) &
! * (exp(pfrz2*supcol)-1.)*rslope3(i,k,1)**2 &
! * rslope(i,k,1)*dtcld,qr(i,k))
temp = rslope3(i,k,1)
temp = temp*temp*rslope(i,k,1)
supcolt=min(supcol,50.)
pfrzdtr = min(20.*(pi*pi)*pfrz1*n0r*denr/den(i,k) &
*(exp(pfrz2*supcolt)-1.)*temp*dtcld, &
qr(i,k))
qg(i,k) = qg(i,k) + pfrzdtr
t(i,k) = t(i,k) + xlf/cpm(i,k)*pfrzdtr
qr(i,k) = qr(i,k)-pfrzdtr
endif
enddo
enddo
!
!
!----------------------------------------------------------------
! update the slope parameters for microphysics computation
!
do k = kts, kte
do i = its, ite
qrs_tmp(i,k,1) = qr(i,k)
qrs_tmp(i,k,2) = qs(i,k)
qrs_tmp(i,k,3) = qg(i,k)
enddo
enddo
call slope_wsm6(qrs_tmp,den_tmp,denfac,t,rslope,rslopeb,rslope2,rslope3, &
work1,its,ite,kts,kte)
!------------------------------------------------------------------
! work1: the thermodynamic term in the denominator associated with
! heat conduction and vapor diffusion
! (ry88, y93, h85)
! work2: parameter associated with the ventilation effects(y93)
!
do k = kts, kte
do i = its, ite
work1(i,k,1) = diffac(xl(i,k),p(i,k),t(i,k),den(i,k),qsat(i,k,1))
work1(i,k,2) = diffac(xls,p(i,k),t(i,k),den(i,k),qsat(i,k,2))
work2(i,k) = venfac(p(i,k),t(i,k),den(i,k))
enddo
enddo
!
!===============================================================
!
! warm rain processes
!
! - follows the processes in RH83 and LFO except for autoconcersion
!
!===============================================================
!
do k = kts, kte
do i = its, ite
supsat = max(q(i,k),qmin)-qsat(i,k,1)
satdt = supsat/dtcld
!---------------------------------------------------------------
! praut: auto conversion rate from cloud to rain [HDC 16]
! (C->R)
!---------------------------------------------------------------
if(qc(i,k).gt.qc0) then
praut(i,k) = qck1*qc(i,k)**(7./3.)
praut(i,k) = min(praut(i,k),qc(i,k)/dtcld)
endif
!---------------------------------------------------------------
! pracw: accretion of cloud water by rain [HL A40] [LFO 51]
! (C->R)
!---------------------------------------------------------------
if(qr(i,k).gt.qcrmin.and.qc(i,k).gt.qmin) then
pracw(i,k) = min(pacrr*rslope3(i,k,1)*rslopeb(i,k,1) &
* qc(i,k)*denfac(i,k),qc(i,k)/dtcld)
endif
!---------------------------------------------------------------
! prevp: evaporation/condensation rate of rain [HDC 14]
! (V->R or R->V)
!---------------------------------------------------------------
if(qr(i,k).gt.0.) then
coeres = rslope2(i,k,1)*sqrt(rslope(i,k,1)*rslopeb(i,k,1))
prevp(i,k) = (rh(i,k,1)-1.)*(precr1*rslope2(i,k,1) &
+ precr2*work2(i,k)*coeres)/work1(i,k,1)
if(prevp(i,k).lt.0.) then
prevp(i,k) = max(prevp(i,k),-qr(i,k)/dtcld)
prevp(i,k) = max(prevp(i,k),satdt/2)
else
prevp(i,k) = min(prevp(i,k),satdt/2)
endif
endif
enddo
enddo
!
!===============================================================
!
! cold rain processes
!
! - follows the revised ice microphysics processes in HDC
! - the processes same as in RH83 and RH84 and LFO behave
! following ice crystal hapits defined in HDC, inclduing
! intercept parameter for snow (n0s), ice crystal number
! concentration (ni), ice nuclei number concentration
! (n0i), ice diameter (d)
!
!===============================================================
!
do k = kts, kte
do i = its, ite
supcol = t0c-t(i,k)
n0sfac(i,k) = max(min(exp(alpha*supcol),n0smax/n0s),1.)
supsat = max(q(i,k),qmin)-qsat(i,k,2)
satdt = supsat/dtcld
ifsat = 0
!-------------------------------------------------------------
! Ni: ice crystal number concentraiton [HDC 5c]
!-------------------------------------------------------------
! xni(i,k) = min(max(5.38e7*(den(i,k) &
! * max(qi(i,k),qmin))**0.75,1.e3),1.e6)
temp = (den(i,k)*max(qi(i,k),qmin))
temp = sqrt(sqrt(temp*temp*temp))
xni(i,k) = min(max(5.38e7*temp,1.e3),1.e6)
eacrs = exp(0.07*(-supcol))
!
xmi = den(i,k)*qi(i,k)/xni(i,k)
diameter = min(dicon * sqrt(xmi),dimax)
vt2i = 1.49e4*diameter**1.31
vt2r=pvtr*rslopeb(i,k,1)*denfac(i,k)
vt2s=pvts*rslopeb(i,k,2)*denfac(i,k)
vt2g=pvtg*rslopeb(i,k,3)*denfac(i,k)
qsum(i,k) = max( (qs(i,k)+qg(i,k)), 1.E-15)
if(qsum(i,k) .gt. 1.e-15) then
vt2ave=(vt2s*qs(i,k)+vt2g*qg(i,k))/(qsum(i,k))
else
vt2ave=0.
endif
if(supcol.gt.0.and.qi(i,k).gt.qmin) then
if(qr(i,k).gt.qcrmin) then
!-------------------------------------------------------------
! praci: accretion of cloud ice by rain [HL A15] [LFO 25]
! (T<T0: I->R)
!-------------------------------------------------------------
acrfac = 2.*rslope3(i,k,1)+2.*diameter*rslope2(i,k,1) &
+ diameter**2*rslope(i,k,1)
praci(i,k) = pi*qi(i,k)*n0r*abs(vt2r-vt2i)*acrfac/4.
! reduce collection efficiency (suggested by B. Wilt)
praci(i,k) = praci(i,k)*min(max(0.0,qr(i,k)/qi(i,k)),1.)**2
praci(i,k) = min(praci(i,k),qi(i,k)/dtcld)
!-------------------------------------------------------------
! piacr: accretion of rain by cloud ice [HL A19] [LFO 26]
! (T<T0: R->S or R->G)
!-------------------------------------------------------------
piacr(i,k) = pi**2*avtr*n0r*denr*xni(i,k)*denfac(i,k) &
* g6pbr*rslope3(i,k,1)*rslope3(i,k,1) &
* rslopeb(i,k,1)/24./den(i,k)
! reduce collection efficiency (suggested by B. Wilt)
piacr(i,k) = piacr(i,k)*min(max(0.0,qi(i,k)/qr(i,k)),1.)**2
piacr(i,k) = min(piacr(i,k),qr(i,k)/dtcld)
endif
!-------------------------------------------------------------
! psaci: accretion of cloud ice by snow [HDC 10]
! (T<T0: I->S)
!-------------------------------------------------------------
if(qs(i,k).gt.qcrmin) then
acrfac = 2.*rslope3(i,k,2)+2.*diameter*rslope2(i,k,2) &
+ diameter**2*rslope(i,k,2)
psaci(i,k) = pi*qi(i,k)*eacrs*n0s*n0sfac(i,k) &
* abs(vt2ave-vt2i)*acrfac/4.
psaci(i,k) = min(psaci(i,k),qi(i,k)/dtcld)
endif
!-------------------------------------------------------------
! pgaci: accretion of cloud ice by graupel [HL A17] [LFO 41]
! (T<T0: I->G)
!-------------------------------------------------------------
if(qg(i,k).gt.qcrmin) then
egi = exp(0.07*(-supcol))
acrfac = 2.*rslope3(i,k,3)+2.*diameter*rslope2(i,k,3) &
+ diameter**2*rslope(i,k,3)
pgaci(i,k) = pi*egi*qi(i,k)*n0g*abs(vt2ave-vt2i)*acrfac/4.
pgaci(i,k) = min(pgaci(i,k),qi(i,k)/dtcld)
endif
endif
!-------------------------------------------------------------
! psacw: accretion of cloud water by snow [HL A7] [LFO 24]
! (T<T0: C->S, and T>=T0: C->R)
!-------------------------------------------------------------
if(qs(i,k).gt.qcrmin.and.qc(i,k).gt.qmin) then
psacw(i,k) = min(pacrc*n0sfac(i,k)*rslope3(i,k,2)*rslopeb(i,k,2) &
! reduce collection efficiency (suggested by B. Wilt)
* min(max(0.0,qs(i,k)/qc(i,k)),1.)**2 &
* qc(i,k)*denfac(i,k),qc(i,k)/dtcld)