-
Notifications
You must be signed in to change notification settings - Fork 0
/
aero_model.F90
2392 lines (2007 loc) · 92.6 KB
/
aero_model.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
!===============================================================================
! Modal Aerosol Model
!===============================================================================
module aero_model
use shr_kind_mod, only: r8 => shr_kind_r8
use constituents, only: pcnst, cnst_name, cnst_get_ind
use ppgrid, only: pcols, pver, pverp
use cam_abortutils, only: endrun
use cam_logfile, only: iulog
use perf_mod, only: t_startf, t_stopf
use camsrfexch, only: cam_in_t, cam_out_t
use aerodep_flx, only: aerodep_flx_prescribed
use physics_types, only: physics_state, physics_ptend, physics_ptend_init
use physics_buffer, only: physics_buffer_desc
use physics_buffer, only: pbuf_get_field, pbuf_get_index, pbuf_set_field
use physconst, only: gravit, rair, rhoh2o
use spmd_utils, only: masterproc
use cam_history, only: outfld, fieldname_len
use chem_mods, only: gas_pcnst, adv_mass
use mo_tracname, only: solsym
use modal_aero_data,only: cnst_name_cw
use modal_aero_data,only: ntot_amode, modename_amode
implicit none
private
public :: aero_model_readnl
public :: aero_model_register
public :: aero_model_init
public :: aero_model_gasaerexch ! create, grow, change, and shrink aerosols.
public :: aero_model_drydep ! aerosol dry deposition and sediment
public :: aero_model_wetdep ! aerosol wet removal
public :: aero_model_emissions ! aerosol emissions
public :: aero_model_surfarea ! aerosol surface area for chemistry
! Misc private data
! number of modes
integer :: nmodes
integer :: pblh_idx = 0
integer :: dgnum_idx = 0
integer :: dgnumwet_idx = 0
integer :: rate1_cw2pr_st_idx = 0
integer :: wetdens_ap_idx = 0
integer :: qaerwat_idx = 0
integer :: fracis_idx = 0
integer :: prain_idx = 0
integer :: nevapr_idx = 0
integer :: rprddp_idx = 0
integer :: rprdsh_idx = 0
! variables for table lookup of aerosol impaction/interception scavenging rates
integer, parameter :: nimptblgrow_mind=-7, nimptblgrow_maxd=12
real(r8) :: dlndg_nimptblgrow
real(r8) :: scavimptblnum(nimptblgrow_mind:nimptblgrow_maxd, ntot_amode)
real(r8) :: scavimptblvol(nimptblgrow_mind:nimptblgrow_maxd, ntot_amode)
! for aero_model_surfarea called from mo_usrrxt
integer :: aitken_idx = -1
integer, dimension(ntot_amode) :: num_idx = -1
integer :: index_tot_mass(ntot_amode,10) = -1
integer :: index_chm_mass(ntot_amode,10) = -1
integer :: ndx_h2so4
character(len=fieldname_len) :: dgnum_name(ntot_amode)
! Namelist variables
character(len=16) :: wetdep_list(pcnst) = ' '
character(len=16) :: drydep_list(pcnst) = ' '
real(r8) :: sol_facti_cloud_borne = 1._r8
real(r8) :: sol_factb_interstitial = 0.1_r8
real(r8) :: sol_factic_interstitial = 0.4_r8
integer :: ndrydep = 0
integer,allocatable :: drydep_indices(:)
integer :: nwetdep = 0
integer,allocatable :: wetdep_indices(:)
logical :: drydep_lq(pcnst)
logical :: wetdep_lq(pcnst)
contains
!=============================================================================
! reads aerosol namelist options
!=============================================================================
subroutine aero_model_readnl(nlfile)
use namelist_utils, only: find_group_name
use units, only: getunit, freeunit
use mpishorthand
character(len=*), intent(in) :: nlfile ! filepath for file containing namelist input
! Local variables
integer :: unitn, ierr
character(len=*), parameter :: subname = 'aero_model_readnl'
! Namelist variables
character(len=16) :: aer_wetdep_list(pcnst) = ' '
character(len=16) :: aer_drydep_list(pcnst) = ' '
namelist /aerosol_nl/ aer_wetdep_list, aer_drydep_list, sol_facti_cloud_borne, &
sol_factb_interstitial, sol_factic_interstitial
!-----------------------------------------------------------------------------
! Read namelist
if (masterproc) then
unitn = getunit()
open( unitn, file=trim(nlfile), status='old' )
call find_group_name(unitn, 'aerosol_nl', status=ierr)
if (ierr == 0) then
read(unitn, aerosol_nl, iostat=ierr)
if (ierr /= 0) then
call endrun(subname // ':: ERROR reading namelist')
end if
end if
close(unitn)
call freeunit(unitn)
end if
#ifdef SPMD
! Broadcast namelist variables
call mpibcast(aer_wetdep_list, len(aer_wetdep_list(1))*pcnst, mpichar, 0, mpicom)
call mpibcast(aer_drydep_list, len(aer_drydep_list(1))*pcnst, mpichar, 0, mpicom)
call mpibcast(sol_facti_cloud_borne, 1, mpir8, 0, mpicom)
call mpibcast(sol_factb_interstitial, 1, mpir8, 0, mpicom)
call mpibcast(sol_factic_interstitial, 1, mpir8, 0, mpicom)
#endif
wetdep_list = aer_wetdep_list
drydep_list = aer_drydep_list
end subroutine aero_model_readnl
!=============================================================================
!=============================================================================
subroutine aero_model_register
use modal_aero_initialize_data, only : modal_aero_register
call modal_aero_register()
end subroutine aero_model_register
!=============================================================================
!=============================================================================
subroutine aero_model_init( pbuf2d )
use mo_chem_utls, only: get_inv_ndx
use cam_history, only: addfld, add_default, phys_decomp
use phys_control, only: phys_getopts
use mo_chem_utls, only: get_rxt_ndx, get_spc_ndx
use modal_aero_data, only: cnst_name_cw
use modal_aero_initialize_data, only: modal_aero_initialize
use rad_constituents, only: rad_cnst_get_info
use dust_model, only: dust_init, dust_names, dust_active, dust_nbin, dust_nnum
use seasalt_model, only: seasalt_init, seasalt_names, seasalt_active,seasalt_nbin
use drydep_mod, only: inidrydep
use wetdep, only: wetdep_init
! args
type(physics_buffer_desc), pointer :: pbuf2d(:,:)
! local vars
character(len=*), parameter :: subrname = 'aero_model_init'
integer :: m, n, id
character(len=20) :: dummy
logical :: history_aerosol ! Output MAM or SECT aerosol tendencies
integer :: l
character(len=6) :: test_name
character(len=64) :: errmes
character(len=2) :: unit_basename ! Units 'kg' or '1'
dgnum_idx = pbuf_get_index('DGNUM')
dgnumwet_idx = pbuf_get_index('DGNUMWET')
call phys_getopts( history_aerosol_out=history_aerosol )
call rad_cnst_get_info(0, nmodes=nmodes)
call modal_aero_initialize(pbuf2d)
call modal_aero_bcscavcoef_init()
call dust_init()
call seasalt_init()
call wetdep_init()
fracis_idx = pbuf_get_index('FRACIS')
prain_idx = pbuf_get_index('PRAIN')
nevapr_idx = pbuf_get_index('NEVAPR')
rprddp_idx = pbuf_get_index('RPRDDP')
rprdsh_idx = pbuf_get_index('RPRDSH')
nwetdep = 0
ndrydep = 0
count_species: do m = 1,pcnst
if ( len_trim(wetdep_list(m)) /= 0 ) then
nwetdep = nwetdep+1
endif
if ( len_trim(drydep_list(m)) /= 0 ) then
ndrydep = ndrydep+1
endif
enddo count_species
if (nwetdep>0) &
allocate(wetdep_indices(nwetdep))
if (ndrydep>0) &
allocate(drydep_indices(ndrydep))
do m = 1,ndrydep
call cnst_get_ind ( drydep_list(m), id, abort=.false. )
if (id>0) then
drydep_indices(m) = id
else
call endrun(subrname//': invalid drydep species: '//trim(drydep_list(m)) )
endif
if (masterproc) then
write(iulog,*) subrname//': '//drydep_list(m)//' will have drydep applied'
endif
enddo
do m = 1,nwetdep
call cnst_get_ind ( wetdep_list(m), id, abort=.false. )
if (id>0) then
wetdep_indices(m) = id
else
call endrun(subrname//': invalid wetdep species: '//trim(wetdep_list(m)) )
endif
if (masterproc) then
write(iulog,*) subrname//': '//wetdep_list(m)//' will have wet removal'
endif
enddo
if (ndrydep>0) then
call inidrydep(rair, gravit)
dummy = 'RAM1'
call addfld (dummy,'frac ',1, 'A','RAM1',phys_decomp)
if ( history_aerosol ) then
call add_default (dummy, 1, ' ')
endif
dummy = 'airFV'
call addfld (dummy,'frac ',1, 'A','FV',phys_decomp)
if ( history_aerosol ) then
call add_default (dummy, 1, ' ')
endif
endif
if (dust_active) then
! emissions diagnostics ....
do m = 1, dust_nbin+dust_nnum
dummy = trim(dust_names(m)) // 'SF'
call addfld (dummy,'kg/m2/s ',1, 'A',trim(dust_names(m))//' dust surface emission',phys_decomp)
if (history_aerosol) then
call add_default (dummy, 1, ' ')
endif
enddo
dummy = 'DSTSFMBL'
call addfld (dummy,'kg/m2/s',1, 'A','Mobilization flux at surface',phys_decomp)
if (history_aerosol) then
call add_default (dummy, 1, ' ')
endif
dummy = 'LND_MBL'
call addfld (dummy,'frac ',1, 'A','Soil erodibility factor',phys_decomp)
if (history_aerosol) then
call add_default (dummy, 1, ' ')
endif
endif
if (seasalt_active) then
dummy = 'SSTSFMBL'
call addfld (dummy,'kg/m2/s',1, 'A','Mobilization flux at surface',phys_decomp)
if (history_aerosol) then
call add_default (dummy, 1, ' ')
endif
do m = 1, seasalt_nbin
dummy = trim(seasalt_names(m)) // 'SF'
call addfld (dummy,'kg/m2/s ',1, 'A',trim(seasalt_names(m))//' seasalt surface emission',phys_decomp)
if (history_aerosol) then
call add_default (dummy, 1, ' ')
endif
enddo
endif
! set flags for drydep tendencies
drydep_lq(:) = .false.
do m=1,ndrydep
id = drydep_indices(m)
drydep_lq(id) = .true.
enddo
! set flags for wetdep tendencies
wetdep_lq(:) = .false.
do m=1,nwetdep
id = wetdep_indices(m)
wetdep_lq(id) = .true.
enddo
wetdens_ap_idx = pbuf_get_index('WETDENS_AP')
qaerwat_idx = pbuf_get_index('QAERWAT')
pblh_idx = pbuf_get_index('pblh')
rate1_cw2pr_st_idx = pbuf_get_index('RATE1_CW2PR_ST')
call pbuf_set_field(pbuf2d, rate1_cw2pr_st_idx, 0.0_r8)
do m = 1,ndrydep
! units
if (drydep_list(m)(1:3) == 'num') then
unit_basename = ' 1'
else
unit_basename = 'kg'
endif
call addfld (trim(drydep_list(m))//'DDF',unit_basename//'/m2/s ', 1, 'A', &
trim(drydep_list(m))//' dry deposition flux at bottom (grav + turb)',phys_decomp)
call addfld (trim(drydep_list(m))//'TBF',unit_basename//'/m2/s', 1, 'A', &
trim(drydep_list(m))//' turbulent dry deposition flux',phys_decomp)
call addfld (trim(drydep_list(m))//'GVF',unit_basename//'/m2/s ', 1, 'A', &
trim(drydep_list(m))//' gravitational dry deposition flux',phys_decomp)
call addfld (trim(drydep_list(m))//'DTQ',unit_basename//'/kg/s ',pver, 'A', &
trim(drydep_list(m))//' dry deposition',phys_decomp)
call addfld (trim(drydep_list(m))//'DDV','m/s ',pver, 'A', &
trim(drydep_list(m))//' deposition velocity',phys_decomp)
if ( history_aerosol ) then
call add_default (trim(drydep_list(m))//'DDF', 1, ' ')
call add_default (trim(drydep_list(m))//'TBF', 1, ' ')
call add_default (trim(drydep_list(m))//'GVF', 1, ' ')
endif
enddo
do m = 1,nwetdep
! units
if (wetdep_list(m)(1:3) == 'num') then
unit_basename = ' 1'
else
unit_basename = 'kg'
endif
call addfld (trim(wetdep_list(m))//'SFWET',unit_basename//'/m2/s ', &
1, 'A','Wet deposition flux at surface',phys_decomp)
call addfld (trim(wetdep_list(m))//'SFSIC',unit_basename//'/m2/s ', &
1, 'A','Wet deposition flux (incloud, convective) at surface',phys_decomp)
call addfld (trim(wetdep_list(m))//'SFSIS',unit_basename//'/m2/s ', &
1, 'A','Wet deposition flux (incloud, stratiform) at surface',phys_decomp)
call addfld (trim(wetdep_list(m))//'SFSBC',unit_basename//'/m2/s ', &
1, 'A','Wet deposition flux (belowcloud, convective) at surface',phys_decomp)
call addfld (trim(wetdep_list(m))//'SFSBS',unit_basename//'/m2/s ', &
1, 'A','Wet deposition flux (belowcloud, stratiform) at surface',phys_decomp)
call addfld (trim(wetdep_list(m))//'WET',unit_basename//'/kg/s ',pver, 'A','wet deposition tendency',phys_decomp)
call addfld (trim(wetdep_list(m))//'SIC',unit_basename//'/kg/s ',pver, 'A', &
trim(wetdep_list(m))//' ic wet deposition',phys_decomp)
call addfld (trim(wetdep_list(m))//'SIS',unit_basename//'/kg/s ',pver, 'A', &
trim(wetdep_list(m))//' is wet deposition',phys_decomp)
call addfld (trim(wetdep_list(m))//'SBC',unit_basename//'/kg/s ',pver, 'A', &
trim(wetdep_list(m))//' bc wet deposition',phys_decomp)
call addfld (trim(wetdep_list(m))//'SBS',unit_basename//'/kg/s ',pver, 'A', &
trim(wetdep_list(m))//' bs wet deposition',phys_decomp)
if ( history_aerosol ) then
call add_default (trim(wetdep_list(m))//'SFWET', 1, ' ')
call add_default (trim(wetdep_list(m))//'SFSIC', 1, ' ')
call add_default (trim(wetdep_list(m))//'SFSIS', 1, ' ')
call add_default (trim(wetdep_list(m))//'SFSBC', 1, ' ')
call add_default (trim(wetdep_list(m))//'SFSBS', 1, ' ')
endif
enddo
do m = 1,gas_pcnst
if ( solsym(m)(1:3) == 'num') then
unit_basename = ' 1' ! Units 'kg' or '1'
else
unit_basename = 'kg' ! Units 'kg' or '1'
end if
call addfld( 'GS_'//trim(solsym(m)), unit_basename//'/m2/s ',1, 'A', &
trim(solsym(m))//' gas chemistry/wet removal (for gas species)', phys_decomp)
call addfld( 'AQ_'//trim(solsym(m)), unit_basename//'/m2/s ',1, 'A', &
trim(solsym(m))//' aqueous chemistry (for gas species)', phys_decomp)
if ( history_aerosol ) then
call add_default( 'GS_'//trim(solsym(m)), 1, ' ')
call add_default( 'AQ_'//trim(solsym(m)), 1, ' ')
endif
enddo
do n = 1,pcnst
if( .not. (cnst_name_cw(n) == ' ') ) then
if (cnst_name_cw(n)(1:3) == 'num') then
unit_basename = ' 1'
else
unit_basename = 'kg'
endif
call addfld( cnst_name_cw(n), unit_basename//'/kg ', pver, 'A', &
trim(cnst_name_cw(n))//' in cloud water',phys_decomp)
call addfld (trim(cnst_name_cw(n))//'SFWET', unit_basename//'/m2/s ',1, 'A', &
trim(cnst_name_cw(n))//' wet deposition flux at surface',phys_decomp)
call addfld (trim(cnst_name_cw(n))//'SFSIC', unit_basename//'/m2/s ',1, 'A', &
trim(cnst_name_cw(n))//' wet deposition flux (incloud, convective) at surface',phys_decomp)
call addfld (trim(cnst_name_cw(n))//'SFSIS', unit_basename//'/m2/s ',1, 'A', &
trim(cnst_name_cw(n))//' wet deposition flux (incloud, stratiform) at surface',phys_decomp)
call addfld (trim(cnst_name_cw(n))//'SFSBC', unit_basename//'/m2/s ',1, 'A', &
trim(cnst_name_cw(n))//' wet deposition flux (belowcloud, convective) at surface',phys_decomp)
call addfld (trim(cnst_name_cw(n))//'SFSBS', unit_basename//'/m2/s ',1, 'A', &
trim(cnst_name_cw(n))//' wet deposition flux (belowcloud, stratiform) at surface',phys_decomp)
call addfld (trim(cnst_name_cw(n))//'DDF', unit_basename//'/m2/s ', 1, 'A', &
trim(cnst_name_cw(n))//' dry deposition flux at bottom (grav + turb)',phys_decomp)
call addfld (trim(cnst_name_cw(n))//'TBF', unit_basename//'/m2/s ', 1, 'A', &
trim(cnst_name_cw(n))//' turbulent dry deposition flux',phys_decomp)
call addfld (trim(cnst_name_cw(n))//'GVF', unit_basename//'/m2/s ', 1, 'A', &
trim(cnst_name_cw(n))//' gravitational dry deposition flux',phys_decomp)
if ( history_aerosol ) then
call add_default( cnst_name_cw(n), 1, ' ' )
call add_default (trim(cnst_name_cw(n))//'GVF', 1, ' ')
call add_default (trim(cnst_name_cw(n))//'SFWET', 1, ' ')
call add_default (trim(cnst_name_cw(n))//'TBF', 1, ' ')
call add_default (trim(cnst_name_cw(n))//'DDF', 1, ' ')
call add_default (trim(cnst_name_cw(n))//'SFSBS', 1, ' ')
call add_default (trim(cnst_name_cw(n))//'SFSIC', 1, ' ')
call add_default (trim(cnst_name_cw(n))//'SFSBC', 1, ' ')
call add_default (trim(cnst_name_cw(n))//'SFSIS', 1, ' ')
endif
endif
enddo
do n=1,ntot_amode
dgnum_name(n) = ' '
write(dgnum_name(n),fmt='(a,i1)') 'dgnumwet',n
call addfld( dgnum_name(n), 'm', pver, 'I', 'Aerosol mode wet diameter', phys_decomp )
if ( history_aerosol ) then
call add_default( dgnum_name(n), 1, ' ' )
endif
end do
ndx_h2so4 = get_spc_ndx('H2SO4')
! for aero_model_surfarea called from mo_usrrxt
do l=1,ntot_amode
if ( trim(modename_amode(l)) == 'aitken' ) then
aitken_idx = l
end if
test_name = ' '
write(test_name,fmt='(a5,i1)') 'num_a',l
num_idx(l) = get_spc_ndx( trim(test_name) )
if (num_idx(l) < 0) then
write(errmes,fmt='(a,i1)') 'usrrxt_inti: cannot find MAM num_idx ',l
write(iulog,*) errmes
call endrun(errmes)
endif
end do
dgnumwet_idx = pbuf_get_index('DGNUMWET')
if ( aitken_idx < 0 ) then
errmes = 'usrrxt_inti: cannot find aitken_idx'
call endrun(errmes)
end if
!
! define indeces associated with the various names (defined in
! chemistry/modal_aero/modal_aero_initialize_data.F90)
!
#if ( defined MODAL_AERO_3MODE )
!
! accumulation mode #1
!
index_tot_mass(1,1) = get_spc_ndx('so4_a1')
index_tot_mass(1,2) = get_spc_ndx('pom_a1')
index_tot_mass(1,3) = get_spc_ndx('soa_a1')
index_tot_mass(1,4) = get_spc_ndx('bc_a1' )
index_tot_mass(1,5) = get_spc_ndx('dst_a1')
index_tot_mass(1,6) = get_spc_ndx('ncl_a1')
index_chm_mass(1,1) = get_spc_ndx('so4_a1')
index_chm_mass(1,2) = get_spc_ndx('soa_a1')
index_chm_mass(1,3) = get_spc_ndx('bc_a1' )
!
! aitken mode
!
index_tot_mass(2,1) = get_spc_ndx('so4_a2')
index_tot_mass(2,2) = get_spc_ndx('soa_a2')
index_tot_mass(2,3) = get_spc_ndx('ncl_a2')
index_chm_mass(2,1) = get_spc_ndx('so4_a2')
index_chm_mass(2,2) = get_spc_ndx('soa_a2')
!
! coarse mode
!
index_tot_mass(3,1) = get_spc_ndx('dst_a3')
index_tot_mass(3,2) = get_spc_ndx('ncl_a3')
index_tot_mass(3,3) = get_spc_ndx('so4_a3')
index_chm_mass(3,1) = get_spc_ndx('so4_a3')
!
#endif
#if ( defined MODAL_AERO_7MODE )
!
! accumulation mode #1
!
index_tot_mass(1,1) = get_spc_ndx('so4_a1')
index_tot_mass(1,2) = get_spc_ndx('nh4_a1')
index_tot_mass(1,3) = get_spc_ndx('pom_a1')
index_tot_mass(1,4) = get_spc_ndx('soa_a1')
index_tot_mass(1,5) = get_spc_ndx('bc_a1' )
index_tot_mass(1,6) = get_spc_ndx('ncl_a1')
index_chm_mass(1,1) = get_spc_ndx('so4_a1')
index_chm_mass(1,2) = get_spc_ndx('nh4_a1')
index_chm_mass(1,3) = get_spc_ndx('soa_a1')
index_chm_mass(1,4) = get_spc_ndx('bc_a1' )
!
! aitken mode
!
index_tot_mass(2,1) = get_spc_ndx('so4_a2')
index_tot_mass(2,2) = get_spc_ndx('nh4_a2')
index_tot_mass(2,3) = get_spc_ndx('soa_a2')
index_tot_mass(2,4) = get_spc_ndx('ncl_a2')
index_chm_mass(2,1) = get_spc_ndx('so4_a2')
index_chm_mass(2,2) = get_spc_ndx('nh4_a2')
index_chm_mass(2,3) = get_spc_ndx('soa_a2')
!
! primary carbon mode not added
!
! fine sea salt
!
index_tot_mass(4,1) = get_spc_ndx('so4_a4')
index_tot_mass(4,2) = get_spc_ndx('nh4_a4')
index_tot_mass(4,3) = get_spc_ndx('ncl_a4')
index_chm_mass(4,1) = get_spc_ndx('so4_a4')
index_chm_mass(4,2) = get_spc_ndx('nh4_a4')
!
! fine soil dust
!
index_tot_mass(5,1) = get_spc_ndx('so4_a5')
index_tot_mass(5,2) = get_spc_ndx('nh4_a5')
index_tot_mass(5,3) = get_spc_ndx('dst_a5')
index_chm_mass(5,1) = get_spc_ndx('so4_a5')
index_chm_mass(5,2) = get_spc_ndx('nh4_a5')
!
! coarse sea salt
!
index_tot_mass(6,1) = get_spc_ndx('so4_a6')
index_tot_mass(6,2) = get_spc_ndx('nh4_a6')
index_tot_mass(6,3) = get_spc_ndx('ncl_a6')
index_chm_mass(6,1) = get_spc_ndx('so4_a6')
index_chm_mass(6,2) = get_spc_ndx('nh4_a6')
!
! coarse soil dust
!
index_tot_mass(7,1) = get_spc_ndx('so4_a7')
index_tot_mass(7,2) = get_spc_ndx('nh4_a7')
index_tot_mass(7,3) = get_spc_ndx('dst_a7')
index_chm_mass(7,1) = get_spc_ndx('so4_a7')
index_chm_mass(7,2) = get_spc_ndx('nh4_a7')
!
#endif
#if ( defined MODAL_AERO_9MODE )
! accumulation mode #1
!
index_tot_mass(1,1) = get_spc_ndx('so4_a1')
index_tot_mass(1,2) = get_spc_ndx('nh4_a1')
index_tot_mass(1,3) = get_spc_ndx('pom_a1')
index_tot_mass(1,4) = get_spc_ndx('soa_a1')
index_tot_mass(1,5) = get_spc_ndx('bc_a1' )
index_tot_mass(1,6) = get_spc_ndx('ncl_a1')
index_chm_mass(1,1) = get_spc_ndx('so4_a1')
index_chm_mass(1,2) = get_spc_ndx('nh4_a1')
index_chm_mass(1,3) = get_spc_ndx('soa_a1')
index_chm_mass(1,4) = get_spc_ndx('bc_a1' )
!
! aitken mode
!
index_tot_mass(2,1) = get_spc_ndx('so4_a2')
index_tot_mass(2,2) = get_spc_ndx('nh4_a2')
index_tot_mass(2,3) = get_spc_ndx('soa_a2')
index_tot_mass(2,4) = get_spc_ndx('ncl_a2')
index_chm_mass(2,1) = get_spc_ndx('so4_a2')
index_chm_mass(2,2) = get_spc_ndx('nh4_a2')
index_chm_mass(2,3) = get_spc_ndx('soa_a2')
!
! primary carbon mode not added
!
! fine sea salt
!
index_tot_mass(4,1) = get_spc_ndx('so4_a4')
index_tot_mass(4,2) = get_spc_ndx('nh4_a4')
index_tot_mass(4,3) = get_spc_ndx('ncl_a4')
index_chm_mass(4,1) = get_spc_ndx('so4_a4')
index_chm_mass(4,2) = get_spc_ndx('nh4_a4')
!
! fine soil dust
!
index_tot_mass(5,1) = get_spc_ndx('so4_a5')
index_tot_mass(5,2) = get_spc_ndx('nh4_a5')
index_tot_mass(5,3) = get_spc_ndx('dst_a5')
index_chm_mass(5,1) = get_spc_ndx('so4_a5')
index_chm_mass(5,2) = get_spc_ndx('nh4_a5')
!
! coarse sea salt
!
index_tot_mass(6,1) = get_spc_ndx('so4_a6')
index_tot_mass(6,2) = get_spc_ndx('nh4_a6')
index_tot_mass(6,3) = get_spc_ndx('ncl_a6')
index_chm_mass(6,1) = get_spc_ndx('so4_a6')
index_chm_mass(6,2) = get_spc_ndx('nh4_a6')
!
! coarse soil dust1 kzm
!
index_tot_mass(7,1) = get_spc_ndx('so4_a7')
index_tot_mass(7,2) = get_spc_ndx('nh4_a7')
index_tot_mass(7,3) = get_spc_ndx('dst_a7')
index_chm_mass(7,1) = get_spc_ndx('so4_a7')
index_chm_mass(7,2) = get_spc_ndx('nh4_a7')
!
! coarse soil dust2 kzm
!
index_tot_mass(8,1) = get_spc_ndx('so4_a8')
index_tot_mass(8,2) = get_spc_ndx('nh4_a8')
index_tot_mass(8,3) = get_spc_ndx('dst_a8')
index_chm_mass(8,1) = get_spc_ndx('so4_a8')
index_chm_mass(8,2) = get_spc_ndx('nh4_a8')
!
! coarse soil dust3 kzm
!
index_tot_mass(9,1) = get_spc_ndx('so4_a9')
index_tot_mass(9,2) = get_spc_ndx('nh4_a9')
index_tot_mass(9,3) = get_spc_ndx('dst_a9')
index_chm_mass(9,1) = get_spc_ndx('so4_a9')
index_chm_mass(9,2) = get_spc_ndx('nh4_a9')
#endif
end subroutine aero_model_init
!=============================================================================
!=============================================================================
subroutine aero_model_drydep ( state, pbuf, obklen, ustar, cam_in, dt, cam_out, ptend )
use dust_sediment_mod, only: dust_sediment_tend
use drydep_mod, only: d3ddflux, calcram
use modal_aero_data, only: qqcw_get_field
use modal_aero_data, only: cnst_name_cw
use modal_aero_data, only: alnsg_amode
use modal_aero_data, only: sigmag_amode
use modal_aero_data, only: nspec_amode
use modal_aero_data, only: numptr_amode
use modal_aero_data, only: numptrcw_amode
use modal_aero_data, only: lmassptr_amode
use modal_aero_data, only: lmassptrcw_amode
use modal_aero_deposition, only: set_srf_drydep
! args
type(physics_state), intent(in) :: state ! Physics state variables
real(r8), intent(in) :: obklen(:)
real(r8), intent(in) :: ustar(:) ! sfc fric vel
type(cam_in_t), target, intent(in) :: cam_in ! import state
real(r8), intent(in) :: dt ! time step
type(cam_out_t), intent(inout) :: cam_out ! export state
type(physics_ptend), intent(out) :: ptend ! indivdual parameterization tendencies
type(physics_buffer_desc), pointer :: pbuf(:)
! local vars
real(r8), pointer :: landfrac(:) ! land fraction
real(r8), pointer :: icefrac(:) ! ice fraction
real(r8), pointer :: ocnfrac(:) ! ocean fraction
real(r8), pointer :: fvin(:) !
real(r8), pointer :: ram1in(:) ! for dry dep velocities from land model for progseasalts
real(r8) :: fv(pcols) ! for dry dep velocities, from land modified over ocean & ice
real(r8) :: ram1(pcols) ! for dry dep velocities, from land modified over ocean & ice
integer :: lchnk ! chunk identifier
integer :: ncol ! number of atmospheric columns
integer :: jvlc ! index for last dimension of vlc_xxx arrays
integer :: lphase ! index for interstitial / cloudborne aerosol
integer :: lspec ! index for aerosol number / chem-mass / water-mass
integer :: m ! aerosol mode index
integer :: mm ! tracer index
integer :: i
real(r8) :: tvs(pcols,pver)
real(r8) :: rho(pcols,pver) ! air density in kg/m3
real(r8) :: sflx(pcols) ! deposition flux
real(r8) :: dep_trb(pcols) !kg/m2/s
real(r8) :: dep_grv(pcols) !kg/m2/s (total of grav and trb)
real(r8) :: pvmzaer(pcols,pverp) ! sedimentation velocity in Pa
real(r8) :: dqdt_tmp(pcols,pver) ! temporary array to hold tendency for 1 species
real(r8) :: rad_drop(pcols,pver)
real(r8) :: dens_drop(pcols,pver)
real(r8) :: sg_drop(pcols,pver)
real(r8) :: rad_aer(pcols,pver)
real(r8) :: dens_aer(pcols,pver)
real(r8) :: sg_aer(pcols,pver)
real(r8) :: vlc_dry(pcols,pver,4) ! dep velocity
real(r8) :: vlc_grv(pcols,pver,4) ! dep velocity
real(r8):: vlc_trb(pcols,4) ! dep velocity
real(r8) :: aerdepdryis(pcols,pcnst) ! aerosol dry deposition (interstitial)
real(r8) :: aerdepdrycw(pcols,pcnst) ! aerosol dry deposition (cloud water)
real(r8), pointer :: fldcw(:,:)
real(r8), pointer :: dgncur_awet(:,:,:)
real(r8), pointer :: wetdens(:,:,:)
real(r8), pointer :: qaerwat(:,:,:)
landfrac => cam_in%landfrac(:)
icefrac => cam_in%icefrac(:)
ocnfrac => cam_in%ocnfrac(:)
fvin => cam_in%fv(:)
ram1in => cam_in%ram1(:)
lchnk = state%lchnk
ncol = state%ncol
! calc ram and fv over ocean and sea ice ...
call calcram( ncol,landfrac,icefrac,ocnfrac,obklen,&
ustar,ram1in,ram1,state%t(:,pver),state%pmid(:,pver),&
state%pdel(:,pver),fvin,fv)
call outfld( 'airFV', fv(:), pcols, lchnk )
call outfld( 'RAM1', ram1(:), pcols, lchnk )
! note that tendencies are not only in sfc layer (because of sedimentation)
! and that ptend is updated within each subroutine for different species
call physics_ptend_init(ptend, state%psetcols, 'aero_model_drydep', lq=drydep_lq)
call pbuf_get_field(pbuf, dgnumwet_idx, dgncur_awet, start=(/1,1,1/), kount=(/pcols,pver,nmodes/) )
call pbuf_get_field(pbuf, wetdens_ap_idx, wetdens, start=(/1,1,1/), kount=(/pcols,pver,nmodes/) )
call pbuf_get_field(pbuf, qaerwat_idx, qaerwat, start=(/1,1,1/), kount=(/pcols,pver,nmodes/) )
tvs(:ncol,:) = state%t(:ncol,:)!*(1+state%q(:ncol,k)
rho(:ncol,:)= state%pmid(:ncol,:)/(rair*state%t(:ncol,:))
!
! calc settling/deposition velocities for cloud droplets (and cloud-borne aerosols)
!
! *** mean drop radius should eventually be computed from ndrop and qcldwtr
rad_drop(:,:) = 5.0e-6_r8
dens_drop(:,:) = rhoh2o
sg_drop(:,:) = 1.46_r8
jvlc = 3
call modal_aero_depvel_part( ncol,state%t(:,:), state%pmid(:,:), ram1, fv, &
vlc_dry(:,:,jvlc), vlc_trb(:,jvlc), vlc_grv(:,:,jvlc), &
rad_drop(:,:), dens_drop(:,:), sg_drop(:,:), 0, lchnk)
jvlc = 4
call modal_aero_depvel_part( ncol,state%t(:,:), state%pmid(:,:), ram1, fv, &
vlc_dry(:,:,jvlc), vlc_trb(:,jvlc), vlc_grv(:,:,jvlc), &
rad_drop(:,:), dens_drop(:,:), sg_drop(:,:), 3, lchnk)
do m = 1, ntot_amode ! main loop over aerosol modes
do lphase = 1, 2 ! loop over interstitial / cloud-borne forms
if (lphase == 1) then ! interstial aerosol - calc settling/dep velocities of mode
! rad_aer = volume mean wet radius (m)
! dgncur_awet = geometric mean wet diameter for number distribution (m)
rad_aer(1:ncol,:) = 0.5_r8*dgncur_awet(1:ncol,:,m) &
*exp(1.5_r8*(alnsg_amode(m)**2))
! dens_aer(1:ncol,:) = wet density (kg/m3)
dens_aer(1:ncol,:) = wetdens(1:ncol,:,m)
sg_aer(1:ncol,:) = sigmag_amode(m)
jvlc = 1
call modal_aero_depvel_part( ncol, state%t(:,:), state%pmid(:,:), ram1, fv, &
vlc_dry(:,:,jvlc), vlc_trb(:,jvlc), vlc_grv(:,:,jvlc), &
rad_aer(:,:), dens_aer(:,:), sg_aer(:,:), 0, lchnk)
jvlc = 2
call modal_aero_depvel_part( ncol, state%t(:,:), state%pmid(:,:), ram1, fv, &
vlc_dry(:,:,jvlc), vlc_trb(:,jvlc), vlc_grv(:,:,jvlc), &
rad_aer(:,:), dens_aer(:,:), sg_aer(:,:), 3, lchnk)
end if
do lspec = 0, nspec_amode(m)+1 ! loop over number + constituents + water
if (lspec == 0) then ! number
if (lphase == 1) then
mm = numptr_amode(m)
jvlc = 1
else
mm = numptrcw_amode(m)
jvlc = 3
endif
else if (lspec <= nspec_amode(m)) then ! non-water mass
if (lphase == 1) then
mm = lmassptr_amode(lspec,m)
jvlc = 2
else
mm = lmassptrcw_amode(lspec,m)
jvlc = 4
endif
else ! water mass
! bypass dry deposition of aerosol water
cycle
if (lphase == 1) then
mm = 0
! mm = lwaterptr_amode(m)
jvlc = 2
else
mm = 0
jvlc = 4
endif
endif
if (mm <= 0) cycle
! if (lphase == 1) then
if ((lphase == 1) .and. (lspec <= nspec_amode(m))) then
ptend%lq(mm) = .TRUE.
! use pvprogseasalts instead (means making the top level 0)
pvmzaer(:ncol,1)=0._r8
pvmzaer(:ncol,2:pverp) = vlc_dry(:ncol,:,jvlc)
call outfld( trim(cnst_name(mm))//'DDV', pvmzaer(:,2:pverp), pcols, lchnk )
if(.true.) then ! use phil's method
! convert from meters/sec to pascals/sec
! pvprogseasalts(:,1) is assumed zero, use density from layer above in conversion
pvmzaer(:ncol,2:pverp) = pvmzaer(:ncol,2:pverp) * rho(:ncol,:)*gravit
! calculate the tendencies and sfc fluxes from the above velocities
call dust_sediment_tend( &
ncol, dt, state%pint(:,:), state%pmid, state%pdel, state%t , &
state%q(:,:,mm), pvmzaer, ptend%q(:,:,mm), sflx )
else !use charlie's method
call d3ddflux( ncol, vlc_dry(:,:,jvlc), state%q(:,:,mm), state%pmid, &
state%pdel, tvs, sflx, ptend%q(:,:,mm), dt )
endif
! apportion dry deposition into turb and gravitational settling for tapes
do i=1,ncol
dep_trb(i)=sflx(i)*vlc_trb(i,jvlc)/vlc_dry(i,pver,jvlc)
dep_grv(i)=sflx(i)*vlc_grv(i,pver,jvlc)/vlc_dry(i,pver,jvlc)
enddo
call outfld( trim(cnst_name(mm))//'DDF', sflx, pcols, lchnk)
call outfld( trim(cnst_name(mm))//'TBF', dep_trb, pcols, lchnk )
call outfld( trim(cnst_name(mm))//'GVF', dep_grv, pcols, lchnk )
call outfld( trim(cnst_name(mm))//'DTQ', ptend%q(:,:,mm), pcols, lchnk)
aerdepdryis(:ncol,mm) = sflx(:ncol)
else if ((lphase == 1) .and. (lspec == nspec_amode(m)+1)) then ! aerosol water
! use pvprogseasalts instead (means making the top level 0)
pvmzaer(:ncol,1)=0._r8
pvmzaer(:ncol,2:pverp) = vlc_dry(:ncol,:,jvlc)
if(.true.) then ! use phil's method
! convert from meters/sec to pascals/sec
! pvprogseasalts(:,1) is assumed zero, use density from layer above in conversion
pvmzaer(:ncol,2:pverp) = pvmzaer(:ncol,2:pverp) * rho(:ncol,:)*gravit
! calculate the tendencies and sfc fluxes from the above velocities
call dust_sediment_tend( &
ncol, dt, state%pint(:,:), state%pmid, state%pdel, state%t , &
qaerwat(:,:,mm), pvmzaer, dqdt_tmp(:,:), sflx )
else !use charlie's method
call d3ddflux( ncol, vlc_dry(:,:,jvlc), qaerwat(:,:,mm), state%pmid, &
state%pdel, tvs, sflx, dqdt_tmp(:,:), dt )
endif
! apportion dry deposition into turb and gravitational settling for tapes
do i=1,ncol
dep_trb(i)=sflx(i)*vlc_trb(i,jvlc)/vlc_dry(i,pver,jvlc)
dep_grv(i)=sflx(i)*vlc_grv(i,pver,jvlc)/vlc_dry(i,pver,jvlc)
enddo
qaerwat(1:ncol,:,mm) = qaerwat(1:ncol,:,mm) + dqdt_tmp(1:ncol,:) * dt
else ! lphase == 2
! use pvprogseasalts instead (means making the top level 0)
pvmzaer(:ncol,1)=0._r8
pvmzaer(:ncol,2:pverp) = vlc_dry(:ncol,:,jvlc)
fldcw => qqcw_get_field(pbuf, mm,lchnk)
if(.true.) then ! use phil's method
! convert from meters/sec to pascals/sec
! pvprogseasalts(:,1) is assumed zero, use density from layer above in conversion
pvmzaer(:ncol,2:pverp) = pvmzaer(:ncol,2:pverp) * rho(:ncol,:)*gravit
! calculate the tendencies and sfc fluxes from the above velocities
call dust_sediment_tend( &
ncol, dt, state%pint(:,:), state%pmid, state%pdel, state%t , &
fldcw(:,:), pvmzaer, dqdt_tmp(:,:), sflx )
else !use charlie's method
call d3ddflux( ncol, vlc_dry(:,:,jvlc), fldcw(:,:), state%pmid, &
state%pdel, tvs, sflx, dqdt_tmp(:,:), dt )
endif
! apportion dry deposition into turb and gravitational settling for tapes
do i=1,ncol
dep_trb(i)=sflx(i)*vlc_trb(i,jvlc)/vlc_dry(i,pver,jvlc)
dep_grv(i)=sflx(i)*vlc_grv(i,pver,jvlc)/vlc_dry(i,pver,jvlc)
enddo
fldcw(1:ncol,:) = fldcw(1:ncol,:) + dqdt_tmp(1:ncol,:) * dt
call outfld( trim(cnst_name_cw(mm))//'DDF', sflx, pcols, lchnk)
call outfld( trim(cnst_name_cw(mm))//'TBF', dep_trb, pcols, lchnk )
call outfld( trim(cnst_name_cw(mm))//'GVF', dep_grv, pcols, lchnk )
aerdepdrycw(:ncol,mm) = sflx(:ncol)
endif
enddo ! lspec = 0, nspec_amode(m)+1
enddo ! lphase = 1, 2
enddo ! m = 1, ntot_amode
! if the user has specified prescribed aerosol dep fluxes then
! do not set cam_out dep fluxes according to the prognostic aerosols
if (.not.aerodep_flx_prescribed()) then
call set_srf_drydep(aerdepdryis, aerdepdrycw, cam_out)
endif
endsubroutine aero_model_drydep
!=============================================================================
!=============================================================================
subroutine aero_model_wetdep( state, dt, dlf, cam_out, ptend, pbuf)
use modal_aero_deposition, only: set_srf_wetdep
use wetdep, only: wetdepa_v2, wetdep_inputs_set, wetdep_inputs_t
use modal_aero_data
use modal_aero_calcsize, only: modal_aero_calcsize_sub
use modal_aero_wateruptake,only: modal_aero_wateruptake_dr
! args
type(physics_state), intent(in) :: state ! Physics state variables
real(r8), intent(in) :: dt ! time step
real(r8), intent(in) :: dlf(:,:) ! shallow+deep convective detrainment [kg/kg/s]
type(cam_out_t), intent(inout) :: cam_out ! export state
type(physics_ptend), intent(out) :: ptend ! indivdual parameterization tendencies
type(physics_buffer_desc), pointer :: pbuf(:)
! local vars
integer :: m ! tracer index
integer :: lchnk ! chunk identifier
integer :: ncol ! number of atmospheric columns
real(r8) :: iscavt(pcols, pver)
integer :: mm
integer :: i,k
real(r8) :: icscavt(pcols, pver)
real(r8) :: isscavt(pcols, pver)
real(r8) :: bcscavt(pcols, pver)
real(r8) :: bsscavt(pcols, pver)
real(r8) :: sol_factb, sol_facti
real(r8) :: sol_factic(pcols,pver)
real(r8) :: sflx(pcols) ! deposition flux
integer :: jnv ! index for scavcoefnv 3rd dimension
integer :: lphase ! index for interstitial / cloudborne aerosol
integer :: lspec ! index for aerosol number / chem-mass / water-mass
integer :: lcoardust, lcoarnacl ! indices for coarse mode dust and seasalt masses
real(r8) :: dqdt_tmp(pcols,pver) ! temporary array to hold tendency for 1 species
real(r8) :: f_act_conv(pcols,pver) ! prescribed aerosol activation fraction for convective cloud ! rce 2010/05/01
real(r8) :: f_act_conv_coarse(pcols,pver) ! similar but for coarse mode ! rce 2010/05/02
real(r8) :: f_act_conv_coarse_dust, f_act_conv_coarse_nacl ! rce 2010/05/02
real(r8) :: fracis_cw(pcols,pver)
real(r8) :: hygro_sum_old(pcols,pver) ! before removal [sum of (mass*hydro/dens)]
real(r8) :: hygro_sum_del(pcols,pver) ! removal change to [sum of (mass*hydro/dens)]
real(r8) :: hygro_sum_old_ik, hygro_sum_new_ik
real(r8) :: prec(pcols) ! precipitation rate
real(r8) :: q_tmp(pcols,pver) ! temporary array to hold "most current" mixing ratio for 1 species
real(r8) :: scavcoefnv(pcols,pver,0:2) ! Dana and Hales coefficient (/mm) for
! cloud-borne num & vol (0),
! interstitial num (1), interstitial vol (2)
real(r8) :: tmpa, tmpb
real(r8) :: tmpdust, tmpnacl
real(r8) :: water_old, water_new ! temporary old/new aerosol water mix-rat
logical :: isprx(pcols,pver) ! true if precipation