forked from NordicESMhub/OsloCTM3
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathscavenging_largescale_uci.f90
2392 lines (2103 loc) · 95.3 KB
/
scavenging_largescale_uci.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
!//=========================================================================
!// Oslo CTM3
!//=========================================================================
!// Ole Amund Sovde, June 2015
!//=========================================================================
!// Large scale scavenging routines.
!//=========================================================================
module scavenging_largescale_uci
!// ----------------------------------------------------------------------
!// MODULE: scavenging_largescale_uci
!// DESCRIPTION: Large scale wet scavenging (Neu ...)
!// Converted from UCI p-scav.f and modified for CTM3.
!//
!// Contains
!// subroutine WASHO
!// subroutine WASH1
!// subroutine WASH2
!// subroutine getHstar
!// subroutine DISGAS2
!// subroutine WASHGAS2
!// subroutine DISGAS
!// subroutine HENRYS
!// subroutine HENRYSallT
!// subroutine RAINGAS
!// subroutine WASHGAS
!// subroutine DIAMEMP
!// subroutine GAMMAX
!// subroutine WETSET_CTM3
!//
!// ----------------------------------------------------------------------
use cmn_precision, only: r8
!// ----------------------------------------------------------------------
implicit none
!// ----------------------------------------------------------------------
!// Parameters
real(r8), parameter :: TICE = 273._r8
real(r8), parameter :: TSOLID = 258._r8 ! Also referred to as TMIX
!// ----------------------------------------------------------------------
character(len=*), parameter :: f90file = 'scavenging_largescale_uci.f90'
!// ----------------------------------------------------------------------
!//---(p-scav.f) p-code 5.6 (Neu, Tang, Prather Aug 2008)
!//Version: qcode_56d; 20090318
!// CTM3: modified VDEP
!// SCAV_FMT is longer
!// TCCNVHENRY is read
!// ----------------------------------------------------------------------
private
public WASHO, WETSET_CTM3, getHstar
!// ----------------------------------------------------------------------
contains
!// ----------------------------------------------------------------------
subroutine WASHO(BTT,BXT,BXX,BYT,BYY,BZT,BZZ,BXY,BXZ,BYZ, &
AIRB,BTEM,DELT,MP)
!// --------------------------------------------------------------------
!// ---p-scav 5.6 (2008 J. Neu, MJP)
!// ---called from pmain to calculate rainout and washout of tracers
!// --------------------------------------------------------------------
use cmn_precision, only: rMom
use cmn_size, only: LPAR, NPAR, IDBLK, JDBLK
use cmn_ctm, only: NSCX
!// --------------------------------------------------------------------
implicit none
!// --------------------------------------------------------------------
!// Input
integer, intent(in) :: MP
real(r8), intent(in), dimension(LPAR,IDBLK,JDBLK) :: AIRB, BTEM
real(r8), intent(in) :: DELT
!// Input/Output
real(r8), intent(inout), dimension(LPAR,NPAR,IDBLK,JDBLK) :: BTT
real(rMom), intent(inout), dimension(LPAR,NPAR,IDBLK,JDBLK) :: &
BXT, BXX, BYT, BYY, BZT, BZZ, BXY, BXZ, BYZ
!// --------------------------------------------------------------------
character(len=*), parameter :: subr = 'WASHO'
!//---------------------------------------------------------------------
!// NSCX = 0 = simpel e-fold vs. ht per tracer
if (NSCX .eq. 0) then
write(6,'(a)') f90file//':'//subr// &
': DO NOT USE THIS ROUTINE in Oslo CTM3'
write(6,'(a)') ' : NSCX should be 1'
stop
call WASH1 (BTT,BXT,BXX,BYT,BYY,BZT,BZZ,BXY,BXZ,BYZ,DELT,MP)
end if
!//NSCX = 1 = new water/ice scavenging (J.Neu 8/2008)
if (NSCX .eq. 1) then
call WASH2 (BTT,BXT,BXX,BYT,BYY,BZT,BZZ,BXY,BXZ,BYZ, &
AIRB,BTEM,DELT,MP)
end if
!// --------------------------------------------------------------------
end subroutine WASHO
!// ----------------------------------------------------------------------
!// ----------------------------------------------------------------------
subroutine WASH1(BTT,BXT,BXX,BYT,BYY,BZT,BZZ,BXY,BXZ,BYZ,DTSCAV,MP)
!// --------------------------------------------------------------------
!//---p-scav 5.6 (2008)
!// --------------------------------------------------------------------
use cmn_precision, only: rMom
use cmn_size, only: LPAR, NPAR, IDBLK, JDBLK, LCONVM, LWEPAR, LPAR
use cmn_ctm, only: MPBLKIB, MPBLKIE, MPBLKJB, MPBLKJE, NTM
use cmn_chem, only: TRWETL
use cmn_met, only: ZOFLE
!// --------------------------------------------------------------------
implicit none
!// --------------------------------------------------------------------
!// Input
integer, intent(in) :: MP
real(r8), intent(in) :: DTSCAV
!// Input/Output
real(r8), intent(inout), dimension(LPAR,NPAR,IDBLK,JDBLK) :: BTT
real(rMom), intent(inout), dimension(LPAR,NPAR,IDBLK,JDBLK) :: &
BXT, BXX, BYT, BYY, BZT, BZZ, BXY, BXZ, BYZ
!// Locals
integer :: I,II,J,JJ,L,N,LTOP,LZ
real(r8) :: FXYW1
!// --------------------------------------------------------------------
!// NSCX = 0 = simpel e-fold vs. ht per tracer
!// --------------------------------------------------------------------
do J = MPBLKJB(MP),MPBLKJE(MP)
JJ = J - MPBLKJB(MP) + 1
do I = MPBLKIB(MP),MPBLKIE(MP)
II = I - MPBLKIB(MP) + 1
!// ---globally uniform "rainout" loss rate (/s) = e-fold vs. height
!// per tracer
!// ---TRWETL(LZ=1:9) stores rate from 1 to 8 km, and then above
!// 8 km (LZ=9)
!// ---applied to both aerosols and gases when NSCX = 0
!// ---top level for rain/wash out scavenging
!// (possibly include LPAUZ?, hence inside IJ)
LTOP = min(LCONVM, LWEPAR, LPAR-1)
do L = 1, LTOP
LZ = max(1,min(9,&
(int(0.5e-3_r8 * (ZOFLE(L,I,J) + ZOFLE(L+1,I,J))))))
do N = 1, NTM
if (TRWETL(LZ,N) .gt. 0._r8) then
FXYW1 = exp(-DTSCAV * TRWETL(LZ,N))
BTT(L,N,II,JJ) = BTT(L,N,II,JJ) * FXYW1
BZT(L,N,II,JJ) = BZT(L,N,II,JJ) * FXYW1
BZZ(L,N,II,JJ) = BZZ(L,N,II,JJ) * FXYW1
BXZ(L,N,II,JJ) = BXZ(L,N,II,JJ) * FXYW1
BYZ(L,N,II,JJ) = BYZ(L,N,II,JJ) * FXYW1
BXT(L,N,II,JJ) = BXT(L,N,II,JJ) * FXYW1
BXX(L,N,II,JJ) = BXX(L,N,II,JJ) * FXYW1
BYT(L,N,II,JJ) = BYT(L,N,II,JJ) * FXYW1
BYY(L,N,II,JJ) = BYY(L,N,II,JJ) * FXYW1
BXY(L,N,II,JJ) = BXY(L,N,II,JJ) * FXYW1
end if
end do
end do
end do
end do
!// --------------------------------------------------------------------
end subroutine WASH1
!// ----------------------------------------------------------------------
!// ----------------------------------------------------------------------
subroutine WASH2(BTT,BXT,BXX,BYT,BYY,BZT,BZZ,BXY,BXZ,BYZ, &
AIRB,BTEM,DTSCAV,MP)
!// --------------------------------------------------------------------
!// ---p-scav 5.6 (2008 J. Neu, MJP, Q. Tang)
!//---called from WASHO (from MAIN) to calc rainout and washout of tracers
!// --------------------------------------------------------------------
use cmn_precision, only: rMom
use cmn_size, only: LPAR, NPAR, IDBLK, JDBLK, LCONVM, LWEPAR, LPAR
use cmn_ctm, only: AREAXY, ETAA, ETAB, &
MPBLKIB, MPBLKIE, MPBLKJB, MPBLKJE, NTM
use cmn_chem, only: TRWETL, TCWETL, TCHENA, TCHENB, TCKAQA, &
TCKAQB, SCV_RETEFF, SCViceFR, SCV_T258, TMASS
use cmn_met, only: ZOFLE, P, PRECLS, CLDLWC, CLDIWC, CLDFR
use utilities, only: CIWMIN
use cmn_oslo, only: TCCNVHENRY, chem_idx
!// --------------------------------------------------------------------
implicit none
!// --------------------------------------------------------------------
!// Input
integer, intent(in) :: MP
real(r8), intent(in) :: DTSCAV
real(r8), intent(in), dimension(LPAR,IDBLK,JDBLK) :: AIRB, BTEM
!// Input/Output
real(r8), intent(inout), dimension(LPAR,NPAR,IDBLK,JDBLK) :: BTT
real(rMom), intent(inout), dimension(LPAR,NPAR,IDBLK,JDBLK) :: &
BXT, BXX, BYT, BYY, BZT, BZZ, BXY, BXZ, BYZ
!// Locals
integer :: I,J,L,N, LTOP, II, JJ
real(r8), dimension(LPAR+1) :: POFLE
real(r8), dimension(LPAR) :: DELZ,POFL
real(r8), dimension(LPAR) :: RLS, RMC, CLWC,CIWC,CFR,TEM
real(r8), dimension(LPAR) :: QM
real(r8), dimension(LPAR,NPAR) :: QTT, QTTNEW
real(r8), dimension(LPAR,NPAR) :: QTTnoSCAV
real(r8) :: GAREA, FXYW1
integer :: LWSHTYP(NPAR),LICETYP(NPAR)
real(r8) :: CLWX,CFXX
real(r8) :: RNEW, RPRECIP, DELTARIMEMASS, DELTARIME
real(r8) :: MASSLOSS
real(r8) :: DOR, DNEW, DEMP, COLEFFSNOW, RHOSNOW
real(r8) :: WEMP, REMP, RRAIN, RWASH
real(r8) :: QTPRECIP, QTCXA, QTAX
real(r8) :: QTDISCF, QTDISRIME, QTDISCXA
real(r8) :: QTDISSTAR
!// --------------------------------------------------------------------
character(len=*), parameter :: subr = 'WASH2'
!//---------------------------------------------------------------------
!// Cloud overlap nomenclature
!// --------------------------------------------------------------------
! TOP=above | 1 = rain under cloud | 2 = rain from clear | (3) = no rain |
!-------------------------------------------------------------------------------
! lyr L | 4=cld core(under rain) | 5=cld(not/rain) | 6=rain thru clr | (7) |
!-------------------------------------------------------------------------------
! BOT=below | 8 = rain from under cloud in lyr L | 9 = rain (ambient)| (10)|
!-------------------------------------------------------------------------------
!
! notes:
! 4 = cloud under rainfall (including cloud above and ambient above)
! 5 = cloud with no rain coming in, can initiate rain.
! 2 & 6 & 9 precip falling thru clear sky (aka ambient)
! in 6, the area and precip shrink (per km)
!
! rain = all forms of precipitation (kg/s) from each sub-box
!
! FR1:FR9 = FRaction of grid box area (0->1)
! PR1:PR9 = PRecipitation thru bottom of region (kg/s)
! QT____ = tracer quantity in sub-box, or precip over time step (kg)
!// Local tracer arrays: QT__ = kg of tracer in region
!// NET = net loss from region (1,2,3) in layer L
!// TOP = amount coming in top region (1,3)
!// PRC = amount in precip regions (0,1,2) 0=1+2
!// RIM = amount picked up riming, all cloud regions
!// EVP = evaporative losses/gains in cloud
!// (1 = all cloud) & ambient (3)
real(r8), dimension(NPAR) :: QTTOP1, QTTOP2, QTTOP3
real(r8), dimension(NPAR) :: QTTOP4, QTTOP6
real(r8), dimension(NPAR) :: QTBOT7, QTBOT8, QTBOT9
real(r8), dimension(NPAR) :: QTT4, QTT5, QTT6
real(r8), dimension(NPAR) :: QTNET4, QTNET5, QTNET6
real(r8), dimension(NPAR) :: QTPRC0, QTPRC4, QTPRC5
real(r8), dimension(NPAR) :: QTWSH4, QTWSH6, QTRIM4, QTEVP4, QTEVP6
real(r8), dimension(NPAR) :: QTEVP4X, QTEVP4Y, QTEVP6X, QTEVP6Y
real(r8), dimension(LPAR) :: CEVAP6, DRLS
!// Cloud fractions
real(r8) :: FR1,FR2,FR3, FR4,FR5,FR6, FR7,FR8,FR9
!// Precipitation rates
real(r8) :: PR1,PR2,PR3, PR4,PR5,PR6, PR7,PR8,PR9
real(r8) :: W1TO4,W1TO6, W2TO4,W2TO6, W3TO4,W3TO6
real(r8) :: D1,D2,D3, D4,D5,D6, D7,D8,D9
real(r8) :: PR4T, D4T, PR6T, D6T, FR6T
real(r8) :: PR6TMP
!// For Henry expressions
real(r8) :: HSTAR
integer, parameter :: HCALLER = 1 !// caller flag for getHstar
!// Parameters
real(r8), parameter :: &
TMIX = TSOLID, &
TFROZ = 240._r8, &
EVAPRATE = 0.25_r8, & !fraction/km
CFMIN = 0.025_r8, &
CFCOREMIN = 0.1_r8, &
DMIN = 0.1_r8, & !mm
VOLPOW = 1._r8/3._r8, &
RHORAIN = 1000._r8, & !kg/m3
RHOSNOWFIX = 100._r8, & !kg/m3
COLEFFRAIN = 0.7_r8, &
COLEFFAER = 0.05_r8, &
QTTMIN = 1.e-20_r8 ! Min QTT to scavenge (tiny QTT may -> NAN)
!// --------------------------------------------------------------------
!// --------------------------------------------------------------------
!// NSCX = 1 = new water/ice scavenging (J.Neu 8/2008)
!// --------------------------------------------------------------------
!// top level for rain/wash out scavenging (possibly include LPAUZ(I,J)?)
LTOP = min(LCONVM, LWEPAR, LPAR-1) ! possibly min (..., LPAUZB(I,J)) ?
QTT(:,:) = 0._r8
QTTNEW(:,:) = 0._r8
QTTnoSCAV(:,:) = 0._r8
!//---set up flags for each species:
do N = 1, NTM
LWSHTYP(N) = 0
LICETYP(N) = 0
!//-----check if soluble
if(TCHENA(N) .gt. 0._r8) then
!// check whether mass-limited or henry's law
if(TCKAQB(N) .gt. 0._r8) then
LWSHTYP(N) = 1
else
LWSHTYP(N) = 2
end if
!// check whether soluble in ice
if(SCViceFR(N) .gt. 0._r8) then
!// Could eventually check on SCV_RETEFF(N) also, but not
!// until it is listed in input file.
LICETYP(N) = 1
else
LICETYP(N) = 2
end if
end if
!// CTM3: possibility of overriding large scale scavenging, so that
!// only convective scavenging is calculated
if (TCCNVHENRY(N) .eq. 5) then
LWSHTYP(N) = 0
LICETYP(N) = 0
else if (TCCNVHENRY(N) .eq. 6) then
LWSHTYP(N) = 3 !// Skip removal on rain
LICETYP(N) = 1 !// Include removal on ice
else if (TCCNVHENRY(N) .eq. 7) then
LWSHTYP(N) = 3 !// Skip removal on rain
LICETYP(N) = 1 !// Include removal on ice
end if
end do
!// Loop through each column
!// --------------------------------------------------------------------
do J = MPBLKJB(MP),MPBLKJE(MP)
JJ = J - MPBLKJB(MP) + 1
do I = MPBLKIB(MP),MPBLKIE(MP)
II = I - MPBLKIB(MP) + 1
!// ------------------------------------------------------------------
!// set ups from (I,J) grid:
RLS(:) = 0._r8
RMC(:) = 0._r8
CFXX = 0._r8
GAREA = AREAXY(I,J)
do L = 1, LTOP
RLS(L) = PRECLS(I,J,L) ! kg/s thru bottom of layer L
!RMC(L) = PRECCNV(I,J,L) !// Convective rain not used
CLWC(L) = CLDLWC(I,J,L) ! kg/kg in box; LWC/CF=in-cloud
CIWC(L) = CLDIWC(I,J,L)
CFR(L) = CLDFR(I,J,L)
TEM(L) = BTEM(L,II,JJ)
QM(L) = AIRB(L,II,JJ)
end do
RLS(LTOP+1) = 0._r8
!// approx height (m) of layer edges (ignore Q)
POFLE(1) = P(I,J)
do L = 2, LPAR+1
POFLE(L) = ETAA(L) + ETAB(L) * POFLE(1)
!// R = 287.*(1-Q) + 461.5*Q -- assume 0.5% bl w.v.==> R = 288.
!// delta-z (m) = dln(P) * R * T / g where R/g = 288/9.81 = 29.36
DELZ(L-1)= ZOFLE(L,I,J) - ZOFLE(L-1,I,J)
end do
do L = 1, LPAR
POFL(L) = 0.5_r8 * (POFLE(L) + POFLE(L+1))
!fract loss of amb.precip.
CEVAP6(L) = EVAPRATE * DELZ(L) * 1.e-3_r8
end do
do N = 1, NTM
!// check if soluble
if(LWSHTYP(N) .gt. 0) then
do L = 1,LPAR
!// Allow scavening of a fraction of grid box
if (TCWETL(N) .eq. 1._r8) then
QTTnoSCAV(L,N) = 0._r8
QTT(L,N) = BTT(L,N,II,JJ)
QTTNEW(L,N) = BTT(L,N,II,JJ)
else
!// Keep out fraction not subject for scavenging
QTTnoSCAV(L,N) = BTT(L,N,II,JJ) * (1._r8 - TCWETL(N))
QTT(L,N) = BTT(L,N,II,JJ) - QTTnoSCAV(L,N)
QTTNEW(L,N) = BTT(L,N,II,JJ) - QTTnoSCAV(L,N)
end if
end do
end if
end do !// do N = 1, NTM
!// initialize input to zero above top layer (LTOP)
QTBOT7(:) = 0._r8
QTBOT8(:) = 0._r8
QTBOT9(:) = 0._r8
D7 = 0._r8
D8 = 0._r8
D9 = 0._r8
PR7 = 0._r8
PR8 = 0._r8
PR9 = 0._r8
FR7 = 0._r8
FR8 = 0._r8
FR9 = 0._r8
W1TO4 = 0._r8
W1TO6 = 0._r8
W2TO4 = 0._r8
W2TO6 = 0._r8
W3TO4 = 0._r8
W3TO6 = 0._r8
do L = LTOP, 1, -1
QTTOP1(:) = QTBOT7(:)
QTTOP2(:) = QTBOT8(:)
QTTOP3(:) = QTBOT9(:)
D1 = D7
D2 = D8
D3 = D9
PR1 = PR7
PR2 = PR8
PR3 = PR9
FR1 = FR7
FR2 = FR8
FR3 = FR9
PR4 = 0._r8
FR4 = 0._r8
PR5 = 0._r8
FR5 = 0._r8
PR6 = 0._r8
FR6 = 0._r8
D4 = 0._r8
D5 = 0._r8
D6 = 0._r8
QTTOP4(:) = 0._r8
QTTOP6(:) = 0._r8
QTPRC0(:) = 0._r8
QTPRC4(:) = 0._r8
QTPRC5(:) = 0._r8
QTRIM4(:) = 0._r8
QTWSH4(:) = 0._r8
QTEVP4(:) = 0._r8
QTEVP4X(:) = 0._r8
QTEVP4Y(:) = 0._r8
QTWSH6(:) = 0._r8
QTEVP6(:) = 0._r8
QTEVP6X(:) = 0._r8
QTEVP6Y(:) = 0._r8
QTDISCF = 0._r8
QTDISRIME = 0._r8
QTDISCXA = 0._r8
RPRECIP = 0._r8
DELTARIMEMASS = 0._r8
DELTARIME = 0._r8
DOR = 0._r8
DNEW = 0._r8
if (RLS(L).le.0._r8 .and. RLS(L+1).le.0._r8) then
!// No rain out at bottom and no rain in at top
PR7 = 0._r8
PR8 = 0._r8
PR9 = 0._r8
FR7 = 0._r8
FR8 = 0._r8
FR9 = 0._r8
D7 = 0._r8
D8 = 0._r8
D9 = 0._r8
QTBOT7(:) = 0._r8
QTBOT8(:) = 0._r8
QTBOT9(:) = 0._r8
else if (RLS(L).le.0._r8 .and. RLS(L+1).gt.0._r8) then
!// No rain out at bottom, but rain in at top
PR7 = 0._r8
PR8 = 0._r8
PR9 = 0._r8
FR7 = 0._r8
FR8 = 0._r8
FR9 = 0._r8
D7 = 0._r8
D8 = 0._r8
D9 = 0._r8
QTBOT7(:) = 0._r8
QTBOT8(:) = 0._r8
QTBOT9(:) = 0._r8
do N = 1, NTM
if (LWSHTYP(N) .gt. 0) then
if (QTTOP1(N) + QTTOP2(N) + QTTOP3(N) .lt. 0._r8) then
QTTNEW(L,N) = QTT(L,N)
else
QTTNEW(L,N) = QTT(L,N) + QTTOP1(N) + QTTOP2(N) + QTTOP3(N)
end if
end if
end do
else
!// Rain out at bottom and rain in at top
! Determine if CFMIN is needed
if (CFR(L).gt.0._r8 .and. CFR(L).lt.CFMIN) then
CFR(L) = CFMIN
CLWC(L) = CLWC(L) / CFR(L)
CIWC(L) = CIWC(L) / CFR(L)
call CIWMIN(TEM(L),CLWC(L),CIWC(L))
CLWC(L) = CLWC(L)*CFR(L)
CIWC(L) = CIWC(L)*CFR(L)
end if
! Determine if CFCOREMIN is needed
PR6TMP = (PR1 + PR2 + PR3) * (1._r8 - CEVAP6(L))
if (RLS(L) .ge. PR6TMP) then
CFXX = max(CFR(L), CFMIN)
CLWC(L) = CLWC(L) / CFXX
CIWC(L) = CIWC(L) / CFXX
call CIWMIN(TEM(L),CLWC(L),CIWC(L))
CLWC(L) = CLWC(L) * CFXX
CIWC(L) = CIWC(L) * CFXX
else
CFXX = CFR(L)
end if
! if (DRLS(L) .gt. 0._r8) then
! CFXX = max(CFR(L),CFCOREMIN)
! CLWC(L) = CLWC(L)/CFXX
! CIWC(L) = CIWC(L)/CFXX
! call CIWMIN(TEM(L),CLWC(L),CIWC(L))
! CLWC(L) = CLWC(L)*CFXX
! CIWC(L) = CIWC(L)*CFXX
! else
! CFXX = CFR(L)
! endif
!// Determine the weighting coeffients and fractions of
!// sub-boxes and D4,D6
if (FR1+FR2+FR3 .gt. 0._r8) then
if (FR1+FR2+FR3 .gt. CFXX) then
FR4 = CFXX
FR6 = FR1 + FR2 + FR3 - CFXX
FR5 = 0._r8
if(FR1.eq.0._r8 .and. FR2.ne.0._r8 .and. FR3.ne.0._r8) then
W1TO4 = 0._r8
W1TO6 = 0._r8
if (FR4 .le. FR2) then
W2TO4 = FR4 / FR2
W2TO6 = 1._r8 - W2TO4
W3TO4 = 0._r8
W3TO6 = 1._r8
else
W2TO4 = 1._r8
W2TO6 = 0._r8
W3TO6 = FR6 / FR3
W3TO4 = 1._r8 - W3TO6
end if
else if(FR2.eq.0._r8.and.FR1.ne.0._r8.and.FR3.ne.0._r8) then
W2TO4 = 0._r8
W2TO6 = 0._r8
if (FR4 .le. FR1) then
W1TO4 = FR4 / FR1
W1TO6 = 1._r8 - W1TO4
W3TO4 = 0._r8
W3TO6 = 1._r8
else
W1TO4 = 1._r8
W1TO6 = 0._r8
W3TO6 = FR6 / FR3
W3TO4 = 1._r8 - W3TO6
end if
else if(FR3.eq.0._r8.and.FR1.ne.0._r8.and.FR2.ne.0._r8) then
W3TO4 = 0._r8
W3TO6 = 0._r8
if (FR4 .le. FR1) then
W1TO4 = FR4 / FR1
W1TO6 = 1._r8 - W1TO4
W2TO4 = 0._r8
W2TO6 = 1._r8
else
W1TO4 = 1._r8
W1TO6 = 0._r8
W2TO6 = FR6 / FR2
W2TO4 = 1._r8 - W2TO6
end if
else if (FR1.eq.0._r8 .and. FR2.eq.0._r8) then
W1TO4 = 0._r8
W1TO6 = 0._r8
W2TO4 = 0._r8
W2TO6 = 0._r8
W3TO4 = FR4 / FR3
W3TO6 = 1._r8 - W3TO4
else if (FR1.eq.0._r8 .and. FR3.eq.0._r8) then
W1TO4 = 0._r8
W1TO6 = 0._r8
W2TO4 = FR4 / FR2
W2TO6 = 1._r8 - W2TO4
W3TO4 = 0._r8
W3TO6 = 0._r8
else if (FR2.eq.0._r8 .and. FR3.eq.0._r8) then
W1TO4 = FR4 / FR1
W1TO6 = 1._r8 - W1TO4
W2TO4 = 0._r8
W2TO6 = 0._r8
W3TO4 = 0._r8
W3TO6 = 0._r8
else if (FR4 .le. FR1) then
W1TO4 = FR4 / FR1
W1TO6 = 1._r8 - W1TO4
W2TO4 = 0._r8
W2TO6 = 1._r8
W3TO4 = 0._r8
W3TO6 = 1._r8
else if (FR4 .le. FR1+FR2) then
W1TO4 = 1._r8
W1TO6 = 0._r8
W2TO4 = (FR4 - FR1) / FR2
W2TO6 = 1._r8 - W2TO4
W3TO4 = 0._r8
W3TO6 = 1._r8
else
W1TO4 = 1._r8
W1TO6 = 0._r8
W2TO4 = 1._r8
W2TO6 = 0._r8
W3TO6 = FR6 / FR3
W3TO4 = 1._r8 - W3TO6
end if
else if ( (FR1 + FR2 + FR3) .eq. CFXX) then
FR4 = FR1 + FR2 + FR3
FR5 = 0._r8 !aviod CFXX-(FR1+FR2+FR3) when they are equal
FR6 = 0._r8
W1TO4 = 1._r8
W1TO6 = 0._r8
W2TO4 = 1._r8
W2TO6 = 0._r8
W3TO4 = 1._r8
W3TO6 = 0._r8
else
FR4 = FR1 + FR2 + FR3
FR5 = CFXX - (FR1 + FR2 + FR3)
FR6 = 0._r8
W1TO4 = 1._r8
W1TO6 = 0._r8
W2TO4 = 1._r8
W2TO6 = 0._r8
W3TO4 = 1._r8
W3TO6 = 0._r8
end if
PR4 = W1TO4*PR1 + W2TO4*PR2 + W3TO4*PR3
PR6 = W1TO6*PR1 + W2TO6*PR2 + W3TO6*PR3
if (PR4 .gt. 0._r8) then
D4 = (W1TO4*PR1*D1 + W2TO4*PR2*D2 + W3TO4*PR3*D3) / PR4
else
D4 = 0._r8
end if
if (PR6 .gt. 0._r8) then
D6 = (W1TO6*PR1*D1 + W2TO6*PR2*D2 + W3TO6*PR3*D3) / PR6
else
D6 = 0._r8
end if
else
FR4 = 0._r8
FR5 = CFXX
FR6 = 0._r8
D4 = 0._r8
D6 = 0._r8
PR4 = 0._r8
PR6 = 0._r8
W1TO4 = 0._r8
W1TO6 = 0._r8
W2TO4 = 0._r8
W2TO6 = 0._r8
W3TO4 = 0._r8
W3TO6 = 0._r8
end if
do N = 1, NTM
QTTOP4(N) = W1TO4*QTTOP1(N) + W2TO4*QTTOP2(N) + W3TO4*QTTOP3(N)
QTTOP6(N) = W1TO6*QTTOP1(N) + W2TO6*QTTOP2(N) + W3TO6*QTTOP3(N)
end do
!// Get TOP values
PR4T = PR4
D4T = D4
PR6T = PR6
D6T = D6
FR6T = FR6
if (CFXX .gt. 0._r8) then
CLWX = CLWC(L) + CIWC(L)
else
CLWX = 0._r8
end if
!-----Evaporate ambient precip and decrease area-------------------------
!-----If ice, diam=diam falling from above If rain, diam=4mm (not used)
!-----Evaporate tracer contained in evaporated precip
!-----Can't evaporate more than we start with-----------------------------
!-----Don't do washout until we adjust ambient precip to match Rbot if needed
!------(after RNEW if statements)
FR6 = max(0._r8, FR6T - CEVAP6(L)*FR6T)
!// Do not use: PR6 = max(0._r8,PR6T*(1._r8-CEVAP6(L))**2) !kg/s
PR6 = max(0._r8, PR6T - CEVAP6(L)*PR6T) !kg/s
if (FR6 .gt. 0._r8) then
if (TEM(L) .lt. TICE) then
D6 = D6T !mm
else
D6 = 4.0_r8 !mm - not necessary
end if
else
D6 = 0._r8
end if
if (PR6T .gt. 0._r8) then
do N = 1, NTM
QTEVP6X(N) = (PR6T-PR6)/PR6T*QTTOP6(N)
end do
end if
RNEW = RLS(L) - PR4T - PR6 !kg/s
!// if RNEW>0, there is growth and/or new precip formation
!if(RNEW.gt.0._r8 .and. CFXX.gt.0._r8) then
if (RNEW .gt. 0._r8) then
!// Min cloudwater requirement for cloud with new precip
!// CWMIN is only needed for new precip formation - do not
!// need for RNEW<0
!// CLWX = max((CLWC(L)+CIWC(L)),CWMIN)
!-------------------------ICE-----------------------------------------
!-----For ice and mixed phase, grow precip in old cloud by riming
!-----Use only portion of cloudwater in old cloud fraction
!-----and rain above old cloud fraction
!-----COLEFF from Lohmann and Roeckner (1996), Loss rate from Rotstayn (1997)
if (TEM(L) .lt. TICE) then
COLEFFSNOW = exp(2.5e-2_r8 * (TEM(L) - 273._r8))
if(TEM(L) .le. TFROZ) then
RHOSNOW = RHOSNOWFIX
else
RHOSNOW = 0.303_r8 * (TEM(L) - TFROZ) * RHOSNOWFIX
end if
if(FR4 .gt. 0._r8) then
if(D4T .gt. 0._r8) then
DELTARIMEMASS = &
(CLWX * QM(L)) * (FR4 / CFXX) &
* (1._r8 - exp((-COLEFFSNOW / (D4T * 1.e-3_r8)) &
* (PR4T / (2._r8 * RHOSNOW * GAREA * FR4)) * DTSCAV))
else
DELTARIMEMASS = 0._r8
end if
else
DELTARIMEMASS = 0._r8
end if
!// Increase in precip rate due to riming (kg/s):
!// Limit to total increase in R in cloud
if (FR4 .gt. 0._r8) then
DELTARIME = min(RNEW, DELTARIMEMASS / DTSCAV)
else
DELTARIME = 0._r8
end if
!// Find diameter of rimed precip, must be at least .1mm
if (PR4T .gt. 0._r8) then
DOR = max(DMIN,(((PR4T + DELTARIME) / PR4T)**VOLPOW) * D4T)
else
DOR = 0._r8
end if
!-----If there is some in-cloud precip left, we have new precip formation
!-----Will be spread over whole cloud fraction
!------Calculate precip rate in old and new cloud fractions
RPRECIP = RNEW - DELTARIME !kg/s
!// Calculate precip rate in old and new cloud fractions
PR4 = PR4T + DELTARIME + RPRECIP * FR4 / CFXX !kg/s
PR5 = RPRECIP * FR5 / CFXX !kg/s
!-----Find diameter of new precip from empirical relation using Rprecip
!----in given area of box- use density of water, not snow, to convert kg/s
!----to mm/s -> as given in Field and Heymsfield, J. Atm. Sci., 60, 544-560,
!---- 2003, doi:10.1175/1520-0469(2003)060%3C0544:AASOIC%3E2.0.CO;2
!-----Also calculate diameter of mixed precip,D4, from empirical relation
!-----using total R in FR4 - this will give larger particles than averaging DOR and
!-----DNEW in the next level
!-----DNEW and D4 must be at least .1mm
if (RPRECIP .gt. 0._r8) then
!// New rain is formed
WEMP = CLWX * QM(L) / (GAREA * CFXX * DELZ(L))
!mm/s local
REMP = RPRECIP / (GAREA * CFXX * RHORAIN * 1.e-3_r8)
call DIAMEMP(WEMP,REMP,DNEW)
DNEW = max(DMIN, DNEW)
!// in sub-box 5
if (FR5 .gt. 0._r8) then
D5 = DNEW
else
D5 = 0._r8
end if
else
D5 = 0._r8
end if
!// in sub-box 4
if (FR4 .gt. 0._r8) then
WEMP = CLWX * QM(L) / (GAREA * CFXX * DELZ(L)) !kg/m3
!mm/s local
REMP = PR4 / (GAREA * FR4 * RHORAIN * 1.e-3_r8)
call DIAMEMP(WEMP,REMP,DEMP)
D4 = ((PR4T + DELTARIME) / PR4) * DOR &
+ (RPRECIP * FR4 / CFXX / PR4) * DNEW
D4 = max(DEMP, D4)
D4 = max(DMIN, D4)
else
D4 = 0._r8
end if
!// ICE SCAVENGING
!// --------------------------------------------------------
!------For ice, rainout only hno3/aerosols using new precip
!------Tracer dissolved given by Kaercher and Voigt (2006) for T<258K
!-----For T>258K, use Henry's Law with Retention coefficient
!------Rain out in whole CF
if (RPRECIP .gt. 0._r8) then
do N = 1, NTM
if (QTT(L,N) .gt. QTTMIN) then
if (LICETYP(N) .eq. 1) then
RRAIN = RPRECIP !kg/s
!// Calculate effective Henry expression
call getHstar(chem_idx(N), TCHENA(N), &
TCHENB(N), TCKAQA(N), TEM(L), HCALLER, HSTAR)
!// Calculate dissolved tracer
call DISGAS2(chem_idx(N),CLWX, CFXX, TMASS(N), HSTAR, &
TEM(L), POFL(L), QM(L), &
QTT(L,N)*CFXX, SCV_RETEFF(N), SCV_T258(N), &
SCViceFR(N), QTDISCF)
!// Rain out
call RAINGAS(RRAIN, DTSCAV, CLWX, CFXX, &
QM(L), QTT(L,N), QTDISCF, QTPRC0(N))
!// Prec.rate in old cloud
QTPRC4(N) = FR4 * QTPRC0(N) / CFXX
!// Prec.rate of new prec.
QTPRC5(N) = FR5 * QTPRC0(N) / CFXX
else if(LICETYP(N) .eq. 2) then
QTPRC4(N) = 0._r8 !// No tracer in old ice prec.
QTPRC5(N) = 0._r8 !// No tracer in new ice prec.
end if
end if
end do
end if !// if (RPRECIP .gt. 0._r8) then
!// For ice, accretion removal for hno3 and aerosols is
!// propotional to riming, no accretion removal for gases
!// remove only in mixed portion of cloud
!// Limit DELTARIMEMASS to RNEW*DTSCAV for ice - evaporation
!// of rimed ice to match
!// RNEW precip rate would result in HNO3 escaping from
!// ice (no trapping)
if (DELTARIME .gt. 0._r8) then
do N = 1, NTM
if (QTT(L,N) .gt. QTTMIN) then
if (LICETYP(N) .eq. 1) then
!!!! this looks odd below, discontinuous at TFROZ
if(TEM(L) .le. TFROZ) then
RHOSNOW = RHOSNOWFIX
else
RHOSNOW = 0.303_r8 * (TEM(L) - TFROZ) * RHOSNOWFIX
end if
!// Calculate effective Henry expression
call getHstar(chem_idx(N), TCHENA(N), &
TCHENB(N), TCKAQA(N), TEM(L), HCALLER, HSTAR)
!// This process is riming, but since this is ice,
!// we apply SCViceFR, so that not all tracer in
!// cloud may rime.
call DISGAS2(chem_idx(N), CLWX * FR4 / CFXX, FR4, TMASS(N), &
HSTAR, TEM(L), POFL(L), &
QM(L), QTT(L,N)*FR4 - QTPRC4(N), &
SCV_RETEFF(N), SCV_T258(N), SCViceFR(N), &
QTDISRIME)
QTDISSTAR = QTDISRIME
QTRIM4(N) = (QTT(L,N)*FR4 - QTPRC4(N)) &
* (1._r8 - exp(-COLEFFSNOW / (D4T * 1.e-3_r8) &
* PR4T / (2._r8 * RHOSNOW * GAREA * FR4) &
* QTDISSTAR / (QTT(L,N)*FR4 - QTPRC4(N))*DTSCAV))
QTRIM4(N) = min(QTRIM4(N), &
DELTARIME * DTSCAV / &
(CLWX * QM(L) * FR4 / CFXX) * QTDISSTAR)
else if (LICETYP(N) .eq. 2) then
QTRIM4(N) = 0._r8 !// No tracer to be removed by ice
end if
end if
end do
end if !// if (DELTARIME .gt. 0._r8) then
!// For ice, no washout in interstitial cloud air
QTWSH4(:) = 0._r8 !! not necessary?
QTEVP4(:) = 0._r8 !! not necessary?
else !// if(TEM(L) .lt. TICE) then
!// RAIN
!// --------------------------------------------------------
!// For rain, accretion increases rain rate but diameter
!// remains constant.
!// Diameter is 4mm (not used)
if (FR4 .gt. 0._r8) then
DELTARIMEMASS = &
CLWX * QM(L) * FR4 / CFXX &
* (1._r8 - exp(-0.24_r8 * COLEFFRAIN &
* ((PR4T / (GAREA * FR4))**0.75_r8) * DTSCAV))
else
DELTARIMEMASS = 0._r8
end if
!// Increase in precip rate due to riming (kg/s):
!// Limit to total increase in R in cloud
if (FR4 .gt. 0._r8) then
DELTARIME = min(RNEW, DELTARIMEMASS / DTSCAV)
else
DELTARIME = 0._r8
end if
!// If there is some in-cloud precip left, we have new
!// precip formation
RPRECIP = RNEW - DELTARIME
PR4 = PR4T + DELTARIME + RPRECIP * FR4 / CFXX !kg/s
PR5 = RPRECIP * FR5 / CFXX !kg/s
if (FR4 .gt. 0._r8) then
D4 = 4._r8
else
D4 = 0._r8
end if
if (FR5 .gt. 0._r8) then
D5 = 4._r8
else
D5 = 0._r8
end if
!// RAIN SCAVENGING
!// --------------------------------------------------------
!// For rain, rainout both hno3/aerosols and gases using
!// new precip
if (RPRECIP .gt. 0._r8) then
!// RPRECIP is new precipitation (RNEW-DELTARIME)
do N = 1, NTM
if (LWSHTYP(N) .eq. 3) then
!// No new precip for type 3
QTPRC4(N) = 0._r8
QTPRC5(N) = 0._r8
else if (TCHENA(N) .gt. 0._r8 .and. &
QTT(L,N).gt.QTTMIN) then
!// There is tracer which can be dissolved
RRAIN = RPRECIP !kg/s
!// Calculate effective Henry expression
call getHstar(chem_idx(N), TCHENA(N), &
TCHENB(N), TCKAQA(N), TEM(L), HCALLER, HSTAR)
!// T>273: Instead of SCViceFR, send in 1._r8
call DISGAS2(chem_idx(N), CLWX, CFXX, TMASS(N), HSTAR, &
TEM(L), POFL(L), QM(L), &
QTT(L,N)*CFXX, SCV_RETEFF(N), SCV_T258(N), &
1._r8, QTDISCF)
call RAINGAS(RRAIN, DTSCAV, CLWX, CFXX, &
QM(L), QTT(L,N), QTDISCF, QTPRC0(N))