updated top-level README and version_decl for V4.5 (#1847)
[WRF.git] / chem / module_mosaic_sect_intr.F
blob55ca9bfca5e5dcdfd04edaa5c1bcecf5f863d0d5
1 module module_mosaic_sect_intr
2 !-----------------------------------------------------------------------
4 ! 2015-06-09 r.c.easter - changes to the MOSAIC box model version for use in WRF-Chem  
6 ! deleted subr sect_iface_dens_size_testbb and the calls to it
7 ! deleted subr mbox_sectional_diagnostics and the call to it
9 !-----------------------------------------------------------------------
12 ! *** NOTE ***
14 ! this routine currently will not work with openmp parallelization
16 ! some of the variables taken from module_data_mosaic_main (cnn, pr_atm, rh, te)
17 ! needs to be passed as subroutine arguments
20   use module_data_mosaic_kind, only:  r8
23   implicit none
26   integer, parameter :: iunits_flagaa = 4
27   integer, parameter :: iunits_flagbb = 1
28   integer, parameter :: iunits_diagout_flagaa = 1
30   !BSINGH - 05/28/2013(RCE updates)
31   real(r8), save, allocatable ::   &
32           dens_aer_tmp(:,:),  &
33           volumcut_sect_tmp(:,:),  &
34           volumcen_sect_tmp(:,:),  &
35           volumlo_sect_tmp(:,:),  &
36           volumhi_sect_tmp(:,:),  &
37           dcut_sect_tmp(:,:),  &
38           dcen_sect_tmp(:,:),  &
39           dlo_sect_tmp(:,:),  &
40           dhi_sect_tmp(:,:)
41   !BSINGH - 05/28/2013(RCE updates ENDS)
43   contains
46   !BSINGH - 05/28/2013(RCE updates)
47 !-----------------------------------------------------------------------
48         subroutine sect_intr_allocate_memory
50         use module_data_mosaic_asecthp, only:  maxd_acomp, maxd_asize, maxd_atype
52         allocate( dens_aer_tmp( maxd_acomp, maxd_atype ) )
53         allocate( volumcut_sect_tmp( 0:maxd_asize, maxd_atype ) )
54         allocate( volumcen_sect_tmp(   maxd_asize, maxd_atype ) )
55         allocate( volumlo_sect_tmp(    maxd_asize, maxd_atype ) )
56         allocate( volumhi_sect_tmp(    maxd_asize, maxd_atype ) )
57         allocate( dcut_sect_tmp( 0:maxd_asize, maxd_atype ) )
58         allocate( dcen_sect_tmp(   maxd_asize, maxd_atype ) )
59         allocate( dlo_sect_tmp(    maxd_asize, maxd_atype ) )
60         allocate( dhi_sect_tmp(    maxd_asize, maxd_atype ) )
62         return
63         end subroutine sect_intr_allocate_memory
64         !BSINGH - 05/28/2013(RCE updates ENDS)
66   !-----------------------------------------------------------------------
67   subroutine sectional_interface_1( &
68        idiagbb_host, &
69        id, ii, jj, kk, ktau, ktauc, dtchem, &
70        rbox_sv1, rbox, &
71        it_mosaic, jaerosolstate_bgn, jaerosolstate, &
72        jhyst_leg, jhyst_leg_sv1, &
73        dp_dry_a, dp_dry_a_sv1, &
74        sigmag_a, sigmag_a_sv1, &
75        mass_dry_a_bgn, mass_dry_a, dens_dry_a_bgn, dens_dry_a, &
76        fact_apmassmr, fact_apnumbmr, fact_apdens_in, fact_apdiam_in, fact_gasmr, &
77        pr_atm, rh, te, cair_mol_m3, cair_mol_cc )
79     use module_data_mosaic_boxmod, only: idiag_sect_coag, idiag_sect_movesect, idiag_sect_newnuc
80     use module_data_mosaic_main, only: mw_air, ntot_max, ntot_used
81     use module_data_mosaic_aero, only: &
82        inh4_a, iso4_a, jh2o, &
83        dens_aer_mac, mw_aer_mac, mw_comp_a, &
84        mcoag_flag1, mmovesect_flag1, mnewnuc_flag1, &
85        naer, naercomp, nbin_a, nbin_a_max
86     use module_data_mosaic_asecthp, only:  ai_phase, &
87        hygro_so4_aer, itype_of_itype_md1md2, &
88        maxd_asize, maxd_atype, nsize_aer, ntype_aer, ntype_md2_aer, &
89        xcut_atype_md2
90     use module_mosaic_movesect1d, only:  move_sections_x3
91     use module_mosaic_movesect3d, only:  move_sect_3d_x1
92     use module_mosaic_coag1d,     only:  mosaic_coag_1box
93     use module_mosaic_coag3d,     only:  mosaic_coag_3d_1box
94     use module_mosaic_newnucb,    only:  mosaic_newnuc_1box
96     implicit none
98     !   subr arguments
99     integer, intent(in) :: idiagbb_host
100     integer, intent(in) :: id, ii, jj, kk, ktau, ktauc, it_mosaic
101     integer, intent(in),    dimension(nbin_a_max) :: jaerosolstate, jaerosolstate_bgn
102     integer, intent(in),    dimension(nbin_a_max) :: jhyst_leg_sv1
103     integer, intent(inout), dimension(nbin_a_max) :: jhyst_leg
105     real(r8), intent(in)    :: dtchem  ! time step [s]
106     real(r8), intent(in)    :: pr_atm  ! pressure [atmospheres]
107     real(r8), intent(in)    :: rh      ! relative humidity [percent]
108     real(r8), intent(in)    :: te      ! temperature [K]
109     real(r8), intent(in)    :: cair_mol_m3, cair_mol_cc  ! air molar densities [mol-air/m3-air] and [mol-air/cm3-air]
110     real(r8), intent(in)    :: fact_apmassmr, fact_apnumbmr, fact_apdens_in, fact_apdiam_in, fact_gasmr
111     real(r8), intent(in)    :: rbox_sv1(ntot_used)
112     real(r8), intent(in),    dimension(nbin_a_max) :: dp_dry_a_sv1, sigmag_a_sv1
113     real(r8), intent(inout) :: rbox(ntot_used)  ! gases = [ppmv],  aero mass/water = [ug-aero/kg-air],  aero numb = [#/kg-air]
114     real(r8), intent(inout), dimension(nbin_a_max) :: dens_dry_a_bgn, dens_dry_a  ! [g-aero/cm3-aero]
115     real(r8), intent(inout), dimension(nbin_a_max) :: mass_dry_a_bgn, mass_dry_a  ! [g-aero/cm3-air]
116     real(r8), intent(inout), dimension(nbin_a_max) :: dp_dry_a, sigmag_a          ! [cm] and [--]
118     !   local variables
119     integer  :: idiagbb_coag, istat_coag
120     integer  :: idiagbb_newnuc, istat_newnuc, itype_newnuc, it2
121     integer  :: idiagbb_movesect
122     integer  :: jhyst_leg_sv2(nbin_a_max)
123     integer  :: jsv
124     integer  :: method_coag, method_movesect
126     logical :: ldiag1, ldiag2
127     logical :: ldiag_movesect, ldiag_newnuc, ldiag_coag
129     real(r8) :: dens_nh4so4a_newnuc
130     real(r8) :: fact_apdiam, fact_apdens
131     real(r8) :: rh_box, rhoair_g_cc
132     real(r8) :: tmpa
134     real(r8) :: cnn(ntot_used), cnn_sv1(ntot_used), cnn_sv2(ntot_used)
135     real(r8) :: rbox_sv2(ntot_used)
136     real(r8) :: rbox_svx(ntot_used,4)
138     real(r8) :: dp_dry_a_sv2(nbin_a_max)
139     real(r8) :: sigmag_a_sv2(nbin_a_max)
141     real(r8) :: drydens_pregrow(maxd_asize,maxd_atype)
142     real(r8) :: drydens_aftgrow(maxd_asize,maxd_atype)
143     real(r8) :: drymass_pregrow(maxd_asize,maxd_atype)
144     real(r8) :: drymass_aftgrow(maxd_asize,maxd_atype)
145     ! drydens_xxxgrow = dry density (g/cm3) of bin
146     ! drymass_xxxgrow = all-component dry-mass mixing ratio (g-AP/m3-air) of bin
147     !    pregrow (aftgrow) means before (after) the last call to the
148     !    gas-aerosol mass-transfer solver
149     ! these arrays are used by the movesect routine to calculate transfer
150     !    between bins due to particle growth/shrinkage
152     real(r8) :: adrydens_box(maxd_asize,maxd_atype)
153     real(r8) :: adrydpav_box(maxd_asize,maxd_atype)
154     real(r8) :: adryqmas_box(maxd_asize,maxd_atype)
155     real(r8) :: awetdens_box(maxd_asize,maxd_atype)
156     real(r8) :: awetdpav_box(maxd_asize,maxd_atype)
157     ! adrydens_box = current dry density (g/cm3) of bin
158     ! adrydpav_box = current mean dry diameter (cm) of bin
159     ! adryqmas_box = current all-component dry-mass mixing ratio (g-AP/m3-air) of bin
160     ! awetdens_box = current wet/ambient density (g/cm3) of bin
161     ! awetdpav_box = current mean wet/ambient diameter (cm) of bin
162     !
163     real(r8) :: adrydens_box_sv(maxd_asize,maxd_atype,4)
164     real(r8) :: adrydpav_box_sv(maxd_asize,maxd_atype,4)
165     real(r8) :: adryqmas_box_sv(maxd_asize,maxd_atype,4)
166     real(r8) :: awetdens_box_sv(maxd_asize,maxd_atype,4)
167     real(r8) :: awetdpav_box_sv(maxd_asize,maxd_atype,4)
170     if (idiagbb_host >= 100) then
171        ldiag1 = .true. ; ldiag2 = .true.
172        ldiag_movesect = .true. ; ldiag_newnuc = .true. ; ldiag_coag = .true. 
173     else
174        ldiag1 = .false. ; ldiag2 = .false.
175        ldiag_movesect = .false. ; ldiag_newnuc = .false. ; ldiag_coag = .false. 
176     end if
178     jhyst_leg_sv2(1:nbin_a) = jhyst_leg(1:nbin_a)
179     dp_dry_a_sv2( 1:nbin_a) = dp_dry_a( 1:nbin_a)
180     sigmag_a_sv2( 1:nbin_a) = sigmag_a( 1:nbin_a)
181     rbox_sv2(1:ntot_used) = rbox(1:ntot_used)
182     fact_apdiam = fact_apdiam_in
183     fact_apdens = fact_apdens_in
186     ! the movesect routine puts "initial" values into the axxxxxxx_box arrays
187     !    so need to have movesect turned on (mmovesect_flag1 > 0)
188     !    for any of the sectional routines to function properly
189     if (mmovesect_flag1 <= 0) then
190        write(*,*)   &
191             '*** skipping sectional_interface_1 -- mmovesect_flag1 <= 0'
192        return
193     end if
195     ! only iunits_flagaa = 4 has been tested, so do not allow other values
196     if (iunits_flagaa /= 4) then
197        call wrf_error_fatal('*** sectional_interface_1 error -- iunits_flagaa /= 4')       
198     end if
199     if (iunits_flagbb /= 1) then
200        call wrf_error_fatal('*** sectional_interface_1 error -- iunits_flagbb /= 1')       
201     end if
202     if (iunits_diagout_flagaa /= 1) then
203        call wrf_error_fatal('*** sectional_interface_1 error -- iunits_diagout_flagaa /= 1')
204     end if
207     ! rbox units factors
208     !
209     ! note that the movesect and coag routines expect the aerosol mass units
210     !    for rbox to be mass mixing ratios or mass concentrations
211     ! they do not apply any molecular weight factors, so they cannot
212     !    work with molar mixing ratios or molar concentrations
213     !
214     ! for trace gases,   rbox = [umol/mol-air = ppmv],     rbox*fact_gasmr = [molecules/cm^3-air]
215     !    fact_gasmr = 1.0e-6_r8
216     ! for aerosol mass,  rbox = [ug-AP/kg-air], rbox*fact_apmassmr = [g-AP/g-air]
217     !    fact_apmassmr = 1.0e-9_r8
218     ! for aerosol numb,  rbox = [#/kg-air],     rbox*fact_apnumbmr = [#/g-air]
219     !    fact_apnumbmr = 1.0e-3_r8
220     ! aerosol densities in mosaic aerchem and sectional routines = [g/cm^3]
221     !     (mosaic densities)*fact_apdens = [g/cm^3]
222     !    fact_apdens = 1.0_r8
223     ! aerosol diameters in mosaic aerchem and sectional routines = [cm]
224     !     (mosaic diameters)*fact_apdiam = [cm]
225     !    fact_apdiam = 1.0_r8
228     dens_nh4so4a_newnuc = dens_aer_mac(iso4_a)
231     ! do some set-up tasks
232     if ( ldiag2 ) call dump_secti( 'dym', &
233        id, ii, jj, kk, ktau, ktauc, it_mosaic, dtchem, &
234        rbox, cair_mol_m3 )
235     if ( ldiag1 ) write(*,*) 'sectiface call sect_iface_bgn_end_calcs( 0 )'
236     call sect_iface_bgn_end_calcs( 0, rbox, &
237          dp_dry_a, sigmag_a, &
238          drydens_aftgrow, drydens_pregrow,   &
239          drymass_aftgrow, drymass_pregrow,   &
240          adrydens_box, awetdens_box, adrydpav_box, awetdpav_box,   &
241          adryqmas_box,   &
242          fact_apmassmr, fact_apnumbmr, fact_apdens, fact_apdiam, fact_gasmr, &
243          it_mosaic, jaerosolstate_bgn, jaerosolstate, jhyst_leg, &
244          mass_dry_a_bgn, mass_dry_a, dens_dry_a_bgn, dens_dry_a, &
245          mw_aer_mac, cair_mol_cc, cair_mol_m3)
247     jsv = 1
248     rbox_svx(:,jsv) = rbox(:)
249     adrydens_box_sv(:,:,jsv) = adrydens_box(:,:)
250     awetdens_box_sv(:,:,jsv) = awetdens_box(:,:)
251     adrydpav_box_sv(:,:,jsv) = adrydpav_box(:,:)
252     awetdpav_box_sv(:,:,jsv) = awetdpav_box(:,:)
253     adryqmas_box_sv(:,:,jsv) = adryqmas_box(:,:)
256     ! move_sections transfers particles between sections
257     ! following condensational growth / evaporative shrinking
258     !
259     method_movesect = mod( max(0,mmovesect_flag1), 100 )
260     idiagbb_movesect = 0
261     if ( ldiag_movesect .or. idiag_sect_movesect > 0) idiagbb_movesect = 70
262     idiagbb_movesect = 70 ! *** temporary for testing
263     if ( ii*jj*kk /= 1 ) idiagbb_movesect = 0
264     if (method_movesect < 50) then
265        if ( ldiag2 ) call dump_secti( 'dyn', &
266           id, ii, jj, kk, ktau, ktauc, it_mosaic, dtchem, &
267           rbox, cair_mol_m3 )
268        if ( ldiag1 ) write(*,*) 'sectiface call move_sections_x3'
269 !      call move_sections_x3( iflag, iclm, jclm, k, m, rbox,   &
270 !           fact_apmassmr, fact_apnumbmr,   &
271 !           fact_apdens, fact_apdiam,   &
272 !           drydens_aftgrow, drydens_pregrow,   &
273 !           drymass_aftgrow, drymass_pregrow,   &
274 !           adrydens_box, awetdens_box, adrydpav_box, awetdpav_box,   &
275 !           adryqmas_box, it_mosaic, mmovesect_flag1, idiag_sect_movesect )
276        idiagbb_movesect = 0
277        if ( ldiag_movesect .and. ii*jj*kk == 1 ) idiagbb_movesect = 70
278        call move_sections_x3( ai_phase, ii, jj, kk, 1, rbox,   &
279             fact_apmassmr, fact_apnumbmr,   &
280             fact_apdens, fact_apdiam,   &
281             drydens_aftgrow, drydens_pregrow,   &
282             drymass_aftgrow, drymass_pregrow,   &
283             adrydens_box, awetdens_box, adrydpav_box, awetdpav_box,   &
284             adryqmas_box, it_mosaic, mmovesect_flag1, idiagbb_movesect )
285     else
286        if ( ldiag2 ) call dump_secti( 'dyo', &
287           id, ii, jj, kk, ktau, ktauc, it_mosaic, dtchem, &
288           rbox, cair_mol_m3 )
289        if ( ldiag1 ) write(*,*) 'sectiface call move_sect_3d_x1'
290        call move_sect_3d_x1( 1, 1, 1, 1, 1, rbox,   &
291             fact_apmassmr, fact_apnumbmr,   &
292             fact_apdens, fact_apdiam,   &
293             drydens_aftgrow, drydens_pregrow,   &
294             drymass_aftgrow, drymass_pregrow,   &
295             adrydens_box, awetdens_box, adrydpav_box, awetdpav_box,   &
296             adryqmas_box, it_mosaic, mmovesect_flag1, idiagbb_movesect )
297     end if
299     jsv = 2
300     rbox_svx(:,jsv) = rbox(:)
301     adrydens_box_sv(:,:,jsv) = adrydens_box(:,:)
302     awetdens_box_sv(:,:,jsv) = awetdens_box(:,:)
303     adrydpav_box_sv(:,:,jsv) = adrydpav_box(:,:)
304     awetdpav_box_sv(:,:,jsv) = awetdpav_box(:,:)
305     adryqmas_box_sv(:,:,jsv) = adryqmas_box(:,:)
308     ! do new particle nucleation
309     !
310     if (mnewnuc_flag1 > 0) then
311        idiagbb_newnuc = 0
312        if ( ldiag_newnuc .or. idiag_sect_newnuc > 0) idiagbb_newnuc = 200
313        if ( ii*jj*kk /= 1 ) idiagbb_newnuc = 0
315        if (ntype_aer == 1) then
316           itype_newnuc = 1
317        else
318           ! newnuc particles go into first wbc bin and kappa bin holding so4 kappa
319           itype_newnuc = 0
320           do it2 = 1, ntype_md2_aer
321              if (hygro_so4_aer <= xcut_atype_md2(it2)) then
322                 itype_newnuc = itype_of_itype_md1md2(1,it2)
323                 exit
324              end if
325           end do
326           if (itype_newnuc == 0) itype_newnuc = itype_of_itype_md1md2(1,ntype_md2_aer)
327        end if
329        rh_box = rh*0.01
330        if ( ldiag2 ) &
331        call dump_secti( 'dyp', id, ii, jj, kk, ktau, ktauc, it_mosaic, dtchem, &
332           rbox, cair_mol_m3 )
333        if ( ldiag1 ) write(*,*) 'sectiface call mosaic_newnuc_1box'
334        call mosaic_newnuc_1box(   &
335             istat_newnuc, idiagbb_newnuc,   &
336             it_mosaic, itype_newnuc, dtchem,   &
337             te, pr_atm, cair_mol_cc, rh_box,   &
338             dens_nh4so4a_newnuc, mw_aer_mac(iso4_a), mw_aer_mac(inh4_a),   &
339             mw_comp_a(jh2o), mw_air,  &
340             fact_apmassmr, fact_apnumbmr, fact_apdens, fact_apdiam, fact_gasmr,   &
341             rbox_sv1, rbox,   &
342             adrydens_box, awetdens_box,   &
343             adrydpav_box, awetdpav_box, adryqmas_box )
344     end if
346     jsv = 3
347     rbox_svx(:,jsv) = rbox(:)
348     adrydens_box_sv(:,:,jsv) = adrydens_box(:,:)
349     awetdens_box_sv(:,:,jsv) = awetdens_box(:,:)
350     adrydpav_box_sv(:,:,jsv) = adrydpav_box(:,:)
351     awetdpav_box_sv(:,:,jsv) = awetdpav_box(:,:)
352     adryqmas_box_sv(:,:,jsv) = adryqmas_box(:,:)
355     ! do particle coagulation
356     !
357     if (mcoag_flag1 > 0) then
358        idiagbb_coag = 0
359        if ( ldiag_coag .or. idiag_sect_coag > 0) idiagbb_coag = 200
360        if ( ii*jj*kk /= 1 ) idiagbb_coag = 0
362        rhoair_g_cc = cair_mol_cc*mw_air
363        method_coag = mod( max(0,mcoag_flag1), 100 )
364        if (method_coag < 50) then
365           if ( ldiag2 ) call dump_secti( 'dyq', &
366              id, ii, jj, kk, ktau, ktauc, it_mosaic, dtchem, &
367              rbox, cair_mol_m3 )
368           if ( ldiag1 ) write(*,*) 'sectiface call mosaic_coag_1box'
369           call mosaic_coag_1box( istat_coag,   &
370                idiagbb_coag, it_mosaic,   &
371                dtchem, te, pr_atm, rhoair_g_cc, rbox,   &
372                fact_apmassmr, fact_apnumbmr, fact_apdens, fact_apdiam,   &
373                adrydens_box, awetdens_box,   &
374                adrydpav_box, awetdpav_box, adryqmas_box, mcoag_flag1 )
375        else
376           if ( ldiag2 ) call dump_secti( 'dyr', &
377              id, ii, jj, kk, ktau, ktauc, it_mosaic, dtchem, &
378              rbox, cair_mol_m3 )
379           if ( ldiag1 ) write(*,*) 'sectiface call mosaic_coag_3d_1box'
380           call mosaic_coag_3d_1box( istat_coag,   &
381                idiagbb_coag, it_mosaic,   &
382                dtchem, te, pr_atm, rhoair_g_cc, rbox,   &
383                fact_apmassmr, fact_apnumbmr, fact_apdens, fact_apdiam,   &
384                adrydens_box, awetdens_box,   &
385                adrydpav_box, awetdpav_box, adryqmas_box, mcoag_flag1 )
386        end if
387     end if
389     jsv = 4
390     rbox_svx(:,jsv) = rbox(:)
391     adrydens_box_sv(:,:,jsv) = adrydens_box(:,:)
392     awetdens_box_sv(:,:,jsv) = awetdens_box(:,:)
393     adrydpav_box_sv(:,:,jsv) = adrydpav_box(:,:)
394     awetdpav_box_sv(:,:,jsv) = awetdpav_box(:,:)
395     adryqmas_box_sv(:,:,jsv) = adryqmas_box(:,:)
398     ! do some clean-up tasks
399     if ( ldiag2 ) call dump_secti( 'dys', &
400        id, ii, jj, kk, ktau, ktauc, it_mosaic, dtchem, &
401        rbox, cair_mol_m3 )
402     if ( ldiag1 ) write(*,*) 'sectiface call sect_iface_bgn_end_calcs( 1 )'
403     call sect_iface_bgn_end_calcs( 1, rbox, &
404          dp_dry_a, sigmag_a, &
405          drydens_aftgrow, drydens_pregrow,   &
406          drymass_aftgrow, drymass_pregrow,   &
407          adrydens_box, awetdens_box, adrydpav_box, awetdpav_box,   &
408          adryqmas_box,   &
409          fact_apmassmr, fact_apnumbmr, fact_apdens, fact_apdiam, fact_gasmr, &
410          it_mosaic, jaerosolstate_bgn, jaerosolstate, jhyst_leg, &
411          mass_dry_a_bgn, mass_dry_a, dens_dry_a_bgn, dens_dry_a, &
412          mw_aer_mac, cair_mol_cc, cair_mol_m3)
415     if ( ldiag2 ) call dump_secti( 'dyt', &
416        id, ii, jj, kk, ktau, ktauc, it_mosaic, dtchem, &
417        rbox, cair_mol_m3 )
418     if ( ldiag1 ) write(*,*) 'sectiface done'
419     return
420   end subroutine sectional_interface_1
423   !-----------------------------------------------------------------------
424   subroutine sect_iface_bgn_end_calcs( ibgn_end, rbox, &
425        dp_dry_a, sigmag_a, &
426        drydens_aftgrow, drydens_pregrow,   &
427        drymass_aftgrow, drymass_pregrow,   &
428        adrydens_box, awetdens_box, adrydpav_box,   &
429        awetdpav_box, adryqmas_box,   &
430        fact_apmassmr, fact_apnumbmr, fact_apdens, fact_apdiam, fact_gasmr, &
431        it_mosaic, jaerosolstate_bgn, jaerosolstate, jhyst_leg, &
432        mass_dry_a_bgn, mass_dry_a, dens_dry_a_bgn, dens_dry_a, &
433        mw_aer_mac, cair_mol_cc, cair_mol_m3)
434     !
435     !   does some set-up and clean-up tasks
436     !
437     use module_data_mosaic_boxmod, only: idiag_sect_movesect
438     use module_data_mosaic_main, only: &
439        avogad, &
440        mw_air, naer_tot, ngas_max, ntot_max, ntot_used
441     use module_data_mosaic_aero, only: &
442        jhyst_lo, jhyst_up, &
443        mhyst_method, mhyst_force_lo, mhyst_force_up, &
444        mhyst_uporlo_jhyst, mhyst_uporlo_waterhyst, &
445        mmovesect_flag1, naer, naercomp, nbin_a, nbin_a_max, &
446        all_solid, all_liquid, mixed, no_aerosol
447     use module_data_mosaic_asecthp, only: &
448        ai_phase, isize_of_ibin, itype_of_ibin, lunerr, &
449        maxd_acomp, maxd_asize, maxd_atype, ncomp_aer, &
450        hyswptr_aer, massptr_aer, &
451        dens_aer, dens_water_aer, hygro_aer 
453     use module_mosaic_movesect1d, only:  test_move_sections
455     use module_peg_util, only:   peg_error_fatal
457     implicit none
459     ! subr parameters
460     integer, intent(in) :: ibgn_end, it_mosaic
461     integer, intent(in), dimension(nbin_a_max) :: jaerosolstate, jaerosolstate_bgn
462     integer, intent(inout), dimension(nbin_a_max) :: jhyst_leg
464     real(r8), intent(in) :: cair_mol_cc, cair_mol_m3
465     real(r8), intent(inout), dimension(nbin_a_max)  :: mass_dry_a_bgn, mass_dry_a, dens_dry_a_bgn, dens_dry_a
466     real(r8), intent(inout), dimension(nbin_a_max)  :: dp_dry_a, sigmag_a
467     real(r8), intent(inout), dimension(naer) :: mw_aer_mac
468     real(r8), intent(inout) :: rbox(ntot_used)
470     real(r8), intent(inout) :: drydens_aftgrow(maxd_asize,maxd_atype)
471     real(r8), intent(inout) :: drydens_pregrow(maxd_asize,maxd_atype)
472     real(r8), intent(inout) :: drymass_aftgrow(maxd_asize,maxd_atype)
473     real(r8), intent(inout) :: drymass_pregrow(maxd_asize,maxd_atype)
475     real(r8), intent(inout) :: adrydens_box(maxd_asize,maxd_atype)
476     real(r8), intent(inout) :: awetdens_box(maxd_asize,maxd_atype)
477     real(r8), intent(inout) :: adrydpav_box(maxd_asize,maxd_atype)
478     real(r8), intent(inout) :: awetdpav_box(maxd_asize,maxd_atype)
479     real(r8), intent(inout) :: adryqmas_box(maxd_asize,maxd_atype)
481     real(r8), intent(in) :: fact_apmassmr, fact_apnumbmr, &
482                             fact_apdens, fact_apdiam, fact_gasmr
485     ! local variables
486     integer :: ibin, iphase, isize, itype, jhyst_tmp
487     integer :: l, ll, mtmp
488     real(r8) :: tmpa, tmph, tmpj, tmpr
489     character(len=200) :: msg
493     ! on first time step, call test_move_sections
494     ! note that test_move_sections only executes when mmovesect_flag1 = 80nn
495     if ((it_mosaic == 1) .and. (ibgn_end == 0)) then
496        call test_move_sections( 1, 1, 1, 1, 1, it_mosaic, mmovesect_flag1, idiag_sect_movesect )
497     end if
500     if ((ibgn_end < 0) .or. (ibgn_end > 1)) return
502     if (ibgn_end == 1) goto 20000
505     !
506     ! ibgn_end = 0
507     !     initialize the a-------_box, dry----_pregrow, and dry----aftgrow 
508     !     for (mhyst_method == mhyst_uporlo_jhyst),
509     !         calc/load rbox(hyswptr_aer) with "hysteresis surrogate" info
510     !
511     adrydens_box(:,:) = 0.0
512     awetdens_box(:,:) = 0.0
513     adrydpav_box(:,:) = 0.0
514     awetdpav_box(:,:) = 0.0
515     adryqmas_box(:,:) = 0.0
517     do ibin = 1, nbin_a
518        isize = isize_of_ibin(ibin)
519        itype = itype_of_ibin(ibin)
521        if (jaerosolstate_bgn(ibin) /= no_aerosol) then
522           drydens_pregrow(isize,itype) = dens_dry_a_bgn(ibin)
523           drymass_pregrow(isize,itype) = mass_dry_a_bgn(ibin)
524        else
525           drydens_pregrow(isize,itype) = -1.0
526           drymass_pregrow(isize,itype) =  0.0
527        end if
528        if (jaerosolstate(ibin) /= no_aerosol) then
529           drydens_aftgrow(isize,itype) = dens_dry_a(ibin)
530           drymass_aftgrow(isize,itype) = mass_dry_a(ibin)
531        else
532           drydens_aftgrow(isize,itype) = -1.0
533           drymass_aftgrow(isize,itype) =  0.0
534        end if
536        ! convert drymass_xxxgrow from mosaic aerchemistry units [g-aero/cm^3-air] 
537        !    to rbox mass mixing ratio units [ug-aero/kg-air]
538        tmpa = 1.0e6_r8/(cair_mol_m3*mw_air*fact_apmassmr)
539        drymass_pregrow(isize,itype) = drymass_pregrow(isize,itype) * tmpa
540        drymass_aftgrow(isize,itype) = drymass_aftgrow(isize,itype) * tmpa
542        ! convert drydens_xxxgrow from [g-aero/cm3-aero] to [???] 
543        !    (currently this is not really needed because fact_apdens = 1.0)
544        tmpa = 1.0_r8/fact_apdens
545        drydens_pregrow(isize,itype) = drydens_pregrow(isize,itype) * tmpa
546        drydens_aftgrow(isize,itype) = drydens_aftgrow(isize,itype) * tmpa
548        if (mhyst_method == mhyst_uporlo_jhyst) then
549           ! this not allowed in wrf-chem 
550           msg = '*** sect_iface_bgn_end_calcs - mhyst_uporlo_jhyst not allowed in wrf-chem'
551           call peg_error_fatal( lunerr, msg )
552        end if
554     end do
556     return
559 20000 continue
560     !
561     ! ibgn_end = 1
562     !     load dp_dry_a with current value
563     !     for (mhyst_method == mhyst_uporlo_jhyst),
564     !         calc/set jhyst_leg using the info in rbox(hyswptr_aer)
565     !     for (mhyst_method == other), set jhyst_leg as needed
566     !
567     do ibin = 1, nbin_a
568        isize = isize_of_ibin(ibin)
569        itype = itype_of_ibin(ibin)
571        ! need to put final value into dp_dry_a here
572        dp_dry_a(ibin) = adrydpav_box(isize,itype)
574        if (mhyst_method == mhyst_uporlo_waterhyst) then
575           continue
576        else if (mhyst_method == mhyst_force_up) then
577           jhyst_leg(ibin) = jhyst_up
578        else if (mhyst_method == mhyst_force_lo) then
579           jhyst_leg(ibin) = jhyst_lo
580        else if (mhyst_method == mhyst_uporlo_jhyst) then
581           ! this not allowed in wrf-chem 
582           msg = '*** sect_iface_bgn_end_calcs - mhyst_uporlo_jhyst not allowed in wrf-chem'
583           call peg_error_fatal( lunerr, msg )
584        else
585           write(msg,'(a,i10)') &
586                '*** sect_iface_bgn_end_calcs - bad mhyst_method =', mhyst_method
587           call peg_error_fatal( lunerr, msg )
588        end if
590     end do
592     return
593   end subroutine sect_iface_bgn_end_calcs
596   !-----------------------------------------------------------------------
597   subroutine dump_secti( dumpchaa, &
598        id, ii, jj, kk, ktau, ktauc, it_mosaic, dtchem, &
599        rbox, cair_mol_m3 )
601     use module_data_mosaic_main, only: mw_air, ntot_max, ntot_used
602     use module_data_mosaic_aero, only: &
603        inh4_a, iso4_a, jh2o, &
604        dens_aer_mac, mw_aer_mac, mw_comp_a, &
605        mcoag_flag1, mmovesect_flag1, mnewnuc_flag1, &
606        naer, naercomp, nbin_a, nbin_a_max
607     use module_data_mosaic_asecthp, only:  maxd_asize, maxd_atype, &
608        nsize_aer, ntype_aer, &
609        numptr_aer, lptr_bc_aer, lptr_oc_aer, lptr_oin_aer
611     implicit none
613     !   subr arguments
614     character(len=*), intent(in) :: dumpchaa
616     integer, intent(in) :: id, ii, jj, kk, ktau, ktauc, it_mosaic
618     real(r8), intent(in)    :: dtchem  ! time step [s]
619     real(r8), intent(in)    :: cair_mol_m3  ! air molar densities [mol-air/m3-air]
620     real(r8), intent(inout) :: rbox(ntot_used)  ! gases = [ppmv],  aero mass/water = [ug-aero/kg-air],  aero numb = [#/kg-air]
622     !   local variables
623     integer :: isize, itype
624     real(r8) :: tmpa
625     character(len=80) :: fmtaa
627     if (ii*jj*kk /= 1) return
629     fmtaa = '(2a,1p,e10.3,7e11.3)'
630     if (nsize_aer(1)*ntype_aer > 8) fmtaa = '(2a,1p,e10.3,7e11.3/(7x,1p,e10.3,7e11.3))'
632 ! number output units = #/cm3 = (rbox #/kg-air)*(airden kg-air/m3)*(1e-6 m3/cm3)
633     tmpa = (cair_mol_m3*mw_air*1.0e-3_r8)*1.0e-6_r8
634     write(131,fmtaa) dumpchaa, ' num', &
635        ( ( rbox(numptr_aer(isize,itype,1))*tmpa, &
636            isize=1,nsize_aer(itype) ), itype=1,ntype_aer )
638 ! mass output units = nmol/m3 = (rbox ug/kg-air)*(airden kg-air/m3)*(1e3 ng/ug)/(molecwght ng/nmol)
639 ! and mosaic molecwght = 1.0 for oin, oc, and bc
640     tmpa = (cair_mol_m3*mw_air*1.0e-3_r8)*1.0e3_r8
641     write(131,fmtaa) dumpchaa, ' oin', &
642        ( ( rbox(lptr_oin_aer(isize,itype,1))*tmpa, &
643            isize=1,nsize_aer(itype) ), itype=1,ntype_aer )
644     write(131,fmtaa) dumpchaa, ' oc ', &
645        ( ( rbox(lptr_oc_aer(isize,itype,1))*tmpa, &
646            isize=1,nsize_aer(itype) ), itype=1,ntype_aer )
647     write(131,fmtaa) dumpchaa, ' bc ', &
648        ( ( rbox(lptr_bc_aer(isize,itype,1))*tmpa, &
649            isize=1,nsize_aer(itype) ), itype=1,ntype_aer )
651     return
652   end subroutine dump_secti
655   !-----------------------------------------------------------------------
658   end module module_mosaic_sect_intr