updated top-level README and version_decl for V4.5 (#1847)
[WRF.git] / chem / module_cbmz_initmixrats.F
blobc97251b485aea2675697f9fcc1bddda59c82b58b
1 !**********************************************************************************  
2 ! This computer software was prepared by Battelle Memorial Institute, hereinafter
3 ! the Contractor, under Contract No. DE-AC05-76RL0 1830 with the Department of 
4 ! Energy (DOE). NEITHER THE GOVERNMENT NOR THE CONTRACTOR MAKES ANY WARRANTY,
5 ! EXPRESS OR IMPLIED, OR ASSUMES ANY LIABILITY FOR THE USE OF THIS SOFTWARE.
7 ! CBMZ module: see module_cbmz.F for references and terms of use
8 !**********************************************************************************  
10 #define CASENAME 4
11       module module_cbmz_initmixrats
14       use module_peg_util
17       integer, parameter :: cbmz_init_wrf_mixrats_flagaa = 1
18                             ! turns subr cbmz_init_wrf_mixrats on/off
21       contains
24 !-----------------------------------------------------------------------
25       subroutine cbmz_init_wrf_mixrats(   &
26                config_flags,              &
27                z_at_w, g,                 &
28                chem, numgas,              &
29                ids,ide, jds,jde, kds,kde, &
30                ims,ime, jms,jme, kms,kme, &
31                its,ite, jts,jte, kts,kte  )
33 ! provides initial values for cbmz gas species
34 !    for gas species that are common to both cbmz and radm2, the initial
35 !       values are provided via the run initialization file
36 !       (The radm2 gas species are initialized from this file
37 !       when chem_in_opt==anything.  This ought to be changed!)
38 !    for gas species that are in cbmz but not in radm2, the initial values
39 !       are provided here
40 !    currently only hcl and "par" have non-zero initial values,
41 !       and other species are near-zero
43 ! when (gas_ic_opt == gas_ic_pnnl) AND (chem_in_opt == 0),
44 !    ozone is set to "Texas August 2000" values
46 ! setting cbmz_init_wrf_mixrats_flagaa = 1/0 turns this subr on/off.
49    USE module_configure, only:  grid_config_rec_type, num_chem, &
50         p_o3, p_ald, p_hc3, p_hc5, p_hc8, p_ket, p_oli, p_olt, p_ora2, &
51     p_hcl, p_par
52    USE module_state_description, only:  param_first_scalar,   &
53                         gas_ic_pnnl
54    USE module_input_chem_data, only:  bdy_chem_value
56    IMPLICIT NONE
59 !-----------------------------------------------------------------------
60 ! subr arguments
62    INTEGER, INTENT(IN) :: numgas, &
63                           ids,ide, jds,jde, kds,kde, &
64                           ims,ime, jms,jme, kms,kme, &
65                           its,ite, jts,jte, kts,kte
67    REAL, INTENT(IN) :: g
69 ! perturbation and base geopotential at layer boundaries
70    REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), &
71          INTENT(IN) :: z_at_w
73 ! advected chemical tracers
74    REAL, DIMENSION( ims:ime, kms:kme, jms:jme, num_chem ), &
75          INTENT(INOUT) :: chem
77    TYPE(grid_config_rec_type), INTENT(IN) :: config_flags
78 !-----------------------------------------------------------------------
80 !   local variables
81         integer i, j, k, kp1
82     real, dimension( its:ite, kts:kte, jts:jte ) :: z
84         if (cbmz_init_wrf_mixrats_flagaa <= 0) return
87 ! calculate the mid-level heights
89     do j = jts, min(jte,jde-1)
90        do k = kts, kte
91           kp1 = min(k+1, kte)
92           do i = its, min(ite,ide-1)
93              z(i,k,j) = (z_at_w(i,k,j)+z_at_w(i,kp1,j))*0.5
94           end do
95        end do
96     end do
98 #if CASENAME !=4
100 !   when (gas_ic_opt == gas_ic_pnnl) AND (chem_in_opt == 0),
101 !   set ozone (and other gas?) to "Texas August 2000" values
103         if ( (config_flags%chem_in_opt == 0) .and.   &
104              (config_flags%gas_ic_opt == gas_ic_pnnl) ) then
105             do j = jts, min(jte,jde-1)
106             do k = kts, kte
107             do i = its, min(ite,ide-1)
108            call bdy_chem_value( chem(i,k,j,p_o3),z(i,k,j), p_o3, numgas )
109             end do
110             end do
111             end do
112         end if
113 #endif
116 !   compute hcl initial mixing ratio based on the article:
117 !   Graedel TE and WC Keene, 1995: Troposhperic budget of reactive chlorine.
118 !       Global Biogeochemical Cycles. 9, (1), 47-77.
119 !   This calculation should mimic the hcl profile in bdy_chem_value_cbmz,
120 !   below.
122 !    Height(m)     HCl concentration
123 !    ---------     -----------------
124 !    <=1000        0.4 ppbv
125 !    1000<z<2500   linear transision zone
126 !    >=2500        0.1 ppbv
128         do j = jts, min(jte,jde-1)
129         do k = kts, kte
130         do i = its, min(ite,ide-1)
131        if( z(i,k,j) <= 1000. ) then
132           chem(i,k,j,p_hcl) = 0.4*1e-3
133        elseif( z(i,k,j) > 1000. &
134             .and. z(i,k,j) <= 2500. ) then
135           chem(i,k,j,p_hcl) = (0.4*1e-3) + (z(i,k,j)-1000.)* &
136                ((0.1*1e-3)-(0.4*1e-3)) / (2500.-1000.)
137        else
138           chem(i,k,j,p_hcl) = 0.1*1e-3
139        end if
140     end do
141     end do
142     end do
144 !wig, 2-May-2007: Moved this logic to make_chem_profile so this species
145 !                 will now come in via the real.exe generated wrfinput.
146 !!$!
147 !!$!   compute par initial mixing ratio from radm2 hydrocarbon species
148 !!$!   using same formula as for par emissions
149 !!$!
150 !!$     do j = jts, min(jte,jde-1)
151 !!$     do k = kts, kte
152 !!$     do i = its, min(ite,ide-1)
153 !!$         chem(i,k,j,p_par) =                               &
154 !!$               0.4*chem(i,k,j,p_ald) + 2.9*chem(i,k,j,p_hc3)   &
155 !!$             + 4.8*chem(i,k,j,p_hc5) + 7.9*chem(i,k,j,p_hc8)   &
156 !!$             + 0.9*chem(i,k,j,p_ket) + 2.8*chem(i,k,j,p_oli)   &
157 !!$             + 1.8*chem(i,k,j,p_olt) + 1.0*chem(i,k,j,p_ora2) 
158 !!$     end do
159 !!$     end do
160 !!$     end do
162         return
163         end subroutine cbmz_init_wrf_mixrats  
167 !-----------------------------------------------------------------------
168         end module module_cbmz_initmixrats
172 !-----------------------------------------------------------------------
173         subroutine bdy_chem_value_cbmz ( chem_bv, z, nch, numgas )
175 ! provides boundary values for cbmz gas species
176 !    for gas species that are common to both cbmz and radm2, the boundary
177 !       values are provided by subr bdy_chem_value
178 !    for gas species that are in cbmz but not in radm2, the boundary values
179 !       are provided here, except for par
180 !    currently only "hcl" has a non-zero boundary value,
181 !       and other species are near-zero
183 ! this is outside of the module declaration because of potential
184 !    module1 --> module2 --> module1 use conflicts
186         use module_configure, only:  grid_config_rec_type,   &
187             p_o3, p_ald, p_hc3, p_hc5, p_hc8, p_ket, p_oli,  &
188             p_olt, p_ora2, p_hcl, p_par
189         use module_input_chem_data, only:  bdy_chem_value
191         implicit none
193 ! arguments
194         REAL,    INTENT(OUT)  :: chem_bv    ! boundary value for chem(-,-,-,nch)
195         REAL,    INTENT(IN)   :: z          ! height
196         INTEGER, INTENT(IN)   :: nch        ! index number of chemical species
197     INTEGER, INTENT(IN)   :: numgas     ! index number of last gas species
198 ! local variables
199         real chem_bv_ald, chem_bv_hc3, chem_bv_hc5,   &
200              chem_bv_hc8, chem_bv_ket, chem_bv_oli,   &
201              chem_bv_olt, chem_bv_ora2
202         real, parameter :: chem_bv_def = 1.0e-20
205     if( nch == p_hcl ) then
206        !This calculation should mimic the hcl profile in
207        !cbmz_init_wrf_mixrats, above.
208        if( z <= 1000. ) then
209           chem_bv = 0.4*1e-3
210        elseif( z > 1000. &
211             .and. z <= 2500. ) then
212           chem_bv = (0.4*1e-3) + (z-1000.)* &
213                ((0.1*1e-3)-(0.4*1e-3)) / (2500.-1000.)
214        else
215           chem_bv = 0.1*1e-3
216        end if
218     else if( nch == p_par ) then
219        call bdy_chem_value( chem_bv_ald,  z, p_ald, numgas )
220 !wig, 2-May-2007: begin
221        ! The extra +1 offsets the blank first index for p_XXX
222        call bdy_chem_value( chem_bv_hc3,  z, numgas+1+1, numgas )
223        call bdy_chem_value( chem_bv_hc5,  z, numgas+1+2, numgas )
224        call bdy_chem_value( chem_bv_hc8,  z, numgas+1+3, numgas )
225 !wig, 2-May-2007: end
226        call bdy_chem_value( chem_bv_ket,  z, p_ket, numgas )
227        call bdy_chem_value( chem_bv_oli,  z, p_oli, numgas )
228        call bdy_chem_value( chem_bv_olt,  z, p_olt, numgas )
229        call bdy_chem_value( chem_bv_ora2, z, p_ora2, numgas )
231        chem_bv = 0.4*chem_bv_ald + 2.9*chem_bv_hc3    &
232             + 4.8*chem_bv_hc5 + 7.9*chem_bv_hc8       &
233             + 0.9*chem_bv_ket + 2.8*chem_bv_oli       &
234             + 1.8*chem_bv_olt + 1.0*chem_bv_ora2
236     else
237        ! chem_bv=0 for all other species
238        chem_bv = chem_bv_def
240     end if
242         return
243         end subroutine bdy_chem_value_cbmz