-
Notifications
You must be signed in to change notification settings - Fork 6
/
Copy pathmodule_mp_tempo_utils.F90
2105 lines (1840 loc) · 88.2 KB
/
module_mp_tempo_utils.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
! Utilities for TEMPO Microphysics
!=================================================================================================================
module module_mp_tempo_utils
use module_mp_tempo_params
#if defined(mpas)
use mpas_kind_types, only: wp => RKIND, sp => R4KIND, dp => R8KIND
use mp_radar
#elif defined(standalone)
use machine, only: wp => kind_phys, sp => kind_sngl_prec, dp => kind_dbl_prec
use module_mp_radar
#else
use machine, only: wp => kind_phys, sp => kind_sngl_prec, dp => kind_dbl_prec
use module_mp_radar
#define ccpp_default 1
#endif
#if defined(ccpp_default) && defined(MPI)
use mpi_f08
#endif
implicit none
contains
!=================================================================================================================
! Normalized lower gamma function calculated either with a series expansion or continued fraction method
! Input:
! a = gamma function argument, x = upper limit of integration
! Output:
! gamma_p = gamma(a, x) / Gamma(a)
real function gamma_p(a, x)
real(wp), intent(in) :: a, x
if ((x < 0.0) .or. (a <= 0.0)) stop "Invalid arguments for function gamma_p"
if (x < (a+1.0)) then
gamma_p = gamma_series(a, x)
else
! gammma_cf computes the upper series
gamma_p = 1.0 - gamma_cf(a, x)
endif
end function gamma_p
!=================================================================================================================
! https://dlmf.nist.gov/8.7 (Equation 8.7.1)
! Solves the normalized lower gamma function: gamma(a,x) / Gamma(a) = x**a * gamma*(a,x)
! Where gamma*(a,x) = exp(-x) * sum (x**k / Gamma(a+k+1))
! Input:
! a = gamma function argument, x = upper limit of integration
! Output:
! Normalized LOWER gamma function: gamma(a, x) / Gamma(a)
real function gamma_series(a, x)
real(wp), intent(in) :: a, x
integer :: k
integer, parameter :: it_max = 100
real(wp), parameter :: smallvalue = 1.e-7
real(wp) :: ap1, sum_term, sum
if (x <= 0.0) stop "Invalid arguments for function gamma_series"
! k = 0 summation term is 1 / Gamma(a+1)
ap1 = a
sum_term = 1.0 / gamma(ap1+1.0)
sum = sum_term
do k = 1, it_max
ap1 = ap1 + 1.0
sum_term = sum_term * x / ap1
sum = sum + sum_term
if (abs(sum_term) < (abs(sum) * smallvalue)) exit
enddo
if (k == it_max) stop "gamma_series solution did not converge"
gamma_series = sum * x**a * exp(-x)
end function gamma_series
!=================================================================================================================
! http://functions.wolfram.com/06.06.10.0003.01
! Solves the normalized upper gamma function: gamma(a,x) / Gamma(a)
! Using a continued fractions method (modifed Lentz Algorithm)
! Input:
! a = gamma function argument, x = lower limit of integration
! Output:
! Normalized UPPER gamma function: gamma(a, x) / Gamma(a)
real function gamma_cf(a, x)
real(wp), intent(in) :: a, x
integer :: k
integer, parameter :: it_max = 100
real(wp), parameter :: smallvalue = 1.e-7
real(wp), parameter :: offset = 1.e-30
real(wp) :: b, d, h0, c, delta, h, aj
b = 1.0 - a + x
d = 1.0 / b
h0 = offset
c = b + (1.0/offset)
delta = c * d
h = h0 * delta
do k = 1, it_max
aj = k * (a-k)
b = b + 2.0
d = b + aj*d
if(abs(d) < offset) d = offset
c = b + aj/c
if(abs(c) < offset) c = offset
d = 1.0 / d
delta = c * d
h = h * delta
if (abs(delta-1.0) < smallvalue) exit
enddo
if (k == it_max) stop "gamma_cf solution did not converge"
gamma_cf = exp(-x+a*log(x)) * h / gamma(a)
end function gamma_cf
!=================================================================================================================
! Calculates log-spaced bins for hydrometer sizes to simplify calculations later
! Input:
! numbins, lowbin, highbin
! Output:
! bins, deltabins
subroutine create_bins(numbins, lowbin, highbin, bins, deltabins)
integer, intent(in) :: numbins
real(dp), intent(in) :: lowbin, highbin
integer :: n
real(dp), dimension(numbins+1) :: xDx
real(dp), dimension(:), intent(out) :: bins
real(dp), dimension(:), intent(out), optional :: deltabins
xDx(1) = lowbin
xDx(numbins+1) = highbin
do n = 2, numbins
xDx(n) = exp(real(n-1, kind=dp)/real(numbins, kind=dp) * log(xDx(numbins+1)/xDx(1)) + log(xDx(1)))
enddo
do n = 1, numbins
bins(n) = sqrt(xDx(n)*xDx(n+1))
enddo
if (present(deltabins)) then
do n = 1, numbins
deltabins(n) = xDx(n+1) - xDx(n)
enddo
endif
end subroutine create_bins
!=================================================================================================================
! Variable collision efficiency for rain collecting cloud water using method of Beard and Grover, 1974
! if a/A less than 0.25; otherwise uses polynomials to get close match of Pruppacher & Klett Fig 14-9.
! Output:
! t_Efrw
subroutine table_Efrw
! Local variables
real(dp) :: vtr, stokes, reynolds, Ef_rw
real(dp) :: p, yc0, F, G, H, z, K0, X
integer :: i, j
do j = 1, nbc
do i = 1, nbr
Ef_rw = 0.0
p = Dc(j) / Dr(i)
if (Dr(i) < 50.e-6 .or. Dc(j) < 3.e-6) then
t_Efrw(i,j) = 0.0
elseif (p > 0.25) then
X = Dc(j) * 1.e6_dp
if (Dr(i) < 75.e-6) then
Ef_rw = 0.026794*X - 0.20604
elseif (Dr(i) < 125.e-6) then
Ef_rw = -0.00066842*X*X + 0.061542*X - 0.37089
elseif (Dr(i) < 175.e-6) then
Ef_rw = 4.091e-06*X*X*X*X - 0.00030908*X*X*X + 0.0066237*X*X - 0.0013687*X - 0.073022
elseif (Dr(i) < 250.e-6) then
Ef_rw = 9.6719e-5*X*X*X - 0.0068901*X*X + 0.17305*X - 0.65988
elseif (Dr(i) < 350.e-6) then
Ef_rw = 9.0488e-5*X*X*X - 0.006585*X*X + 0.16606*X - 0.56125
else
Ef_rw = 0.00010721*X*X*X - 0.0072962*X*X + 0.1704*X - 0.46929
endif
else
vtr = -0.1021 + 4.932e3*Dr(i) - 0.9551e6*Dr(i)*Dr(i) + 0.07934e9*Dr(i)*Dr(i)*Dr(i) &
- 0.002362e12*Dr(i)*Dr(i)*Dr(i)*Dr(i)
stokes = Dc(j) * Dc(j) * vtr * rho_w2 / (9.*1.718e-5*Dr(i))
reynolds = 9. * stokes / (p*p*rho_w2)
F = log(reynolds)
G = -0.1007_dp - 0.358_dp*F + 0.0261_dp*F*F
K0 = exp(G)
z = log(stokes / (K0+1.e-15_dp))
H = 0.1465_dp + 1.302_dp*z - 0.607_dp*z*z + 0.293_dp*z*z*z
yc0 = 2.0_dp / PI * atan(H)
Ef_rw = (yc0+p)*(yc0+p) / ((1.+p)*(1.+p))
endif
t_Efrw(i,j) = max(0.0, min(real(Ef_rw, kind=wp), 0.95))
enddo
enddo
end subroutine table_Efrw
!=================================================================================================================
! Variable collision efficiency for snow collecting cloud water using method of Wang and Ji, 2000 except
! equate melted snow diameter to their "effective collision cross-section."
! Output:
! t_Efsw
subroutine table_Efsw
! Local variables
real(dp) :: Ds_m, vts, vtc, stokes, reynolds, Ef_sw
real(dp) :: p, yc0, F, G, H, z, K0
integer :: i, j
do j = 1, nbc
vtc = 1.19e4_dp * (1.0e4_dp*Dc(j)*Dc(j)*0.25_dp)
do i = 1, nbs
vts = av_s*Ds(i)**bv_s * exp(-fv_s*Ds(i)) - vtc
Ds_m = (am_s*Ds(i)**bm_s / am_r)**obmr
p = Dc(j) / Ds_m
if (p > 0.25 .or. Ds(i) < D0s .or. Dc(j) < 6.e-6 .or. vts < 1.e-3) then
t_Efsw(i,j) = 0.0
else
stokes = Dc(j) * Dc(j) * vts * rho_w2 / (9.*1.718e-5*Ds_m)
reynolds = 9. * stokes / (p*p*rho_w2)
F = log(reynolds)
G = -0.1007_dp - 0.358_dp*F + 0.0261_dp*F*F
K0 = exp(G)
z = log(stokes / (K0+1.e-15_dp))
H = 0.1465_dp + 1.302_dp*z - 0.607_dp*z*z + 0.293_dp*z*z*z
yc0 = 2.0_dp / PI * atan(H)
Ef_sw = (yc0+p)*(yc0+p) / ((1.+p)*(1.+p))
t_Efsw(i,j) = max(0.0, min(real(Ef_sw, kind=wp), 0.95))
endif
enddo
enddo
end subroutine table_Efsw
!=================================================================================================================
! Droplet evaporation
! Output:
! tpc_wev, tnc_wev
subroutine table_dropEvap
! Local variables
integer :: i, j, k, n
real(dp), dimension(nbc) :: N_c, massc
real(dp) :: summ, summ2, lamc, N0_c
integer :: nu_c
do n = 1, nbc
massc(n) = am_r*Dc(n)**bm_r
enddo
do k = 1, nbc
nu_c = min(nu_c_max, nint(nu_c_scale/t_Nc(k)) + nu_c_min)
do j = 1, ntb_c
lamc = (t_Nc(k)*am_r* ccg(2,nu_c)*ocg1(nu_c) / r_c(j))**obmr
N0_c = t_Nc(k)*ocg1(nu_c) * lamc**cce(1,nu_c)
do i = 1, nbc
N_c(i) = N0_c* Dc(i)**nu_c*exp(-lamc*Dc(i))*dtc(i)
summ = 0.
summ2 = 0.
do n = 1, i
summ = summ + massc(n)*N_c(n)
summ2 = summ2 + N_c(n)
enddo
tpc_wev(i,j,k) = summ
tnc_wev(i,j,k) = summ2
enddo
enddo
enddo
end subroutine table_dropEvap
!=================================================================================================================
! Rain collecting graupel (and inverse). Explicit CE integration.
subroutine qr_acr_qg(NRHGtable)
implicit none
INTEGER, INTENT(IN) ::NRHGtable
!..Local variables
INTEGER:: i, j, k, m, n, n2, n3, idx_bg
INTEGER:: km, km_s, km_e
DOUBLE PRECISION, DIMENSION(nbg):: N_g
DOUBLE PRECISION, DIMENSION(nbg,NRHGtable):: vg
DOUBLE PRECISION, DIMENSION(nbr):: vr, N_r
DOUBLE PRECISION:: N0_r, N0_g, lam_exp, lamg, lamr
DOUBLE PRECISION:: massg, massr, dvg, dvr, t1, t2, z1, z2, y1, y2
!+---+
do n2 = 1, nbr
! vr(n2) = av_r*Dr(n2)**bv_r * DEXP(-fv_r*Dr(n2))
vr(n2) = -0.1021 + 4.932E3*Dr(n2) - 0.9551E6*Dr(n2)*Dr(n2) &
+ 0.07934E9*Dr(n2)*Dr(n2)*Dr(n2) &
- 0.002362E12*Dr(n2)*Dr(n2)*Dr(n2)*Dr(n2)
enddo
do n3 = 1, NRHGtable
do n = 1, nbg
if (NRHGtable == NRHG) idx_bg = n3
if (NRHGtable == NRHG1) idx_bg = idx_bg1
vg(n,n3) = av_g(idx_bg)*Dg(n)**bv_g(idx_bg)
enddo
enddo
km_s = 0
km_e = ntb_r*ntb_r1 - 1
do km = km_s, km_e
m = km / ntb_r1 + 1
k = mod( km , ntb_r1 ) + 1
lam_exp = (N0r_exp(k)*am_r*crg(1)/r_r(m))**ore1
lamr = lam_exp * (crg(3)*org2*org1)**obmr
N0_r = N0r_exp(k)/(crg(2)*lam_exp) * lamr**cre(2)
do n2 = 1, nbr
N_r(n2) = N0_r*Dr(n2)**mu_r *DEXP(-lamr*Dr(n2))*dtr(n2)
enddo
do n3 = 1, NRHGtable
if (NRHGtable == NRHG) idx_bg = n3
if (NRHGtable == NRHG1) idx_bg = idx_bg1
do j = 1, ntb_g
do i = 1, ntb_g1
lam_exp = (N0g_exp(i)*am_g(idx_bg)*cgg(1,1)/r_g(j))**oge1
lamg = lam_exp * (cgg(3,1)*ogg2*ogg1)**obmg
N0_g = N0g_exp(i)/(cgg(2,1)*lam_exp) * lamg**cge(2,1)
do n = 1, nbg
N_g(n) = N0_g*Dg(n)**mu_g * DEXP(-lamg*Dg(n))*dtg(n)
enddo
t1 = 0.0d0
t2 = 0.0d0
z1 = 0.0d0
z2 = 0.0d0
y1 = 0.0d0
y2 = 0.0d0
do n2 = 1, nbr
massr = am_r * Dr(n2)**bm_r
do n = 1, nbg
massg = am_g(idx_bg) * Dg(n)**bm_g
dvg = 0.5d0*((vr(n2) - vg(n,n3)) + DABS(vr(n2)-vg(n,n3)))
dvr = 0.5d0*((vg(n,n3) - vr(n2)) + DABS(vg(n,n3)-vr(n2)))
t1 = t1+ PI*.25*Ef_rg*(Dg(n)+Dr(n2))*(Dg(n)+Dr(n2)) &
*dvg*massg * N_g(n)* N_r(n2)
z1 = z1+ PI*.25*Ef_rg*(Dg(n)+Dr(n2))*(Dg(n)+Dr(n2)) &
*dvg*massr * N_g(n)* N_r(n2)
y1 = y1+ PI*.25*Ef_rg*(Dg(n)+Dr(n2))*(Dg(n)+Dr(n2)) &
*dvg * N_g(n)* N_r(n2)
t2 = t2+ PI*.25*Ef_rg*(Dg(n)+Dr(n2))*(Dg(n)+Dr(n2)) &
*dvr*massr * N_g(n)* N_r(n2)
y2 = y2+ PI*.25*Ef_rg*(Dg(n)+Dr(n2))*(Dg(n)+Dr(n2)) &
*dvr * N_g(n)* N_r(n2)
z2 = z2+ PI*.25*Ef_rg*(Dg(n)+Dr(n2))*(Dg(n)+Dr(n2)) &
*dvr*massg * N_g(n)* N_r(n2)
enddo
97 continue
enddo
tcg_racg(i,j,n3,k,m) = t1
tmr_racg(i,j,n3,k,m) = DMIN1(z1, r_r(m)*1.0d0)
tcr_gacr(i,j,n3,k,m) = t2
tnr_racg(i,j,n3,k,m) = y1
tnr_gacr(i,j,n3,k,m) = y2
enddo
enddo
enddo
enddo
end subroutine qr_acr_qg
!=================================================================================================================
! Rain collecting snow (and inverse). Explicit CE integration.
subroutine qr_acr_qs
implicit none
!..Local variables
INTEGER:: i, j, k, m, n, n2
INTEGER:: km, km_s, km_e
DOUBLE PRECISION, DIMENSION(nbr):: vr, D1, N_r
DOUBLE PRECISION, DIMENSION(nbs):: vs, N_s
DOUBLE PRECISION:: loga_, a_, b_, second, M0, M2, M3, Mrat, oM3
DOUBLE PRECISION:: N0_r, lam_exp, lamr, slam1, slam2
DOUBLE PRECISION:: dvs, dvr, masss, massr
DOUBLE PRECISION:: t1, t2, t3, t4, z1, z2, z3, z4
DOUBLE PRECISION:: y1, y2, y3, y4
!+---+
do n2 = 1, nbr
! vr(n2) = av_r*Dr(n2)**bv_r * DEXP(-fv_r*Dr(n2))
vr(n2) = -0.1021 + 4.932E3*Dr(n2) - 0.9551E6*Dr(n2)*Dr(n2) &
+ 0.07934E9*Dr(n2)*Dr(n2)*Dr(n2) &
- 0.002362E12*Dr(n2)*Dr(n2)*Dr(n2)*Dr(n2)
D1(n2) = (vr(n2)/av_s)**(1./bv_s)
enddo
do n = 1, nbs
vs(n) = 1.5*av_s*Ds(n)**bv_s * DEXP(-fv_s*Ds(n))
enddo
km_s = 0
km_e = ntb_r*ntb_r1 - 1
do km = km_s, km_e
m = km / ntb_r1 + 1
k = mod( km , ntb_r1 ) + 1
lam_exp = (N0r_exp(k)*am_r*crg(1)/r_r(m))**ore1
lamr = lam_exp * (crg(3)*org2*org1)**obmr
N0_r = N0r_exp(k)/(crg(2)*lam_exp) * lamr**cre(2)
do n2 = 1, nbr
N_r(n2) = N0_r*Dr(n2)**mu_r * DEXP(-lamr*Dr(n2))*dtr(n2)
enddo
do j = 1, ntb_t
do i = 1, ntb_s
!..From the bm_s moment, compute plus one moment. If we are not
!.. using bm_s=2, then we must transform to the pure 2nd moment
!.. (variable called "second") and then to the bm_s+1 moment.
M2 = r_s(i)*oams *1.0d0
if (bm_s.gt.2.0-1.E-3 .and. bm_s.lt.2.0+1.E-3) then
loga_ = sa(1) + sa(2)*Tc(j) + sa(3)*bm_s &
+ sa(4)*Tc(j)*bm_s + sa(5)*Tc(j)*Tc(j) &
+ sa(6)*bm_s*bm_s + sa(7)*Tc(j)*Tc(j)*bm_s &
+ sa(8)*Tc(j)*bm_s*bm_s + sa(9)*Tc(j)*Tc(j)*Tc(j) &
+ sa(10)*bm_s*bm_s*bm_s
a_ = 10.0**loga_
b_ = sb(1) + sb(2)*Tc(j) + sb(3)*bm_s &
+ sb(4)*Tc(j)*bm_s + sb(5)*Tc(j)*Tc(j) &
+ sb(6)*bm_s*bm_s + sb(7)*Tc(j)*Tc(j)*bm_s &
+ sb(8)*Tc(j)*bm_s*bm_s + sb(9)*Tc(j)*Tc(j)*Tc(j) &
+ sb(10)*bm_s*bm_s*bm_s
second = (M2/a_)**(1./b_)
else
second = M2
endif
loga_ = sa(1) + sa(2)*Tc(j) + sa(3)*cse(1) &
+ sa(4)*Tc(j)*cse(1) + sa(5)*Tc(j)*Tc(j) &
+ sa(6)*cse(1)*cse(1) + sa(7)*Tc(j)*Tc(j)*cse(1) &
+ sa(8)*Tc(j)*cse(1)*cse(1) + sa(9)*Tc(j)*Tc(j)*Tc(j) &
+ sa(10)*cse(1)*cse(1)*cse(1)
a_ = 10.0**loga_
b_ = sb(1)+sb(2)*Tc(j)+sb(3)*cse(1) + sb(4)*Tc(j)*cse(1) &
+ sb(5)*Tc(j)*Tc(j) + sb(6)*cse(1)*cse(1) &
+ sb(7)*Tc(j)*Tc(j)*cse(1) + sb(8)*Tc(j)*cse(1)*cse(1) &
+ sb(9)*Tc(j)*Tc(j)*Tc(j)+sb(10)*cse(1)*cse(1)*cse(1)
M3 = a_ * second**b_
oM3 = 1./M3
Mrat = M2*(M2*oM3)*(M2*oM3)*(M2*oM3)
M0 = (M2*oM3)**mu_s
slam1 = M2 * oM3 * Lam0
slam2 = M2 * oM3 * Lam1
do n = 1, nbs
N_s(n) = Mrat*(Kap0*DEXP(-slam1*Ds(n)) &
+ Kap1*M0*Ds(n)**mu_s * DEXP(-slam2*Ds(n)))*dts(n)
enddo
t1 = 0.0d0
t2 = 0.0d0
t3 = 0.0d0
t4 = 0.0d0
z1 = 0.0d0
z2 = 0.0d0
z3 = 0.0d0
z4 = 0.0d0
y1 = 0.0d0
y2 = 0.0d0
y3 = 0.0d0
y4 = 0.0d0
do n2 = 1, nbr
massr = am_r * Dr(n2)**bm_r
do n = 1, nbs
masss = am_s * Ds(n)**bm_s
dvs = 0.5d0*((vr(n2) - vs(n)) + DABS(vr(n2)-vs(n)))
dvr = 0.5d0*((vs(n) - vr(n2)) + DABS(vs(n)-vr(n2)))
if (massr .gt. 1.5*masss) then
t1 = t1+ PI*.25*Ef_rs*(Ds(n)+Dr(n2))*(Ds(n)+Dr(n2)) &
*dvs*masss * N_s(n)* N_r(n2)
z1 = z1+ PI*.25*Ef_rs*(Ds(n)+Dr(n2))*(Ds(n)+Dr(n2)) &
*dvs*massr * N_s(n)* N_r(n2)
y1 = y1+ PI*.25*Ef_rs*(Ds(n)+Dr(n2))*(Ds(n)+Dr(n2)) &
*dvs * N_s(n)* N_r(n2)
else
t3 = t3+ PI*.25*Ef_rs*(Ds(n)+Dr(n2))*(Ds(n)+Dr(n2)) &
*dvs*masss * N_s(n)* N_r(n2)
z3 = z3+ PI*.25*Ef_rs*(Ds(n)+Dr(n2))*(Ds(n)+Dr(n2)) &
*dvs*massr * N_s(n)* N_r(n2)
y3 = y3+ PI*.25*Ef_rs*(Ds(n)+Dr(n2))*(Ds(n)+Dr(n2)) &
*dvs * N_s(n)* N_r(n2)
endif
if (massr .gt. 1.5*masss) then
t2 = t2+ PI*.25*Ef_rs*(Ds(n)+Dr(n2))*(Ds(n)+Dr(n2)) &
*dvr*massr * N_s(n)* N_r(n2)
y2 = y2+ PI*.25*Ef_rs*(Ds(n)+Dr(n2))*(Ds(n)+Dr(n2)) &
*dvr * N_s(n)* N_r(n2)
z2 = z2+ PI*.25*Ef_rs*(Ds(n)+Dr(n2))*(Ds(n)+Dr(n2)) &
*dvr*masss * N_s(n)* N_r(n2)
else
t4 = t4+ PI*.25*Ef_rs*(Ds(n)+Dr(n2))*(Ds(n)+Dr(n2)) &
*dvr*massr * N_s(n)* N_r(n2)
y4 = y4+ PI*.25*Ef_rs*(Ds(n)+Dr(n2))*(Ds(n)+Dr(n2)) &
*dvr * N_s(n)* N_r(n2)
z4 = z4+ PI*.25*Ef_rs*(Ds(n)+Dr(n2))*(Ds(n)+Dr(n2)) &
*dvr*masss * N_s(n)* N_r(n2)
endif
enddo
enddo
tcs_racs1(i,j,k,m) = t1
tmr_racs1(i,j,k,m) = DMIN1(z1, r_r(m)*1.0d0)
tcs_racs2(i,j,k,m) = t3
tmr_racs2(i,j,k,m) = z3
tcr_sacr1(i,j,k,m) = t2
tms_sacr1(i,j,k,m) = z2
tcr_sacr2(i,j,k,m) = t4
tms_sacr2(i,j,k,m) = z4
tnr_racs1(i,j,k,m) = y1
tnr_racs2(i,j,k,m) = y3
tnr_sacr1(i,j,k,m) = y2
tnr_sacr2(i,j,k,m) = y4
enddo
enddo
enddo
end subroutine qr_acr_qs
!=================================================================================================================
! This is a literal adaptation of Bigg (1954) probability of drops of a particular volume freezing.
! Given this probability, simply freeze the proportion of drops summing their masses.
subroutine freezeH2O
implicit none
!..Local variables
INTEGER:: i, j, k, m, n, n2
DOUBLE PRECISION :: N_r, N_c
DOUBLE PRECISION, DIMENSION(nbr):: massr
DOUBLE PRECISION, DIMENSION(nbc):: massc
DOUBLE PRECISION:: sum1, sum2, sumn1, sumn2, &
prob, vol, Texp, orho_w, &
lam_exp, lamr, N0_r, lamc, N0_c, y
INTEGER:: nu_c
REAL:: T_adjust
orho_w = 1./rho_w2
do n2 = 1, nbr
massr(n2) = am_r*Dr(n2)**bm_r
enddo
do n = 1, nbc
massc(n) = am_r*Dc(n)**bm_r
enddo
!..Freeze water (smallest drops become cloud ice, otherwise graupel).
do m = 1, ntb_IN
T_adjust = MAX(-3.0, MIN(3.0 - ALOG10(Nt_IN(m)), 3.0))
do k = 1, 45
! print*, ' Freezing water for temp = ', -k
Texp = DEXP( REAL(k,KIND=dp) - T_adjust*1.0D0 ) - 1.0D0
do j = 1, ntb_r1
do i = 1, ntb_r
lam_exp = (N0r_exp(j)*am_r*crg(1)/r_r(i))**ore1
lamr = lam_exp * (crg(3)*org2*org1)**obmr
N0_r = N0r_exp(j)/(crg(2)*lam_exp) * lamr**cre(2)
sum1 = 0.0d0
sum2 = 0.0d0
sumn1 = 0.0d0
sumn2 = 0.0d0
do n2 = nbr, 1, -1
N_r = N0_r*Dr(n2)**mu_r*DEXP(-lamr*Dr(n2))*dtr(n2)
vol = massr(n2)*orho_w
prob = MAX(0.0D0, 1.0D0 - DEXP(-120.0D0*vol*5.2D-4 * Texp))
if (massr(n2) .lt. xm0g) then
sumn1 = sumn1 + prob*N_r
sum1 = sum1 + prob*N_r*massr(n2)
else
sumn2 = sumn2 + prob*N_r
sum2 = sum2 + prob*N_r*massr(n2)
endif
if ((sum1+sum2).ge.r_r(i)) EXIT
enddo
tpi_qrfz(i,j,k,m) = sum1
tni_qrfz(i,j,k,m) = sumn1
tpg_qrfz(i,j,k,m) = sum2
tnr_qrfz(i,j,k,m) = sumn2
enddo
enddo
do j = 1, nbc
nu_c = MIN(15, NINT(1000.E6/t_Nc(j)) + 2)
do i = 1, ntb_c
lamc = (t_Nc(j)*am_r* ccg(2,nu_c) * ocg1(nu_c) / r_c(i))**obmr
N0_c = t_Nc(j)*ocg1(nu_c) * lamc**cce(1,nu_c)
sum1 = 0.0d0
sumn2 = 0.0d0
do n = nbc, 1, -1
vol = massc(n)*orho_w
prob = MAX(0.0D0, 1.0D0 - DEXP(-120.0D0*vol*5.2D-4 * Texp))
N_c = N0_c*Dc(n)**nu_c*EXP(-lamc*Dc(n))*dtc(n)
sumn2 = MIN(t_Nc(j), sumn2 + prob*N_c)
sum1 = sum1 + prob*N_c*massc(n)
if (sum1 .ge. r_c(i)) EXIT
enddo
tpi_qcfz(i,j,k,m) = sum1
tni_qcfz(i,j,k,m) = sumn2
enddo
enddo
enddo
enddo
end subroutine freezeH2O
!=================================================================================================================
! Cloud ice converting to snow since portion greater than min snow size. Given cloud ice content (kg/m**3),
! number concentration (#/m**3) and gamma shape parameter, mu_i, break the distrib into bins and figure out
! the mass/number of ice with sizes larger than D0s. Also, compute incomplete gamma function for the
! integration of ice depositional growth from diameter=0 to D0s. Amount of ice depositional growth is this
! portion of distrib while larger diameters contribute to snow growth (as in Harrington et al. 1995).
subroutine qi_aut_qs
implicit none
!..Local variables
INTEGER:: i, j, n2
DOUBLE PRECISION, DIMENSION(nbi):: N_i
DOUBLE PRECISION:: N0_i, lami, Di_mean, t1, t2
REAL:: xlimit_intg
do j = 1, ntb_i1
do i = 1, ntb_i
lami = (am_i*cig(2)*oig1*Nt_i(j)/r_i(i))**obmi
Di_mean = (bm_i + mu_i + 1.) / lami
N0_i = Nt_i(j)*oig1 * lami**cie(1)
t1 = 0.0d0
t2 = 0.0d0
if (SNGL(Di_mean) .gt. 5.*D0s) then
t1 = r_i(i)
t2 = Nt_i(j)
tpi_ide(i,j) = 0.0D0
elseif (SNGL(Di_mean) .lt. D0i) then
t1 = 0.0D0
t2 = 0.0D0
tpi_ide(i,j) = 1.0D0
else
xlimit_intg = lami*D0s
tpi_ide(i,j) = gamma_p(mu_i+2.0, xlimit_intg) * 1.0D0
do n2 = 1, nbi
N_i(n2) = N0_i*Di(n2)**mu_i * DEXP(-lami*Di(n2))*dti(n2)
if (Di(n2).ge.D0s) then
t1 = t1 + N_i(n2) * am_i*Di(n2)**bm_i
t2 = t2 + N_i(n2)
endif
enddo
endif
tps_iaus(i,j) = t1
tni_iaus(i,j) = t2
enddo
enddo
end subroutine qi_aut_qs
#if defined (ccpp_default)
!=================================================================================================================
!..Fill the table of CCN activation data created from parcel model run
!.. by Trude Eidhammer with inputs of aerosol number concentration,
!.. vertical velocity, temperature, lognormal mean aerosol radius, and
!.. hygroscopicity, kappa. The data are read from external file and
!.. contain activated fraction of CCN for given conditions.
!+---+-----------------------------------------------------------------+
!+---+-----------------------------------------------------------------+
!>\ingroup aathompson
!! Fill the table of CCN activation data created from parcel model run
!! by Trude Eidhammer with inputs of aerosol number concentration,
!! vertical velocity, temperature, lognormal mean aerosol radius, and
!! hygroscopicity, kappa. The data are read from external file and
!! contain activated fraction of CCN for given conditions.
subroutine table_ccnAct(errmess,errflag)
implicit none
!..Error handling variables
CHARACTER(len=*), INTENT(INOUT) :: errmess
INTEGER, INTENT(INOUT) :: errflag
!..Local variables
INTEGER:: iunit_mp_th1, i
LOGICAL:: opened
iunit_mp_th1 = -1
DO i = 20,99
INQUIRE ( i , OPENED = opened )
IF ( .NOT. opened ) THEN
iunit_mp_th1 = i
GOTO 2010
ENDIF
ENDDO
2010 CONTINUE
IF ( iunit_mp_th1 < 0 ) THEN
write(0,*) 'module_mp_tempo: table_ccnAct: '// &
'Can not find unused fortran unit to read in lookup table.'
return
ENDIF
!WRITE(*, '(A,I2)') 'module_mp_tempo: opening CCN_ACTIVATE.BIN on unit ',iunit_mp_th1
OPEN(iunit_mp_th1,FILE='CCN_ACTIVATE.BIN', &
FORM='UNFORMATTED',STATUS='OLD',CONVERT='BIG_ENDIAN',ERR=9009)
!sms$serial begin
READ(iunit_mp_th1,ERR=9010) tnccn_act
!sms$serial end
RETURN
9009 CONTINUE
WRITE( errmess , '(A,I2)' ) 'module_mp_tempo: error opening CCN_ACTIVATE.BIN on unit ',iunit_mp_th1
errflag = 1
RETURN
9010 CONTINUE
WRITE( errmess , '(A,I2)' ) 'module_mp_tempo: error reading CCN_ACTIVATE.BIN on unit ',iunit_mp_th1
errflag = 1
RETURN
end subroutine table_ccnAct
!=================================================================================================================
! Rain collecting graupel (and inverse). Explicit CE integration.
subroutine qr_acr_qg_par(NRHGtable, qrqg_file)
implicit none
INTEGER, INTENT(IN) ::NRHGtable
!..Local variables
INTEGER:: i, j, k, m, n, n2, n3, idx_bg
INTEGER:: km, km_s, km_e
DOUBLE PRECISION, DIMENSION(nbg):: N_g
DOUBLE PRECISION, DIMENSION(nbg,NRHGtable):: vg
DOUBLE PRECISION, DIMENSION(nbr):: vr, N_r
DOUBLE PRECISION:: N0_r, N0_g, lam_exp, lamg, lamr
DOUBLE PRECISION:: massg, massr, dvg, dvr, t1, t2, z1, z2, y1, y2
LOGICAL force_read_thompson, write_thompson_tables
LOGICAL lexist,lopen
INTEGER good,ierr
character(len=*), intent(in) :: qrqg_file
force_read_thompson = .false.
write_thompson_tables = .false.
!+---+
good = 0
INQUIRE(FILE=qrqg_file, EXIST=lexist)
#ifdef MPI
call MPI_BARRIER(mpi_communicator,ierr)
#endif
IF ( lexist ) THEN
OPEN(63,file=qrqg_file,form="unformatted",err=1234)
!sms$serial begin
READ(63,err=1234) tcg_racg
READ(63,err=1234) tmr_racg
READ(63,err=1234) tcr_gacr
!! READ(63,err=1234) tmg_gacr
READ(63,err=1234) tnr_racg
READ(63,err=1234) tnr_gacr
!sms$serial end
good = 1
1234 CONTINUE
IF ( good .NE. 1 ) THEN
INQUIRE(63,opened=lopen)
IF (lopen) THEN
IF( force_read_thompson ) THEN
write(0,*) "Error reading "//qrqg_file//" Aborting because force_read_thompson is .true."
return
ENDIF
CLOSE(63)
ELSE
IF( force_read_thompson ) THEN
write(0,*) "Error opening "//qrqg_file//" Aborting because force_read_thompson is .true."
return
ENDIF
ENDIF
ELSE
INQUIRE(63,opened=lopen)
IF (lopen) THEN
CLOSE(63)
ENDIF
ENDIF
ELSE
IF( force_read_thompson ) THEN
write(0,*) "Non-existent "//qrqg_file//" Aborting because force_read_thompson is .true."
return
ENDIF
ENDIF
IF (.NOT. good .EQ. 1 ) THEN
if (thompson_table_writer) then
write_thompson_tables = .true.
write(0,*) "ThompMP: computing qr_acr_qg"
endif
do n2 = 1, nbr
! vr(n2) = av_r*Dr(n2)**bv_r * DEXP(-fv_r*Dr(n2))
vr(n2) = -0.1021 + 4.932E3*Dr(n2) - 0.9551E6*Dr(n2)*Dr(n2) &
+ 0.07934E9*Dr(n2)*Dr(n2)*Dr(n2) &
- 0.002362E12*Dr(n2)*Dr(n2)*Dr(n2)*Dr(n2)
enddo
! do n = 1, nbg
! vg(n) = av_g*Dg(n)**bv_g
! enddo
do n3 = 1, NRHGtable
do n = 1, nbg
if (NRHGtable == NRHG) idx_bg = n3
if (NRHGtable == NRHG1) idx_bg = idx_bg1
vg(n,n3) = av_g(idx_bg)*Dg(n)**bv_g(idx_bg)
enddo
enddo
!..Note values returned from wrf_dm_decomp1d are zero-based, add 1 for
!.. fortran indices. J. Michalakes, 2009Oct30.
#if ( defined( DM_PARALLEL ) && ( ! defined( STUBMPI ) ) )
CALL wrf_dm_decomp1d ( ntb_r*ntb_r1, km_s, km_e )
#else
km_s = 0
km_e = ntb_r*ntb_r1 - 1
#endif
do km = km_s, km_e
m = km / ntb_r1 + 1
k = mod( km , ntb_r1 ) + 1
lam_exp = (N0r_exp(k)*am_r*crg(1)/r_r(m))**ore1
lamr = lam_exp * (crg(3)*org2*org1)**obmr
N0_r = N0r_exp(k)/(crg(2)*lam_exp) * lamr**cre(2)
do n2 = 1, nbr
N_r(n2) = N0_r*Dr(n2)**mu_r *DEXP(-lamr*Dr(n2))*dtr(n2)
enddo
do n3 = 1, NRHGtable
if (NRHGtable == NRHG) idx_bg = n3
if (NRHGtable == NRHG1) idx_bg = idx_bg1
do j = 1, ntb_g
do i = 1, ntb_g1
lam_exp = (N0g_exp(i)*am_g(idx_bg)*cgg(1,1)/r_g(j))**oge1
lamg = lam_exp * (cgg(3,1)*ogg2*ogg1)**obmg
N0_g = N0g_exp(i)/(cgg(2,1)*lam_exp) * lamg**cge(2,1)
do n = 1, nbg
N_g(n) = N0_g*Dg(n)**mu_g * DEXP(-lamg*Dg(n))*dtg(n)
enddo
t1 = 0.0d0
t2 = 0.0d0
z1 = 0.0d0
z2 = 0.0d0
y1 = 0.0d0
y2 = 0.0d0
do n2 = 1, nbr
massr = am_r * Dr(n2)**bm_r
do n = 1, nbg
massg = am_g(idx_bg) * Dg(n)**bm_g
dvg = 0.5d0*((vr(n2) - vg(n,n3)) + DABS(vr(n2)-vg(n,n3)))
dvr = 0.5d0*((vg(n,n3) - vr(n2)) + DABS(vg(n,n3)-vr(n2)))
t1 = t1+ PI*.25*Ef_rg*(Dg(n)+Dr(n2))*(Dg(n)+Dr(n2)) &
*dvg*massg * N_g(n)* N_r(n2)
z1 = z1+ PI*.25*Ef_rg*(Dg(n)+Dr(n2))*(Dg(n)+Dr(n2)) &
*dvg*massr * N_g(n)* N_r(n2)
y1 = y1+ PI*.25*Ef_rg*(Dg(n)+Dr(n2))*(Dg(n)+Dr(n2)) &
*dvg * N_g(n)* N_r(n2)
t2 = t2+ PI*.25*Ef_rg*(Dg(n)+Dr(n2))*(Dg(n)+Dr(n2)) &
*dvr*massr * N_g(n)* N_r(n2)
y2 = y2+ PI*.25*Ef_rg*(Dg(n)+Dr(n2))*(Dg(n)+Dr(n2)) &
*dvr * N_g(n)* N_r(n2)
z2 = z2+ PI*.25*Ef_rg*(Dg(n)+Dr(n2))*(Dg(n)+Dr(n2)) &
*dvr*massg * N_g(n)* N_r(n2)
enddo
97 continue
enddo
tcg_racg(i,j,n3,k,m) = t1
tmr_racg(i,j,n3,k,m) = DMIN1(z1, r_r(m)*1.0d0)
tcr_gacr(i,j,n3,k,m) = t2
!! tmg_gacr(i,j,k,m) = DMIN1(z2, r_g(j)*1.0d0)
tnr_racg(i,j,n3,k,m) = y1
tnr_gacr(i,j,n3,k,m) = y2
enddo
enddo
enddo
enddo
IF ( write_thompson_tables ) THEN
write(0,*) "Writing "//qrqg_file//" in Tempo MP init"
OPEN(63,file=qrqg_file,form="unformatted",err=9234)
WRITE(63,err=9234) tcg_racg
WRITE(63,err=9234) tmr_racg
WRITE(63,err=9234) tcr_gacr
WRITE(63,err=9234) tnr_racg
WRITE(63,err=9234) tnr_gacr
CLOSE(63)
RETURN ! ----- RETURN
9234 CONTINUE
write(0,*) "Error writing "//qrqg_file
return
ENDIF
ENDIF
end subroutine qr_acr_qg_par
!=================================================================================================================
! Rain collecting snow (and inverse). Explicit CE integration.
subroutine qr_acr_qs_par
implicit none
!..Local variables
INTEGER:: i, j, k, m, n, n2
INTEGER:: km, km_s, km_e
DOUBLE PRECISION, DIMENSION(nbr):: vr, D1, N_r
DOUBLE PRECISION, DIMENSION(nbs):: vs, N_s
DOUBLE PRECISION:: loga_, a_, b_, second, M0, M2, M3, Mrat, oM3
DOUBLE PRECISION:: N0_r, lam_exp, lamr, slam1, slam2
DOUBLE PRECISION:: dvs, dvr, masss, massr
DOUBLE PRECISION:: t1, t2, t3, t4, z1, z2, z3, z4
DOUBLE PRECISION:: y1, y2, y3, y4
LOGICAL force_read_thompson, write_thompson_tables
LOGICAL lexist,lopen
INTEGER good,ierr
!+---+
force_read_thompson = .false.
write_thompson_tables = .false.
good = 0
INQUIRE(FILE=qr_acr_qs_file, EXIST=lexist)
#ifdef MPI
call MPI_BARRIER(mpi_communicator,ierr)
#endif
IF ( lexist ) THEN
!write(0,*) "ThompMP: read "//qr_acr_qs_file//" instead of computing"
OPEN(63,file=qr_acr_qs_file,form="unformatted",err=1234)
!sms$serial begin
READ(63,err=1234)tcs_racs1
READ(63,err=1234)tmr_racs1
READ(63,err=1234)tcs_racs2
READ(63,err=1234)tmr_racs2
READ(63,err=1234)tcr_sacr1
READ(63,err=1234)tms_sacr1
READ(63,err=1234)tcr_sacr2
READ(63,err=1234)tms_sacr2
READ(63,err=1234)tnr_racs1
READ(63,err=1234)tnr_racs2
READ(63,err=1234)tnr_sacr1
READ(63,err=1234)tnr_sacr2
!sms$serial end
good = 1
1234 CONTINUE
IF ( good .NE. 1 ) THEN
INQUIRE(63,opened=lopen)
IF (lopen) THEN
IF( force_read_thompson ) THEN
write(0,*) "Error reading "//qr_acr_qs_file//" Aborting because force_read_thompson is .true."