Line data Source code
1 : !--------------------------------------------------------------------------------
2 : ! Manages the absorber concentrations in the layers RRTMG operates
3 : ! including an extra layer over the model if needed.
4 : !
5 : ! Creator: Francis Vitt
6 : ! 9 May 2011
7 : !--------------------------------------------------------------------------------
8 : module rrtmg_state
9 :
10 : use shr_kind_mod, only: r8 => shr_kind_r8
11 : use ppgrid, only: pcols, pver, pverp
12 : use cam_history, only: outfld
13 :
14 : implicit none
15 : private
16 : save
17 :
18 : public :: rrtmg_state_t
19 : public :: rrtmg_state_init
20 : public :: rrtmg_state_create
21 : public :: rrtmg_state_update
22 : public :: rrtmg_state_destroy
23 : public :: num_rrtmg_levs
24 :
25 : type rrtmg_state_t
26 :
27 : real(r8), allocatable :: h2ovmr(:,:) ! h2o volume mixing ratio
28 : real(r8), allocatable :: o3vmr(:,:) ! o3 volume mixing ratio
29 : real(r8), allocatable :: co2vmr(:,:) ! co2 volume mixing ratio
30 : real(r8), allocatable :: ch4vmr(:,:) ! ch4 volume mixing ratio
31 : real(r8), allocatable :: o2vmr(:,:) ! o2 volume mixing ratio
32 : real(r8), allocatable :: n2ovmr(:,:) ! n2o volume mixing ratio
33 : real(r8), allocatable :: cfc11vmr(:,:) ! cfc11 volume mixing ratio
34 : real(r8), allocatable :: cfc12vmr(:,:) ! cfc12 volume mixing ratio
35 : real(r8), allocatable :: cfc22vmr(:,:) ! cfc22 volume mixing ratio
36 : real(r8), allocatable :: ccl4vmr(:,:) ! ccl4 volume mixing ratio
37 :
38 : real(r8), allocatable :: pmidmb(:,:) ! Level pressure (hPa)
39 : real(r8), allocatable :: pintmb(:,:) ! Model interface pressure (hPa)
40 : real(r8), allocatable :: tlay(:,:) ! mid point temperature
41 : real(r8), allocatable :: tlev(:,:) ! interface temperature
42 :
43 : end type rrtmg_state_t
44 :
45 : integer :: num_rrtmg_levs ! number of pressure levels greate than 1.e-4_r8 mbar
46 :
47 : real(r8), parameter :: amdw = 1.607793_r8 ! Molecular weight of dry air / water vapor
48 : real(r8), parameter :: amdc = 0.658114_r8 ! Molecular weight of dry air / carbon dioxide
49 : real(r8), parameter :: amdo = 0.603428_r8 ! Molecular weight of dry air / ozone
50 : real(r8), parameter :: amdm = 1.805423_r8 ! Molecular weight of dry air / methane
51 : real(r8), parameter :: amdn = 0.658090_r8 ! Molecular weight of dry air / nitrous oxide
52 : real(r8), parameter :: amdo2 = 0.905140_r8 ! Molecular weight of dry air / oxygen
53 : real(r8), parameter :: amdc1 = 0.210852_r8 ! Molecular weight of dry air / CFC11
54 : real(r8), parameter :: amdc2 = 0.239546_r8 ! Molecular weight of dry air / CFC12
55 :
56 : contains
57 :
58 : !--------------------------------------------------------------------------------
59 : ! sets the number of model levels RRTMG operates
60 : !--------------------------------------------------------------------------------
61 1536 : subroutine rrtmg_state_init
62 :
63 : use ref_pres,only : pref_edge
64 : implicit none
65 :
66 : ! The following cuts off RRTMG at roughly the point where it becomes
67 : ! invalid due to low pressure.
68 52224 : num_rrtmg_levs = count( pref_edge(:) > 1._r8 ) ! pascals (1.e-2 mbar)
69 :
70 1536 : end subroutine rrtmg_state_init
71 :
72 : !--------------------------------------------------------------------------------
73 : ! creates (alloacates) an rrtmg_state object
74 : !--------------------------------------------------------------------------------
75 :
76 38400 : function rrtmg_state_create( pstate, cam_in ) result( rstate )
77 1536 : use physics_types, only: physics_state
78 : use camsrfexch, only: cam_in_t
79 : use physconst, only: stebol
80 :
81 : implicit none
82 :
83 : type(physics_state), intent(in) :: pstate
84 : type(cam_in_t), intent(in) :: cam_in
85 :
86 : type(rrtmg_state_t) :: rstate
87 :
88 : real(r8) dy ! Temporary layer pressure thickness
89 : real(r8) :: tint(pcols,pverp) ! Model interface temperature
90 : integer :: ncol, i, kk, k
91 :
92 115200 : allocate( rstate%h2ovmr(pcols,num_rrtmg_levs) )
93 76800 : allocate( rstate%o3vmr(pcols,num_rrtmg_levs) )
94 76800 : allocate( rstate%co2vmr(pcols,num_rrtmg_levs) )
95 76800 : allocate( rstate%ch4vmr(pcols,num_rrtmg_levs) )
96 76800 : allocate( rstate%o2vmr(pcols,num_rrtmg_levs) )
97 76800 : allocate( rstate%n2ovmr(pcols,num_rrtmg_levs) )
98 76800 : allocate( rstate%cfc11vmr(pcols,num_rrtmg_levs) )
99 76800 : allocate( rstate%cfc12vmr(pcols,num_rrtmg_levs) )
100 76800 : allocate( rstate%cfc22vmr(pcols,num_rrtmg_levs) )
101 76800 : allocate( rstate%ccl4vmr(pcols,num_rrtmg_levs) )
102 :
103 76800 : allocate( rstate%pmidmb(pcols,num_rrtmg_levs) )
104 115200 : allocate( rstate%pintmb(pcols,num_rrtmg_levs+1) )
105 76800 : allocate( rstate%tlay(pcols,num_rrtmg_levs) )
106 76800 : allocate( rstate%tlev(pcols,num_rrtmg_levs+1) )
107 :
108 38400 : ncol = pstate%ncol
109 :
110 : ! Calculate interface temperatures (following method
111 : ! used in radtpl for the longwave), using surface upward flux and
112 : ! stebol constant in mks units
113 591360 : do i = 1,ncol
114 552960 : tint(i,1) = pstate%t(i,1)
115 552960 : tint(i,pverp) = sqrt(sqrt(cam_in%lwup(i)/stebol))
116 17733120 : do k = 2,pver
117 17141760 : dy = (pstate%lnpint(i,k) - pstate%lnpmid(i,k)) / (pstate%lnpmid(i,k-1) - pstate%lnpmid(i,k))
118 17694720 : tint(i,k) = pstate%t(i,k) - dy * (pstate%t(i,k) - pstate%t(i,k-1))
119 : end do
120 : end do
121 :
122 1305600 : do k = 1, num_rrtmg_levs
123 :
124 1267200 : kk = max(k + (pverp-num_rrtmg_levs)-1,1)
125 :
126 19514880 : rstate%pmidmb(:ncol,k) = pstate%pmid(:ncol,kk) * 1.e-2_r8
127 19514880 : rstate%pintmb(:ncol,k) = pstate%pint(:ncol,kk) * 1.e-2_r8
128 :
129 19514880 : rstate%tlay(:ncol,k) = pstate%t(:ncol,kk)
130 19553280 : rstate%tlev(:ncol,k) = tint(:ncol,kk)
131 :
132 : enddo
133 :
134 : ! bottom interface
135 591360 : rstate%pintmb(:ncol,num_rrtmg_levs+1) = pstate%pint(:ncol,pverp) * 1.e-2_r8 ! mbar
136 591360 : rstate%tlev(:ncol,num_rrtmg_levs+1) = tint(:ncol,pverp)
137 :
138 : ! top layer thickness
139 38400 : if (num_rrtmg_levs==pverp) then
140 591360 : rstate%pmidmb(:ncol,1) = 0.5_r8 * rstate%pintmb(:ncol,2)
141 591360 : rstate%pintmb(:ncol,1) = 1.e-4_r8 ! mbar
142 : endif
143 :
144 76800 : endfunction rrtmg_state_create
145 :
146 : !--------------------------------------------------------------------------------
147 : ! updates the concentration fields
148 : !--------------------------------------------------------------------------------
149 76800 : subroutine rrtmg_state_update(pstate,pbuf,icall,rstate)
150 38400 : use physics_types, only: physics_state
151 : use physics_buffer, only: physics_buffer_desc
152 : use rad_constituents, only: rad_cnst_get_gas
153 :
154 : implicit none
155 :
156 : type(physics_state), intent(in), target :: pstate
157 : type(physics_buffer_desc), pointer :: pbuf(:)
158 : integer, intent(in) :: icall ! index through climate/diagnostic radiation calls
159 : type(rrtmg_state_t) :: rstate
160 :
161 76800 : real(r8), pointer, dimension(:,:) :: sp_hum ! specific humidity
162 76800 : real(r8), pointer, dimension(:,:) :: n2o ! nitrous oxide mass mixing ratio
163 76800 : real(r8), pointer, dimension(:,:) :: ch4 ! methane mass mixing ratio
164 76800 : real(r8), pointer, dimension(:,:) :: o2 ! O2 mass mixing ratio
165 76800 : real(r8), pointer, dimension(:,:) :: cfc11 ! cfc11 mass mixing ratio
166 76800 : real(r8), pointer, dimension(:,:) :: cfc12 ! cfc12 mass mixing ratio
167 76800 : real(r8), pointer, dimension(:,:) :: o3 ! Ozone mass mixing ratio
168 76800 : real(r8), pointer, dimension(:,:) :: co2 ! co2 mass mixing ratio
169 :
170 : integer :: ncol, i, kk, k, lchnk
171 : real(r8) :: H, P_top, P_surface
172 : real(r8), dimension(pcols) :: P_int, P_mid, alpha, beta, a, b, chi_mid, chi_0, chi_eff
173 :
174 76800 : ncol = pstate%ncol
175 76800 : lchnk = pstate%lchnk
176 :
177 : ! Get specific humidity
178 76800 : call rad_cnst_get_gas(icall,'H2O', pstate, pbuf, sp_hum)
179 : ! Get oxygen mass mixing ratio.
180 76800 : call rad_cnst_get_gas(icall,'O2', pstate, pbuf, o2)
181 : ! Get ozone mass mixing ratio.
182 76800 : call rad_cnst_get_gas(icall,'O3', pstate, pbuf, o3)
183 : ! Get CO2 mass mixing ratio
184 76800 : call rad_cnst_get_gas(icall,'CO2', pstate, pbuf, co2)
185 : ! Get N2O mass mixing ratio
186 76800 : call rad_cnst_get_gas(icall,'N2O', pstate, pbuf, n2o)
187 : ! Get CH4 mass mixing ratio
188 76800 : call rad_cnst_get_gas(icall,'CH4', pstate, pbuf, ch4)
189 : ! Get CFC mass mixing ratios
190 76800 : call rad_cnst_get_gas(icall,'CFC11', pstate, pbuf, cfc11)
191 76800 : call rad_cnst_get_gas(icall,'CFC12', pstate, pbuf, cfc12)
192 :
193 2611200 : do k = 1, num_rrtmg_levs
194 :
195 2534400 : kk = max(k + (pverp-num_rrtmg_levs)-1,1)
196 :
197 39029760 : rstate%ch4vmr(:ncol,k) = ch4(:ncol,kk) * amdm
198 39029760 : rstate%h2ovmr(:ncol,k) = (sp_hum(:ncol,kk) / (1._r8 - sp_hum(:ncol,kk))) * amdw
199 39029760 : rstate%o3vmr(:ncol,k) = o3(:ncol,kk) * amdo
200 39029760 : rstate%co2vmr(:ncol,k) = co2(:ncol,kk) * amdc
201 39029760 : rstate%ch4vmr(:ncol,k) = ch4(:ncol,kk) * amdm
202 39029760 : rstate%o2vmr(:ncol,k) = o2(:ncol,kk) * amdo2
203 39029760 : rstate%n2ovmr(:ncol,k) = n2o(:ncol,kk) * amdn
204 39029760 : rstate%cfc11vmr(:ncol,k) = cfc11(:ncol,kk) * amdc1
205 39029760 : rstate%cfc12vmr(:ncol,k) = cfc12(:ncol,kk) * amdc2
206 39029760 : rstate%cfc22vmr(:ncol,k) = 0._r8
207 39106560 : rstate%ccl4vmr(:ncol,k) = 0._r8
208 :
209 : enddo
210 :
211 : ! For the purpose of attenuating solar fluxes above the CAM model top, we assume that ozone
212 : ! mixing decreases linearly in each column from the value in the top layer of CAM to zero at
213 : ! the pressure level set by P_top. P_top has been set to 50 Pa (0.5 hPa) based on model tuning
214 : ! to produce temperatures at the top of CAM that are most consistent with WACCM at similar pressure levels.
215 :
216 76800 : P_top = 50.0_r8 ! pressure (Pa) at which we assume O3 = 0 in linear decay from CAM top
217 1182720 : P_int(:ncol) = pstate%pint(:ncol,1) ! pressure (Pa) at upper interface of CAM
218 1182720 : P_mid(:ncol) = pstate%pmid(:ncol,1) ! pressure (Pa) at midpoint of top layer of CAM
219 76800 : alpha(:) = 0.0_r8
220 76800 : beta(:) = 0.0_r8
221 1182720 : alpha(:ncol) = log(P_int(:ncol)/P_top)
222 1182720 : beta(:ncol) = log(P_mid(:ncol)/P_int(:ncol))/log(P_mid(:ncol)/P_top)
223 :
224 1182720 : a(:ncol) = ( (1._r8 + alpha(:ncol)) * exp(-alpha(:ncol)) - 1._r8 ) / alpha(:ncol)
225 1182720 : b(:ncol) = 1_r8 - exp(-alpha(:ncol))
226 :
227 11289600 : where(alpha .gt. 0) ! only apply where top level is below 80 km
228 76800 : chi_mid(:) = o3(:,1)*amdo ! molar mixing ratio of O3 at midpoint of top layer
229 : chi_0(:) = chi_mid(:) / (1._r8 + beta(:))
230 : chi_eff(:) = chi_0(:) * (a(:) + b(:))
231 76800 : rstate%o3vmr(:,1) = chi_eff(:)
232 : chi_eff(:) = chi_eff(:) * P_int(:) / amdo / 9.8_r8 ! O3 column above in kg m-2
233 : chi_eff(:) = chi_eff(:) / 2.1415e-5_r8 ! O3 column above in DU
234 : endwhere
235 :
236 76800 : call outfld('O3colAbove', chi_eff(:ncol), ncol, lchnk)
237 :
238 153600 : end subroutine rrtmg_state_update
239 :
240 : !--------------------------------------------------------------------------------
241 : ! de-allocates an rrtmg_state object
242 : !--------------------------------------------------------------------------------
243 38400 : subroutine rrtmg_state_destroy(rstate)
244 :
245 : implicit none
246 :
247 : type(rrtmg_state_t) :: rstate
248 :
249 38400 : deallocate(rstate%h2ovmr)
250 38400 : deallocate(rstate%o3vmr)
251 38400 : deallocate(rstate%co2vmr)
252 38400 : deallocate(rstate%ch4vmr)
253 38400 : deallocate(rstate%o2vmr)
254 38400 : deallocate(rstate%n2ovmr)
255 38400 : deallocate(rstate%cfc11vmr)
256 38400 : deallocate(rstate%cfc12vmr)
257 38400 : deallocate(rstate%cfc22vmr)
258 38400 : deallocate(rstate%ccl4vmr)
259 :
260 38400 : deallocate(rstate%pmidmb)
261 38400 : deallocate(rstate%pintmb)
262 38400 : deallocate(rstate%tlay)
263 38400 : deallocate(rstate%tlev)
264 :
265 76800 : endsubroutine rrtmg_state_destroy
266 :
267 0 : end module rrtmg_state
|