Update version info for release v4.6.1 (#2122)
[WRF.git] / chem / module_cbm4_initmixrats.F
blob62b35a67072528eef27217b53ccfaa71392412ec
1 module module_cbm4_initmixrats
2   
3   use module_peg_util
4   
5   
6   integer, parameter :: cbm4_init_wrf_mixrats_flagaa = 1
7   ! turns subr cbm4_init_wrf_mixrats on/off
8   
9   
10 contains
11   
12   
13   !-----------------------------------------------------------------------
14   subroutine cbm4_init_wrf_mixrats(   &
15        config_flags,              &
16        z_at_w, g,                 &
17        chem, numgas,              &
18        ids,ide, jds,jde, kds,kde, &
19        ims,ime, jms,jme, kms,kme, &
20        its,ite, jts,jte, kts,kte  )
22 ! provides initial values for cbm4 gas species
23 !    for gas species that are common to both cbm4 and radm2, the initial
24 !       values are provided via the run initialization file
25 !       (The radm2 gas species are initialized from this file
26 !       when chem_in_opt==anything.  This ought to be changed!)
27 !    for gas species that are in cbm4 but not in radm2, the initial values
28 !       are provided here
29 !    currently only hcl and "par" have non-zero initial values,
30 !       and other species are near-zero
32 ! when (gas_ic_opt == gas_ic_pnnl) AND (chem_in_opt == 0),
33 !    ozone is set to "Texas August 2000" values
35 ! setting cbm4_init_wrf_mixrats_flagaa = 1/0 turns this subr on/off.
38     USE module_configure, only:  grid_config_rec_type, num_chem, &
39          p_o3, p_ald2, p_par
40     USE module_state_description, only:  param_first_scalar,   &
41          gas_ic_pnnl
42     USE module_input_chem_data, only:  bdy_chem_value
43     
44     IMPLICIT NONE
47 !-----------------------------------------------------------------------
48 ! subr arguments
50    INTEGER, INTENT(IN) :: numgas, &
51                           ids,ide, jds,jde, kds,kde, &
52                           ims,ime, jms,jme, kms,kme, &
53                           its,ite, jts,jte, kts,kte
55    REAL, INTENT(IN) :: g
57 ! perturbation and base geopotential at layer boundaries
58    REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), &
59          INTENT(IN) :: z_at_w
61 ! advected chemical tracers
62    REAL, DIMENSION( ims:ime, kms:kme, jms:jme, num_chem ), &
63          INTENT(INOUT) :: chem
65    TYPE(grid_config_rec_type), INTENT(IN) :: config_flags
66 !-----------------------------------------------------------------------
68 !   local variables
69         integer i, j, k, kp1
70     real, dimension( its:ite, kts:kte, jts:jte ) :: z
72         if (cbm4_init_wrf_mixrats_flagaa <= 0) return
75 ! calculate the mid-level heights
77     do j = jts, min(jte,jde-1)
78        do k = kts, kte
79           kp1 = min(k+1, kte)
80           do i = its, min(ite,ide-1)
81              z(i,k,j) = (z_at_w(i,k,j)+z_at_w(i,kp1,j))*0.5
82           end do
83        end do
84     end do
87 !   when (gas_ic_opt == gas_ic_cbm4) AND (chem_in_opt == 0),
88 !   set ozone (and other gas?) to "Texas August 2000" values
90         if ( (config_flags%chem_in_opt == 0) .and.   &
91              (config_flags%gas_ic_opt == gas_ic_pnnl) ) then
92             do j = jts, min(jte,jde-1)
93             do k = kts, kte
94             do i = its, min(ite,ide-1)
95            call bdy_chem_value( chem(i,k,j,p_o3),z(i,k,j), p_o3, numgas )
96             end do
97             end do
98             end do
99         end if
103 !wig, 2-May-2007: Moved this logic to make_chem_profile so this species
104 !                 will now come in via the real.exe generated wrfinput.
105 !!$!
106 !!$!   compute par initial mixing ratio from radm2 hydrocarbon species
107 !!$!   using same formula as for par emissions
108 !!$!
109 !!$     do j = jts, min(jte,jde-1)
110 !!$     do k = kts, kte
111 !!$     do i = its, min(ite,ide-1)
112 !!$         chem(i,k,j,p_par) =                               &
113 !!$               0.4*chem(i,k,j,p_ald) + 2.9*chem(i,k,j,p_hc3)   &
114 !!$             + 4.8*chem(i,k,j,p_hc5) + 7.9*chem(i,k,j,p_hc8)   &
115 !!$             + 0.9*chem(i,k,j,p_ket) + 2.8*chem(i,k,j,p_oli)   &
116 !!$             + 1.8*chem(i,k,j,p_olt) + 1.0*chem(i,k,j,p_ora2) 
117 !!$     end do
118 !!$     end do
119 !!$     end do
121         return
122       end subroutine cbm4_init_wrf_mixrats
129 !-----------------------------------------------------------------------
130         subroutine bdy_chem_value_cbm4 ( chem_bv, z, nch, numgas )
132 ! provides boundary values for cbm4 gas species
133 !    for gas species that are common to both cbm4 and radm2, the boundary
134 !       values are provided by subr bdy_chem_value
135 !    for gas species that are in cbm4 but not in radm2, the boundary values
136 !       are provided here, except for par
137 !    currently only "hcl" has a non-zero boundary value,
138 !       and other species are near-zero
140 ! this is outside of the module declaration because of potential
141 !    module1 --> module2 --> module1 use conflicts
143         use module_configure, only:  grid_config_rec_type,   &
144             p_par, p_ald2
145         use module_input_chem_data, only:  bdy_chem_value
147         implicit none
149 ! arguments
150         REAL,    INTENT(OUT)  :: chem_bv    ! boundary value for chem(-,-,-,nch)
151         REAL,    INTENT(IN)   :: z          ! height
152         INTEGER, INTENT(IN)   :: nch        ! index number of chemical species
153     INTEGER, INTENT(IN)   :: numgas     ! index number of last gas species
154 ! local variables
155         real chem_bv_ald, chem_bv_hc3, chem_bv_hc5,   &
156              chem_bv_hc8, chem_bv_ket, chem_bv_oli,   &
157              chem_bv_olt, chem_bv_ora2
158         real, parameter :: chem_bv_def = 1.0e-20
161 !    if( nch == p_hcl ) then
162 !       !This calculation should mimic the hcl profile in
163 !       !cbm4_init_wrf_mixrats, above.
164 !       if( z <= 1000. ) then
165 !          chem_bv = 0.4*1e-3
166 !       elseif( z > 1000. &
167 !            .and. z <= 2500. ) then
168 !          chem_bv = (0.4*1e-3) + (z-1000.)* &
169 !               ((0.1*1e-3)-(0.4*1e-3)) / (2500.-1000.)
170 !       else
171 !          chem_bv = 0.1*1e-3
172 !       end if
174     if( nch == p_par ) then
175        call bdy_chem_value( chem_bv_ald,  z, p_ald2, numgas )
176 !wig, 2-May-2007: begin
177        ! The extra +1 offsets the blank first index for p_XXX
178        call bdy_chem_value( chem_bv_hc3,  z, numgas+1+1, numgas )
179        call bdy_chem_value( chem_bv_hc5,  z, numgas+1+2, numgas )
180        call bdy_chem_value( chem_bv_hc8,  z, numgas+1+3, numgas )
181 !wig, 2-May-2007: end
182        call bdy_chem_value( chem_bv_ket,  z, numgas+1+4, numgas )
183        call bdy_chem_value( chem_bv_oli,  z, numgas+1+5, numgas )
184        call bdy_chem_value( chem_bv_olt,  z, numgas+1+6, numgas )
185        call bdy_chem_value( chem_bv_ora2, z, numgas+1+7, numgas )
187        chem_bv = 0.4*chem_bv_ald + 2.9*chem_bv_hc3    &
188             + 4.8*chem_bv_hc5 + 7.9*chem_bv_hc8       &
189             + 0.9*chem_bv_ket + 2.8*chem_bv_oli       &
190             + 1.8*chem_bv_olt + 1.0*chem_bv_ora2
192     else
193        ! chem_bv=0 for all other species
194        chem_bv = chem_bv_def
196     end if
198         return
199       end subroutine bdy_chem_value_cbm4
200    
202 !-----------------------------------------------------------------------
203     end module module_cbm4_initmixrats