Update version info for release v4.6.1 (#2122)
[WRF.git] / chem / module_mosaic_soa_vbs.F
blob4f5ffe93b28bbfcc42e13caf9368767d77109e28
1         module module_mosaic_soa_vbs
4         use module_data_mosaic_kind, only: r8
6         implicit none
9         contains
12 !-------------------------------------------------------------------------------
13         subroutine mosaic_soa_vbs_intr( &
14            dtchem, p_atm, t_k, swdown_cell, &
15            jaerosolstate, &
16            aer, gas, water_a, area_wet_a, dp_wet_a, &
17            kg, sat_soa, total_species, &
18            ma, mc, mosaic_vars_aa )
20         use module_data_mosaic_aero, only: &
21            iasoaX_g, ipcg1_b_c_g, ismpa_g, &
22            msoa_vbs_info, &
23            naer, nanion, nbin_a_max, ncation, &
24            ngas_aerchtot, ngas_ioa, ngas_soa, ngas_volatile, &
25            mosaic_vars_aa_type
27 ! arguments
28         integer,  intent(inout), dimension(nbin_a_max)               :: jaerosolstate
30         real(r8) :: dtchem
31         real(r8) :: p_atm
32         real(r8) :: t_k
33         real(r8) :: swdown_cell
35         real(r8), intent(inout), dimension(naer,3,nbin_a_max)        :: aer
36         real(r8), intent(inout), dimension(ngas_aerchtot)            :: gas
37         real(r8), intent(inout), dimension(nbin_a_max)               :: water_a
38         real(r8), intent(inout), dimension(nbin_a_max)               :: area_wet_a
39         real(r8), intent(inout), dimension(nbin_a_max)               :: dp_wet_a
40         real(r8), intent(inout), dimension(ngas_aerchtot,nbin_a_max) :: kg
41         real(r8), intent(inout), dimension(ngas_volatile)            :: sat_soa
42         real(r8), intent(inout), dimension(ngas_volatile)            :: total_species
43         real(r8), intent(inout), dimension(ncation,nbin_a_max)       :: mc
44         real(r8), intent(inout), dimension(nanion,nbin_a_max)        :: ma
46         type (mosaic_vars_aa_type), intent(inout) :: mosaic_vars_aa
48 ! local
49         integer :: ii, jj, kk
50         integer :: start_svoc, nsoa_tmp
51         integer :: vbs_nbin, vbs_uq_aqsoa, vbs_uq_par
54         vbs_nbin = msoa_vbs_info(1)
55         vbs_uq_aqsoa = msoa_vbs_info(2)
56         vbs_uq_par = msoa_vbs_info(3)
57        
58         ii = mosaic_vars_aa%hostgridinfo(2)
59         jj = mosaic_vars_aa%hostgridinfo(3)
60         kk = mosaic_vars_aa%hostgridinfo(4)
61         if (ii*jj*kk == 1) write(*,'(/a,3i10)') &
62            'mosaic_soa_vbs_intr - vbs_nbin, uq_aqsoa, uq+par =', &
63            vbs_nbin, vbs_uq_aqsoa, vbs_uq_par
67 ! *** eventually need to do something about these glysoa subrs
69 !       if (glysoa_param == glysoa_param_simple)  call glysoa_simple(dtchem)
70 !       if (glysoa_param == glysoa_param_complex) call glysoa_complex(dtchem)
72         start_svoc = 1
73         nsoa_tmp       = 0
74         if (vbs_nbin .eq. 0) then
75            ! simple version, Hodzic and Jimenez, GMD, 2011
76            start_svoc = ismpa_g
77            ! 4-bin version, Knote et al., ACPD, 2014
78         else if (vbs_nbin .eq. 4) then
79            start_svoc = iasoaX_g
80         else
81            ! 9-bin version
82            start_svoc = ipcg1_b_c_g
83 !          nsoa_tmp = ngas_volatile-start_svoc
84         end if
85         nsoa_tmp = ngas_ioa + ngas_soa - start_svoc + 1
87 !       call equilibrium(start_svoc,nsoa_tmp)
88         call  equilibrium( start_svoc, nsoa_tmp, &
89            aer, gas, kg, sat_soa, total_species )
92         return
93         end subroutine mosaic_soa_vbs_intr
96 !-------------------------------------------------------------------------------
97 ! Calculates the equilibrium gas-particle partitioning for SOA species
98         subroutine  equilibrium( start_ind, N, &
99            aer, gas, kg, sat_soa, total_species )
100 !       subroutine  equilibrium( start_ind, N )
101 ! This routine was implemented by Manish Shrivastava on 12/24/2009 to do gas-particle partitioning of SOA assuming thermodynamic equilibrium.
102 ! Modified by Alma Hodzic 12/2012 to implement the partitioning for mozart-mosaic species (based on the initial code implemented by Manish Shrivastava and originated from CAMx) 
103 ! This would give MOSAIC cpabilities of running both dynamic and equilibrium gas-particle partitioning
104 ! Calls the subr soap. Subroutine soap calls subr spfcn
105 !        use module_data_mosaic_main
106 !        use module_data_mosaic_aero
108         use module_data_mosaic_aero, only: &
109            ioc_a, jtotal, &
110            naer, nbin_a, nbin_a_max, ngas_aerchtot, ngas_volatile
112         implicit none
114 ! arguments
115         integer, intent(in) :: start_ind, N
116         real(r8), intent(inout), dimension(naer,3,nbin_a_max)        :: aer
117         real(r8), intent(inout), dimension(ngas_aerchtot)            :: gas
118         real(r8), intent(inout), dimension(ngas_aerchtot,nbin_a_max) :: kg
119         real(r8), intent(inout), dimension(ngas_volatile)            :: sat_soa
120         real(r8), intent(inout), dimension(ngas_volatile)            :: total_species
122 ! local variables
123 !       integer, parameter :: N=ngas_soa !Total number of soa species
124         integer, parameter :: itermax=2000
125         integer :: idxfresh(N),idxaged(N)   !counter for fresh and aged soa species
126         integer :: flagsoap(N) ! flagsoap determines if the species 'i' is fresh (flagsoap(i)=2) or aged(flagsoap(i)=1
127         integer :: nsolfresh,nsolaged,ntrack,icontfresh,icontaged ! counters corresponding to fresh and aged species for mapping
128         integer :: ibin,iter ! Bin nos.
129         integer :: iv, jp
130         integer :: i
132         real(r8), parameter :: tinys=1.0d-15
133         real(r8) :: dq,frqfresh(nbin_a),frqaged(nbin_a)
134         real(r8) :: frqtotfresh,frqtotaged,frt
135         real(r8) :: xsumfresh(nbin_a),xsumaged(nbin_a)
136         real(r8) :: mnkfresh,mxkfresh,mnkaged,mxkaged
137 !       real betak
138         real(r8) ::  Csatfresh(N), Ctotfresh(N)
139         real(r8) ::  Cgasfresh(N),Caerfresh(N) ! Csat: Saturation conc., Ctot: Total organic mass
140 !       in gas+aerosol phase, Cgas:gas phase, Caer: Particle
141         real(r8) ::    Csataged(N), Ctotaged(N)
142         real(r8) ::  Cgasaged(N),Caeraged(N)
143         real(r8) :: cpxfresh,cpxaged !Moles of pre-existing fresh and aged particle phase organic mass
144         real(r8) :: dum, sum_dum, sum_soa, small_oc
146 !       real, parameter :: tolmin = 1.E-12
147 !       real, parameter :: conmin = 1.E-20
148 !       real totOA,minitw !total OA in particle phase
149         real(r8) :: cpx !pre-existing OA umol/m3
150         real(r8) :: Ctot(N),Caer(N),Cgas(N),Csat(N)
151         real(r8) :: Paer(ngas_volatile)
154 !       LOGICAL check
155         jp=jtotal
156         iter=0
157          cpxaged=0.0
158         cpxfresh=0.0 ! Assume no pres-existing OA forms a solution
159         nsolfresh=0
160          nsolaged=0
161          icontfresh=0
162          icontaged=0
163          dq=0.0
165 ! Paer holds the organic aerosol values in each volatility bin (sum of all size bins)
166           do iv=1,ngas_volatile
167            Paer(iv)=0.0
168           enddo
170 ! Initialize flagsoap
171           do i=1,N
172              flagsoap(i)=1
173              Ctot(i) = 0.0
174              Ctotaged(i) = 0.0
175              Ctotfresh(i) = 0.0
176              Caer(i) = 0.0
177              Caeraged(i) = 0.0
178              Caerfresh(i) = 0.0
179              Cgas(i) = 0.0
180              Cgasaged(i) = 0.0
181              Cgasfresh(i) = 0.0
182              Csat(i) = 0.0
183              Csataged(i) = 0.0
184              Csatfresh(i) = 0.0
185           enddo
187 ! Calculate Ctot and Paer
188 !              do iv = ipcg1_b_c_g, ngas_volatile
189 !              do iv = start_ind, ngas_ioa + ngas_soa
190               do iv = start_ind, (start_ind + N - 1)
191         total_species(iv) = gas(iv)
192         do ibin = 1, nbin_a
193           total_species(iv) = total_species(iv) + aer(iv,jtotal,ibin)
194            Paer(iv)=Paer(iv)+aer(iv,jtotal,ibin)
195         enddo
196       enddo
198 ! Calculate pre-existing moles of OA (cpx) as sum of all size bins
199         do ibin=1,nbin_a
200         cpxaged= cpxaged+aer(ioc_a,jp,ibin)
201          enddo
203 !  Maps arrays starting from start_ind or ipcg1_b_c_g on to corresponding arrays starting from 1 for just soa species
204         do i=1,N
205            Ctot(i)=total_species(start_ind+i-1)
206            Caer(i)=Paer(start_ind+i-1)
207            Csat(i)=sat_soa(start_ind+i-1)
208            Cgas(i)=gas(start_ind+i-1)
209          enddo
211 ! Initialize mapping array indices
212           do i=1,N
213             idxfresh(i)=0
214             idxaged(i)=0
215           enddo
217 !     Seperate the fresh and aged species and treat them as 2 different solutions. Note this approach differes from PMCAMx
218 !     In PMCAMx if flagsoap(i) was set to zero those species were not considered solution forming.
219          do i=1,N
220             flagsoap(i)=1 ! Biomass burning(carbon and oxygen species) +traditional soa species
221          enddo
223 !         do i=1,9
224 !            flagsoap(i)=1 ! Biomass burning(carbon and oxygen species) +traditional soa species
225 !          enddo
226 !         do i=10,18
227 !            flagsoap(i)=1 ! Biomass burning(carbon and oxygen species) +traditional soa species
228 !          enddo
229 !        do i=19,26
230 !            flagsoap(i)=1 ! Biomass burning(carbon and oxygen species) +traditional soa species
231 !          enddo
232 !        do i=27,34
233 !            flagsoap(i)=1 ! Biomass burning(carbon and oxygen species) +traditional soa species
234 !          enddo
235 !        do i=35,43
236 !            flagsoap(i)=1 ! Biomass burning(carbon and oxygen species) +traditional soa species
237 !          enddo
238 !         do i=44,52
239 !            flagsoap(i)=1 ! Biomass burning(carbon and oxygen species) +traditional soa species
240 !          enddo
241 !        do i=53,60
242 !            flagsoap(i)=1 ! Biomass burning(carbon and oxygen species) +traditional soa species
243 !          enddo
244 !      do i=61,68
245 !            flagsoap(i)=1 ! Biomass burning(carbon and oxygen species) +traditional soa species
246 !          enddo
247 !        do i=69,84
248 !            flagsoap(i)=1 !Oxidized fossil oxygen
249 !        enddo
251       do i=1,N
252          if (flagsoap(i).eq.2) then ! fresh primary species forming 1 solution
253            icontfresh=icontfresh+1  ! count the number of fresh species
254             idxfresh(icontfresh) = i  !Map the species
255             Csatfresh(icontfresh)=Csat(i)
256             Ctotfresh(icontfresh)=Ctot(i)
257             Caerfresh(icontfresh)=Caer(i)
258             Cgasfresh(icontfresh)=Cgas(i)
259             nsolfresh=nsolfresh+1
260          elseif (flagsoap(i).eq.1) then                       ! Aged SOA species forming another solution
261             icontaged=icontaged+1
262             idxaged(icontaged) = i
263             Csataged(icontaged)=Csat(i)
264             Ctotaged(icontaged)=Ctot(i)
265             Caeraged(icontaged)=Caer(i)
266             Cgasaged(icontaged)=Cgas(i)
267             nsolaged=nsolaged+1
268          endif
269       enddo
271 !      Caluclate the initial equilibrium partitioning by the bisection method (CMU PMCAMx approach)
272 !       If all fresh abd aged species form a solution
273 !         call soap(ngas_soa,Ctot,Csat,Caer,Cgas,cpx)
275 !       if fresh and aged species form seperate solutions
276       if (nsolfresh.gt.0)  call soap(nsolfresh,Ctotfresh, &
277                     Csatfresh,Caerfresh,Cgasfresh,cpxfresh)
278       if (nsolaged.gt.0)  call soap(nsolaged,Ctotaged, &
279                   Csataged,Caeraged,Cgasaged,cpxaged)
281 !     Map the fresh and aged species back into original arrays
282 !     Now assign the equilibrium gas-particle partitioning arrays
283         ntrack=0
284        do i=1,N ! Map the fresh and aged species back into array from 1 to N after calculating equilibrium
285          if (idxfresh(i).gt.0) then
286          Caer(idxfresh(i))= Caerfresh(i)
287          Cgas(idxfresh(i))= Cgasfresh(i)
288          Ctot(idxfresh(i))= Ctotfresh(i)
289          ntrack=ntrack+1
290          endif
291          if (idxaged(i).gt.0) then
292          Caer(idxaged(i))= Caeraged(i)
293          Cgas(idxaged(i))= Cgasaged(i)
294          Ctot(idxaged(i))= Ctotaged(i)
295          ntrack=ntrack+1
296          endif
297        enddo
299 !       Check for total number of species
300         if (ntrack.ne.N) then
301            call wrf_error_fatal('Error in mapping fresh and primary species arrays')
302         endif
304 ! From here on distribute the organic aerosol in size bins following Koo et al. 2003 " Integrated approaches to modeling
305 ! the organic and inorganic atmospheric aerosol components"
306 ! The original code from PMCAMx was modified to include 2 solutions for fresh and primary species
307 ! by Manish Shrivastava on 01/11/2010
308 ! Calculate total organic aerosol OA(in nmoles/m3) in each bin for either of fresh and aged aerosols
310          do ibin=1,nbin_a
311            xsumfresh(ibin)=0.0
312            xsumaged(ibin)=0.0
313               xsumaged(ibin)= xsumaged(ibin)+aer(ioc_a,jp,ibin)!Caluclate pre-existing primary in each bin for aged aerosol
314 !         do iv = start_ind, ngas_ioa + ngas_soa
315          do iv = start_ind, (start_ind + N - 1)
316            if (flagsoap(iv-start_ind+1).eq.2) then
317                xsumfresh(ibin)= xsumfresh(ibin)+aer(iv,jtotal,ibin)
318            elseif (flagsoap(iv-start_ind+1).eq.1) then
319               xsumaged(ibin)= xsumaged(ibin)+aer(iv,jtotal,ibin)
320                 elseif (flagsoap(iv-start_ind+1).eq.0) then
321                  write(*,*) 'Error in mapping flagsoap to start_ind'
322            endif
323          enddo
324 !         do iv = ipcg1_b_c_g, ngas_volatile
325 !           if (flagsoap(iv-ipcg1_b_c_g+1).eq.2) then
326 !               xsumfresh(ibin)= xsumfresh(ibin)+aer(iv,jtotal,ibin)
327 !           elseif (flagsoap(iv-ipcg1_b_c_g+1).eq.1) then
328 !              xsumaged(ibin)= xsumaged(ibin)+aer(iv,jtotal,ibin)
329 !                elseif (flagsoap(iv-ipcg1_b_c_g+1).eq.0) then
330 !                 write(*,*) 'Error in mapping flagsoap to ipcg1_b_c_g'
331 !           endif
332 !         enddo
334 ! Give a small non-zero value to xsum if it is zero in the section
335           if (xsumfresh(ibin).eq.0.0) xsumfresh(ibin)=tinys
336           if (xsumaged(ibin).eq.0.0) xsumaged(ibin)=tinys
337         enddo
340 ! Calculate dq as (gas concentration) G(t)-G(t+h):
341 ! Caluclate driving force at previous time step (Cgas,i-XiCsati) for both fresh and aged solutions
342 !          do iv = start_ind, ngas_ioa + ngas_soa
343           do iv = start_ind, (start_ind + N - 1)
344            if (Ctot(iv-start_ind+1).lt.1d-10) goto 120 ! If a given species concentration is too low skip
345             dq=gas(iv)-Cgas(iv-start_ind+1) !Since both fresh and aged species have been remapped to an array going from 1 to N
346 !          do iv = ipcg1_b_c_g, ngas_volatile
347 !           if (Ctot(iv-ipcg1_b_c_g+1).lt.1d-10) goto 120 ! If a given species concentration is too low skip
348 !            dq=gas(iv)-Cgas(iv-ipcg1_b_c_g+1) !Since both fresh and aged species have been remapped to an array going from 1 to N
349            frqtotfresh=0.0d0
350            frqtotaged=0.0d0
351            mnkfresh=0.0d0
352            mnkaged=0.0d0
353            mxkfresh=0.0d0
354            mxkaged=0.0d0
356              do ibin=1,nbin_a
357 ! fraceq(iv,ibin) is calculated as the rate of mass transfer
358 ! The weighting fractions frqfresh(ibin) amd frqaged(ibin) are caluclated assuming mole fractions from previous time step
359 ! This assumtion could be relaxed by iterativetely solving this equation
360            if (flagsoap(iv-start_ind+1).eq.2) then
361               frqfresh(ibin)= kg(iv,ibin)*(gas(iv) & ! replaced fraceq(iv,ibin) by kg(iv,ibin) on 01/19/10
362               -(aer(iv,jtotal,ibin))/xsumfresh(ibin) &
363               *Csat(iv-start_ind+1))
364           endif
366            if (flagsoap(iv-start_ind+1).eq.1) then
367               frqaged(ibin)= kg(iv,ibin)*(gas(iv) & ! replaced fraceq(iv,ibin) by kg(iv,ibin) on 01/19/10
368              -(aer(iv,jtotal,ibin))/xsumaged(ibin) &
369               *Csat(iv-start_ind+1))
370           endif
372 !           if (flagsoap(iv-ipcg1_b_c_g+1).eq.2) then
373 !              frqfresh(ibin)= kg(iv,ibin)*(gas(iv) & ! replaced fraceq(iv,ibin) by kg(iv,ibin) on 01/19/10
374 !              -(aer(iv,jtotal,ibin))/xsumfresh(ibin) &
375 !              *Csat(iv-ipcg1_b_c_g+1))
376 !          endif
378 !           if (flagsoap(iv-ipcg1_b_c_g+1).eq.1) then
379 !              frqaged(ibin)= kg(iv,ibin)*(gas(iv) & ! replaced fraceq(iv,ibin) by kg(iv,ibin) on 01/19/10
380 !             -(aer(iv,jtotal,ibin))/xsumaged(ibin) &
381 !              *Csat(iv-ipcg1_b_c_g+1))
382 !          endif
383             mnkfresh=min(mnkfresh,frqfresh(ibin))
384             mnkaged=min(mnkaged,frqaged(ibin))
386             mxkfresh=max(mxkfresh,frqfresh(ibin))
387             mxkaged=max(mxkaged,frqaged(ibin))
388           enddo ! for ibin
390 !          Repeat code from this point on for aged aerosol species
391             if (flagsoap(iv-start_ind+1).eq.2) then
392 !            if (flagsoap(iv-ipcg1_b_c_g+1).eq.2) then
393 ! Condensation is favored in the next time step in this bin
394           if(dq.gt.0.and.mnkfresh.lt.0.and.mxkfresh.gt.0) then
395              do ibin=1,nbin_a
396                frqfresh(ibin)=max(frqfresh(ibin)-mnkfresh,0.0d0)
397               enddo
398 ! evaporation is favored in the next time step in this bin
399           elseif(dq.lt.0.and.mxkfresh.gt.0.and.mnkfresh.lt.0) then
400               do ibin=1,nbin_a
401               frqfresh(ibin)=min(frqfresh(ibin)-mxkfresh,0.0d0)
402               enddo
403            endif
404            do ibin=1,nbin_a
405             frqtotfresh=frqtotfresh+frqfresh(ibin)
406            enddo
408 ! Re-normalize frqfresh(ibin)
409 ! Additional code to check for frqtotfresh and frqtotaged
410 ! Added by Manish Shrivastava on 02/19/2010
412           do ibin=1,nbin_a
413            frqfresh(ibin)=frqfresh(ibin)/frqtotfresh
414            enddo
416             elseif(flagsoap(iv-start_ind+1).eq.1) then
417 !            elseif(flagsoap(iv-ipcg1_b_c_g+1).eq.1) then
418           if(dq.gt.0.and.mnkaged.lt.0.and.mxkaged.gt.0) then
419              do ibin=1,nbin_a
420                frqaged(ibin)=max(frqaged(ibin)-mnkaged,0.0d0)
421               enddo
422           elseif(dq.lt.0.and.mxkaged.gt.0.and.mnkaged.lt.0) then
423               do ibin=1,nbin_a
424               frqaged(ibin)=min(frqaged(ibin)-mxkaged,0.0d0)
425               enddo
426            endif
428            do ibin=1,nbin_a
429             frqtotaged=frqtotaged+frqaged(ibin)
430            enddo
432            do ibin=1,nbin_a
433            frqaged(ibin)=frqaged(ibin)/frqtotaged
434            enddo
436            endif ! for flagsoap
437 !     Condense all condensing species
438            if(dq.gt.0.0d0) then
440             !  Map the species back into the original MOSAIC arrays
441              do ibin=1,nbin_a
442                  if (flagsoap(iv-start_ind+1).eq.2) then
443                  aer(iv,jtotal,ibin)= aer(iv,jtotal,ibin)+dq*frqfresh(ibin)
444                  endif
445                 if (flagsoap(iv-start_ind+1).eq.1) then
446                 aer(iv,jtotal,ibin)= aer(iv,jtotal,ibin)+dq*frqaged(ibin)
447                 endif
448              enddo
449 ! Set the gas phase species to equilibrium value
450                 gas(iv)=Cgas(iv-start_ind+1)
452 !             do ibin=1,nbin_a
453 !                 if (flagsoap(iv-ipcg1_b_c_g+1).eq.2) then
454 !                 aer(iv,jtotal,ibin)= aer(iv,jtotal,ibin)+dq*frqfresh(ibin)
455 !                 endif
456 !                if (flagsoap(iv-ipcg1_b_c_g+1).eq.1) then
457 !                aer(iv,jtotal,ibin)= aer(iv,jtotal,ibin)+dq*frqaged(ibin)
458 !                endif
459 !             enddo
460 !! Set the gas phase species to equilibrium value
461 !                gas(iv)=Cgas(iv-ipcg1_b_c_g+1)
463 !     Evaporate all evaporating species
464          elseif(dq.lt.0.0d0) then
465             iter=0
466 100         frt=1.0d0
467                do ibin=1,nbin_a
468                    if (flagsoap(iv-start_ind+1).eq.2) then
469 !                   if (flagsoap(iv-ipcg1_b_c_g+1).eq.2) then
470 ! Cannot evaporate more than whats in the bin ie ratio (aer(iv,jtotal,ibin)/dq*frqfresh(ibin)) should be less than equal to 1
471                  if(frqfresh(ibin).gt.0.0d0) &
472          frt=MAX(MIN(aer(iv,jtotal,ibin)/abs(-dq*frqfresh(ibin)),frt),0.0d0)
473 !                  elseif(flagsoap(iv-ipcg1_b_c_g+1).eq.1) then
474                   elseif(flagsoap(iv-start_ind+1).eq.1) then
475                if(frqaged(ibin).gt.0.0d0) &
476          frt=MAX(MIN(aer(iv,jtotal,ibin)/abs(-dq*frqaged(ibin)),frt),0.0d0)
477                   endif ! for flagsoap
478                enddo ! for ibin
482          frqtotfresh=0.0d0
483          frqtotaged=0.0d0
485              do ibin=1,nbin_a
486         if (flagsoap(iv-start_ind+1).eq.2) then
487 !        if (flagsoap(iv-ipcg1_b_c_g+1).eq.2) then
488                aer(iv,jtotal,ibin)= &
489 ! Since dq is negative this is evaporating aerosols
490                MAX(aer(iv,jtotal,ibin)+frt*dq*frqfresh(ibin),0.0d0)
491          if(aer(iv,jtotal,ibin).lt.tinys) frqfresh(ibin)=0.0d0
492               frqtotfresh=frqtotfresh+frqfresh(ibin)
493 !        elseif (flagsoap(iv-ipcg1_b_c_g+1).eq.1) then
494         elseif (flagsoap(iv-start_ind+1).eq.1) then
495                aer(iv,jtotal,ibin)= &
496                MAX(aer(iv,jtotal,ibin)+frt*dq*frqaged(ibin),0.0d0)
497          if(aer(iv,jtotal,ibin).lt.tinys) frqaged(ibin)=0.0d0
498               frqtotaged=frqtotaged+frqaged(ibin)
499          endif ! for flagsoap
500              enddo ! for ibin
502 ! Check if we should evaporate more
503           dq=(1.0d0-frt)*dq
504 !         if (flagsoap(iv-ipcg1_b_c_g+1).eq.2) then
505          if (flagsoap(iv-start_ind+1).eq.2) then
506            if(dq.lt.-1.d-8) then ! check if d-8 is better
507              if(frqtotfresh.gt.tinys) then ! we have sections which are not empty
508               if(iter.le.itermax) then ! check infinite loop
509                 iter = iter + 1
510                 do ibin = 1,nbin_a
511                   frqfresh(ibin) = frqfresh(ibin) / frqtotfresh
512                 enddo ! for ibin
513              goto 100
514             endif ! for iter
515           endif ! frqtotfresh.gt.tinys
516            endif ! dq.lt.-1.d-7
517 !          elseif (flagsoap(iv-ipcg1_b_c_g+1).eq.1) then
518           elseif (flagsoap(iv-start_ind+1).eq.1) then
519            if(dq.lt.-1.d-8) then
520              if(frqtotaged.gt.tinys) then ! we have sections which are not empty
521               if(iter.le.itermax) then ! check infinite loop
522                 iter = iter + 1
523                 do ibin = 1,nbin_a
524                   frqaged(ibin) = frqaged(ibin) / frqtotaged
525                 enddo
526                goto 100
527           endif
528             endif
529             endif
531             ! we need to evaporate more to achieve equilibrium
532             ! but we completely evaporated the species in all sections
533             ! or exceeded itermax
534            endif ! for flagsoap
536 ! now set the gas species concentration conservatively
537 !           gas(iv)=Ctot(iv-ipcg1_b_c_g+1)
538            gas(iv)=Ctot(iv-start_ind+1)
539              do ibin=1,nbin_a
540                gas(iv)=gas(iv)-aer(iv,jtotal,ibin)
541              enddo
542         endif ! if dq.gt.0
544 120       continue
545            enddo ! for iv=start_ind
547        end subroutine equilibrium
550 !--------------------------------------------------------------------------------
551 !---------------------------------
552 ! Calculates the equilibrium gas-particle partitioning for SOA species when MOZART_MOSAIC_4BIN_KPP is used
553 ! This routine was modified by Alma Hodzic based on the initial code implemented by Manish Shrivastava and originated from CAMx 
555 !++ alma - removed the subr  equilibrium_smp
556 !        subroutine  equilibrium_smp
558 !       end subroutine equilibrium_smp
562 !--------------------------------------------------------------------------------
563 !    This subr spfcn calculates the objective function fval to solve gas-particle partitioning of SOA
564 !    subr spfcn is called from within the subr soap
565         subroutine spfcn(N,Ctot,Csat,Ca,cpx,tom,fval)
566 !        use module_data_mosaic_main
567 !        use module_data_mosaic_aero
568       implicit none
569        real(r8):: Ctot(N),Csat(N),Ca(N),tom,fval,cpx
571          integer i,N
572         fval=0.0
573          do i=1, N
574          Ca(i)=Ctot(i)*tom/(tom+Csat(i)/1)! Replace the divisor 1 by Molecular Weights if the units for Csat(i) are ug/m3 or ng/m3
575         fval=fval+Ca(i)/1 ! The divisor is set to 1 as the species are in nmol/m3
576         enddo
577           fval=fval+cpx-tom
578         return
580        end subroutine spfcn
583 !--------------------------------------------------------------------------------
584         subroutine soap(N,Ctot,Csat,Ca,Cgas,cpx)
585 ! SOAP calculates the gas-partitioning of SOA. Adapted from PMCAMx and uses the bisection approach.
586 ! SOAP calls subr spfcn which calculates the objective function for solving gas-particle partitioning
587 !        use module_data_mosaic_main
588 !        use module_data_mosaic_aero
590         real(r8),  parameter :: xtol = 5.0e-5
591           real(r8):: Ctot(N),Csat(N),cpx,Ca(N),Cgas(N)
592           real(r8):: xend,dx,xmid,fend,fmid,sun
593          integer i,N,znum
594         
595          sun=0.0
596           do i=1,N
597             if (Csat(i).gt.0) then
598             sun=sun+Ctot(i)/Csat(i) !If a species does not exist its Csat is zero
599            else
600            endif
601           enddo
602          if(cpx.lt.1e-9.and.sun.le.1.0) then !if ctots for all species are less than corr. csats and cpre is negligible
603            do i=1,N
604              Cgas(i)=Ctot(i)
605              Ca(i)=0.0
606            enddo
607          goto 900
608         endif
610        xend=0.0
611        do i=1, N
612          xend=xend+Ctot(i)/1 ! Replace the divisor 1 by molecular weight if the units of Ctot(i) are in ng/m3 or ug/m3
613          enddo
614         xend=xend+cpx ! total number of moles
615        if (xend.gt.1e-10) then 
616            call spfcn(N,Ctot,Csat,Ca,cpx,xend,fend) ! Calculates the objective function
617         else
618 !          write(*,*) "Total no of moles less than 1e-10 bypassing soap" 
619               goto 100
620       endif
621           if(abs(fend).le.xtol*xend) goto 99 ! Check for tolerance
622           if (fend.gt.0.0) then ! The objective function is supposed to be less than equal to zero
623          write(*,*) "Error in SOAP"
624          goto 50
625         endif
626            dx=xend-cpx
627         do znum=1,200
628         dx=0.5*dx
629          xmid=xend-dx ! Find the midpoint following the bisection approach
630            call spfcn (N,Ctot,Csat,Ca,cpx,xmid,fmid) ! Re-calculate the objective function
631           if(abs(fmid).le.xtol*xmid.or.dx.le.xtol*xmid) goto 100 ! converged
632            if (fmid.lt.0.0) xend=xmid
633          enddo
634 50         call wrf_message("Error in SOAP")
635          call wrf_error_fatal("Error: max number of iterations reached")
636       
638 99     xmid=xend
639 100    continue
640         do i=1, N
641         Ca(i)=min(Ctot(i), Ca(i))
642         Cgas(i)=Ctot(i)-Ca(i)
643        enddo
644 900   continue
645         
647 !     write(*,*) xmid
648      return
650        end subroutine soap
654 !-------------------------------------------------------------------------------
655       subroutine soa_vbs_load_params( ipass )
657 ! sets species pointers (gas, aerosol, electrolyte) and 
658 !    physical parameters (molecular weight, density, ...)
659 !    that are different from those in the mosaic offline box-model code
661       use module_data_mosaic_aero, only: &
662          aer_name, gas_name, &
663          dens_aer_mac, dens_comp_a, &
664          mw_aer_mac, mw_comp_a, mw_gas, &
665          ngas_aerchtot, partial_molar_vol, v_molar_gas
667       use module_data_mosaic_aero, only: &
668          ih2so4_g,     ihno3_g,      ihcl_g,      inh3_g,         &
669          imsa_g,                                                  &
670          iaro1_g,      iaro2_g,      ialk1_g,     iole1_g,        &
671          iapi1_g,      iapi2_g,      ilim1_g,     ilim2_g,        &
672          iso4_a,     ino3_a,     icl_a,     inh4_a,     ico3_a,   &
673          imsa_a,     ina_a,      ica_a,     ioc_a,      ibc_a,    &
674          ioin_a,     iaro1_a,    iaro2_a,   ialk1_a,    iole1_a,  &
675          iapi1_a,    iapi2_a,    ilim1_a,   ilim2_a,              &
676          jnh4so4,    jlvcite,    jnh4hso4,   jnh4no3,    jnh4cl,  &
677          jna2so4,    jna3hso4,   jnahso4,    jnano3,     jnacl,   &
678          jcaso4,     jcano3,     jcacl2,     jcaco3,     jh2so4,  &
679          jhno3,      jhcl,       jhhso4,                          &
680          jnh4msa,    jnamsa,     jcamsa2,    jmsa,                &
681          joc,        jbc,        join,       jaro1,      jaro2,   &
682          jalk1,      jole1,      japi1,      japi2,      jlim1,   &
683          jlim2,      jh2o
684   
685       use module_data_mosaic_aero, only: &
686          in2o5_g,      iclno2_g,                                 &
687          ipcg1_b_c_g,  ipcg2_b_c_g,  ipcg3_b_c_g,  ipcg4_b_c_g,  &
688          ipcg5_b_c_g,  ipcg6_b_c_g,  ipcg7_b_c_g,  ipcg8_b_c_g,  &
689          ipcg9_b_c_g,                                            &
690          ipcg1_b_o_g,  ipcg2_b_o_g,  ipcg3_b_o_g,  ipcg4_b_o_g,  &
691          ipcg5_b_o_g,  ipcg6_b_o_g,  ipcg7_b_o_g,  ipcg8_b_o_g,  &
692          ipcg9_b_o_g,                                            &
693          iopcg1_b_c_g, iopcg2_b_c_g, iopcg3_b_c_g, iopcg4_b_c_g, &
694          iopcg5_b_c_g, iopcg6_b_c_g, iopcg7_b_c_g, iopcg8_b_c_g, &
695          iopcg1_b_o_g, iopcg2_b_o_g, iopcg3_b_o_g, iopcg4_b_o_g, &
696          iopcg5_b_o_g, iopcg6_b_o_g, iopcg7_b_o_g, iopcg8_b_o_g, &
697          ipcg1_f_c_g,  ipcg2_f_c_g,  ipcg3_f_c_g,  ipcg4_f_c_g,  &
698          ipcg5_f_c_g,  ipcg6_f_c_g,  ipcg7_f_c_g,  ipcg8_f_c_g,  &
699          ipcg9_f_c_g,                                            &
700          ipcg1_f_o_g,  ipcg2_f_o_g,  ipcg3_f_o_g,  ipcg4_f_o_g,  &
701          ipcg5_f_o_g,  ipcg6_f_o_g,  ipcg7_f_o_g,  ipcg8_f_o_g,  &
702          ipcg9_f_o_g,                                            &
703          iopcg1_f_c_g, iopcg2_f_c_g, iopcg3_f_c_g, iopcg4_f_c_g, &
704          iopcg5_f_c_g, iopcg6_f_c_g, iopcg7_f_c_g, iopcg8_f_c_g, &
705          iopcg1_f_o_g, iopcg2_f_o_g, iopcg3_f_o_g, iopcg4_f_o_g, &
706          iopcg5_f_o_g, iopcg6_f_o_g, iopcg7_f_o_g, iopcg8_f_o_g, &
707          iant1_c_g,    iant2_c_g,    iant3_c_g,    iant4_c_g,    &
708          iant1_o_g,    iant2_o_g,    iant3_o_g,    iant4_o_g,    &
709          ibiog1_c_g,   ibiog2_c_g,   ibiog3_c_g,   ibiog4_c_g,   &
710          ibiog1_o_g,   ibiog2_o_g,   ibiog3_o_g,   ibiog4_o_g,   &
711          ismpa_g,      ismpbb_g,                                 &
712          iasoa1_g,     iasoa2_g,     iasoa3_g,     iasoa4_g,     &
713          iasoaX_g,                                               &
714          ibsoa1_g,     ibsoa2_g,     ibsoa3_g,     ibsoa4_g,     &
715          ibsoaX_g,                                               &
716          igly,         iho
717 !        iiepox_g,     igly_g
719       use module_data_mosaic_aero, only: &
720 !        idust_a,                                                &
721 !        itr1r1_a,     itr1r2_a,     itr1r3_a,     itr1r4_a,     &
722          ipcg1_b_c_a,  ipcg2_b_c_a,  ipcg3_b_c_a,  ipcg4_b_c_a,  &
723          ipcg5_b_c_a,  ipcg6_b_c_a,  ipcg7_b_c_a,  ipcg8_b_c_a,  &
724          ipcg9_b_c_a,                                            &
725          ipcg1_b_o_a,  ipcg2_b_o_a,  ipcg3_b_o_a,  ipcg4_b_o_a,  &
726          ipcg5_b_o_a,  ipcg6_b_o_a,  ipcg7_b_o_a,  ipcg8_b_o_a,  &
727          ipcg9_b_o_a,                                            &
728          iopcg1_b_c_a, iopcg2_b_c_a, iopcg3_b_c_a, iopcg4_b_c_a, &
729          iopcg5_b_c_a, iopcg6_b_c_a, iopcg7_b_c_a, iopcg8_b_c_a, &
730          iopcg1_b_o_a, iopcg2_b_o_a, iopcg3_b_o_a, iopcg4_b_o_a, &
731          iopcg5_b_o_a, iopcg6_b_o_a, iopcg7_b_o_a, iopcg8_b_o_a, &
732          ipcg1_f_c_a,  ipcg2_f_c_a,  ipcg3_f_c_a,  ipcg4_f_c_a,  &
733          ipcg5_f_c_a,  ipcg6_f_c_a,  ipcg7_f_c_a,  ipcg8_f_c_a,  &
734          ipcg9_f_c_a,                                            &
735          ipcg1_f_o_a,  ipcg2_f_o_a,  ipcg3_f_o_a,  ipcg4_f_o_a,  &
736          ipcg5_f_o_a,  ipcg6_f_o_a,  ipcg7_f_o_a,  ipcg8_f_o_a,  &
737          ipcg9_f_o_a,                                            &
738          iopcg1_f_c_a, iopcg2_f_c_a, iopcg3_f_c_a, iopcg4_f_c_a, &
739          iopcg5_f_c_a, iopcg6_f_c_a, iopcg7_f_c_a, iopcg8_f_c_a, &
740          iopcg1_f_o_a, iopcg2_f_o_a, iopcg3_f_o_a, iopcg4_f_o_a, &
741          iopcg5_f_o_a, iopcg6_f_o_a, iopcg7_f_o_a, iopcg8_f_o_a, &
742          iant1_c_a,    iant2_c_a,    iant3_c_a,    iant4_c_a,    &
743          iant1_o_a,    iant2_o_a,    iant3_o_a,    iant4_o_a,    &
744          ibiog1_c_a,   ibiog2_c_a,   ibiog3_c_a,   ibiog4_c_a,   &
745          ibiog1_o_a,   ibiog2_o_a,   ibiog3_o_a,   ibiog4_o_a,   &
746          ismpa_a,      ismpbb_a,                                 &
747          iasoa1_a,     iasoa2_a,     iasoa3_a,     iasoa4_a,     &
748          iasoaX_a,                                               &
749          ibsoa1_a,     ibsoa2_a,     ibsoa3_a,     ibsoa4_a,     &
750          ibsoaX_a,                                               &
751          iglysoa_r1_a, iglysoa_r2_a, iglysoa_sfc_a, iglysoa_nh4_a, &
752          iglysoa_oh_a
753 !        iiepox_a,     igly_a,       iiepoxos_a,   itetrol_a,    &
754 !        itanv_a,      isopnv_a,     iternv_a,     iseqnv_a,     &
755 !        isianv_a
757       use module_data_mosaic_aero, only: &
758 !        jdust,                                          &
759 !        jtr1r1,     jtr1r2,     jtr1r3,     jtr1r4,     &
760          jpcg1_b_c,  jpcg2_b_c,  jpcg3_b_c,  jpcg4_b_c,  &
761          jpcg5_b_c,  jpcg6_b_c,  jpcg7_b_c,  jpcg8_b_c,  &
762          jpcg9_b_c,                                      &
763          jpcg1_b_o,  jpcg2_b_o,  jpcg3_b_o,  jpcg4_b_o,  &
764          jpcg5_b_o,  jpcg6_b_o,  jpcg7_b_o,  jpcg8_b_o,  &
765          jpcg9_b_o,                                      &
766          jopcg1_b_c, jopcg2_b_c, jopcg3_b_c, jopcg4_b_c, &
767          jopcg5_b_c, jopcg6_b_c, jopcg7_b_c, jopcg8_b_c, &
768          jopcg1_b_o, jopcg2_b_o, jopcg3_b_o, jopcg4_b_o, &
769          jopcg5_b_o, jopcg6_b_o, jopcg7_b_o, jopcg8_b_o, &
770          jpcg1_f_c,  jpcg2_f_c,  jpcg3_f_c,  jpcg4_f_c,  &
771          jpcg5_f_c,  jpcg6_f_c,  jpcg7_f_c,  jpcg8_f_c,  &
772          jpcg9_f_c,                                      &
773          jpcg1_f_o,  jpcg2_f_o,  jpcg3_f_o,  jpcg4_f_o,  &
774          jpcg5_f_o,  jpcg6_f_o,  jpcg7_f_o,  jpcg8_f_o,  &
775          jpcg9_f_o,                                      &
776          jopcg1_f_c, jopcg2_f_c, jopcg3_f_c, jopcg4_f_c, &
777          jopcg5_f_c, jopcg6_f_c, jopcg7_f_c, jopcg8_f_c, &
778          jopcg1_f_o, jopcg2_f_o, jopcg3_f_o, jopcg4_f_o, &
779          jopcg5_f_o, jopcg6_f_o, jopcg7_f_o, jopcg8_f_o, &
780          jant1_c,    jant2_c,    jant3_c,    jant4_c,    &
781          jant1_o,    jant2_o,    jant3_o,    jant4_o,    &
782          jbiog1_c,   jbiog2_c,   jbiog3_c,   jbiog4_c,   &
783          jbiog1_o,   jbiog2_o,   jbiog3_o,   jbiog4_o,   &
784          jsmpa,      jsmpbb,                             &
785          jasoa1,     jasoa2,     jasoa3,     jasoa4,     &
786          jasoaX,                                         &
787          jbsoa1,     jbsoa2,     jbsoa3,     jbsoa4,     &
788          jbsoaX,                                         &
789          jglysoa_r1, jglysoa_r2, jglysoa_sfc, jglysoa_nh4, &
790          jglysoa_oh
791 !        jiepox,     jgly,       jiepoxos,   jtetrol,    &
792 !        jtanv,      jsopnv,     jternv,     jseqnv,     &
793 !        jsianv
795 ! arguments
796       integer, intent(in) :: ipass
799       if (ipass == 1) then
801 !                                                 jdust      = 26 ! insoluble - part of naercomp
802       ipcg1_b_c_a   =   6 ;  ipcg1_b_c_g  =  6 ;  jpcg1_b_c  = 26
803       ipcg2_b_c_a   =   7 ;  ipcg2_b_c_g  =  7 ;  jpcg2_b_c  = 27
804       ipcg3_b_c_a   =   8 ;  ipcg3_b_c_g  =  8 ;  jpcg3_b_c  = 28
805       ipcg4_b_c_a   =   9 ;  ipcg4_b_c_g  =  9 ;  jpcg4_b_c  = 29
806       ipcg5_b_c_a   =  10 ;  ipcg5_b_c_g  = 10 ;  jpcg5_b_c  = 30
807       ipcg6_b_c_a   =  11 ;  ipcg6_b_c_g  = 11 ;  jpcg6_b_c  = 31
808       ipcg7_b_c_a   =  12 ;  ipcg7_b_c_g  = 12 ;  jpcg7_b_c  = 32
809       ipcg8_b_c_a   =  13 ;  ipcg8_b_c_g  = 13 ;  jpcg8_b_c  = 33
810       ipcg9_b_c_a   =  14 ;  ipcg9_b_c_g  = 14 ;  jpcg9_b_c  = 34
811       ipcg1_b_o_a   =  15 ;  ipcg1_b_o_g  = 15 ;  jpcg1_b_o  = 35
812       ipcg2_b_o_a   =  16 ;  ipcg2_b_o_g  = 16 ;  jpcg2_b_o  = 36
813       ipcg3_b_o_a   =  17 ;  ipcg3_b_o_g  = 17 ;  jpcg3_b_o  = 37
814       ipcg4_b_o_a   =  18 ;  ipcg4_b_o_g  = 18 ;  jpcg4_b_o  = 38
815       ipcg5_b_o_a   =  19 ;  ipcg5_b_o_g  = 19 ;  jpcg5_b_o  = 39
816       ipcg6_b_o_a   =  20 ;  ipcg6_b_o_g  = 20 ;  jpcg6_b_o  = 40
817       ipcg7_b_o_a   =  21 ;  ipcg7_b_o_g  = 21 ;  jpcg7_b_o  = 41
818       ipcg8_b_o_a   =  22 ;  ipcg8_b_o_g  = 22 ;  jpcg8_b_o  = 42
819       ipcg9_b_o_a   =  23 ;  ipcg9_b_o_g  = 23 ;  jpcg9_b_o  = 43
821       iopcg1_b_c_a  =  24 ;  iopcg1_b_c_g = 24 ;  jopcg1_b_c = 44
822       iopcg2_b_c_a  =  25 ;  iopcg2_b_c_g = 25 ;  jopcg2_b_c = 45
823       iopcg3_b_c_a  =  26 ;  iopcg3_b_c_g = 26 ;  jopcg3_b_c = 46
824       iopcg4_b_c_a  =  27 ;  iopcg4_b_c_g = 27 ;  jopcg4_b_c = 47
825       iopcg5_b_c_a  =  28 ;  iopcg5_b_c_g = 28 ;  jopcg5_b_c = 48
826       iopcg6_b_c_a  =  29 ;  iopcg6_b_c_g = 29 ;  jopcg6_b_c = 49
827       iopcg7_b_c_a  =  30 ;  iopcg7_b_c_g = 30 ;  jopcg7_b_c = 50
828       iopcg8_b_c_a  =  31 ;  iopcg8_b_c_g = 31 ;  jopcg8_b_c = 51
829       iopcg1_b_o_a  =  32 ;  iopcg1_b_o_g = 32 ;  jopcg1_b_o = 52
830       iopcg2_b_o_a  =  33 ;  iopcg2_b_o_g = 33 ;  jopcg2_b_o = 53
831       iopcg3_b_o_a  =  34 ;  iopcg3_b_o_g = 34 ;  jopcg3_b_o = 54
832       iopcg4_b_o_a  =  35 ;  iopcg4_b_o_g = 35 ;  jopcg4_b_o = 55
833       iopcg5_b_o_a  =  36 ;  iopcg5_b_o_g = 36 ;  jopcg5_b_o = 56
834       iopcg6_b_o_a  =  37 ;  iopcg6_b_o_g = 37 ;  jopcg6_b_o = 57
835       iopcg7_b_o_a  =  38 ;  iopcg7_b_o_g = 38 ;  jopcg7_b_o = 58
836       iopcg8_b_o_a  =  39 ;  iopcg8_b_o_g = 39 ;  jopcg8_b_o = 59
838       ipcg1_f_c_a   =  40 ;  ipcg1_f_c_g  = 40 ;  jpcg1_f_c  = 60
839       ipcg2_f_c_a   =  41 ;  ipcg2_f_c_g  = 41 ;  jpcg2_f_c  = 61
840       ipcg3_f_c_a   =  42 ;  ipcg3_f_c_g  = 42 ;  jpcg3_f_c  = 62
841       ipcg4_f_c_a   =  43 ;  ipcg4_f_c_g  = 43 ;  jpcg4_f_c  = 63
842       ipcg5_f_c_a   =  44 ;  ipcg5_f_c_g  = 44 ;  jpcg5_f_c  = 64
843       ipcg6_f_c_a   =  45 ;  ipcg6_f_c_g  = 45 ;  jpcg6_f_c  = 65
844       ipcg7_f_c_a   =  46 ;  ipcg7_f_c_g  = 46 ;  jpcg7_f_c  = 66
845       ipcg8_f_c_a   =  47 ;  ipcg8_f_c_g  = 47 ;  jpcg8_f_c  = 67
846       ipcg9_f_c_a   =  48 ;  ipcg9_f_c_g  = 48 ;  jpcg9_f_c  = 68
847       ipcg1_f_o_a   =  49 ;  ipcg1_f_o_g  = 49 ;  jpcg1_f_o  = 69
848       ipcg2_f_o_a   =  50 ;  ipcg2_f_o_g  = 50 ;  jpcg2_f_o  = 70
849       ipcg3_f_o_a   =  51 ;  ipcg3_f_o_g  = 51 ;  jpcg3_f_o  = 71
850       ipcg4_f_o_a   =  52 ;  ipcg4_f_o_g  = 52 ;  jpcg4_f_o  = 72
851       ipcg5_f_o_a   =  53 ;  ipcg5_f_o_g  = 53 ;  jpcg5_f_o  = 73
852       ipcg6_f_o_a   =  54 ;  ipcg6_f_o_g  = 54 ;  jpcg6_f_o  = 74
853       ipcg7_f_o_a   =  55 ;  ipcg7_f_o_g  = 55 ;  jpcg7_f_o  = 75
854       ipcg8_f_o_a   =  56 ;  ipcg8_f_o_g  = 56 ;  jpcg8_f_o  = 76
855       ipcg9_f_o_a   =  57 ;  ipcg9_f_o_g  = 57 ;  jpcg9_f_o  = 77
857       iopcg1_f_c_a  =  58 ;  iopcg1_f_c_g = 58 ;  jopcg1_f_c = 78
858       iopcg2_f_c_a  =  59 ;  iopcg2_f_c_g = 59 ;  jopcg2_f_c = 79
859       iopcg3_f_c_a  =  60 ;  iopcg3_f_c_g = 60 ;  jopcg3_f_c = 80
860       iopcg4_f_c_a  =  61 ;  iopcg4_f_c_g = 61 ;  jopcg4_f_c = 81
861       iopcg5_f_c_a  =  62 ;  iopcg5_f_c_g = 62 ;  jopcg5_f_c = 82
862       iopcg6_f_c_a  =  63 ;  iopcg6_f_c_g = 63 ;  jopcg6_f_c = 83
863       iopcg7_f_c_a  =  64 ;  iopcg7_f_c_g = 64 ;  jopcg7_f_c = 84
864       iopcg8_f_c_a  =  65 ;  iopcg8_f_c_g = 65 ;  jopcg8_f_c = 85
865       iopcg1_f_o_a  =  66 ;  iopcg1_f_o_g = 66 ;  jopcg1_f_o = 86
866       iopcg2_f_o_a  =  67 ;  iopcg2_f_o_g = 67 ;  jopcg2_f_o = 87
867       iopcg3_f_o_a  =  68 ;  iopcg3_f_o_g = 68 ;  jopcg3_f_o = 88
868       iopcg4_f_o_a  =  69 ;  iopcg4_f_o_g = 69 ;  jopcg4_f_o = 89
869       iopcg5_f_o_a  =  70 ;  iopcg5_f_o_g = 70 ;  jopcg5_f_o = 90
870       iopcg6_f_o_a  =  71 ;  iopcg6_f_o_g = 71 ;  jopcg6_f_o = 91
871       iopcg7_f_o_a  =  72 ;  iopcg7_f_o_g = 72 ;  jopcg7_f_o = 92
872       iopcg8_f_o_a  =  73 ;  iopcg8_f_o_g = 73 ;  jopcg8_f_o = 93
874       ismpa_a       =  74 ;  ismpa_g      = 74 ;  jsmpa      = 94
875       ismpbb_a      =  75 ;  ismpbb_g     = 75 ;  jsmpbb     = 95
877       iant1_c_a     =  76 ;  iant1_c_g    = 76 ;  jant1_c    = 96
878       iant2_c_a     =  77 ;  iant2_c_g    = 77 ;  jant2_c    = 97
879       iant3_c_a     =  78 ;  iant3_c_g    = 78 ;  jant3_c    = 98
880       iant4_c_a     =  79 ;  iant4_c_g    = 79 ;  jant4_c    = 99
881       iant1_o_a     =  80 ;  iant1_o_g    = 80 ;  jant1_o    =100
882       iant2_o_a     =  81 ;  iant2_o_g    = 81 ;  jant2_o    =101
883       iant3_o_a     =  82 ;  iant3_o_g    = 82 ;  jant3_o    =102
884       iant4_o_a     =  83 ;  iant4_o_g    = 83 ;  jant4_o    =103
885       ibiog1_c_a    =  84 ;  ibiog1_c_g   = 84 ;  jbiog1_c   =104
886       ibiog2_c_a    =  85 ;  ibiog2_c_g   = 85 ;  jbiog2_c   =105
887       ibiog3_c_a    =  86 ;  ibiog3_c_g   = 86 ;  jbiog3_c   =106
888       ibiog4_c_a    =  87 ;  ibiog4_c_g   = 87 ;  jbiog4_c   =107
889       ibiog1_o_a    =  88 ;  ibiog1_o_g   = 88 ;  jbiog1_o   =108
890       ibiog2_o_a    =  89 ;  ibiog2_o_g   = 89 ;  jbiog2_o   =109
891       ibiog3_o_a    =  90 ;  ibiog3_o_g   = 90 ;  jbiog3_o   =110
892       ibiog4_o_a    =  91 ;  ibiog4_o_g   = 91 ;  jbiog4_o   =111
894       iasoaX_a      =  92 ;  iasoaX_g     = 92 ;  jasoaX     =112
895       iasoa1_a      =  93 ;  iasoa1_g     = 93 ;  jasoa1     =113
896       iasoa2_a      =  94 ;  iasoa2_g     = 94 ;  jasoa2     =114
897       iasoa3_a      =  95 ;  iasoa3_g     = 95 ;  jasoa3     =115
898       iasoa4_a      =  96 ;  iasoa4_g     = 96 ;  jasoa4     =116
899       ibsoaX_a      =  97 ;  ibsoaX_g     = 97 ;  jbsoaX     =117
900       ibsoa1_a      =  98 ;  ibsoa1_g     = 98 ;  jbsoa1     =118
901       ibsoa2_a      =  99 ;  ibsoa2_g     = 99 ;  jbsoa2     =119
902       ibsoa3_a      = 100 ;  ibsoa3_g     =100 ;  jbsoa3     =120
903       ibsoa4_a      = 101 ;  ibsoa4_g     =101 ;  jbsoa4     =121
905       iglysoa_r1_a  = 102 ;                       jglysoa_r1 =122
906       iglysoa_r2_a  = 103 ;                       jglysoa_r2 =123
907       iglysoa_sfc_a = 104 ;                       jglysoa_sfc=124
908       iglysoa_nh4_a = 105 ;                       jglysoa_nh4=125
909       iglysoa_oh_a  = 106 ;                       jglysoa_oh =126
911                                                   jh2o       =127 ! water - part of naercomp
913                              in2o5_g      =102  ! ioa --> NO3-
914                              iclno2_g     =103  ! ioa N2O5+Cl- -->
915                              igly         =104
916                              iho          =105
917 !                            ico2_g       = 14  ! currently not used
919       ico3_a        = 107  
920       ina_a         = 108  
921       ica_a         = 109  
922       ioin_a        = 110  
923       ioc_a         = 111  
924       ibc_a         = 112  
926 !     iiepox_a      =  92 ;  iiepox_g     = 92 ;  jiepox     =113
927 !     igly_a        =  93 ;  igly_g       = 93 ;  jgly       =114
928 !                            in2o5_g      = 94 
929 !                            iclno2_g     = 95 
930 !     iiepoxos_a    =  94                      ;  jiepoxos   =115
931 !     itetrol_a     =  95                      ;  jtetrol    =116
932 !     itanv_a       =  96                      ;  jtanv      =117
933 !     isopnv_a      =  97                      ;  jsopnv     =118
934 !     iternv_a      =  98                      ;  jternv     =119
935 !     iseqnv_a      =  99                      ;  jseqnv     =120
936 !     isianv_a      = 100                      ;  jsianv     =121
937 !                                                 jh2o       =122 ! water - part of naercomp
938 !     idust_a       = 105  
939 !     itr1r1_a      = 108                      ;  jtr1r1     =123
940 !     itr1r2_a      = 109                      ;  jtr1r2     =124
941 !     itr1r3_a      = 110                      ;  jtr1r3     =125
942 !     itr1r4_a      = 111                      ;  jtr1r4     =126
944       return
945       end if ! (ipass == 1) then
948       if (ipass == 2) then
950 ! names of aer species
951 !     aer_name(idust_a)      = "dust"
952 !     aer_name(itr1r1_a)     = "tr1r1"
953 !     aer_name(itr1r2_a)     = "tr1r2"
954 !     aer_name(itr1r3_a)     = "tr1r3"
955 !     aer_name(itr1r4_a)     = "tr1r4"
956 !     aer_name(itanv_a)      = "tanv"
957 !     aer_name(isopnv_a)     = "sopnv"
958 !     aer_name(iternv_a)     = "ternv"
959 !     aer_name(iseqnv_a)     = "seqnv"
960 !     aer_name(isianv_a)     = "sianv"
961       aer_name(ipcg1_b_c_a)  = "pcg1_b_c"
962       aer_name(ipcg2_b_c_a)  = "pcg2_b_c"
963       aer_name(ipcg3_b_c_a)  = "pcg3_b_c"
964       aer_name(ipcg4_b_c_a)  = "pcg4_b_c"
965       aer_name(ipcg5_b_c_a)  = "pcg5_b_c"
966       aer_name(ipcg6_b_c_a)  = "pcg6_b_c"
967       aer_name(ipcg7_b_c_a)  = "pcg7_b_c"
968       aer_name(ipcg8_b_c_a)  = "pcg8_b_c"
969       aer_name(ipcg9_b_c_a)  = "pcg9_b_c"
970       aer_name(ipcg1_b_o_a)  = "pcg1_b_o"
971       aer_name(ipcg2_b_o_a)  = "pcg2_b_o"
972       aer_name(ipcg3_b_o_a)  = "pcg3_b_o"
973       aer_name(ipcg4_b_o_a)  = "pcg4_b_o"
974       aer_name(ipcg5_b_o_a)  = "pcg5_b_o"
975       aer_name(ipcg6_b_o_a)  = "pcg6_b_o"
976       aer_name(ipcg7_b_o_a)  = "pcg7_b_o"
977       aer_name(ipcg8_b_o_a)  = "pcg8_b_o"
978       aer_name(ipcg9_b_o_a)  = "pcg9_b_o"
979       aer_name(iopcg1_b_c_a) = "opcg1_b_c"
980       aer_name(iopcg2_b_c_a) = "opcg2_b_c"
981       aer_name(iopcg3_b_c_a) = "opcg3_b_c"
982       aer_name(iopcg4_b_c_a) = "opcg4_b_c"
983       aer_name(iopcg5_b_c_a) = "opcg5_b_c"
984       aer_name(iopcg6_b_c_a) = "opcg6_b_c"
985       aer_name(iopcg7_b_c_a) = "opcg7_b_c"
986       aer_name(iopcg8_b_c_a) = "opcg8_b_c"
987       aer_name(iopcg1_b_o_a) = "opcg1_b_o"
988       aer_name(iopcg2_b_o_a) = "opcg2_b_o"
989       aer_name(iopcg3_b_o_a) = "opcg3_b_o"
990       aer_name(iopcg4_b_o_a) = "opcg4_b_o"
991       aer_name(iopcg5_b_o_a) = "opcg5_b_o"
992       aer_name(iopcg6_b_o_a) = "opcg6_b_o"
993       aer_name(iopcg7_b_o_a) = "opcg7_b_o"
994       aer_name(iopcg8_b_o_a) = "opcg8_b_o"
995       aer_name(ipcg1_f_c_a)  = "pcg1_f_c"
996       aer_name(ipcg2_f_c_a)  = "pcg2_f_c"
997       aer_name(ipcg3_f_c_a)  = "pcg3_f_c"
998       aer_name(ipcg4_f_c_a)  = "pcg4_f_c"
999       aer_name(ipcg5_f_c_a)  = "pcg5_f_c"
1000       aer_name(ipcg6_f_c_a)  = "pcg6_f_c"
1001       aer_name(ipcg7_f_c_a)  = "pcg7_f_c"
1002       aer_name(ipcg8_f_c_a)  = "pcg8_f_c"
1003       aer_name(ipcg9_f_c_a)  = "pcg9_f_c"
1004       aer_name(ipcg1_f_o_a)  = "pcg1_f_o"
1005       aer_name(ipcg2_f_o_a)  = "pcg2_f_o"
1006       aer_name(ipcg3_f_o_a)  = "pcg3_f_o"
1007       aer_name(ipcg4_f_o_a)  = "pcg4_f_o"
1008       aer_name(ipcg5_f_o_a)  = "pcg5_f_o"
1009       aer_name(ipcg6_f_o_a)  = "pcg6_f_o"
1010       aer_name(ipcg7_f_o_a)  = "pcg7_f_o"
1011       aer_name(ipcg8_f_o_a)  = "pcg8_f_o"
1012       aer_name(ipcg9_f_o_a)  = "pcg9_f_o"
1013       aer_name(iopcg1_f_c_a) = "opcg1_f_c"
1014       aer_name(iopcg2_f_c_a) = "opcg2_f_c"
1015       aer_name(iopcg3_f_c_a) = "opcg3_f_c"
1016       aer_name(iopcg4_f_c_a) = "opcg4_f_c"
1017       aer_name(iopcg5_f_c_a) = "opcg5_f_c"
1018       aer_name(iopcg6_f_c_a) = "opcg6_f_c"
1019       aer_name(iopcg7_f_c_a) = "opcg7_f_c"
1020       aer_name(iopcg8_f_c_a) = "opcg8_f_c"
1021       aer_name(iopcg1_f_o_a) = "opcg1_f_o"
1022       aer_name(iopcg2_f_o_a) = "opcg2_f_o"
1023       aer_name(iopcg3_f_o_a) = "opcg3_f_o"
1024       aer_name(iopcg4_f_o_a) = "opcg4_f_o"
1025       aer_name(iopcg5_f_o_a) = "opcg5_f_o"
1026       aer_name(iopcg6_f_o_a) = "opcg6_f_o"
1027       aer_name(iopcg7_f_o_a) = "opcg7_f_o"
1028       aer_name(iopcg8_f_o_a) = "opcg8_f_o"
1029       aer_name(ismpa_a)      = "smpa"
1030       aer_name(ismpbb_a)     = "smpbb"
1031 !     aer_name(igly_a)       = "gly"
1032 !     aer_name(iiepox_a)     = "iepox"
1033 !     aer_name(iiepoxos_a)   = "iepoxos"
1034 !     aer_name(itetrol_a)    = "tetrol"
1035       aer_name(iant1_c_a)    = "ant1_c"
1036       aer_name(iant2_c_a)    = "ant2_c"
1037       aer_name(iant3_c_a)    = "ant3_c"
1038       aer_name(iant4_c_a)    = "ant4_c"
1039       aer_name(iant1_o_a)    = "ant1_o"
1040       aer_name(iant2_o_a)    = "ant2_o"
1041       aer_name(iant3_o_a)    = "ant3_o"
1042       aer_name(iant4_o_a)    = "ant4_o"
1043       aer_name(ibiog1_c_a)   = "biog1_c"
1044       aer_name(ibiog2_c_a)   = "biog2_c"
1045       aer_name(ibiog3_c_a)   = "biog3_c"
1046       aer_name(ibiog4_c_a)   = "biog4_c"
1047       aer_name(ibiog1_o_a)   = "biog1_o"
1048       aer_name(ibiog2_o_a)   = "biog2_o"
1049       aer_name(ibiog3_o_a)   = "biog3_o"
1050       aer_name(ibiog4_o_a)   = "biog4_o"
1052       aer_name(iglysoa_r1_a) ="glysoa_r1"
1053       aer_name(iglysoa_r2_a) ="glysoa_r2"
1054       aer_name(iglysoa_sfc_a)="glysoa_sfc"
1055       aer_name(iglysoa_nh4_a)="glysoa_nh4"
1056       aer_name(iglysoa_oh_a) ="glysoa_oh"
1057       aer_name(iasoaX_a)     ="asoaX"
1058       aer_name(iasoa1_a)     ="asoa1"
1059       aer_name(iasoa2_a)     ="asoa2"
1060       aer_name(iasoa3_a)     ="asoa3"
1061       aer_name(iasoa4_a)     ="asoa4"
1062       aer_name(ibsoaX_a)     ="bsoaX"
1063       aer_name(ibsoa1_a)     ="bsoa1"
1064       aer_name(ibsoa2_a)     ="bsoa2"
1065       aer_name(ibsoa3_a)     ="bsoa3"
1066       aer_name(ibsoa4_a)     ="bsoa4"
1068 ! names of gas species
1069       gas_name(ipcg1_b_c_g)  = "pcg1_b_c"
1070       gas_name(ipcg2_b_c_g)  = "pcg2_b_c"
1071       gas_name(ipcg3_b_c_g)  = "pcg3_b_c"
1072       gas_name(ipcg4_b_c_g)  = "pcg4_b_c"
1073       gas_name(ipcg5_b_c_g)  = "pcg5_b_c"
1074       gas_name(ipcg6_b_c_g)  = "pcg6_b_c"
1075       gas_name(ipcg7_b_c_g)  = "pcg7_b_c"
1076       gas_name(ipcg8_b_c_g)  = "pcg8_b_c"
1077       gas_name(ipcg9_b_c_g)  = "pcg9_b_c"
1078       gas_name(ipcg1_b_o_g)  = "pcg1_b_o"
1079       gas_name(ipcg2_b_o_g)  = "pcg2_b_o"
1080       gas_name(ipcg3_b_o_g)  = "pcg3_b_o"
1081       gas_name(ipcg4_b_o_g)  = "pcg4_b_o"
1082       gas_name(ipcg5_b_o_g)  = "pcg5_b_o"
1083       gas_name(ipcg6_b_o_g)  = "pcg6_b_o"
1084       gas_name(ipcg7_b_o_g)  = "pcg7_b_o"
1085       gas_name(ipcg8_b_o_g)  = "pcg8_b_o"
1086       gas_name(ipcg9_b_o_g)  = "pcg9_b_o"
1087       gas_name(iopcg1_b_c_g) = "opcg1_b_c"
1088       gas_name(iopcg2_b_c_g) = "opcg2_b_c"
1089       gas_name(iopcg3_b_c_g) = "opcg3_b_c"
1090       gas_name(iopcg4_b_c_g) = "opcg4_b_c"
1091       gas_name(iopcg5_b_c_g) = "opcg5_b_c"
1092       gas_name(iopcg6_b_c_g) = "opcg6_b_c"
1093       gas_name(iopcg7_b_c_g) = "opcg7_b_c"
1094       gas_name(iopcg8_b_c_g) = "opcg8_b_c"
1095       gas_name(iopcg1_b_o_g) = "opcg1_b_o"
1096       gas_name(iopcg2_b_o_g) = "opcg2_b_o"
1097       gas_name(iopcg3_b_o_g) = "opcg3_b_o"
1098       gas_name(iopcg4_b_o_g) = "opcg4_b_o"
1099       gas_name(iopcg5_b_o_g) = "opcg5_b_o"
1100       gas_name(iopcg6_b_o_g) = "opcg6_b_o"
1101       gas_name(iopcg7_b_o_g) = "opcg7_b_o"
1102       gas_name(iopcg8_b_o_g) = "opcg8_b_o"
1103       gas_name(ipcg1_f_c_g)  = "pcg1_f_c"
1104       gas_name(ipcg2_f_c_g)  = "pcg2_f_c"
1105       gas_name(ipcg3_f_c_g)  = "pcg3_f_c"
1106       gas_name(ipcg4_f_c_g)  = "pcg4_f_c"
1107       gas_name(ipcg5_f_c_g)  = "pcg5_f_c"
1108       gas_name(ipcg6_f_c_g)  = "pcg6_f_c"
1109       gas_name(ipcg7_f_c_g)  = "pcg7_f_c"
1110       gas_name(ipcg8_f_c_g)  = "pcg8_f_c"
1111       gas_name(ipcg9_f_c_g)  = "pcg9_f_c"
1112       gas_name(ipcg1_f_o_g)  = "pcg1_f_o"
1113       gas_name(ipcg2_f_o_g)  = "pcg2_f_o"
1114       gas_name(ipcg3_f_o_g)  = "pcg3_f_o"
1115       gas_name(ipcg4_f_o_g)  = "pcg4_f_o"
1116       gas_name(ipcg5_f_o_g)  = "pcg5_f_o"
1117       gas_name(ipcg6_f_o_g)  = "pcg6_f_o"
1118       gas_name(ipcg7_f_o_g)  = "pcg7_f_o"
1119       gas_name(ipcg8_f_o_g)  = "pcg8_f_o"
1120       gas_name(ipcg9_f_o_g)  = "pcg9_f_o"
1121       gas_name(iopcg1_f_c_g) = "opcg1_f_c"
1122       gas_name(iopcg2_f_c_g) = "opcg2_f_c"
1123       gas_name(iopcg3_f_c_g) = "opcg3_f_c"
1124       gas_name(iopcg4_f_c_g) = "opcg4_f_c"
1125       gas_name(iopcg5_f_c_g) = "opcg5_f_c"
1126       gas_name(iopcg6_f_c_g) = "opcg6_f_c"
1127       gas_name(iopcg7_f_c_g) = "opcg7_f_c"
1128       gas_name(iopcg8_f_c_g) = "opcg8_f_c"
1129       gas_name(iopcg1_f_o_g) = "opcg1_f_o"
1130       gas_name(iopcg2_f_o_g) = "opcg2_f_o"
1131       gas_name(iopcg3_f_o_g) = "opcg3_f_o"
1132       gas_name(iopcg4_f_o_g) = "opcg4_f_o"
1133       gas_name(iopcg5_f_o_g) = "opcg5_f_o"
1134       gas_name(iopcg6_f_o_g) = "opcg6_f_o"
1135       gas_name(iopcg7_f_o_g) = "opcg7_f_o"
1136       gas_name(iopcg8_f_o_g) = "opcg8_f_o"
1137       gas_name(ismpa_g)      = "smpa"
1138       gas_name(ismpbb_g)     = "smpbb"
1139 !     gas_name(igly_g)       = "gly"
1140 !     gas_name(iiepox_g)     = "iepox"
1141       gas_name(iant1_c_g)    = "ant1_c"
1142       gas_name(iant2_c_g)    = "ant2_c"
1143       gas_name(iant3_c_g)    = "ant3_c"
1144       gas_name(iant4_c_g)    = "ant4_c"
1145       gas_name(iant1_o_g)    = "ant1_o"
1146       gas_name(iant2_o_g)    = "ant2_o"
1147       gas_name(iant3_o_g)    = "ant3_o"
1148       gas_name(iant4_o_g)    = "ant4_o"
1149       gas_name(ibiog1_c_g)   = "biog1_c"
1150       gas_name(ibiog2_c_g)   = "biog2_c"
1151       gas_name(ibiog3_c_g)   = "biog3_c"
1152       gas_name(ibiog4_c_g)   = "biog4_c"
1153       gas_name(ibiog1_o_g)   = "biog1_o"
1154       gas_name(ibiog2_o_g)   = "biog2_o"
1155       gas_name(ibiog3_o_g)   = "biog3_o"
1156       gas_name(ibiog4_o_g)   = "biog4_o"
1157       gas_name(in2o5_g)      = "n2o5 "
1158       gas_name(iclno2_g)     = "clno2"
1159       
1160       gas_name(iasoaX_g)="asoaX"
1161       gas_name(iasoa1_g)="asoa1"
1162       gas_name(iasoa2_g)="asoa2"
1163       gas_name(iasoa3_g)="asoa3"
1164       gas_name(iasoa4_g)="asoa4"
1165       gas_name(ibsoaX_g)="bsoaX"
1166       gas_name(ibsoa1_g)="bsoa1"
1167       gas_name(ibsoa2_g)="bsoa2"
1168       gas_name(ibsoa3_g)="bsoa3"
1169       gas_name(ibsoa4_g)="bsoa4"
1171 ! densities of compounds in g/cc
1172       dens_comp_a(:)        = 1.0       ! default
1174       dens_comp_a(jh2o)     = 1.0
1175       dens_comp_a(jnh4so4)  = 1.8
1176       dens_comp_a(jlvcite)  = 1.8
1177       dens_comp_a(jnh4hso4) = 1.8
1178       dens_comp_a(jnh4msa)  = 1.8       ! assumed same as nh4hso4
1179       dens_comp_a(jnh4no3)  = 1.7
1180       dens_comp_a(jnh4cl)   = 1.5
1181       dens_comp_a(jnacl)    = 2.2
1182       dens_comp_a(jnano3)   = 2.2
1183       dens_comp_a(jna2so4)  = 2.2
1184       dens_comp_a(jna3hso4) = 2.2
1185       dens_comp_a(jnahso4)  = 2.2
1186       dens_comp_a(jnamsa)   = 2.2       ! assumed same as nahso4
1187       dens_comp_a(jcaso4)   = 2.6
1188       dens_comp_a(jcamsa2)  = 2.6       ! assumed same as caso4
1189       dens_comp_a(jcano3)   = 2.6
1190       dens_comp_a(jcacl2)   = 2.6
1191       dens_comp_a(jcaco3)   = 2.6
1192       dens_comp_a(jh2so4)   = 1.8
1193       dens_comp_a(jhhso4)   = 1.8
1194       dens_comp_a(jhno3)    = 1.8
1195       dens_comp_a(jhcl)     = 1.8
1196       dens_comp_a(jmsa)     = 1.8       ! assumed same as h2so4
1197       dens_comp_a(joc)      = 1.0
1198       dens_comp_a(jbc)      = 1.8
1199       dens_comp_a(join)     = 2.6
1201       dens_comp_a(iasoaX_a) = 1.5
1202       dens_comp_a(iasoa1_a) = 1.5
1203       dens_comp_a(iasoa2_a) = 1.5
1204       dens_comp_a(iasoa3_a) = 1.5
1205       dens_comp_a(iasoa4_a) = 1.5
1206       dens_comp_a(ibsoaX_a) = 1.5
1207       dens_comp_a(ibsoa1_a) = 1.5
1208       dens_comp_a(ibsoa2_a) = 1.5
1209       dens_comp_a(ibsoa3_a) = 1.5
1210       dens_comp_a(ibsoa4_a) = 1.5
1212 !     dens_comp_a(jdust)    = 2.6
1213 !     dens_comp_a(jtr1r1)   = 2.6
1214 !     dens_comp_a(jtr1r2)   = 2.6
1215 !     dens_comp_a(jtr1r3)   = 2.6
1216 !     dens_comp_a(jtr1r4)   = 2.6
1217 !     dens_comp_a(jtanv)    = 1.0
1220 ! molecular weights of generic aerosol species
1221       mw_aer_mac(:)        = 250.0 ! default
1223       mw_aer_mac(iso4_a)   =  96.0
1224       mw_aer_mac(ino3_a)   =  62.0
1225       mw_aer_mac(icl_a)    =  35.5
1226       mw_aer_mac(imsa_a)   =  95.0 ! ch3so3
1227       mw_aer_mac(ico3_a)   =  60.0
1228       mw_aer_mac(inh4_a)   =  18.0
1229       mw_aer_mac(ina_a)    =  23.0
1230       mw_aer_mac(ica_a)    =  40.0
1231       mw_aer_mac(ioin_a)   =   1.0  ! not used
1232 !     mw_aer_mac(idust_a)  =   1.0  ! not used
1233 !     mw_aer_mac(itr1r1_a) =   1.0  ! not used
1234 !     mw_aer_mac(itr1r2_a) =   1.0  ! not used
1235 !     mw_aer_mac(itr1r3_a) =   1.0  ! not used
1236 !     mw_aer_mac(itr1r4_a) =   1.0  ! not used
1237       mw_aer_mac(ibc_a)    =   1.0  ! not used
1239       mw_aer_mac(ioc_a)      = 250.0
1240 !     mw_aer_mac(igly_a)     =  58.0
1241 !     mw_aer_mac(iiepox_a)   = 118.0
1242 !     mw_aer_mac(iiepoxos_a) = 216.0
1243 !     mw_aer_mac(itetrol_a)  = 136.0
1246 ! molecular weights of compounds
1247       mw_comp_a(:)        = 250.0 ! default
1249       mw_comp_a(jh2o)     =  18.0
1250       mw_comp_a(jnh4so4)  = 132.0
1251       mw_comp_a(jlvcite)  = 247.0
1252       mw_comp_a(jnh4hso4) = 115.0
1253       mw_comp_a(jnh4msa)  = 113.0
1254       mw_comp_a(jnh4no3)  =  80.0
1255       mw_comp_a(jnh4cl)   =  53.5
1256       mw_comp_a(jnacl)    =  58.5
1257       mw_comp_a(jnano3)   =  85.0
1258       mw_comp_a(jna2so4)  = 142.0
1259       mw_comp_a(jna3hso4) = 262.0
1260       mw_comp_a(jnahso4)  = 120.0
1261       mw_comp_a(jnamsa)   = 118.0
1262       mw_comp_a(jcaso4)   = 136.0
1263       mw_comp_a(jcamsa2)  = 230.0
1264       mw_comp_a(jcano3)   = 164.0
1265       mw_comp_a(jcacl2)   = 111.0
1266       mw_comp_a(jcaco3)   = 100.0
1267       mw_comp_a(jh2so4)   =  98.0
1268       mw_comp_a(jhhso4)   =  98.0
1269       mw_comp_a(jhno3)    =  63.0
1270       mw_comp_a(jhcl)     =  36.5
1271       mw_comp_a(jmsa)     =  96.0
1272       mw_comp_a(joc)      = 250.0
1273       mw_comp_a(jbc)      =   1.0
1274       mw_comp_a(join)     =   1.0
1275 !     mw_comp_a(jdust)    =   1.0
1277 !     mw_comp_a(jtr1r1)   =   1.0
1278 !     mw_comp_a(jtr1r2)   =   1.0
1279 !     mw_comp_a(jtr1r3)   =   1.0
1280 !     mw_comp_a(jtr1r4)   =   1.0
1281 !     mw_comp_a(jgly)     =  58.0
1282 !     mw_comp_a(jiepox)   = 118.0
1283 !     mw_comp_a(jiepoxos) = 216.0
1284 !     mw_comp_a(jtetrol)  = 136.0
1287 ! densities of generic aerosol species
1288       dens_aer_mac(:)        = 1.0        ! default
1290       dens_aer_mac(iso4_a)   = 1.8        ! used
1291       dens_aer_mac(ino3_a)   = 1.8        ! used
1292       dens_aer_mac(icl_a)    = 2.2        ! used
1293       dens_aer_mac(imsa_a)   = 1.8        ! used
1294       dens_aer_mac(ico3_a)   = 2.6        ! used
1295       dens_aer_mac(inh4_a)   = 1.8        ! used
1296       dens_aer_mac(ina_a)    = 2.2        ! used
1297       dens_aer_mac(ica_a)    = 2.6        ! used
1298       dens_aer_mac(ioin_a)   = 2.6        ! used
1299 !     dens_aer_mac(idust_a)  = 2.6       ! used
1300 !     dens_aer_mac(itr1r1_a) = 2.6      ! used
1301 !     dens_aer_mac(itr1r2_a) = 2.6      ! used
1302 !     dens_aer_mac(itr1r3_a) = 2.6      ! used
1303 !     dens_aer_mac(itr1r4_a) = 2.6      ! used
1304       dens_aer_mac(ioc_a)    = 1.0        ! used
1305       dens_aer_mac(ibc_a)    = 1.7        ! used
1307       dens_aer_mac(iasoa1_a) = 1.5
1308       dens_aer_mac(iasoa2_a) = 1.5
1309       dens_aer_mac(iasoa3_a) = 1.5
1310       dens_aer_mac(iasoa4_a) = 1.5
1311       dens_aer_mac(iasoaX_a) = 1.5
1312       dens_aer_mac(ibsoa1_a) = 1.5
1313       dens_aer_mac(ibsoa2_a) = 1.5
1314       dens_aer_mac(ibsoa3_a) = 1.5
1315       dens_aer_mac(ibsoa4_a) = 1.5
1316       dens_aer_mac(ibsoaX_a) = 1.5
1319 ! kappa_aer_mac not used in previous wrfchem mosaic
1322 ! partial molar volumes of condensing species
1323       partial_molar_vol(:)        = 250.0  ! default
1325       partial_molar_vol(ih2so4_g) =  51.83
1326       partial_molar_vol(ihno3_g)  =  31.45
1327       partial_molar_vol(ihcl_g)   =  20.96
1328       partial_molar_vol(inh3_g)   =  24.03
1329       partial_molar_vol(imsa_g)   =  53.33
1331 !     partial_molar_vol(iiepox_g) = 118.0
1332 !     partial_molar_vol(igly_g)   =  58.0
1333       partial_molar_vol(in2o5_g)  = 200.0       ! assumed...
1334       partial_molar_vol(iclno2_g) = 200.0       ! assumed...
1337 ! used to calculate diffusivities of condensing gases
1338 ! (other values set in module_mosaic_init_aerpar are ok)
1339       v_molar_gas(in2o5_g)  = 60.40
1340       v_molar_gas(iclno2_g) = 52.70
1343 ! molecular weights of condensing gases
1344       mw_gas(1:ngas_aerchtot) = 250.0  ! default
1345       mw_gas(ih2so4_g) = 98.0
1346       mw_gas(ihno3_g)  = 63.0
1347       mw_gas(ihcl_g)   = 36.5
1348       mw_gas(inh3_g)   = 17.0
1349       mw_gas(imsa_g)   = 96.0
1350       mw_gas(in2o5_g)  = 108.0
1351       mw_gas(iclno2_g) = 81.5
1352 !     mw_gas(igly)=58.0
1353 !     mw_gas(iho)=17.0
1356       return
1357       end if ! (ipass == 2) then
1360       return
1361       end subroutine soa_vbs_load_params
1364 !-------------------------------------------------------------------------------
1365       subroutine soa_vbs_update_thermcons( t_k, po_soa, sat_soa )
1367       use module_data_mosaic_aero, only: &
1368          msoa_vbs_info, &
1369          ngas_ioa, ngas_soa, ngas_volatile, &
1370          ipcg1_b_c_g,  ipcg2_b_c_g,  ipcg3_b_c_g,  ipcg4_b_c_g,  &
1371          ipcg5_b_c_g,  ipcg6_b_c_g,  ipcg7_b_c_g,  ipcg8_b_c_g,  &
1372          ipcg9_b_c_g,                                            &
1373          ipcg1_b_o_g,  ipcg2_b_o_g,  ipcg3_b_o_g,  ipcg4_b_o_g,  &
1374          ipcg5_b_o_g,  ipcg6_b_o_g,  ipcg7_b_o_g,  ipcg8_b_o_g,  &
1375          ipcg9_b_o_g,                                            &
1376          iopcg1_b_c_g, iopcg2_b_c_g, iopcg3_b_c_g, iopcg4_b_c_g, &
1377          iopcg5_b_c_g, iopcg6_b_c_g, iopcg7_b_c_g, iopcg8_b_c_g, &
1378          iopcg1_b_o_g, iopcg2_b_o_g, iopcg3_b_o_g, iopcg4_b_o_g, &
1379          iopcg5_b_o_g, iopcg6_b_o_g, iopcg7_b_o_g, iopcg8_b_o_g, &
1380          ipcg1_f_c_g,  ipcg2_f_c_g,  ipcg3_f_c_g,  ipcg4_f_c_g,  &
1381          ipcg5_f_c_g,  ipcg6_f_c_g,  ipcg7_f_c_g,  ipcg8_f_c_g,  &
1382          ipcg9_f_c_g,                                            &
1383          ipcg1_f_o_g,  ipcg2_f_o_g,  ipcg3_f_o_g,  ipcg4_f_o_g,  &
1384          ipcg5_f_o_g,  ipcg6_f_o_g,  ipcg7_f_o_g,  ipcg8_f_o_g,  &
1385          ipcg9_f_o_g,                                            &
1386          iopcg1_f_c_g, iopcg2_f_c_g, iopcg3_f_c_g, iopcg4_f_c_g, &
1387          iopcg5_f_c_g, iopcg6_f_c_g, iopcg7_f_c_g, iopcg8_f_c_g, &
1388          iopcg1_f_o_g, iopcg2_f_o_g, iopcg3_f_o_g, iopcg4_f_o_g, &
1389          iopcg5_f_o_g, iopcg6_f_o_g, iopcg7_f_o_g, iopcg8_f_o_g, &
1390          iant1_c_g,    iant2_c_g,    iant3_c_g,    iant4_c_g,    &
1391          iant1_o_g,    iant2_o_g,    iant3_o_g,    iant4_o_g,    &
1392          ibiog1_c_g,   ibiog2_c_g,   ibiog3_c_g,   ibiog4_c_g,   &
1393          ibiog1_o_g,   ibiog2_o_g,   ibiog3_o_g,   ibiog4_o_g,   &
1394          ismpa_g,      ismpbb_g,                                 &
1395          iasoa1_g,     iasoa2_g,     iasoa3_g,     iasoa4_g,     &
1396          iasoaX_g,                                               &
1397          ibsoa1_g,     ibsoa2_g,     ibsoa3_g,     ibsoa4_g,     &
1398          ibsoaX_g
1399 !        iiepox_g,     igly_g
1401 ! arguments
1402       real(r8), intent(in   ) :: t_k
1403       real(r8), intent(inout) :: po_soa(ngas_volatile)
1404       real(r8), intent(inout) :: sat_soa(ngas_volatile)
1406 ! local variables
1407       integer :: itmpa, iv
1408       integer :: vbs_nbin(1)
1411       vbs_nbin(1) = msoa_vbs_info(1)
1414 ! vapor pressures of soa species
1416       po_soa = 1.0e5_r8  ! this default value gives csat=6.1e9 ug/m^3 (~1e9 ppbv) for 298K and MW=150
1418       if (vbs_nbin(1) .eq. 9) then ! default
1420       po_soa(ipcg1_b_c_g)  = fn_po(9.91d-9, 83.0d0, t_k) ! [pascal]
1421       po_soa(ipcg2_b_c_g)  = fn_po(9.91d-7, 106.0d0, t_k) ! [pascal]
1422       po_soa(ipcg3_b_c_g)  = fn_po(9.91d-6, 100.0d0, t_k) ! [pascal]
1423       po_soa(ipcg4_b_c_g)  = fn_po(9.91d-5, 94.0d0, t_k) ! [pascal]
1424       po_soa(ipcg5_b_c_g)  = fn_po(9.91d-4, 88.0d0, t_k) ! [pascal]
1425       po_soa(ipcg6_b_c_g)  = fn_po(9.91d-3, 82.0d0, t_k) ! [pascal]
1426       po_soa(ipcg7_b_c_g)  = fn_po(9.91d-2, 64.0d0, t_k) ! [pascal]
1427       po_soa(ipcg8_b_c_g)  = fn_po(9.91d-1, 70.0d0, t_k) ! [pascal]
1428       po_soa(ipcg9_b_c_g)  = fn_po(9.91d0,  64.0d0, t_k) ! [pascal]
1429       po_soa(iopcg1_b_c_g) = fn_po(9.91d-8, 112.0d0, t_k) ! [pascal]
1430       po_soa(iopcg2_b_c_g) = fn_po(9.91d-7, 106.0d0, t_k) ! [pascal]
1431       po_soa(iopcg3_b_c_g) = fn_po(9.91d-6, 100.0d0, t_k) ! [pascal]
1432       po_soa(iopcg4_b_c_g) = fn_po(9.91d-5, 94.0d0, t_k) ! [pascal]
1433       po_soa(iopcg5_b_c_g) = fn_po(9.91d-4, 88.0d0, t_k) ! [pascal]
1434       po_soa(iopcg6_b_c_g) = fn_po(9.91d-3, 82.0d0, t_k) ! [pascal]
1435       po_soa(iopcg7_b_c_g) = fn_po(9.91d-2, 76.0d0, t_k) ! [pascal]
1436       po_soa(iopcg8_b_c_g) = fn_po(9.91d-1, 70.0d0, t_k) ! [pascal]
1438       po_soa(ipcg1_b_o_g)  = fn_po(9.91d-8, 112.0d0, t_k) ! [pascal]
1439       po_soa(ipcg2_b_o_g)  = fn_po(9.91d-7, 106.0d0, t_k) ! [pascal]
1440       po_soa(ipcg3_b_o_g)  = fn_po(9.91d-6, 100.0d0, t_k) ! [pascal]
1441       po_soa(ipcg4_b_o_g)  = fn_po(9.91d-5, 94.0d0, t_k) ! [pascal]
1442       po_soa(ipcg5_b_o_g)  = fn_po(9.91d-4, 88.0d0, t_k) ! [pascal]
1443       po_soa(ipcg6_b_o_g)  = fn_po(9.91d-3, 82.0d0, t_k) ! [pascal]
1444       po_soa(ipcg7_b_o_g)  = fn_po(9.91d-2, 76.0d0, t_k) ! [pascal]
1445       po_soa(ipcg8_b_o_g)  = fn_po(9.91d-1, 70.0d0, t_k) ! [pascal]
1446       po_soa(ipcg9_b_o_g)  = fn_po(9.91d0,  64.0d0, t_k) ! [pascal]
1447       po_soa(iopcg1_b_o_g) = fn_po(9.91d-8, 112.0d0, t_k) ! [pascal]
1448       po_soa(iopcg2_b_o_g) = fn_po(9.91d-7, 106.0d0, t_k) ! [pascal]
1449       po_soa(iopcg3_b_o_g) = fn_po(9.91d-6, 100.0d0, t_k) ! [pascal]
1450       po_soa(iopcg4_b_o_g) = fn_po(9.91d-5, 94.0d0, t_k) ! [pascal]
1451       po_soa(iopcg5_b_o_g) = fn_po(9.91d-4, 88.0d0, t_k) ! [pascal]
1452       po_soa(iopcg6_b_o_g) = fn_po(9.91d-3, 82.0d0, t_k) ! [pascal]
1453       po_soa(iopcg7_b_o_g) = fn_po(9.91d-2, 76.0d0, t_k) ! [pascal]
1454       po_soa(iopcg8_b_o_g) = fn_po(9.91d-1, 70.0d0, t_k) ! [pascal]
1456       po_soa(ipcg1_f_c_g)  = fn_po(9.91d-8, 112.0d0, t_k) ! [pascal]
1457       po_soa(ipcg2_f_c_g)  = fn_po(9.91d-7, 106.0d0, t_k) ! [pascal]
1458       po_soa(ipcg3_f_c_g)  = fn_po(9.91d-6, 100.0d0, t_k) ! [pascal]
1459       po_soa(ipcg4_f_c_g)  = fn_po(9.91d-5, 94.0d0, t_k) ! [pascal]
1460       po_soa(ipcg5_f_c_g)  = fn_po(9.91d-4, 88.0d0, t_k) ! [pascal]
1461       po_soa(ipcg6_f_c_g)  = fn_po(9.91d-3, 82.0d0, t_k) ! [pascal]
1462       po_soa(ipcg7_f_c_g)  = fn_po(9.91d-2, 76.0d0, t_k) ! [pascal]
1463       po_soa(ipcg8_f_c_g)  = fn_po(9.91d-1, 70.0d0, t_k) ! [pascal]
1464       po_soa(ipcg9_f_c_g)  = fn_po(9.91d0,  64.0d0, t_k) ! [pascal]
1465       po_soa(iopcg1_f_c_g) = fn_po(9.91d-8, 112.0d0, t_k) ! [pascal]
1466       po_soa(iopcg2_f_c_g) = fn_po(9.91d-7, 106.0d0, t_k) ! [pascal]
1467       po_soa(iopcg3_f_c_g) = fn_po(9.91d-6, 100.0d0, t_k) ! [pascal]
1468       po_soa(iopcg4_f_c_g) = fn_po(9.91d-5, 94.0d0, t_k) ! [pascal]
1469       po_soa(iopcg5_f_c_g) = fn_po(9.91d-4, 88.0d0, t_k) ! [pascal]
1470       po_soa(iopcg6_f_c_g) = fn_po(9.91d-3, 82.0d0, t_k) ! [pascal]
1471       po_soa(iopcg7_f_c_g) = fn_po(9.91d-2, 76.0d0, t_k) ! [pascal]
1472       po_soa(iopcg8_f_c_g) = fn_po(9.91d-1, 70.0d0, t_k) ! [pascal]
1474       po_soa(ipcg1_f_o_g)  = fn_po(9.91d-8, 112.0d0, t_k) ! [pascal]
1475       po_soa(ipcg2_f_o_g)  = fn_po(9.91d-7, 106.0d0, t_k) ! [pascal]
1476       po_soa(ipcg3_f_o_g)  = fn_po(9.91d-6, 100.0d0, t_k) ! [pascal]
1477       po_soa(ipcg4_f_o_g)  = fn_po(9.91d-5, 94.0d0, t_k) ! [pascal]
1478       po_soa(ipcg5_f_o_g)  = fn_po(9.91d-4, 88.0d0, t_k) ! [pascal]
1479       po_soa(ipcg6_f_o_g)  = fn_po(9.91d-3, 82.0d0, t_k) ! [pascal]
1480       po_soa(ipcg7_f_o_g)  = fn_po(9.91d-2, 76.0d0, t_k) ! [pascal]
1481       po_soa(ipcg8_f_o_g)  = fn_po(9.91d-1, 70.0d0, t_k) ! [pascal]
1482       po_soa(ipcg9_f_o_g)  = fn_po(9.91d0,  64.0d0, t_k) ! [pascal]
1483       po_soa(iopcg1_f_o_g) = fn_po(9.91d-8, 112.0d0, t_k) ! [pascal]
1484       po_soa(iopcg2_f_o_g) = fn_po(9.91d-7, 106.0d0, t_k) ! [pascal]
1485       po_soa(iopcg3_f_o_g) = fn_po(9.91d-6, 100.0d0, t_k) ! [pascal]
1486       po_soa(iopcg4_f_o_g) = fn_po(9.91d-5, 94.0d0, t_k) ! [pascal]
1487       po_soa(iopcg5_f_o_g) = fn_po(9.91d-4, 88.0d0, t_k) ! [pascal]
1488       po_soa(iopcg6_f_o_g) = fn_po(9.91d-3, 82.0d0, t_k) ! [pascal]
1489       po_soa(iopcg7_f_o_g) = fn_po(9.91d-2, 76.0d0, t_k) ! [pascal]
1490       po_soa(iopcg8_f_o_g) = fn_po(9.91d-1, 70.0d0, t_k) ! [pascal]
1492       po_soa(iant1_o_g) = fn_po(9.91d-6, 100.0d0, t_k) ! [pascal]
1493       po_soa(iant2_o_g) = fn_po(9.91d-5, 94.0d0, t_k) ! [pascal]
1494       po_soa(iant3_o_g) = fn_po(9.91d-4, 88.0d0, t_k) ! [pascal]
1495       po_soa(iant4_o_g) = fn_po(9.91d-3, 82.0d0, t_k) ! [pascal]
1497       po_soa(iant1_c_g) = fn_po(9.91d-6, 106.0d0, t_k) ! [pascal]
1498       po_soa(iant2_c_g) = fn_po(9.91d-5, 100.0d0, t_k) ! [pascal]
1499       po_soa(iant3_c_g) = fn_po(9.91d-4, 94.0d0, t_k) ! [pascal]
1500       po_soa(iant4_c_g) = fn_po(9.91d-3, 88.0d0, t_k) ! [pascal]
1502       po_soa(ibiog1_c_g) = fn_po(9.91d-6, 106.0d0, t_k) ! [pascal]
1503       po_soa(ibiog2_c_g) = fn_po(9.91d-5, 100.0d0, t_k) ! [pascal]
1504       po_soa(ibiog3_c_g) = fn_po(9.91d-4, 94.0d0, t_k) ! [pascal]
1505       po_soa(ibiog4_c_g) = fn_po(9.91d-3, 88.0d0, t_k) ! [pascal]
1507       po_soa(ibiog1_o_g) = fn_po(9.91d-6, 106.0d0, t_k) ! [pascal]
1508       po_soa(ibiog2_o_g) = fn_po(9.91d-5, 100.0d0, t_k) ! [pascal]
1509       po_soa(ibiog3_o_g) = fn_po(9.91d-4, 94.0d0, t_k) ! [pascal]
1510       po_soa(ibiog4_o_g) = fn_po(9.91d-3, 88.0d0, t_k) ! [pascal]
1512       end if ! (vbs_nbin(1) .eq. 9) then ! default
1515       if (vbs_nbin(1).eq.4) then
1516         po_soa(iasoaX_g) = fn_po(9.91d-10, 40.0d0, T_K) ! [Pascal]
1517         po_soa(iasoa1_g) = fn_po(9.91d-6, dhr_approx(0.0d0), T_K) !  [Pascal]
1518         po_soa(iasoa2_g) = fn_po(9.91d-5, dhr_approx(1.0d0), T_K) !  [Pascal]
1519         po_soa(iasoa3_g) = fn_po(9.91d-4, dhr_approx(2.0d0), T_K) !  [Pascal]
1520         po_soa(iasoa4_g) = fn_po(9.91d-3, dhr_approx(3.0d0), T_K) !  [Pascal]
1521         po_soa(ibsoaX_g) = fn_po(9.91d-10, 40.0d0, T_K) ! [Pascal]
1522         po_soa(ibsoa1_g) = fn_po(9.91d-6, dhr_approx(0.0d0), T_K) !  [Pascal]
1523         po_soa(ibsoa2_g) = fn_po(9.91d-5, dhr_approx(1.0d0), T_K) !  [Pascal]
1524         po_soa(ibsoa3_g) = fn_po(9.91d-4, dhr_approx(2.0d0), T_K) !  [Pascal]
1525         po_soa(ibsoa4_g) = fn_po(9.91d-3, dhr_approx(3.0d0), T_K) !  [Pascal]
1526       endif
1529       if (vbs_nbin(1).eq.2) then
1530         po_soa(ipcg1_b_c_g)  = fn_po(9.91d-8, 83.0d0, t_k) ! [pascal]
1531         po_soa(ipcg2_b_c_g)  = fn_po(9.91d-1, 83.0d0, t_k) ! [pascal]
1532         po_soa(iopcg1_b_c_g) = fn_po(9.91d-8, 83.0d0, t_k) ! [pascal]
1533         po_soa(ipcg1_b_o_g)  = fn_po(9.91d-8, 83.0d0, t_k) ! [pascal]
1534         po_soa(ipcg2_b_o_g)  = fn_po(9.91d-1, 83.0d0, t_k) ! [pascal]
1535         po_soa(iopcg1_b_o_g) = fn_po(9.91d-8, 83.0d0, t_k) ! [pascal]
1536         po_soa(ipcg1_f_c_g)  = fn_po(9.91d-8, 83.0d0, t_k) ! [pascal]
1537         po_soa(ipcg2_f_c_g)  = fn_po(9.91d-1, 83.0d0, t_k) ! [pascal]
1538         po_soa(iopcg1_f_c_g) = fn_po(9.91d-8, 83.0d0, t_k) ! [pascal]
1539         po_soa(ipcg1_f_o_g)  = fn_po(9.91d-8, 83.0d0, t_k) ! [pascal]
1540         po_soa(ipcg2_f_o_g)  = fn_po(9.91d-1, 83.0d0, t_k) ! [pascal]
1541         po_soa(iopcg1_f_o_g) = fn_po(9.91d-8, 83.0d0, t_k) ! [pascal]
1542         po_soa(iant1_c_g)    = fn_po(9.91d-6, 83.0d0, t_k) ! [pascal]
1543         po_soa(iant1_o_g)    = fn_po(9.91d-6, 83.0d0, t_k) ! [pascal]
1544         po_soa(ibiog1_c_g)   = fn_po(9.91d-6, 83.0d0, t_k) ! [pascal]
1545         po_soa(ibiog1_o_g)   = fn_po(9.91d-6, 83.0d0, t_k) ! [pascal]
1546       endif
1548       if (vbs_nbin(1).eq.3) then
1549 ! these values for pcg and opcg gases are the same as vbs_nbin(1)==2
1550 ! above
1551         po_soa(ipcg1_b_c_g)  = fn_po(9.91d-8, 83.0d0, t_k) ! [pascal]
1552         po_soa(ipcg2_b_c_g)  = fn_po(9.91d-1, 83.0d0, t_k) ! [pascal]
1553         po_soa(iopcg1_b_c_g) = fn_po(9.91d-8, 83.0d0, t_k) ! [pascal]
1554         po_soa(ipcg1_b_o_g)  = fn_po(9.91d-8, 83.0d0, t_k) ! [pascal]
1555         po_soa(ipcg2_b_o_g)  = fn_po(9.91d-1, 83.0d0, t_k) ! [pascal]
1556         po_soa(iopcg1_b_o_g) = fn_po(9.91d-8, 83.0d0, t_k) ! [pascal]
1557         po_soa(ipcg1_f_c_g)  = fn_po(9.91d-8, 83.0d0, t_k) ! [pascal]
1558         po_soa(ipcg2_f_c_g)  = fn_po(9.91d-1, 83.0d0, t_k) ! [pascal]
1559         po_soa(iopcg1_f_c_g) = fn_po(9.91d-8, 83.0d0, t_k) ! [pascal]
1560         po_soa(ipcg1_f_o_g)  = fn_po(9.91d-8, 83.0d0, t_k) ! [pascal]
1561         po_soa(ipcg2_f_o_g)  = fn_po(9.91d-1, 83.0d0, t_k) ! [pascal]
1562         po_soa(iopcg1_f_o_g) = fn_po(9.91d-8, 83.0d0, t_k) ! [pascal]
1563 ! these values for ant and bio gases are from manish wrfchem 3.5
1564         po_soa(iant1_c_g)    = fn_po(9.91d-7, 106.0d0, t_k) ! [pascal]
1565         po_soa(iant2_c_g)    = fn_po(9.91d-6, 100.0d0, t_k) ! [pascal]
1566         po_soa(iant3_c_g)    = fn_po(9.91d-5,  94.0d0, t_k) ! [pascal]
1567         po_soa(iant4_c_g)    = fn_po(9.91d-4,  88.0d0, t_k) ! [pascal]
1568         po_soa(ibiog1_c_g)   = fn_po(9.91d-7, 106.0d0, t_k) ! [pascal]
1569         po_soa(ibiog2_c_g)   = fn_po(9.91d-6, 100.0d0, t_k) ! [pascal]
1570         po_soa(ibiog3_c_g)   = fn_po(9.91d-5,  94.0d0, t_k) ! [pascal]
1571         po_soa(ibiog1_o_g)   = fn_po(9.91d-7, 106.0d0, t_k) ! [pascal]
1572         po_soa(ibiog2_o_g)   = fn_po(9.91d-6, 100.0d0, t_k) ! [pascal]
1573       endif
1575 ! scm_temp - temporary for scm
1577 ! chem_opt=198 (vbs_nbin=2) a01 soa
1578 ! pcg1_b_c_a01,pcg2_b_c_a01,pcg1_b_o_a01,pcg2_b_o_a01,opcg1_b_c_a01,opcg1_b_o_a01,pcg1_f_c_a01,pcg2_f_c_a01,pcg1_f_o_a01,pcg2_f_o_a01,opcg1_f_c_a01,opcg1_f_o_a01,
1579 ! ant1_c_a01,ant1_o_a01,
1580 ! biog1_c_a01,biog1_o_a01,
1582 ! chem_opt=204 (vbs_nbin=3) a01 soa
1583 ! pcg1_b_c_a01,pcg1_b_o_a01,opcg1_b_c_a01,opcg1_b_o_a01,pcg1_f_c_a01,pcg1_f_o_a01,opcg1_f_c_a01,opcg1_f_o_a01,
1584 ! ant1_c_a01, ant2_c_a01, ant3_c_a01, ant4_c_a01,
1585 ! biog1_c_a01, biog2_c_a01, biog3_c_a01, biog1_o_a01, biog2_o_a01,
1587 ! to get chem_opt=204 to work, do the following in both the mosaic1 and mosaic2 code
1588 !    duplicate the vbs_nbin=2 block above for vbs_nbin=3
1589 !    add lines for ant2-4 and biog2-3 from the vbs_nbin=9 block
1590 ! *** best to not bother about this, as it is a manish problem ***
1593       if (vbs_nbin(1).eq.0) then
1594         po_soa(ismpa_g)    = fn_po(9.91d-8, 83.0d0, t_k) ! [pascal]
1595         po_soa(ismpbb_g)   = fn_po(9.91d-8, 83.0d0, t_k) ! [pascal]
1596         po_soa(ibiog1_c_g) = fn_po(9.91d-6, 83.0d0, t_k) ! [pascal]
1597         po_soa(ibiog1_o_g) = fn_po(9.91d-6, 83.0d0, t_k) ! [pascal]
1598       endif
1602 ! set saturation vapor concentrations
1604       if (vbs_nbin(1).eq.0) then
1605         itmpa = ismpa_g
1606       else if (vbs_nbin(1).eq.4) then
1607         itmpa = iasoaX_g
1608       else
1609         itmpa = ipcg1_b_c_g
1610       end if
1612       sat_soa(:) = 0.0
1613       do iv = itmpa, ngas_ioa + ngas_soa
1614         sat_soa(iv) = 1.e9*po_soa(iv)/(8.314*t_k)       !  [nmol/m^3(air)]
1615       enddo
1618       return
1619       end subroutine soa_vbs_update_thermcons
1622 !-------------------------------------------------------------------------------
1623       ! Function to approximate enthalpy of vaporization for
1624       ! semi-volatile organic aerosols as a function of volatility
1625       ! from Epstein et al., ES&T, 2010
1626       ! http://pubs.acs.org/doi/abs/10.1021/es902497z
1627       real(r8) function dhr_approx(log10_Csat_298)
1629       real(r8), intent(in) :: log10_Csat_298
1631       dhr_approx = -11.0 * log10_Csat_298 + 131.0 ! kJ/mol
1633       end function dhr_approx
1636 !-------------------------------------------------------------------------------
1637       real(r8) function fn_po(po_298, dh, t)        ! touch
1638 ! subr. arguments
1639       real(r8), intent(in) :: po_298, dh, t
1640 ! local variables (none)
1642       fn_po = po_298*exp(-(dh/8.314e-3)*(1./t - 3.354016435e-3))
1644       return
1645       end function fn_po
1649       end module module_mosaic_soa_vbs