forked from ESCOMP/PUMAS
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathmicro_mg_utils.F90
3003 lines (2333 loc) · 101 KB
/
micro_mg_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
module micro_mg_utils
!--------------------------------------------------------------------------
!
! This module contains process rates and utility functions used by the MG
! microphysics.
!
! Original MG authors: Andrew Gettelman, Hugh Morrison
! Contributions from: Peter Caldwell, Xiaohong Liu and Steve Ghan
!
! Separated from MG 1.5 by B. Eaton.
! Separated module switched to MG 2.0 and further changes by S. Santos.
!
! for questions contact Hugh Morrison, Andrew Gettelman
! e-mail: [email protected], [email protected]
!
!--------------------------------------------------------------------------
!
! List of required external functions that must be supplied:
! gamma --> standard mathematical gamma function (if gamma is an
! intrinsic, define HAVE_GAMMA_INTRINSICS)
!
!--------------------------------------------------------------------------
!
! Constants that must be specified in the "init" method (module variables):
!
! kind kind of reals (to verify correct linkage only) -
! gravit acceleration due to gravity m s-2
! rair dry air gas constant for air J kg-1 K-1
! rh2o gas constant for water vapor J kg-1 K-1
! cpair specific heat at constant pressure for dry air J kg-1 K-1
! tmelt temperature of melting point for water K
! latvap latent heat of vaporization J kg-1
! latice latent heat of fusion J kg-1
!
!--------------------------------------------------------------------------
#ifndef HAVE_GAMMA_INTRINSICS
use shr_spfn_mod, only: gamma => shr_spfn_gamma
#endif
!++ag DEBUG: remove later
use cam_logfile, only: iulog
!--ag
implicit none
private
save
public :: &
micro_mg_utils_init, &
size_dist_param_liq, &
size_dist_param_basic, &
avg_diameter, &
rising_factorial, &
ice_deposition_sublimation, &
ice_deposition_sublimation_mg4, & !! ktc
sb2001v2_liq_autoconversion,&
sb2001v2_accre_cld_water_rain,&
kk2000_liq_autoconversion, &
ice_autoconversion, &
immersion_freezing, &
contact_freezing, &
snow_self_aggregation, &
ice_self_aggregation, & !! mg4 added for using instead of snow_self_agg
accrete_cloud_water_snow, &
accrete_cloud_water_ice, & !! mg4 added for using instead of cloud_water_snow
secondary_ice_production, &
accrete_rain_snow, &
accrete_rain_ice, & !! mg4 added for using instead of rain_snow
heterogeneous_rain_freezing, &
accrete_cloud_water_rain, &
self_collection_rain, &
accrete_cloud_ice_snow, &
evaporate_sublimate_precip, &
bergeron_process_snow, &
!++ag
graupel_collecting_snow, &
graupel_collecting_rain, &
graupel_collecting_cld_water, &
graupel_riming_liquid_snow, &
graupel_rain_riming_snow, &
graupel_rime_splintering, &
evaporate_sublimate_precip_graupel, &
! graupel_sublimate_evap
!--ag
evaporate_sublimate_precip_mg4, &
evaporate_sublimate_precip_graupel_mg4, &
access_lookup_table, & !! mg4
access_lookup_table_coll, & !! mg4
init_lookup_table !! mg4
! 8 byte real and integer
integer, parameter, public :: r8 = selected_real_kind(12)
integer, parameter, public :: i8 = selected_int_kind(18)
public :: MGHydrometeorProps
type :: MGHydrometeorProps
! Density (kg/m^3)
real(r8) :: rho
! Information for size calculations.
! Basic calculation of mean size is:
! lambda = (shape_coef*nic/qic)^(1/eff_dim)
! Then lambda is constrained by bounds.
real(r8) :: eff_dim
real(r8) :: shape_coef
real(r8) :: lambda_bounds(2)
! Minimum average particle mass (kg).
! Limit is applied at the beginning of the size distribution calculations.
real(r8) :: min_mean_mass
end type MGHydrometeorProps
interface MGHydrometeorProps
module procedure NewMGHydrometeorProps
end interface
type(MGHydrometeorProps), public :: mg_liq_props
type(MGHydrometeorProps), public :: mg_ice_props
type(MGHydrometeorProps), public :: mg_rain_props
type(MGHydrometeorProps), public :: mg_snow_props
type(MGHydrometeorProps), public :: mg_graupel_props
type(MGHydrometeorProps), public :: mg_hail_props
interface size_dist_param_liq
module procedure size_dist_param_liq_vect
module procedure size_dist_param_liq_line
end interface
interface size_dist_param_basic
module procedure size_dist_param_basic_vect
module procedure size_dist_param_basic_line
end interface
!=================================================
! Public module parameters (mostly for MG itself)
!=================================================
! Pi to 20 digits; more than enough to reach the limit of double precision.
real(r8), parameter, public :: pi = 3.14159265358979323846_r8
! "One minus small number": number near unity for round-off issues.
real(r8), parameter, public :: omsm = 1._r8 - 1.e-5_r8
! Smallest mixing ratio considered in microphysics.
real(r8), parameter, public :: qsmall = 1.e-18_r8
! minimum allowed cloud fraction
real(r8), parameter, public :: mincld = 0.0001_r8
real(r8), parameter, public :: rhosn = 250._r8 ! bulk density snow
real(r8), parameter, public :: rhoi = 500._r8 ! bulk density ice
real(r8), parameter, public :: rhow = 1000._r8 ! bulk density liquid
real(r8), parameter, public :: rhows = 917._r8 ! bulk density water solid
!Hail and Graupel (set in MG3)
real(r8), parameter, public :: rhog = 500._r8
real(r8), parameter, public :: rhoh = 900._r8
! fall speed parameters, V = aD^b (V is in m/s)
! droplets
real(r8), parameter, public :: ac = 3.e7_r8
real(r8), parameter, public :: bc = 2._r8
! snow
real(r8), parameter, public :: as = 11.72_r8
real(r8), parameter, public :: bs = 0.41_r8
! cloud ice
real(r8), parameter, public :: ai = 700._r8
real(r8), parameter, public :: bi = 1._r8
! small cloud ice (r< 10 um) - sphere, bulk density
real(r8), parameter, public :: aj = ac*((rhoi/rhows)**(bc/3._r8))*rhows/rhow
real(r8), parameter, public :: bj = bc
! rain
real(r8), parameter, public :: ar = 841.99667_r8
real(r8), parameter, public :: br = 0.8_r8
! graupel
real(r8), parameter, public :: ag = 19.3_r8
real(r8), parameter, public :: bg = 0.37_r8
! hail
real(r8), parameter, public :: ah = 114.5_r8
real(r8), parameter, public :: bh = 0.5_r8
! mass of new crystal due to aerosol freezing and growth (kg)
! Make this consistent with the lower bound, to support UTLS and
! stratospheric ice, and the smaller ice size limit.
real(r8), parameter, public :: mi0 = 4._r8/3._r8*pi*rhoi*(1.e-6_r8)**3
! mass of new graupel particle (assume same as mi0 for now, may want to make bigger?)
!real(r8), parameter, public :: mg0 = 4._r8/3._r8*pi*rhoi*(1.e-6_r8)**3
!or set based on M2005:
real(r8), parameter, public :: mg0 = 1.6e-10_r8
! radius of contact nuclei
real(r8), parameter, public :: mmult = 4._r8/3._r8*pi*rhoi*(5.e-6_r8)**3
!=================================================
! Private module parameters
!=================================================
! Signaling NaN bit pattern that represents a limiter that's turned off.
integer(i8), parameter :: limiter_off = int(Z'7FF1111111111111', i8)
! alternate threshold used for some in-cloud mmr
real(r8), parameter, public :: icsmall = 1.e-8_r8
! particle mass-diameter relationship
! currently we assume spherical particles for cloud ice/snow
! m = cD^d
! exponent
real(r8), parameter :: dsph = 3._r8
! Bounds for mean diameter for different constituents.
real(r8), parameter :: lam_bnd_rain(2) = 1._r8/[500.e-6_r8, 20.e-6_r8]
real(r8), parameter :: lam_bnd_snow(2) = 1._r8/[2000.e-6_r8, 10.e-6_r8]
! Minimum average mass of particles.
real(r8), parameter :: min_mean_mass_liq = 1.e-20_r8
real(r8), parameter :: min_mean_mass_ice = 1.e-20_r8
! ventilation parameters
! for snow
real(r8), parameter :: f1s = 0.86_r8
real(r8), parameter :: f2s = 0.28_r8
! for rain
real(r8), parameter :: f1r = 0.78_r8
real(r8), parameter :: f2r = 0.308_r8
! collection efficiencies
! aggregation of cloud ice and snow
real(r8), parameter :: eii = 0.5_r8
! collection efficiency, ice-droplet collisions
real(r8), parameter, public :: ecid = 0.7_r8
! collection efficiency between droplets/rain and snow/rain
real(r8), parameter, public :: ecr = 1.0_r8
! immersion freezing parameters, bigg 1953
real(r8), parameter :: bimm = 100._r8
real(r8), parameter :: aimm = 0.66_r8
! Mass of each raindrop created from autoconversion.
real(r8), parameter :: droplet_mass_25um = 4._r8/3._r8*pi*rhow*(25.e-6_r8)**3
real(r8), parameter :: droplet_mass_40um = 4._r8/3._r8*pi*rhow*(40.e-6_r8)**3
!!!!!!
! MG4 look up table parameters
!!!!!!
! ice microphysics lookup table array dimensions
integer, parameter, public :: tsize = 3
integer, parameter, public :: isize = 20
integer, parameter, public :: jsize = 20
integer, parameter, public :: rcollsize = 30
! number of ice microphysical quantities used from lookup table
integer, parameter, public :: tabsize = 14
! number of ice-rain collection microphysical quantities used from lookup table
integer, parameter, public :: colltabsize = 2
!ice lookup table values
real(r8) :: itab(tsize,isize,jsize,tabsize)
!ice lookup table values for ice-rain collision/collection
real(r8) :: itabcoll(tsize,isize,jsize,rcollsize,colltabsize)
private :: itab,itabcoll
!=========================================================
! Constants set in initialization
!=========================================================
! Set using arguments to micro_mg_init
real(r8) :: rv ! water vapor gas constant
real(r8) :: cpp ! specific heat of dry air
real(r8) :: tmelt ! freezing point of water (K)
real(r8) :: ra ! dry air gas constant
! latent heats of:
real(r8) :: xxlv ! vaporization
real(r8) :: xlf ! freezing
real(r8) :: xxls ! sublimation
! additional constants to help speed up code
real(r8) :: gamma_bs_plus3
real(r8) :: gamma_half_br_plus5
real(r8) :: gamma_half_bs_plus5
real(r8) :: gamma_2bs_plus2
!=========================================================
! Utilities that are cheaper if the compiler knows that
! some argument is an integer.
!=========================================================
interface rising_factorial
module procedure rising_factorial_r8
module procedure rising_factorial_integer
end interface rising_factorial
interface var_coef
module procedure var_coef_r8
module procedure var_coef_integer
end interface var_coef
!==========================================================================
contains
!==========================================================================
! Initialize module variables.
!
! "kind" serves no purpose here except to check for unlikely linking
! issues; always pass in the kind for a double precision real.
!
! "errstring" is the only output; it is blank if there is no error, or set
! to a message if there is an error.
!
! Check the list at the top of this module for descriptions of all other
! arguments.
subroutine micro_mg_utils_init( kind, rair, rh2o, cpair, tmelt_in, latvap, &
latice, dcs, errstring)
integer, intent(in) :: kind
real(r8), intent(in) :: rair
real(r8), intent(in) :: rh2o
real(r8), intent(in) :: cpair
real(r8), intent(in) :: tmelt_in
real(r8), intent(in) :: latvap
real(r8), intent(in) :: latice
real(r8), intent(in) :: dcs
character(128), intent(out) :: errstring
! Name this array to workaround an XLF bug (otherwise could just use the
! expression that sets it).
real(r8) :: ice_lambda_bounds(2)
!-----------------------------------------------------------------------
errstring = ' '
if( kind .ne. r8 ) then
errstring = 'micro_mg_init: KIND of reals does not match'
return
endif
! declarations for MG code (transforms variable names)
rv= rh2o ! water vapor gas constant
cpp = cpair ! specific heat of dry air
tmelt = tmelt_in
ra = rair !dry air gas constant
! latent heats
xxlv = latvap ! latent heat vaporization
xlf = latice ! latent heat freezing
xxls = xxlv + xlf ! latent heat of sublimation
! Define constants to help speed up code (this limits calls to gamma function)
gamma_bs_plus3=gamma(3._r8+bs)
gamma_half_br_plus5=gamma(5._r8/2._r8+br/2._r8)
gamma_half_bs_plus5=gamma(5._r8/2._r8+bs/2._r8)
gamma_2bs_plus2=gamma(2._r8*bs+2._r8)
! Don't specify lambda bounds for cloud liquid, as they are determined by
! pgam dynamically.
mg_liq_props = MGHydrometeorProps(rhow, dsph, &
min_mean_mass=min_mean_mass_liq)
! Mean ice diameter can not grow bigger than twice the autoconversion
! threshold for snow.
ice_lambda_bounds = 1._r8/[2._r8*dcs, 1.e-6_r8]
mg_ice_props = MGHydrometeorProps(rhoi, dsph, &
ice_lambda_bounds, min_mean_mass_ice)
mg_rain_props = MGHydrometeorProps(rhow, dsph, lam_bnd_rain)
mg_snow_props = MGHydrometeorProps(rhosn, dsph, lam_bnd_snow)
mg_graupel_props = MGHydrometeorProps(rhog, dsph, lam_bnd_snow)
mg_hail_props = MGHydrometeorProps(rhoh, dsph, lam_bnd_snow)
end subroutine micro_mg_utils_init
! Constructor for a constituent property object.
function NewMGHydrometeorProps(rho, eff_dim, lambda_bounds, min_mean_mass) &
result(res)
real(r8), intent(in) :: rho, eff_dim
real(r8), intent(in), optional :: lambda_bounds(2), min_mean_mass
type(MGHydrometeorProps) :: res
res%rho = rho
res%eff_dim = eff_dim
if (present(lambda_bounds)) then
res%lambda_bounds = lambda_bounds
else
res%lambda_bounds = no_limiter()
end if
if (present(min_mean_mass)) then
res%min_mean_mass = min_mean_mass
else
res%min_mean_mass = no_limiter()
end if
res%shape_coef = rho*pi*gamma(eff_dim+1._r8)/6._r8
end function NewMGHydrometeorProps
!========================================================================
!FORMULAS
!========================================================================
! Use gamma function to implement rising factorial extended to the reals.
pure function rising_factorial_r8(x, n) result(res)
real(r8), intent(in) :: x, n
real(r8) :: res
res = gamma(x+n)/gamma(x)
end function rising_factorial_r8
! Rising factorial can be performed much cheaper if n is a small integer.
pure function rising_factorial_integer(x, n) result(res)
real(r8), intent(in) :: x
integer, intent(in) :: n
real(r8) :: res
integer :: i
real(r8) :: factor
res = 1._r8
factor = x
do i = 1, n
res = res * factor
factor = factor + 1._r8
end do
end function rising_factorial_integer
! Calculate correction due to latent heat for evaporation/sublimation
elemental function calc_ab(t, qv, xxl) result(ab)
real(r8), intent(in) :: t ! Temperature
real(r8), intent(in) :: qv ! Saturation vapor pressure
real(r8), intent(in) :: xxl ! Latent heat
real(r8) :: ab
real(r8) :: dqsdt
dqsdt = xxl*qv / (rv * t**2)
ab = 1._r8 + dqsdt*xxl/cpp
end function calc_ab
! get cloud droplet size distribution parameters
elemental subroutine size_dist_param_liq_line(props, qcic, ncic, rho, pgam, lamc)
type(MGHydrometeorProps), intent(in) :: props
real(r8), intent(in) :: qcic
real(r8), intent(inout) :: ncic
real(r8), intent(in) :: rho
real(r8), intent(out) :: pgam
real(r8), intent(out) :: lamc
type(MGHydrometeorProps) :: props_loc
if (qcic > qsmall) then
! Local copy of properties that can be modified.
! (Elemental routines that operate on arrays can't modify scalar
! arguments.)
props_loc = props
! Get pgam from fit Rotstayn and Liu 2003 (changed from Martin 1994 for CAM6)
pgam = 1.0_r8 - 0.7_r8 * exp(-0.008_r8*1.e-6_r8*ncic*rho)
pgam = 1._r8/(pgam**2) - 1._r8
pgam = max(pgam, 2._r8)
! Set coefficient for use in size_dist_param_basic.
! The 3D case is so common and optimizable that we specialize it:
if (props_loc%eff_dim == 3._r8) then
props_loc%shape_coef = pi / 6._r8 * props_loc%rho * &
rising_factorial(pgam+1._r8, 3)
else
props_loc%shape_coef = pi / 6._r8 * props_loc%rho * &
rising_factorial(pgam+1._r8, props_loc%eff_dim)
end if
! Limit to between 2 and 50 microns mean size.
props_loc%lambda_bounds = (pgam+1._r8)*1._r8/[50.e-6_r8, 2.e-6_r8]
call size_dist_param_basic(props_loc, qcic, ncic, lamc)
else
! pgam not calculated in this case, so set it to a value likely to
! cause an error if it is accidentally used
! (gamma function undefined for negative integers)
pgam = -100._r8
lamc = 0._r8
end if
end subroutine size_dist_param_liq_line
! get cloud droplet size distribution parameters
subroutine size_dist_param_liq_vect(props, qcic, ncic, rho, pgam, lamc, mgncol)
type(mghydrometeorprops), intent(in) :: props
integer, intent(in) :: mgncol
real(r8), dimension(mgncol), intent(in) :: qcic
real(r8), dimension(mgncol), intent(inout) :: ncic
real(r8), dimension(mgncol), intent(in) :: rho
real(r8), dimension(mgncol), intent(out) :: pgam
real(r8), dimension(mgncol), intent(out) :: lamc
type(mghydrometeorprops) :: props_loc
integer :: i
do i=1,mgncol
if (qcic(i) > qsmall) then
! Local copy of properties that can be modified.
! (Elemental routines that operate on arrays can't modify scalar
! arguments.)
props_loc = props
! Get pgam from fit Rotstayn and Liu 2003 (changed from Martin 1994 for CAM6)
pgam(i) = 1.0_r8 - 0.7_r8 * exp(-0.008_r8*1.e-6_r8*ncic(i)*rho(i))
pgam(i) = 1._r8/(pgam(i)**2) - 1._r8
pgam(i) = max(pgam(i), 2._r8)
endif
enddo
do i=1,mgncol
if (qcic(i) > qsmall) then
! Set coefficient for use in size_dist_param_basic.
! The 3D case is so common and optimizable that we specialize
! it:
if (props_loc%eff_dim == 3._r8) then
props_loc%shape_coef = pi / 6._r8 * props_loc%rho * &
rising_factorial(pgam(i)+1._r8, 3)
else
props_loc%shape_coef = pi / 6._r8 * props_loc%rho * &
rising_factorial(pgam(i)+1._r8, props_loc%eff_dim)
end if
! Limit to between 2 and 50 microns mean size.
props_loc%lambda_bounds(1) = (pgam(i)+1._r8)*1._r8/50.e-6_r8
props_loc%lambda_bounds(2) = (pgam(i)+1._r8)*1._r8/2.e-6_r8
call size_dist_param_basic(props_loc, qcic(i), ncic(i), lamc(i))
endif
enddo
do i=1,mgncol
if (qcic(i) <= qsmall) then
! pgam not calculated in this case, so set it to a value likely to
! cause an error if it is accidentally used
! (gamma function undefined for negative integers)
pgam(i) = -100._r8
lamc(i) = 0._r8
end if
enddo
end subroutine size_dist_param_liq_vect
! Basic routine for getting size distribution parameters.
elemental subroutine size_dist_param_basic_line(props, qic, nic, lam, n0)
type(MGHydrometeorProps), intent(in) :: props
real(r8), intent(in) :: qic
real(r8), intent(inout) :: nic
real(r8), intent(out) :: lam
real(r8), intent(out), optional :: n0
if (qic > qsmall) then
! add upper limit to in-cloud number concentration to prevent
! numerical error
if (limiter_is_on(props%min_mean_mass)) then
nic = min(nic, qic / props%min_mean_mass)
end if
! lambda = (c n/q)^(1/d)
lam = (props%shape_coef * nic/qic)**(1._r8/props%eff_dim)
! check for slope
! adjust vars
if (lam < props%lambda_bounds(1)) then
lam = props%lambda_bounds(1)
nic = lam**(props%eff_dim) * qic/props%shape_coef
else if (lam > props%lambda_bounds(2)) then
lam = props%lambda_bounds(2)
nic = lam**(props%eff_dim) * qic/props%shape_coef
end if
else
lam = 0._r8
end if
if (present(n0)) n0 = nic * lam
end subroutine size_dist_param_basic_line
subroutine size_dist_param_basic_vect(props, qic, nic, lam, mgncol, n0)
type (mghydrometeorprops), intent(in) :: props
integer, intent(in) :: mgncol
real(r8), dimension(mgncol), intent(in) :: qic
real(r8), dimension(mgncol), intent(inout) :: nic
real(r8), dimension(mgncol), intent(out) :: lam
real(r8), dimension(mgncol), intent(out), optional :: n0
integer :: i
do i=1,mgncol
if (qic(i) > qsmall) then
! add upper limit to in-cloud number concentration to prevent
! numerical error
if (limiter_is_on(props%min_mean_mass)) then
nic(i) = min(nic(i), qic(i) / props%min_mean_mass)
end if
! lambda = (c n/q)^(1/d)
lam(i) = (props%shape_coef * nic(i)/qic(i))**(1._r8/props%eff_dim)
! check for slope
! adjust vars
if (lam(i) < props%lambda_bounds(1)) then
lam(i) = props%lambda_bounds(1)
nic(i) = lam(i)**(props%eff_dim) * qic(i)/props%shape_coef
else if (lam(i) > props%lambda_bounds(2)) then
lam(i) = props%lambda_bounds(2)
nic(i) = lam(i)**(props%eff_dim) * qic(i)/props%shape_coef
end if
else
lam(i) = 0._r8
end if
enddo
if (present(n0)) n0 = nic * lam
end subroutine size_dist_param_basic_vect
real(r8) elemental function avg_diameter(q, n, rho_air, rho_sub)
! Finds the average diameter of particles given their density, and
! mass/number concentrations in the air.
! Assumes that diameter follows an exponential distribution.
real(r8), intent(in) :: q ! mass mixing ratio
real(r8), intent(in) :: n ! number concentration (per volume)
real(r8), intent(in) :: rho_air ! local density of the air
real(r8), intent(in) :: rho_sub ! density of the particle substance
avg_diameter = (pi * rho_sub * n/(q*rho_air))**(-1._r8/3._r8)
end function avg_diameter
elemental function var_coef_r8(relvar, a) result(res)
! Finds a coefficient for process rates based on the relative variance
! of cloud water.
real(r8), intent(in) :: relvar
real(r8), intent(in) :: a
real(r8) :: res
res = rising_factorial(relvar, a) / relvar**a
end function var_coef_r8
elemental function var_coef_integer(relvar, a) result(res)
! Finds a coefficient for process rates based on the relative variance
! of cloud water.
real(r8), intent(in) :: relvar
integer, intent(in) :: a
real(r8) :: res
res = rising_factorial(relvar, a) / relvar**a
end function var_coef_integer
!========================================================================
!MICROPHYSICAL PROCESS CALCULATIONS
!========================================================================
!========================================================================
! Initial ice deposition and sublimation loop.
! Run before the main loop
! This subroutine written by Peter Caldwell
subroutine ice_deposition_sublimation(t, qv, qi, ni, &
icldm, rho, dv,qvl, qvi, &
berg, vap_dep, ice_sublim, mgncol)
!INPUT VARS:
!===============================================
integer, intent(in) :: mgncol
real(r8), dimension(mgncol), intent(in) :: t
real(r8), dimension(mgncol), intent(in) :: qv
real(r8), dimension(mgncol), intent(in) :: qi
real(r8), dimension(mgncol), intent(in) :: ni
real(r8), dimension(mgncol), intent(in) :: icldm
real(r8), dimension(mgncol), intent(in) :: rho
real(r8), dimension(mgncol), intent(in) :: dv
real(r8), dimension(mgncol), intent(in) :: qvl
real(r8), dimension(mgncol), intent(in) :: qvi
!OUTPUT VARS:
!===============================================
real(r8), dimension(mgncol), intent(out) :: vap_dep !ice deposition (cell-ave value)
real(r8), dimension(mgncol), intent(out) :: ice_sublim !ice sublimation (cell-ave value)
real(r8), dimension(mgncol), intent(out) :: berg !bergeron enhancement (cell-ave value)
!INTERNAL VARS:
!===============================================
real(r8) :: ab
real(r8) :: epsi
real(r8) :: qiic
real(r8) :: niic
real(r8) :: lami
real(r8) :: n0i
integer :: i
do i=1,mgncol
if (qi(i)>=qsmall) then
!GET IN-CLOUD qi, ni
!===============================================
qiic = qi(i)/icldm(i)
niic = ni(i)/icldm(i)
!Compute linearized condensational heating correction
ab=calc_ab(t(i), qvi(i), xxls)
!Get slope and intercept of gamma distn for ice.
call size_dist_param_basic(mg_ice_props, qiic, niic, lami, n0i)
!Get depletion timescale=1/eps
epsi = 2._r8*pi*n0i*rho(i)*Dv(i)/(lami*lami)
!Compute deposition/sublimation
vap_dep(i) = epsi/ab*(qv(i) - qvi(i))
!Make this a grid-averaged quantity
vap_dep(i)=vap_dep(i)*icldm(i)
!Split into deposition or sublimation.
if (t(i) < tmelt .and. vap_dep(i)>0._r8) then
ice_sublim(i)=0._r8
else
! make ice_sublim negative for consistency with other evap/sub processes
ice_sublim(i)=min(vap_dep(i),0._r8)
vap_dep(i)=0._r8
end if
!sublimation occurs @ any T. Not so for berg.
if (t(i) < tmelt) then
!Compute bergeron rate assuming cloud for whole step.
berg(i) = max(epsi/ab*(qvl(i) - qvi(i)), 0._r8)
else !T>frz
berg(i)=0._r8
end if !T<frz
else !where qi<qsmall
berg(i)=0._r8
vap_dep(i)=0._r8
ice_sublim(i)=0._r8
end if !qi>qsmall
enddo
end subroutine ice_deposition_sublimation
subroutine ice_deposition_sublimation_mg4(t, qv, qi, niic, &
icldm, rho, dv,qvl, qvi, &
berg, vap_dep, ice_sublim, &
af1pr5, af1pr14, rhof, mu, sc, &
mgncol)
!INPUT VARS:
!===============================================
integer, intent(in) :: mgncol
real(r8), dimension(mgncol), intent(in) :: t
real(r8), dimension(mgncol), intent(in) :: qv
real(r8), dimension(mgncol), intent(in) :: qi
real(r8), dimension(mgncol), intent(in) :: niic ! trude: input nicc as other routines
real(r8), dimension(mgncol), intent(in) :: icldm
real(r8), dimension(mgncol), intent(in) :: rho
real(r8), dimension(mgncol), intent(in) :: dv
real(r8), dimension(mgncol), intent(in) :: qvl
real(r8), dimension(mgncol), intent(in) :: qvi
real(r8), dimension(mgncol), intent(in) :: af1pr5, af1pr14
real(r8), dimension(mgncol), intent(in) :: rhof
real(r8), dimension(mgncol), intent(in) :: mu
real(r8), dimension(mgncol), intent(in) :: sc
!OUTPUT VARS:
!===============================================
real(r8), dimension(mgncol), intent(out) :: vap_dep !ice deposition (cell-ave value)
real(r8), dimension(mgncol), intent(out) :: ice_sublim !ice sublimation (cell-ave value)
real(r8), dimension(mgncol), intent(out) :: berg !bergeron enhancement (cell-ave value)
!INTERNAL VARS:
!===============================================
real(r8) :: ab
real(r8) :: epsi
real(r8) :: qiic
real(r8) :: lami
real(r8) :: n0i
real(r8), parameter :: thrd = 1._r8/3._r8
integer :: i
do i=1,mgncol
if (qi(i)>=qsmall) then
!GET IN-CLOUD qi, ni
!===============================================
qiic = qi(i)/icldm(i)
!Compute linearized condensational heating correction
ab=calc_ab(t(i), qvi(i), xxls)
if ( t(i).lt.tmelt) then
epsi = (af1pr5(i)+af1pr14(i)*sc(i)**thrd*(rhof(i)*rho(i)/mu(i))**0.5_r8)*2._r8*pi* &
rho(i)*dv(i)
else
epsi = 0._r8
endif
!Compute deposition/sublimation
vap_dep(i) = epsi/ab*(qv(i) - qvi(i))
!Make this a grid-averaged quantity
vap_dep(i)=vap_dep(i)*icldm(i)
!Split into deposition or sublimation.
if (t(i) < tmelt .and. vap_dep(i)>0._r8) then
ice_sublim(i)=0._r8
else
! make ice_sublim negative for consistency with other evap/sub processes
ice_sublim(i)=min(vap_dep(i),0._r8)
vap_dep(i)=0._r8
end if
!sublimation occurs @ any T. Not so for berg.
if (t(i) < tmelt) then
!Compute bergeron rate assuming cloud for whole step.
berg(i) = max(epsi/ab*(qvl(i) - qvi(i)), 0._r8)
else !T>frz
berg(i)=0._r8
end if !T<frz
else !where qi<qsmall
berg(i)=0._r8
vap_dep(i)=0._r8
ice_sublim(i)=0._r8
end if !qi>qsmall
enddo
end subroutine ice_deposition_sublimation_mg4
!========================================================================
! autoconversion of cloud liquid water to rain
! formula from Khrouditnov and Kogan (2000), modified for sub-grid distribution of qc
! minimum qc of 1 x 10^-8 prevents floating point error
subroutine kk2000_liq_autoconversion(microp_uniform, qcic, &
ncic, rho, relvar, prc, nprc, nprc1, mgncol)
integer, intent(in) :: mgncol
logical, intent(in) :: microp_uniform
real(r8), dimension(mgncol), intent(in) :: qcic
real(r8), dimension(mgncol), intent(in) :: ncic
real(r8), dimension(mgncol), intent(in) :: rho
real(r8), dimension(mgncol), intent(in) :: relvar
real(r8), dimension(mgncol), intent(out) :: prc
real(r8), dimension(mgncol), intent(out) :: nprc
real(r8), dimension(mgncol), intent(out) :: nprc1
real(r8), dimension(mgncol) :: prc_coef
integer :: i
! Take variance into account, or use uniform value.
if (.not. microp_uniform) then
prc_coef(:) = var_coef(relvar(:), 2.47_r8)
else
prc_coef(:) = 1._r8
end if
do i=1,mgncol
if (qcic(i) >= icsmall) then
! nprc is increase in rain number conc due to autoconversion
! nprc1 is decrease in cloud droplet conc due to autoconversion
! assume exponential sub-grid distribution of qc, resulting in additional
! factor related to qcvar below
! switch for sub-columns, don't include sub-grid qc
prc(i) = prc_coef(i) * &
0.01_r8 * 1350._r8 * qcic(i)**2.47_r8 * (ncic(i)*1.e-6_r8*rho(i))**(-1.1_r8)
nprc(i) = prc(i) * (1._r8/droplet_mass_25um)
nprc1(i) = prc(i)*ncic(i)/qcic(i)
else
prc(i)=0._r8
nprc(i)=0._r8
nprc1(i)=0._r8
end if
enddo
end subroutine kk2000_liq_autoconversion
!========================================================================
subroutine sb2001v2_liq_autoconversion(pgam,qc,nc,qr,rho,relvar,au,nprc,nprc1,mgncol)
!
! ---------------------------------------------------------------------
! AUTO_SB: calculates the evolution of mass- and number mxg-ratio for
! drizzle drops due to autoconversion. The autoconversion rate assumes
! f(x)=A*x**(nu_c)*exp(-Bx) in drop MASS x.
! Code from Hugh Morrison, Sept 2014
! autoconversion
! use simple lookup table of dnu values to get mass spectral shape parameter
! equivalent to the size spectral shape parameter pgam
integer, intent(in) :: mgncol
real(r8), dimension(mgncol), intent (in) :: pgam
real(r8), dimension(mgncol), intent (in) :: qc ! = qc (cld water mixing ratio)
real(r8), dimension(mgncol), intent (in) :: nc ! = nc (cld water number conc /kg)
real(r8), dimension(mgncol), intent (in) :: qr ! = qr (rain water mixing ratio)
real(r8), dimension(mgncol), intent (in) :: rho ! = rho : density profile
real(r8), dimension(mgncol), intent (in) :: relvar
real(r8), dimension(mgncol), intent (out) :: au ! = prc autoconversion rate
real(r8), dimension(mgncol), intent (out) :: nprc1 ! = number tendency
real(r8), dimension(mgncol), intent (out) :: nprc ! = number tendency fixed size for rain
! parameters for droplet mass spectral shape,
!used by Seifert and Beheng (2001)
! warm rain scheme only (iparam = 1)
real(r8), parameter :: dnu(16) = [0._r8,-0.557_r8,-0.430_r8,-0.307_r8, &
-0.186_r8,-0.067_r8,0.050_r8,0.167_r8,0.282_r8,0.397_r8,0.512_r8, &
0.626_r8,0.739_r8,0.853_r8,0.966_r8,0.966_r8]
! parameters for Seifert and Beheng (2001) autoconversion/accretion
real(r8), parameter :: kc = 9.44e9_r8
real(r8), parameter :: kr = 5.78e3_r8
real(r8) :: dum, dum1, nu, pra_coef
integer :: dumi, i
do i=1,mgncol
pra_coef = var_coef(relvar(i), 2.47_r8)
if (qc(i) > qsmall) then
dumi=int(pgam(i))
nu=dnu(dumi)+(dnu(dumi+1)-dnu(dumi))* &
(pgam(i)-dumi)
dum = 1._r8-qc(i)/(qc(i)+qr(i))
dum1 = 600._r8*dum**0.68_r8*(1._r8-dum**0.68_r8)**3
au(i) = kc/(20._r8*2.6e-7_r8)* &
(nu+2._r8)*(nu+4._r8)/(nu+1._r8)**2._r8* &
(rho(i)*qc(i)/1000._r8)**4._r8/(rho(i)*nc(i)/1.e6_r8)**2._r8* &
(1._r8+dum1/(1._r8-dum)**2)*1000._r8 / rho(i)
nprc1(i) = au(i)*2._r8/2.6e-7_r8*1000._r8
nprc(i) = au(i)/droplet_mass_40um
else
au(i) = 0._r8
nprc1(i) = 0._r8
nprc(i)=0._r8
end if
enddo
end subroutine sb2001v2_liq_autoconversion
!========================================================================
!SB2001 Accretion V2
subroutine sb2001v2_accre_cld_water_rain(qc,nc,qr,rho,relvar,pra,npra,mgncol)
!
! ---------------------------------------------------------------------
! ACCR_SB calculates the evolution of mass mxng-ratio due to accretion
! and self collection following Seifert & Beheng (2001).
!
integer, intent(in) :: mgncol
real(r8), dimension(mgncol), intent (in) :: qc ! = qc (cld water mixing ratio)
real(r8), dimension(mgncol), intent (in) :: nc ! = nc (cld water number conc /kg)
real(r8), dimension(mgncol), intent (in) :: qr ! = qr (rain water mixing ratio)
real(r8), dimension(mgncol), intent (in) :: rho ! = rho : density profile
real(r8), dimension(mgncol), intent (in) :: relvar
! Output tendencies
real(r8), dimension(mgncol), intent(out) :: pra ! MMR
real(r8), dimension(mgncol), intent(out) :: npra ! Number
! parameters for Seifert and Beheng (2001) autoconversion/accretion
real(r8), parameter :: kc = 9.44e9_r8
real(r8), parameter :: kr = 5.78e3_r8