Line data Source code
1 : !================================================================================================
2 : ! output/input data necessary to drive radiation offline
3 : ! Francis Vitt -- Created 15 Dec 2009
4 : !================================================================================================
5 : module radiation_data
6 :
7 : use shr_kind_mod, only: r8=>shr_kind_r8
8 : use ppgrid, only : pcols, pver, pverp, begchunk, endchunk
9 : use cam_history, only: addfld, add_default, horiz_only, outfld
10 : use rad_constituents, only: rad_cnst_get_info, rad_cnst_get_gas, rad_cnst_get_aer_mmr
11 : use radconstants, only: nradgas, gaslist
12 : use cam_history_support, only: fieldname_len, fillvalue
13 : use spmd_utils, only: masterproc
14 : use cam_abortutils, only: endrun
15 :
16 : use drv_input_data,only: drv_input_data_t
17 : use drv_input_data,only: drv_input_data_get, drv_input_2d_t, drv_input_3d_t, drv_input_4d_t, drv_input_2di_t
18 : use physics_types, only: physics_state
19 : use physics_buffer,only: physics_buffer_desc, pbuf_get_chunk
20 :
21 : implicit none
22 : private
23 :
24 : public :: rad_data_readnl
25 : public :: rad_data_register
26 : public :: rad_data_init
27 : public :: rad_data_write
28 : public :: rad_data_read
29 : public :: rad_data_enable
30 :
31 : integer :: volcgeom_ifld, volcgeom1_ifld, volcgeom2_ifld, volcgeom3_ifld
32 : integer :: cld_ifld,rel_ifld,rei_ifld
33 : integer :: dei_ifld,mu_ifld,lambdac_ifld,iciwp_ifld,iclwp_ifld
34 : integer :: des_ifld,icswp_ifld,cldfsnow_ifld
35 : integer :: cldemis_ifld, cldtau_ifld, cicewp_ifld, cliqwp_ifld, nmxrgn_ifld, pmxrgn_ifld
36 : integer :: qrs_ifld, qrl_ifld
37 : integer :: dgnumwet_ifld, qaerwat_ifld
38 :
39 : character(len=fieldname_len), parameter :: &
40 : lndfrc_fldn = 'rad_lndfrc ' , &
41 : icefrc_fldn = 'rad_icefrc ' , &
42 : snowh_fldn = 'rad_snowh ' , &
43 : asdir_fldn = 'rad_asdir ' , &
44 : asdif_fldn = 'rad_asdif ' , &
45 : aldir_fldn = 'rad_aldir ' , &
46 : aldif_fldn = 'rad_aldif ' , &
47 : coszen_fldn = 'rad_coszen ' , &
48 : asdir_pos_fldn = 'rad_asdir_pos ' , &
49 : asdif_pos_fldn = 'rad_asdif_pos ' , &
50 : aldir_pos_fldn = 'rad_aldir_pos ' , &
51 : aldif_pos_fldn = 'rad_aldif_pos ' , &
52 : lwup_fldn = 'rad_lwup ' , &
53 : ts_fldn = 'rad_ts ' , &
54 : temp_fldn = 'rad_temp ' , &
55 : pdel_fldn = 'rad_pdel ' , &
56 : pdeldry_fldn = 'rad_pdeldry ' , &
57 : pmid_fldn = 'rad_pmid ' , &
58 : watice_fldn = 'rad_watice ' , &
59 : watliq_fldn = 'rad_watliq ' , &
60 : watvap_fldn = 'rad_watvap ' , &
61 : zint_fldn = 'rad_zint ' , &
62 : pint_fldn = 'rad_pint ' , &
63 : cld_fldn = 'rad_cld ' , &
64 : cldemis_fldn = 'rad_cldemis ' , &
65 : cldtau_fldn = 'rad_cldtau ' , &
66 : cicewp_fldn = 'rad_cicewp ' , &
67 : cliqwp_fldn = 'rad_cliqwp ' , &
68 : nmxrgn_fldn = 'rad_nmxrgn ' , &
69 : pmxrgn_fldn = 'rad_pmxrgn ' , &
70 : cldfsnow_fldn = 'rad_cldfsnow ' , &
71 : rel_fldn = 'rad_rel ' , &
72 : rei_fldn = 'rad_rei ' , &
73 : dei_fldn = 'rad_dei ' , &
74 : des_fldn = 'rad_des ' , &
75 : mu_fldn = 'rad_mu ' , &
76 : lambdac_fldn = 'rad_lambdac ' , &
77 : iciwp_fldn = 'rad_iciwp ' , &
78 : iclwp_fldn = 'rad_iclwp ' , &
79 : icswp_fldn = 'rad_icswp ' , &
80 : qrs_fldn = 'rad_qrs ' , &
81 : qrl_fldn = 'rad_qrl ' , &
82 : volcgeom_fldn = 'rad_volc_geom ' , &
83 : volcgeom1_fldn = 'rad_volc_geom1 ' , &
84 : volcgeom2_fldn = 'rad_volc_geom2 ' , &
85 : volcgeom3_fldn = 'rad_volc_geom3 '
86 :
87 : ! for modal aerosols
88 : character(len=fieldname_len), allocatable :: dgnumwet_fldn(:)
89 : character(len=fieldname_len), allocatable :: qaerwat_fldn(:)
90 : integer :: nmodes=0
91 :
92 : ! rad constituents mixing ratios
93 : integer :: ngas, naer=0
94 : character(len=64), allocatable :: gasnames(:)
95 : character(len=64), allocatable :: aernames(:)
96 :
97 : ! control options
98 : logical :: rad_data_output = .false.
99 : integer :: rad_data_histfile_num = 2
100 : character(len=1) :: rad_data_avgflag = 'I'
101 :
102 : ! MG microphys check
103 : logical :: mg_microphys
104 :
105 : logical :: fixed_dyn_heating = .false.
106 :
107 : ! for fixed dynamical heating ...
108 : logical :: do_fdh = .false.
109 :
110 : integer :: tcorr_idx = -1
111 : integer :: qrsin_idx = -1
112 : integer :: qrlin_idx = -1
113 : integer :: tropp_idx = -1
114 :
115 : logical :: enabled = .false.
116 : logical :: gmean_3modes = .false.
117 :
118 : contains
119 :
120 :
121 : !================================================================================================
122 : !================================================================================================
123 1536 : subroutine rad_data_register
124 : use physics_buffer, only: pbuf_add_field, dtype_r8
125 :
126 1536 : if ( do_fdh ) then
127 0 : call pbuf_add_field('tcorr', 'global', dtype_r8, (/pcols,pver/), tcorr_idx)
128 0 : call pbuf_add_field('qrsin', 'physpkg', dtype_r8, (/pcols,pver/), qrsin_idx)
129 0 : call pbuf_add_field('qrlin', 'physpkg', dtype_r8, (/pcols,pver/), qrlin_idx)
130 0 : call pbuf_add_field('tropp', 'physpkg', dtype_r8, (/pcols/), tropp_idx)
131 : end if
132 :
133 1536 : end subroutine rad_data_register
134 :
135 : !================================================================================================
136 : !================================================================================================
137 1536 : subroutine rad_data_readnl(nlfile)
138 :
139 : ! Read rad_data_nl namelist group. Parse input.
140 :
141 1536 : use namelist_utils, only: find_group_name
142 : use units, only: getunit, freeunit
143 : use mpishorthand
144 :
145 : character(len=*), intent(in) :: nlfile ! filepath for file containing namelist input
146 : logical :: rad_data_fdh = .false.
147 :
148 : ! Local variables
149 : integer :: unitn, ierr, i
150 : character(len=*), parameter :: subname = 'rad_data_readnl'
151 :
152 : namelist /rad_data_nl/ rad_data_output, rad_data_histfile_num, rad_data_avgflag, rad_data_fdh
153 :
154 : !-----------------------------------------------------------------------------
155 :
156 1536 : if (masterproc) then
157 2 : unitn = getunit()
158 2 : open( unitn, file=trim(nlfile), status='old' )
159 2 : call find_group_name(unitn, 'rad_data_nl', status=ierr)
160 2 : if (ierr == 0) then
161 0 : read(unitn, rad_data_nl, iostat=ierr)
162 0 : if (ierr /= 0) then
163 0 : call endrun(subname // ':: ERROR reading namelist')
164 : end if
165 : end if
166 2 : close(unitn)
167 2 : call freeunit(unitn)
168 : end if
169 :
170 : #ifdef SPMD
171 : ! Broadcast namelist variables
172 1536 : call mpibcast (rad_data_output, 1, mpilog , 0, mpicom)
173 1536 : call mpibcast (rad_data_fdh, 1, mpilog , 0, mpicom)
174 1536 : call mpibcast (rad_data_histfile_num, 1, mpiint , 0, mpicom)
175 1536 : call mpibcast (rad_data_avgflag, 1, mpichar , 0, mpicom)
176 : #endif
177 1536 : do_fdh = rad_data_fdh
178 1536 : enabled = rad_data_output
179 :
180 1536 : end subroutine rad_data_readnl
181 :
182 : !================================================================================================
183 : !================================================================================================
184 0 : subroutine rad_data_enable()
185 0 : enabled = .true.
186 0 : end subroutine rad_data_enable
187 :
188 : !================================================================================================
189 : !================================================================================================
190 1536 : subroutine rad_data_init( pbuf2d )
191 : use phys_control, only: phys_getopts
192 : use physics_buffer, only: pbuf_get_index
193 : use physics_buffer, only: pbuf_set_field
194 : implicit none
195 :
196 : ! args
197 : type(physics_buffer_desc), pointer :: pbuf2d(:,:)
198 :
199 : ! local vars
200 : integer :: i, err
201 : integer :: m,l, nspec
202 : character(len=64) :: name
203 : character(len=32) :: aername = ' '
204 : character(len=128):: long_name
205 : character(len=64) :: long_name_description
206 : character(len=16) :: microp_scheme ! microphysics scheme
207 : character(len=16) :: rad_scheme
208 :
209 1536 : if (.not.enabled) return
210 :
211 0 : call phys_getopts(microp_scheme_out=microp_scheme, radiation_scheme_out=rad_scheme)
212 0 : mg_microphys = (trim(microp_scheme) == 'MG')
213 :
214 0 : volcgeom_ifld = pbuf_get_index('VOLC_RAD_GEOM',errcode=err) ! might need 3 more for 3-mode inputs
215 0 : volcgeom1_ifld = pbuf_get_index('VOLC_RAD_GEOM1',errcode=err) ! might need 3 more for 3-mode inputs
216 0 : volcgeom2_ifld = pbuf_get_index('VOLC_RAD_GEOM2',errcode=err) ! might need 3 more for 3-mode inputs
217 0 : volcgeom3_ifld = pbuf_get_index('VOLC_RAD_GEOM3',errcode=err) ! might need 3 more for 3-mode inputs
218 :
219 0 : gmean_3modes = volcgeom1_ifld > 0 .and. volcgeom2_ifld > 0 .and. volcgeom3_ifld > 0 .and. volcgeom_ifld < 1
220 :
221 0 : if (volcgeom_ifld > 0) then
222 0 : call addfld(volcgeom_fldn, (/ 'lev' /), 'I','m', 'combined volcanic aerosol geometric-mode radius' )
223 : endif
224 0 : if (gmean_3modes) then
225 0 : call addfld(volcgeom1_fldn, (/ 'lev' /), 'I','m', 'mode 1 volcanic aerosol geometric-mode radius' )
226 0 : call addfld(volcgeom2_fldn, (/ 'lev' /), 'I','m', 'mode 2 volcanic aerosol geometric-mode radius' )
227 0 : call addfld(volcgeom3_fldn, (/ 'lev' /), 'I','m', 'mode 3 volcanic aerosol geometric-mode radius' )
228 : endif
229 :
230 0 : cld_ifld = pbuf_get_index('CLD')
231 0 : rel_ifld = pbuf_get_index('REL')
232 0 : rei_ifld = pbuf_get_index('REI')
233 0 : qrs_ifld = pbuf_get_index('QRS')
234 0 : qrl_ifld = pbuf_get_index('QRL')
235 0 : if (mg_microphys) then
236 0 : dei_ifld = pbuf_get_index('DEI')
237 0 : des_ifld = pbuf_get_index('DES')
238 0 : mu_ifld = pbuf_get_index('MU')
239 0 : lambdac_ifld = pbuf_get_index('LAMBDAC')
240 0 : iciwp_ifld = pbuf_get_index('ICIWP')
241 0 : iclwp_ifld = pbuf_get_index('ICLWP')
242 0 : icswp_ifld = pbuf_get_index('ICSWP')
243 0 : cldfsnow_ifld = pbuf_get_index('CLDFSNOW')
244 : else
245 0 : cldemis_ifld= pbuf_get_index('CLDEMIS')
246 0 : cldtau_ifld = pbuf_get_index('CLDTAU')
247 0 : cicewp_ifld = pbuf_get_index('CICEWP')
248 0 : cliqwp_ifld = pbuf_get_index('CLIQWP')
249 0 : nmxrgn_ifld = pbuf_get_index('NMXRGN')
250 0 : pmxrgn_ifld = pbuf_get_index('PMXRGN')
251 : endif
252 :
253 0 : if ( do_fdh ) then
254 0 : call addfld('rad_TROP_P', horiz_only, 'A','Pa','Pressure at which tropopause is defined for radiation' )
255 0 : call addfld('rad_TCORR', (/ 'lev' /), 'I','K', 'Fixed dynamical heating temperature correction ' )
256 0 : call pbuf_set_field(pbuf2d, tcorr_idx, 0.0_r8)
257 0 : call pbuf_set_field(pbuf2d, qrsin_idx, 0.0_r8)
258 0 : call pbuf_set_field(pbuf2d, qrlin_idx, 0.0_r8)
259 0 : call pbuf_set_field(pbuf2d, tropp_idx, -1.0_r8)
260 : endif
261 :
262 0 : call rad_cnst_get_info(0, ngas=ngas, naero=naer, nmodes=nmodes)
263 :
264 : ! The code to output the gases assumes that the rad_constituents module has
265 : ! ordered them in the same way that they are ordered in the "gaslist" array
266 : ! in module radconstants, and that there are nradgas of them. This ordering
267 : ! is performed in the internal init_lists routine in rad_constituents.
268 0 : if (ngas /= nradgas) then
269 0 : call endrun('rad_data_init: ERROR: ngas /= nradgas')
270 : end if
271 :
272 0 : allocate( gasnames(ngas) )
273 0 : call rad_cnst_get_info(0, gasnames=gasnames)
274 :
275 0 : if (naer > 0) then
276 0 : allocate( aernames(naer) )
277 0 : call rad_cnst_get_info(0, aernames=aernames)
278 : endif
279 :
280 0 : if (nmodes>0) then
281 0 : allocate(dgnumwet_fldn(nmodes),qaerwat_fldn(nmodes))
282 0 : dgnumwet_ifld = pbuf_get_index('DGNUMWET')
283 0 : qaerwat_ifld = pbuf_get_index('QAERWAT')
284 0 : do m = 1, nmodes
285 0 : write(dgnumwet_fldn(m), 1001 ) m
286 0 : write(qaerwat_fldn(m), 1003 ) m
287 : enddo
288 : endif
289 :
290 0 : if (.not.rad_data_output) return
291 :
292 : call addfld (lndfrc_fldn, horiz_only, rad_data_avgflag, 'fraction', &
293 0 : 'radiation input: land fraction')
294 : call addfld (icefrc_fldn, horiz_only, rad_data_avgflag, 'fraction', &
295 0 : 'radiation input: ice fraction')
296 : call addfld (snowh_fldn, horiz_only, rad_data_avgflag, 'm', &
297 0 : 'radiation input: water equivalent snow depth')
298 : call addfld (asdir_fldn, horiz_only, rad_data_avgflag, '1', &
299 0 : 'radiation input: s hort wave direct albedo', flag_xyfill=.true.)
300 : call addfld (asdif_fldn, horiz_only, rad_data_avgflag, '1', &
301 0 : 'radiation input: s hort wave difuse albedo', flag_xyfill=.true.)
302 : call addfld (aldir_fldn, horiz_only, rad_data_avgflag, '1', &
303 0 : 'radiation input: long wave direct albedo', flag_xyfill=.true.)
304 : call addfld (aldif_fldn, horiz_only, rad_data_avgflag, '1', &
305 0 : 'radiation input: long wave difuse albedo', flag_xyfill=.true.)
306 :
307 : call addfld (coszen_fldn, horiz_only, rad_data_avgflag, '1', &
308 0 : 'radiation input: cosine solar zenith when positive', flag_xyfill=.true.)
309 : call addfld (asdir_pos_fldn, horiz_only, rad_data_avgflag, '1', &
310 0 : 'radiation input: s hort wave direct albedo weighted by coszen', flag_xyfill=.true.)
311 : call addfld (asdif_pos_fldn, horiz_only, rad_data_avgflag, '1', &
312 0 : 'radiation input: s hort wave difuse albedo weighted by coszen', flag_xyfill=.true.)
313 : call addfld (aldir_pos_fldn, horiz_only, rad_data_avgflag, '1', &
314 0 : 'radiation input: long wave direct albedo weighted by coszen', flag_xyfill=.true.)
315 : call addfld (aldif_pos_fldn, horiz_only, rad_data_avgflag, '1', &
316 0 : 'radiation input: long wave difuse albedo weighted by coszen', flag_xyfill=.true.)
317 :
318 : call addfld (lwup_fldn, horiz_only, rad_data_avgflag, 'W/m2', &
319 0 : 'radiation input: long wave up radiation flux ')
320 : call addfld (ts_fldn, horiz_only, rad_data_avgflag, 'K', &
321 0 : 'radiation input: surface temperature')
322 :
323 : call addfld (temp_fldn, (/ 'lev' /), rad_data_avgflag, 'K', &
324 0 : 'radiation input: midpoint temperature')
325 : call addfld (pdel_fldn, (/ 'lev' /), rad_data_avgflag, 'Pa', &
326 0 : 'radiation input: pressure layer thickness')
327 : call addfld (pdeldry_fldn, (/ 'lev' /), rad_data_avgflag,'Pa', &
328 0 : 'radiation input: dry pressure layer thickness')
329 : call addfld (pmid_fldn, (/ 'lev' /), rad_data_avgflag, 'Pa', &
330 0 : 'radiation input: midpoint pressure')
331 : call addfld (watice_fldn, (/ 'lev' /), rad_data_avgflag, 'kg/kg', &
332 0 : 'radiation input: cloud ice')
333 : call addfld (watliq_fldn, (/ 'lev' /), rad_data_avgflag, 'kg/kg', &
334 0 : 'radiation input: cloud liquid water')
335 : call addfld (watvap_fldn, (/ 'lev' /), rad_data_avgflag, 'kg/kg', &
336 0 : 'radiation input: water vapor')
337 :
338 : call addfld (zint_fldn, (/ 'ilev' /), rad_data_avgflag, 'km', &
339 0 : 'radiation input: interface height')
340 : call addfld (pint_fldn, (/ 'ilev' /), rad_data_avgflag, 'Pa', &
341 0 : 'radiation input: interface pressure')
342 :
343 : call addfld (cld_fldn, (/ 'lev' /), rad_data_avgflag, 'fraction', &
344 0 : 'radiation input: cloud fraction')
345 : call addfld (rel_fldn, (/ 'lev' /), rad_data_avgflag, 'micron', &
346 0 : 'radiation input: effective liquid drop radius')
347 : call addfld (rei_fldn, (/ 'lev' /), rad_data_avgflag, 'micron', &
348 0 : 'radiation input: effective ice partical radius')
349 : call addfld (qrs_fldn, (/ 'lev' /), rad_data_avgflag, 'J/s/kg', &
350 0 : 'radiation input: solar heating rate')
351 : call addfld (qrl_fldn, (/ 'lev' /), rad_data_avgflag, 'J/s/kg', &
352 0 : 'radiation input: longwave heating rate')
353 :
354 0 : if (mg_microphys) then
355 : call addfld (dei_fldn, (/ 'lev' /), rad_data_avgflag, 'micron', &
356 0 : 'radiation input: effective ice partical diameter')
357 : call addfld (des_fldn, (/ 'lev' /), rad_data_avgflag, 'micron', &
358 0 : 'radiation input: effective snow partical diameter')
359 : call addfld (mu_fldn, (/ 'lev' /), rad_data_avgflag, ' ', &
360 0 : 'radiation input: ice gamma parameter for optics (radiation)')
361 : call addfld (lambdac_fldn, (/ 'lev' /), rad_data_avgflag,' ', &
362 0 : 'radiation input: slope of droplet distribution for optics (radiation)')
363 : call addfld (iciwp_fldn, (/ 'lev' /), rad_data_avgflag, 'kg/m2', &
364 0 : 'radiation input: In-cloud ice water path')
365 : call addfld (iclwp_fldn, (/ 'lev' /), rad_data_avgflag, 'kg/m2', &
366 0 : 'radiation input: In-cloud liquid water path')
367 : call addfld (icswp_fldn, (/ 'lev' /), rad_data_avgflag, 'kg/m2', &
368 0 : 'radiation input: In-cloud snow water path')
369 : call addfld (cldfsnow_fldn, (/ 'lev' /), rad_data_avgflag, 'fraction', &
370 0 : 'radiation input: cloud liquid drops + snow')
371 : else
372 :
373 : call addfld (nmxrgn_fldn, horiz_only, rad_data_avgflag, ' ', &
374 0 : 'radiation input: ')
375 : call addfld (pmxrgn_fldn, (/ 'ilev' /), rad_data_avgflag, 'Pa', &
376 0 : 'radiation input: ')
377 : call addfld (cldemis_fldn, (/ 'lev' /), rad_data_avgflag, ' ', &
378 0 : 'radiation input: cloud property ')
379 : call addfld (cldtau_fldn, (/ 'lev' /), rad_data_avgflag, ' ', &
380 0 : 'radiation input: cloud property ')
381 : call addfld (cicewp_fldn, (/ 'lev' /), rad_data_avgflag, ' ', &
382 0 : 'radiation input: cloud property ')
383 : call addfld (cliqwp_fldn, (/ 'lev' /), rad_data_avgflag, ' ', &
384 0 : 'radiation input: cloud property ')
385 : endif
386 :
387 0 : call add_default (lndfrc_fldn, rad_data_histfile_num, ' ')
388 0 : call add_default (icefrc_fldn, rad_data_histfile_num, ' ')
389 0 : call add_default (snowh_fldn, rad_data_histfile_num, ' ')
390 0 : call add_default (asdir_fldn, rad_data_histfile_num, ' ')
391 0 : call add_default (asdif_fldn, rad_data_histfile_num, ' ')
392 0 : call add_default (aldir_fldn, rad_data_histfile_num, ' ')
393 0 : call add_default (aldif_fldn, rad_data_histfile_num, ' ')
394 :
395 0 : call add_default (coszen_fldn, rad_data_histfile_num, ' ')
396 0 : call add_default (asdir_pos_fldn, rad_data_histfile_num, ' ')
397 0 : call add_default (asdif_pos_fldn, rad_data_histfile_num, ' ')
398 0 : call add_default (aldir_pos_fldn, rad_data_histfile_num, ' ')
399 0 : call add_default (aldif_pos_fldn, rad_data_histfile_num, ' ')
400 :
401 0 : call add_default (lwup_fldn, rad_data_histfile_num, ' ')
402 0 : call add_default (ts_fldn, rad_data_histfile_num, ' ')
403 0 : call add_default (temp_fldn, rad_data_histfile_num, ' ')
404 0 : call add_default (pdel_fldn, rad_data_histfile_num, ' ')
405 0 : call add_default (pdeldry_fldn, rad_data_histfile_num, ' ')
406 0 : call add_default (pmid_fldn, rad_data_histfile_num, ' ')
407 0 : call add_default (watice_fldn, rad_data_histfile_num, ' ')
408 0 : call add_default (watliq_fldn, rad_data_histfile_num, ' ')
409 0 : call add_default (watvap_fldn, rad_data_histfile_num, ' ')
410 0 : call add_default (zint_fldn, rad_data_histfile_num, ' ')
411 0 : call add_default (pint_fldn, rad_data_histfile_num, ' ')
412 :
413 0 : call add_default (cld_fldn, rad_data_histfile_num, ' ')
414 0 : call add_default (rel_fldn, rad_data_histfile_num, ' ')
415 0 : call add_default (rei_fldn, rad_data_histfile_num, ' ')
416 0 : call add_default (qrs_fldn, rad_data_histfile_num, ' ')
417 0 : call add_default (qrl_fldn, rad_data_histfile_num, ' ')
418 :
419 0 : if (mg_microphys) then
420 0 : call add_default (dei_fldn, rad_data_histfile_num, ' ')
421 0 : call add_default (des_fldn, rad_data_histfile_num, ' ')
422 0 : call add_default (mu_fldn, rad_data_histfile_num, ' ')
423 0 : call add_default (lambdac_fldn, rad_data_histfile_num, ' ')
424 0 : call add_default (iciwp_fldn, rad_data_histfile_num, ' ')
425 0 : call add_default (iclwp_fldn, rad_data_histfile_num, ' ')
426 0 : call add_default (icswp_fldn, rad_data_histfile_num, ' ')
427 0 : call add_default (cldfsnow_fldn, rad_data_histfile_num, ' ')
428 : else
429 0 : call add_default (cldemis_fldn, rad_data_histfile_num, ' ')
430 0 : call add_default (cldtau_fldn, rad_data_histfile_num, ' ')
431 0 : call add_default (cicewp_fldn, rad_data_histfile_num, ' ')
432 0 : call add_default (cliqwp_fldn, rad_data_histfile_num, ' ')
433 0 : call add_default (nmxrgn_fldn, rad_data_histfile_num, ' ')
434 0 : call add_default (pmxrgn_fldn, rad_data_histfile_num, ' ')
435 : endif
436 :
437 : ! stratospheric aersols geometric mean radius
438 0 : if (volcgeom_ifld > 0) then
439 0 : call add_default (volcgeom_fldn, rad_data_histfile_num, ' ')
440 : endif
441 0 : if (gmean_3modes) then
442 0 : call add_default (volcgeom1_fldn, rad_data_histfile_num, ' ')
443 0 : call add_default (volcgeom2_fldn, rad_data_histfile_num, ' ')
444 0 : call add_default (volcgeom3_fldn, rad_data_histfile_num, ' ')
445 : endif
446 :
447 : ! rad constituents
448 :
449 0 : long_name_description = ' mass mixing ratio used in rad climate calculation'
450 :
451 0 : do i = 1, ngas
452 0 : long_name = trim(gasnames(i))//trim(long_name_description)
453 0 : name = 'rad_'//gasnames(i)
454 0 : call addfld(trim(name), (/ 'lev' /), rad_data_avgflag, 'kg/kg', trim(long_name))
455 0 : call add_default (trim(name), rad_data_histfile_num, ' ')
456 : end do
457 :
458 0 : if (naer > 0) then
459 :
460 0 : do i = 1, naer
461 0 : long_name = trim(aernames(i))//trim(long_name_description)
462 0 : name = 'rad_'//aernames(i)
463 0 : call addfld(trim(name), (/ 'lev' /), rad_data_avgflag, 'kg/kg', trim(long_name))
464 0 : call add_default (trim(name), rad_data_histfile_num, ' ')
465 : end do
466 : endif
467 0 : if (nmodes>0) then
468 :
469 : ! for modal aerosol model
470 0 : do m = 1, nmodes
471 0 : write(long_name, 1002) m
472 0 : call addfld ( dgnumwet_fldn(m), (/ 'lev' /), rad_data_avgflag, '', trim(long_name) )
473 0 : call add_default (dgnumwet_fldn(m), rad_data_histfile_num, ' ')
474 :
475 0 : write(long_name, 1004) m
476 0 : call addfld ( qaerwat_fldn(m), (/ 'lev' /), rad_data_avgflag, '', trim(long_name) )
477 0 : call add_default (qaerwat_fldn(m), rad_data_histfile_num, ' ')
478 :
479 : ! get mode info
480 0 : call rad_cnst_get_info(0, m, nspec=nspec)
481 : ! aerosol species loop
482 0 : do l = 1, nspec
483 0 : call rad_cnst_get_info(0,m,l, spec_name=aername)
484 0 : name = 'rad_'//trim(aername)
485 0 : call addfld(trim(name), (/ 'lev' /), rad_data_avgflag, 'kg/kg', trim(long_name))
486 0 : call add_default (trim(name), rad_data_histfile_num, ' ')
487 : end do
488 : end do
489 : endif
490 :
491 : 1001 format ( 'rad_dgnumwet_', I1 )
492 : 1002 format ( 'radiation input: DGNUMWET for mode ', I1 )
493 : 1003 format ( 'rad_qaerwat_', I1 )
494 : 1004 format ( 'radiation input: QAERWAT for mode ', I1 )
495 :
496 1536 : end subroutine rad_data_init
497 :
498 : !================================================================================================
499 : !================================================================================================
500 1492272 : subroutine rad_data_write( pbuf, state, cam_in, coszen )
501 :
502 1536 : use physics_types, only: physics_state
503 : use camsrfexch, only: cam_in_t
504 : use constituents, only: cnst_get_ind
505 : use physics_buffer, only: pbuf_get_field, pbuf_old_tim_idx
506 : use drv_input_data, only: drv_input_data_freq
507 : use physconst, only: cpair
508 :
509 : implicit none
510 : type(physics_buffer_desc), pointer :: pbuf(:)
511 :
512 : type(physics_state), intent(in), target :: state
513 : type(cam_in_t), intent(in) :: cam_in
514 : real(r8), intent(in) :: coszen(pcols)
515 :
516 : ! Local variables
517 : integer :: i, k
518 : integer :: m,l, nspec
519 : character(len=32) :: name, aername
520 1492272 : real(r8), pointer :: mmr(:,:)
521 :
522 : integer :: lchnk, itim_old, ifld
523 : integer :: ixcldice ! cloud ice water index
524 : integer :: ixcldliq ! cloud liquid water index
525 : integer :: icol
526 : integer :: ncol
527 :
528 : ! surface albedoes weighted by (positive cosine zenith angle)
529 : real(r8):: coszrs_pos(pcols) ! = max(coszrs,0)
530 : real(r8):: asdir_pos (pcols) !
531 : real(r8):: asdif_pos (pcols) !
532 : real(r8):: aldir_pos (pcols) !
533 : real(r8):: aldif_pos (pcols) !
534 :
535 1492272 : real(r8), pointer, dimension(:,:) :: ptr
536 1492272 : integer , pointer, dimension(:) :: iptr1d
537 : real(r8), dimension(pcols) :: rptr1d
538 :
539 : ! for Fixed Dynamical Heating
540 1492272 : real(r8), pointer, dimension(:,:) :: tcorr
541 1492272 : real(r8), pointer, dimension(:,:) :: qrsin
542 1492272 : real(r8), pointer, dimension(:,:) :: qrlin
543 1492272 : real(r8), pointer, dimension(:) :: tropp
544 1492272 : real(r8), pointer, dimension(:,:) :: qrs
545 1492272 : real(r8), pointer, dimension(:,:) :: qrl
546 1492272 : real(r8), pointer, dimension(:,:,:) :: dgnumwet_ptr
547 1492272 : real(r8), pointer, dimension(:,:,:) :: qaerwat_ptr
548 :
549 1492272 : if (.not.enabled) return
550 0 : call pbuf_get_field(pbuf, qrs_ifld, qrs )
551 0 : call pbuf_get_field(pbuf, qrl_ifld, qrl )
552 :
553 0 : lchnk = state%lchnk
554 0 : ncol = state%ncol
555 :
556 : ! store temperature correction above tropopause in pbuf for Fixed Dynamical Heating
557 0 : if(do_fdh) then
558 :
559 0 : call pbuf_get_field(pbuf, tcorr_idx, tcorr)
560 0 : call pbuf_get_field(pbuf, qrsin_idx, qrsin)
561 0 : call pbuf_get_field(pbuf, qrlin_idx, qrlin)
562 0 : call pbuf_get_field(pbuf, tropp_idx, tropp)
563 0 : do k = 1, pver
564 0 : do i = 1, ncol
565 0 : if ( state%pmid(i,k) < tropp(i)) then
566 0 : tcorr(i,k) = tcorr(i,k) + drv_input_data_freq*(qrs(i,k) - qrsin(i,k) + qrl(i,k) - qrlin(i,k)) /cpair
567 : endif
568 : enddo
569 : enddo
570 :
571 0 : call outfld('rad_TROP_P', tropp(:ncol), ncol, lchnk)
572 0 : call outfld('rad_TCORR', tcorr(:ncol,:), ncol, lchnk)
573 :
574 : end if
575 :
576 0 : if (.not.rad_data_output) return
577 :
578 0 : call outfld(qrs_fldn, qrs, pcols, lchnk )
579 0 : call outfld(qrl_fldn, qrl, pcols, lchnk )
580 :
581 : ! get index of (liquid+ice) cloud water
582 0 : call cnst_get_ind('CLDICE', ixcldice)
583 0 : call cnst_get_ind('CLDLIQ', ixcldliq)
584 :
585 0 : do icol = 1, ncol
586 0 : coszrs_pos(icol) = max(coszen(icol),0._r8)
587 : enddo
588 0 : asdir_pos(:ncol) = cam_in%asdir(:ncol) * coszrs_pos(:ncol)
589 0 : asdif_pos(:ncol) = cam_in%asdif(:ncol) * coszrs_pos(:ncol)
590 0 : aldir_pos(:ncol) = cam_in%aldir(:ncol) * coszrs_pos(:ncol)
591 0 : aldif_pos(:ncol) = cam_in%aldif(:ncol) * coszrs_pos(:ncol)
592 :
593 0 : call outfld(lndfrc_fldn, cam_in%landfrac, pcols, lchnk)
594 0 : call outfld(icefrc_fldn, cam_in%icefrac, pcols, lchnk)
595 0 : call outfld(snowh_fldn, cam_in%snowhland, pcols, lchnk)
596 0 : call outfld(temp_fldn, state%t, pcols, lchnk )
597 0 : call outfld(pdel_fldn, state%pdel, pcols, lchnk )
598 0 : call outfld(pdeldry_fldn,state%pdeldry, pcols, lchnk )
599 0 : call outfld(watice_fldn, state%q(:,:,ixcldice), pcols, lchnk )
600 0 : call outfld(watliq_fldn, state%q(:,:,ixcldliq), pcols, lchnk )
601 0 : call outfld(watvap_fldn, state%q(:,:,1), pcols, lchnk )
602 0 : call outfld(zint_fldn, state%zi, pcols, lchnk )
603 0 : call outfld(pint_fldn, state%pint, pcols, lchnk )
604 0 : call outfld(pmid_fldn, state%pmid, pcols, lchnk )
605 :
606 0 : call outfld(asdir_fldn, cam_in%asdir, pcols, lchnk )
607 0 : call outfld(asdif_fldn, cam_in%asdif, pcols, lchnk )
608 0 : call outfld(aldir_fldn, cam_in%aldir, pcols, lchnk )
609 0 : call outfld(aldif_fldn, cam_in%aldif, pcols, lchnk )
610 :
611 0 : call outfld(coszen_fldn, coszrs_pos, pcols, lchnk )
612 0 : call outfld(asdir_pos_fldn, asdir_pos, pcols, lchnk )
613 0 : call outfld(asdif_pos_fldn, asdif_pos, pcols, lchnk )
614 0 : call outfld(aldir_pos_fldn, aldir_pos, pcols, lchnk )
615 0 : call outfld(aldif_pos_fldn, aldif_pos, pcols, lchnk )
616 :
617 0 : call outfld(lwup_fldn, cam_in%lwup, pcols, lchnk )
618 0 : call outfld(ts_fldn, cam_in%ts, pcols, lchnk )
619 :
620 0 : itim_old = pbuf_old_tim_idx()
621 :
622 0 : call pbuf_get_field(pbuf, cld_ifld, ptr, start=(/1,1,itim_old/), kount=(/pcols,pver,1/) )
623 0 : call outfld(cld_fldn, ptr, pcols, lchnk )
624 :
625 0 : call pbuf_get_field(pbuf, rel_ifld, ptr )
626 0 : call outfld(rel_fldn, ptr, pcols, lchnk )
627 :
628 0 : call pbuf_get_field(pbuf, rei_ifld, ptr )
629 0 : call outfld(rei_fldn, ptr, pcols, lchnk )
630 :
631 0 : if (mg_microphys) then
632 :
633 0 : call pbuf_get_field(pbuf, dei_ifld, ptr )
634 0 : call outfld(dei_fldn, ptr, pcols, lchnk )
635 :
636 0 : call pbuf_get_field(pbuf, des_ifld, ptr )
637 0 : call outfld(des_fldn, ptr, pcols, lchnk )
638 :
639 0 : call pbuf_get_field(pbuf, mu_ifld, ptr )
640 0 : call outfld(mu_fldn, ptr, pcols, lchnk )
641 :
642 0 : call pbuf_get_field(pbuf, lambdac_ifld, ptr )
643 0 : call outfld(lambdac_fldn, ptr, pcols, lchnk )
644 :
645 0 : call pbuf_get_field(pbuf, iciwp_ifld, ptr )
646 0 : call outfld(iciwp_fldn, ptr, pcols, lchnk )
647 :
648 0 : call pbuf_get_field(pbuf, iclwp_ifld, ptr )
649 0 : call outfld(iclwp_fldn, ptr, pcols, lchnk )
650 :
651 0 : call pbuf_get_field(pbuf, icswp_ifld, ptr )
652 0 : call outfld(icswp_fldn, ptr, pcols, lchnk )
653 :
654 0 : call pbuf_get_field(pbuf, cldfsnow_ifld, ptr, start=(/1,1,itim_old/), kount=(/pcols,pver,1/) )
655 0 : call outfld(cldfsnow_fldn, ptr, pcols, lchnk )
656 :
657 : else
658 :
659 0 : call pbuf_get_field(pbuf, cldemis_ifld, ptr )
660 0 : call outfld(cldemis_fldn, ptr, pcols, lchnk )
661 :
662 0 : call pbuf_get_field(pbuf, cldtau_ifld, ptr )
663 0 : call outfld(cldtau_fldn, ptr, pcols, lchnk )
664 :
665 0 : call pbuf_get_field(pbuf, cicewp_ifld, ptr )
666 0 : call outfld(cicewp_fldn, ptr, pcols, lchnk )
667 :
668 0 : call pbuf_get_field(pbuf, cliqwp_ifld, ptr )
669 0 : call outfld(cliqwp_fldn, ptr, pcols, lchnk )
670 :
671 0 : call pbuf_get_field(pbuf, nmxrgn_ifld, iptr1d )
672 0 : rptr1d = dble(iptr1d)
673 0 : call outfld(nmxrgn_fldn, rptr1d, pcols, lchnk )
674 :
675 0 : call pbuf_get_field(pbuf, pmxrgn_ifld, ptr )
676 0 : call outfld(pmxrgn_fldn, ptr, pcols, lchnk )
677 :
678 : endif
679 :
680 : ! output mixing ratio of rad constituents
681 :
682 0 : do i = 1, ngas
683 0 : name = 'rad_'//gasnames(i)
684 0 : call rad_cnst_get_gas(0, gaslist(i), state, pbuf, mmr)
685 0 : call outfld(trim(name), mmr, pcols, lchnk)
686 : end do
687 :
688 0 : if ( naer>0 ) then
689 : ! bulk aerosol model
690 0 : do i = 1, naer
691 0 : name = 'rad_'//aernames(i)
692 0 : call rad_cnst_get_aer_mmr(0, i, state, pbuf, mmr)
693 0 : call outfld(trim(name), mmr, pcols, lchnk)
694 : end do
695 : endif
696 0 : if (nmodes>0) then
697 :
698 0 : call pbuf_get_field(pbuf, dgnumwet_ifld, dgnumwet_ptr)
699 0 : call pbuf_get_field(pbuf, qaerwat_ifld, qaerwat_ptr)
700 :
701 : ! for modal aerosol model
702 0 : do m = 1, nmodes
703 0 : ptr => dgnumwet_ptr(:,:,m)
704 0 : call outfld( dgnumwet_fldn(m), ptr, pcols, lchnk )
705 0 : ptr => qaerwat_ptr(:,:,m)
706 0 : call outfld( qaerwat_fldn(m), ptr, pcols, lchnk )
707 :
708 : ! get mode info
709 0 : call rad_cnst_get_info(0, m, nspec=nspec)
710 : ! aerosol species loop
711 0 : do l = 1, nspec
712 0 : call rad_cnst_get_info(0,m,l, spec_name=aername)
713 0 : call rad_cnst_get_aer_mmr(0, m, l, 'a', state, pbuf, mmr)
714 0 : name = 'rad_'//aername
715 0 : call outfld(trim(name), mmr, pcols, lchnk)
716 : enddo
717 : enddo
718 : endif
719 :
720 : ! stratospheric aersols geometric mean radius
721 0 : if (volcgeom_ifld > 0) then
722 0 : call pbuf_get_field(pbuf, volcgeom_ifld, ptr)
723 0 : call outfld(volcgeom_fldn, ptr, pcols, lchnk)
724 : endif
725 0 : if (gmean_3modes) then
726 0 : call pbuf_get_field(pbuf, volcgeom1_ifld, ptr)
727 0 : call outfld(volcgeom1_fldn, ptr, pcols, lchnk)
728 0 : call pbuf_get_field(pbuf, volcgeom2_ifld, ptr)
729 0 : call outfld(volcgeom2_fldn, ptr, pcols, lchnk)
730 0 : call pbuf_get_field(pbuf, volcgeom3_ifld, ptr)
731 0 : call outfld(volcgeom3_fldn, ptr, pcols, lchnk)
732 : endif
733 :
734 1492272 : end subroutine rad_data_write
735 :
736 : !=================================================================================
737 : !=================================================================================
738 0 : subroutine rad_data_read(indata, phys_state, pbuf2d, cam_in, recno )
739 :
740 1492272 : use camsrfexch, only: cam_in_t
741 : use physics_buffer, only: pbuf_get_field, pbuf_old_tim_idx
742 : use constituents, only: cnst_get_ind
743 : use tropopause, only: tropopause_find_cam, TROP_ALG_HYBSTOB, TROP_ALG_CLIMATE
744 :
745 : implicit none
746 :
747 : type(drv_input_data_t), intent(inout) :: indata
748 : type(physics_buffer_desc), pointer :: pbuf2d(:,:)
749 : type(physics_state), target, intent(inout) :: phys_state(begchunk:endchunk)
750 : type(cam_in_t), target, intent(inout) :: cam_in(begchunk:endchunk)
751 : integer, intent(in) :: recno
752 :
753 : ! local vars
754 0 : type(physics_buffer_desc), pointer :: pbuf(:)
755 :
756 0 : type(drv_input_3d_t) :: cld_ptrs(begchunk:endchunk)
757 0 : type(drv_input_3d_t) :: cldfsnow_ptrs(begchunk:endchunk)
758 :
759 0 : type(drv_input_3d_t) :: rel_ptrs(begchunk:endchunk)
760 0 : type(drv_input_3d_t) :: rei_ptrs(begchunk:endchunk)
761 0 : type(drv_input_3d_t) :: dei_ptrs(begchunk:endchunk)
762 0 : type(drv_input_3d_t) :: des_ptrs(begchunk:endchunk)
763 0 : type(drv_input_3d_t) :: mu_ptrs(begchunk:endchunk)
764 0 : type(drv_input_3d_t) :: lambdac_ptrs(begchunk:endchunk)
765 0 : type(drv_input_3d_t) :: iciwp_ptrs(begchunk:endchunk)
766 0 : type(drv_input_3d_t) :: iclwp_ptrs(begchunk:endchunk)
767 0 : type(drv_input_3d_t) :: icswp_ptrs(begchunk:endchunk)
768 :
769 0 : type(drv_input_3d_t) :: qrs_ptrs(begchunk:endchunk)
770 0 : type(drv_input_3d_t) :: qrl_ptrs(begchunk:endchunk)
771 :
772 0 : type(drv_input_3d_t) :: watvap_ptrs(begchunk:endchunk)
773 0 : type(drv_input_3d_t) :: watliq_ptrs(begchunk:endchunk)
774 0 : type(drv_input_3d_t) :: watice_ptrs(begchunk:endchunk)
775 :
776 0 : type(drv_input_3d_t) :: zi_ptrs(begchunk:endchunk)
777 0 : type(drv_input_3d_t) :: pint_ptrs(begchunk:endchunk)
778 0 : type(drv_input_3d_t) :: lnpint_ptrs(begchunk:endchunk)
779 :
780 0 : type(drv_input_3d_t) :: temp_ptrs(begchunk:endchunk)
781 0 : type(drv_input_3d_t) :: pdel_ptrs(begchunk:endchunk)
782 0 : type(drv_input_3d_t) :: pdeldry_ptrs(begchunk:endchunk)
783 0 : type(drv_input_3d_t) :: lnpmid_ptrs(begchunk:endchunk)
784 0 : type(drv_input_3d_t) :: pmid_ptrs(begchunk:endchunk)
785 :
786 0 : type(drv_input_2d_t) :: lndfrac_ptrs(begchunk:endchunk)
787 0 : type(drv_input_2d_t) :: icefrac_ptrs(begchunk:endchunk)
788 0 : type(drv_input_2d_t) :: snowhland_ptrs(begchunk:endchunk)
789 0 : type(drv_input_2d_t) :: asdir_ptrs(begchunk:endchunk)
790 0 : type(drv_input_2d_t) :: asdif_ptrs(begchunk:endchunk)
791 0 : type(drv_input_2d_t) :: aldir_ptrs(begchunk:endchunk)
792 0 : type(drv_input_2d_t) :: aldif_ptrs(begchunk:endchunk)
793 :
794 0 : type(drv_input_2di_t):: nmxrgn_ptrs(begchunk:endchunk)
795 0 : type(drv_input_3d_t) :: pmxrgn_ptrs(begchunk:endchunk)
796 :
797 0 : type(drv_input_3d_t) :: cldemis_ptrs(begchunk:endchunk)
798 0 : type(drv_input_3d_t) :: cldtau_ptrs(begchunk:endchunk)
799 0 : type(drv_input_3d_t) :: cicewp_ptrs(begchunk:endchunk)
800 0 : type(drv_input_3d_t) :: cliqwp_ptrs(begchunk:endchunk)
801 :
802 0 : type(drv_input_2d_t) :: lwup_ptrs(begchunk:endchunk)
803 0 : type(drv_input_2d_t) :: ts_ptrs(begchunk:endchunk)
804 :
805 0 : type(drv_input_3d_t) :: volcgeom_ptrs(begchunk:endchunk)
806 0 : type(drv_input_3d_t) :: volcgeom1_ptrs(begchunk:endchunk)
807 0 : type(drv_input_3d_t) :: volcgeom2_ptrs(begchunk:endchunk)
808 0 : type(drv_input_3d_t) :: volcgeom3_ptrs(begchunk:endchunk)
809 :
810 : integer :: i, k, c, ncol, itim
811 :
812 : integer :: ixcldice ! cloud ice water index
813 : integer :: ixcldliq ! cloud liquid water index
814 :
815 : ! for Fixed Dynamical Heating
816 0 : real(r8), pointer, dimension(:,:) :: tcorr
817 0 : real(r8), pointer, dimension(:,:) :: qrsin
818 0 : real(r8), pointer, dimension(:,:) :: qrlin
819 0 : real(r8), pointer, dimension(:) :: tropp
820 : integer :: troplev(pcols)
821 :
822 : ! for modal aerosols
823 0 : type(drv_input_4d_t) :: dgnumwet_ptrs(begchunk:endchunk)
824 0 : type(drv_input_4d_t) :: qaerwat_ptrs(begchunk:endchunk)
825 :
826 : ! get index of (liquid+ice) cloud water
827 0 : call cnst_get_ind('CLDICE', ixcldice)
828 0 : call cnst_get_ind('CLDLIQ', ixcldliq)
829 :
830 : ! phys buffer time index
831 0 : itim = pbuf_old_tim_idx()
832 :
833 : ! setup the data pointers
834 : !$OMP PARALLEL DO PRIVATE (C,pbuf)
835 0 : do c=begchunk,endchunk
836 0 : pbuf => pbuf_get_chunk(pbuf2d, c)
837 :
838 0 : watvap_ptrs (c)%array => phys_state(c)%q(:,:,1)
839 0 : watliq_ptrs (c)%array => phys_state(c)%q(:,:,ixcldliq)
840 0 : watice_ptrs (c)%array => phys_state(c)%q(:,:,ixcldice)
841 :
842 0 : zi_ptrs (c)%array => phys_state(c)%zi
843 0 : pint_ptrs (c)%array => phys_state(c)%pint
844 0 : lnpint_ptrs (c)%array => phys_state(c)%lnpint
845 :
846 0 : temp_ptrs (c)%array => phys_state(c)%t
847 0 : pdel_ptrs (c)%array => phys_state(c)%pdel
848 0 : pdeldry_ptrs(c)%array => phys_state(c)%pdeldry
849 0 : lnpmid_ptrs (c)%array => phys_state(c)%lnpmid
850 0 : pmid_ptrs (c)%array => phys_state(c)%pmid
851 :
852 0 : lndfrac_ptrs (c)%array => cam_in(c)%landfrac
853 0 : icefrac_ptrs (c)%array => cam_in(c)%icefrac
854 0 : snowhland_ptrs(c)%array => cam_in(c)%snowhland
855 0 : asdir_ptrs (c)%array => cam_in(c)%asdir
856 0 : asdif_ptrs (c)%array => cam_in(c)%asdif
857 0 : aldir_ptrs (c)%array => cam_in(c)%aldir
858 0 : aldif_ptrs (c)%array => cam_in(c)%aldif
859 0 : lwup_ptrs (c)%array => cam_in(c)%lwup
860 0 : ts_ptrs (c)%array => cam_in(c)%ts
861 :
862 0 : call pbuf_get_field(pbuf, cld_ifld, cld_ptrs (c)%array, start=(/1,1,itim/), kount=(/pcols,pver,1/) )
863 0 : call pbuf_get_field(pbuf, rel_ifld, rel_ptrs(c)%array )
864 0 : call pbuf_get_field(pbuf, rei_ifld, rei_ptrs(c)%array )
865 :
866 0 : call pbuf_get_field(pbuf, qrs_ifld, qrs_ptrs(c)%array )
867 0 : call pbuf_get_field(pbuf, qrl_ifld, qrl_ptrs(c)%array )
868 :
869 0 : if (mg_microphys) then
870 0 : call pbuf_get_field(pbuf, dei_ifld, dei_ptrs(c)%array )
871 0 : call pbuf_get_field(pbuf, des_ifld, des_ptrs(c)%array )
872 0 : call pbuf_get_field(pbuf, mu_ifld, mu_ptrs (c)%array )
873 0 : call pbuf_get_field(pbuf, lambdac_ifld, lambdac_ptrs(c)%array )
874 0 : call pbuf_get_field(pbuf, iciwp_ifld, iciwp_ptrs (c)%array )
875 0 : call pbuf_get_field(pbuf, iclwp_ifld, iclwp_ptrs (c)%array )
876 0 : call pbuf_get_field(pbuf, icswp_ifld, icswp_ptrs (c)%array )
877 0 : call pbuf_get_field(pbuf, cldfsnow_ifld, cldfsnow_ptrs(c)%array, start=(/1,1,itim/), kount=(/pcols,pver,1/) )
878 : else
879 0 : call pbuf_get_field(pbuf, nmxrgn_ifld, nmxrgn_ptrs(c)%array )
880 0 : call pbuf_get_field(pbuf, pmxrgn_ifld, pmxrgn_ptrs(c)%array )
881 0 : call pbuf_get_field(pbuf, cldemis_ifld, cldemis_ptrs(c)%array )
882 0 : call pbuf_get_field(pbuf, cldtau_ifld, cldtau_ptrs(c)%array )
883 0 : call pbuf_get_field(pbuf, cicewp_ifld, cicewp_ptrs(c)%array )
884 0 : call pbuf_get_field(pbuf, cliqwp_ifld, cliqwp_ptrs(c)%array )
885 : endif
886 :
887 0 : if (nmodes>0) then
888 0 : call pbuf_get_field(pbuf, dgnumwet_ifld, dgnumwet_ptrs(c)%array )
889 0 : call pbuf_get_field(pbuf, qaerwat_ifld, qaerwat_ptrs(c)%array )
890 : endif
891 :
892 : ! stratospheric aersols geometric mean radius
893 0 : if (volcgeom_ifld > 0) then
894 0 : call pbuf_get_field(pbuf, volcgeom_ifld, volcgeom_ptrs(c)%array )
895 : endif
896 0 : if (gmean_3modes) then
897 0 : call pbuf_get_field(pbuf, volcgeom1_ifld, volcgeom1_ptrs(c)%array )
898 0 : call pbuf_get_field(pbuf, volcgeom2_ifld, volcgeom2_ptrs(c)%array )
899 0 : call pbuf_get_field(pbuf, volcgeom3_ifld, volcgeom3_ptrs(c)%array )
900 : endif
901 :
902 : enddo
903 :
904 :
905 0 : if (nmodes>0) then
906 0 : call get_aeromode_data( indata, dgnumwet_fldn, recno, dgnumwet_ptrs )
907 0 : call get_aeromode_data( indata, qaerwat_fldn, recno, qaerwat_ptrs )
908 : endif
909 :
910 : ! get the 2D data
911 :
912 0 : call drv_input_data_get( indata, lndfrc_fldn, recno, lndfrac_ptrs )
913 0 : call drv_input_data_get( indata, icefrc_fldn, recno, icefrac_ptrs )
914 0 : call drv_input_data_get( indata, snowh_fldn, recno, snowhland_ptrs )
915 0 : call drv_input_data_get( indata, asdir_fldn, recno, asdir_ptrs )
916 0 : call drv_input_data_get( indata, asdif_fldn, recno, asdif_ptrs )
917 0 : call drv_input_data_get( indata, aldir_fldn, recno, aldir_ptrs )
918 0 : call drv_input_data_get( indata, aldif_fldn, recno, aldif_ptrs )
919 0 : call drv_input_data_get( indata, lwup_fldn, recno, lwup_ptrs )
920 0 : call drv_input_data_get( indata, ts_fldn, recno, ts_ptrs )
921 :
922 : ! get the 3D data
923 :
924 0 : call drv_input_data_get( indata, cld_fldn, 'lev', pver, recno, cld_ptrs )
925 0 : call drv_input_data_get( indata, rel_fldn, 'lev', pver, recno, rel_ptrs )
926 0 : call drv_input_data_get( indata, rei_fldn, 'lev', pver, recno, rei_ptrs )
927 0 : call drv_input_data_get( indata, qrs_fldn, 'lev', pver, recno, qrs_ptrs )
928 0 : call drv_input_data_get( indata, qrl_fldn, 'lev', pver, recno, qrl_ptrs )
929 :
930 0 : if (mg_microphys) then
931 0 : call drv_input_data_get( indata, dei_fldn, 'lev', pver, recno, dei_ptrs )
932 0 : call drv_input_data_get( indata, des_fldn, 'lev', pver, recno, des_ptrs )
933 0 : call drv_input_data_get( indata, mu_fldn, 'lev', pver, recno, mu_ptrs )
934 0 : call drv_input_data_get( indata, lambdac_fldn,'lev', pver, recno, lambdac_ptrs )
935 0 : call drv_input_data_get( indata, iciwp_fldn, 'lev', pver, recno, iciwp_ptrs )
936 0 : call drv_input_data_get( indata, iclwp_fldn, 'lev', pver, recno, iclwp_ptrs )
937 0 : call drv_input_data_get( indata, icswp_fldn, 'lev', pver, recno, icswp_ptrs )
938 0 : call drv_input_data_get( indata, cldfsnow_fldn,'lev',pver, recno, cldfsnow_ptrs )
939 : else
940 0 : call drv_input_data_get( indata, nmxrgn_fldn, recno, nmxrgn_ptrs )
941 0 : call drv_input_data_get( indata, pmxrgn_fldn, 'ilev',pverp,recno, pmxrgn_ptrs )
942 0 : call drv_input_data_get( indata, cldemis_fldn,'lev', pver, recno, cldemis_ptrs )
943 0 : call drv_input_data_get( indata, cldtau_fldn, 'lev', pver, recno, cldtau_ptrs )
944 0 : call drv_input_data_get( indata, cicewp_fldn, 'lev', pver, recno, cicewp_ptrs )
945 0 : call drv_input_data_get( indata, cliqwp_fldn, 'lev', pver, recno, cliqwp_ptrs )
946 : endif
947 :
948 : ! stratospheric aersols geometric mean radius
949 0 : if (volcgeom_ifld > 0) then
950 0 : call drv_input_data_get( indata, volcgeom_fldn, 'lev', pver, recno, volcgeom_ptrs )
951 : endif
952 0 : if (gmean_3modes) then
953 0 : call drv_input_data_get( indata, volcgeom1_fldn, 'lev', pver, recno, volcgeom1_ptrs )
954 0 : call drv_input_data_get( indata, volcgeom2_fldn, 'lev', pver, recno, volcgeom2_ptrs )
955 0 : call drv_input_data_get( indata, volcgeom3_fldn, 'lev', pver, recno, volcgeom3_ptrs )
956 : endif
957 :
958 0 : call drv_input_data_get( indata, watvap_fldn, 'lev', pver, recno, watvap_ptrs )
959 0 : call drv_input_data_get( indata, watliq_fldn, 'lev', pver, recno, watliq_ptrs )
960 0 : call drv_input_data_get( indata, watice_fldn, 'lev', pver, recno, watice_ptrs )
961 0 : call drv_input_data_get( indata, temp_fldn, 'lev', pver, recno, temp_ptrs )
962 0 : call drv_input_data_get( indata, pdel_fldn, 'lev', pver, recno, pdel_ptrs )
963 0 : call drv_input_data_get( indata, pdeldry_fldn,'lev', pver, recno, pdeldry_ptrs )
964 0 : call drv_input_data_get( indata, pmid_fldn, 'lev', pver, recno, pmid_ptrs )
965 :
966 0 : call drv_input_data_get( indata, zint_fldn, 'ilev',pverp,recno, zi_ptrs )
967 0 : call drv_input_data_get( indata, pint_fldn, 'ilev',pverp,recno, pint_ptrs )
968 :
969 0 : do c=begchunk,endchunk
970 0 : ncol = phys_state(c)%ncol
971 0 : if (all(pmid_ptrs(c)%array(:ncol,:) > 0._r8 )) then
972 0 : do i=1,ncol
973 0 : lnpmid_ptrs(c)%array(i,:) = log(pmid_ptrs(c)%array(i,:))
974 0 : lnpint_ptrs(c)%array(i,:) = log(pint_ptrs(c)%array(i,:))
975 : enddo
976 : endif
977 :
978 : ! adjust temperatue input above tropopause for Fixed Dynamical Heating ....
979 0 : if (do_fdh) then
980 0 : pbuf => pbuf_get_chunk(pbuf2d, c)
981 :
982 0 : call pbuf_get_field(pbuf, tropp_idx, tropp)
983 0 : call pbuf_get_field(pbuf, tcorr_idx, tcorr)
984 0 : call pbuf_get_field(pbuf, qrsin_idx, qrsin)
985 0 : call pbuf_get_field(pbuf, qrlin_idx, qrlin)
986 :
987 : !REMOVECAM - no longer need this when CAM is retired and pcols no longer exists
988 0 : troplev(:) = 0
989 0 : tropp(:) = 0._r8
990 : !REMOVECAM_END
991 0 : call tropopause_find_cam(phys_state(c), troplev, tropP=tropp, primary=TROP_ALG_CLIMATE, &
992 0 : backup=TROP_ALG_CLIMATE)
993 :
994 0 : qrsin(:,:) = qrs_ptrs(c)%array(:,:)
995 0 : qrlin(:,:) = qrl_ptrs(c)%array(:,:)
996 :
997 0 : do k = 1, pver
998 0 : do i = 1, ncol
999 0 : if ( phys_state(c)%pmid(i,k) < tropp(i) ) then
1000 0 : temp_ptrs(c)%array(i,k) = temp_ptrs(c)%array(i,k) + tcorr(i,k)
1001 : endif
1002 : enddo
1003 : enddo
1004 : endif
1005 :
1006 : enddo
1007 :
1008 0 : call get_rad_cnst_data(indata, phys_state, pbuf2d, recno)
1009 :
1010 0 : end subroutine rad_data_read
1011 :
1012 : !================================================================================================
1013 : ! Private routines
1014 : !================================================================================================
1015 :
1016 :
1017 : !================================================================================================
1018 : !================================================================================================
1019 0 : subroutine get_aeromode_data(indata, infld_names, recno, chunk_ptrs)
1020 0 : use drv_input_data, only: drv_input_data_read
1021 :
1022 : type(drv_input_data_t), intent(inout) :: indata
1023 : character(len=*), intent(in) :: infld_names(nmodes)
1024 : integer, intent(in) :: recno
1025 : type(drv_input_4d_t), intent(inout) :: chunk_ptrs(begchunk:endchunk)
1026 :
1027 0 : real(r8) :: data (pcols, pver, begchunk:endchunk)
1028 :
1029 : integer :: c, n
1030 :
1031 0 : do n = 1,nmodes
1032 :
1033 0 : data = drv_input_data_read( indata, infld_names(n), 'lev', pver, recno )
1034 0 : do c=begchunk,endchunk
1035 0 : chunk_ptrs(c)%array(:,:,n) = data(:,:,c)
1036 : enddo
1037 :
1038 : enddo
1039 :
1040 0 : end subroutine get_aeromode_data
1041 :
1042 : !================================================================================================
1043 : !================================================================================================
1044 0 : subroutine get_rad_cnst_data(indata, state, pbuf2d, recno)
1045 :
1046 0 : use physics_types, only: physics_state
1047 : use physics_buffer, only: physics_buffer_desc
1048 : use rad_constituents, only: rad_cnst_get_info
1049 :
1050 : implicit none
1051 :
1052 : ! Arguments
1053 : type(drv_input_data_t), intent(inout) :: indata
1054 : type(physics_buffer_desc), pointer :: pbuf2d(:,:)
1055 : type(physics_state), intent(inout) :: state(begchunk:endchunk)
1056 : integer, intent(in) :: recno
1057 :
1058 : ! Local variables
1059 : integer :: i, ncol
1060 : integer :: m, l, nspec
1061 : integer :: idx
1062 : character(len=32) :: aername
1063 :
1064 : !-----------------------------------------------------------------------------
1065 :
1066 : ! read in mixing ratio of rad constituents
1067 :
1068 0 : do i = 1,ngas
1069 0 : call read_rad_gas_data(indata, gasnames(i), i, state, pbuf2d, recno )
1070 : enddo
1071 :
1072 0 : do i = 1,naer
1073 0 : call read_rad_aer_data(indata, aernames(i), i, state, pbuf2d, recno )
1074 : enddo
1075 :
1076 : ! aerosol mode loop
1077 0 : do m = 1, nmodes
1078 :
1079 : ! get mode info
1080 0 : call rad_cnst_get_info(0, m, nspec=nspec)
1081 : ! aerosol species loop
1082 0 : do l = 1, nspec
1083 0 : call rad_cnst_get_info(0,m,l, spec_name=aername)
1084 0 : call read_rad_mam_data( indata, aername, m, l, state, pbuf2d, recno )
1085 : end do
1086 : end do
1087 :
1088 0 : end subroutine get_rad_cnst_data
1089 :
1090 : !=================================================================================
1091 : !=================================================================================
1092 0 : subroutine read_rad_gas_data(indata, name, idx, state, pbuf2d, recno )
1093 0 : use rad_constituents, only: rad_cnst_get_gas
1094 : use drv_input_data, only: drv_input_data_read
1095 :
1096 : type(drv_input_data_t), intent(inout) :: indata
1097 : character(len=32), intent(in) :: name
1098 : integer, intent(in) :: idx
1099 : type(physics_state),target, intent(inout) :: state(begchunk:endchunk)
1100 : type(physics_buffer_desc), pointer, dimension(:,:) :: pbuf2d
1101 : integer, intent(in) :: recno
1102 :
1103 0 : type(physics_buffer_desc), pointer, dimension(:) :: phys_buffer_chunk
1104 : character(len=32) :: radname
1105 : integer :: c, ncol
1106 0 : real(r8) :: mass(pcols,pver,begchunk:endchunk)
1107 0 : real(r8), pointer :: mmr_ptr(:,:)
1108 :
1109 0 : radname = 'rad_'//trim(name)
1110 0 : mass = drv_input_data_read( indata, radname, 'lev', pver, recno )
1111 :
1112 0 : do c = begchunk,endchunk
1113 0 : ncol = state(c)%ncol
1114 0 : phys_buffer_chunk => pbuf_get_chunk(pbuf2d, c)
1115 0 : call rad_cnst_get_gas(0, gaslist(idx), state(c), phys_buffer_chunk, mmr_ptr)
1116 0 : mmr_ptr(:ncol,:) = mass(:ncol,:,c)
1117 : enddo
1118 :
1119 0 : end subroutine read_rad_gas_data
1120 :
1121 : !=================================================================================
1122 : !=================================================================================
1123 0 : subroutine read_rad_aer_data(indata, name, idx, state, pbuf2d, recno )
1124 0 : use rad_constituents, only: rad_cnst_get_aer_mmr
1125 : use drv_input_data, only: drv_input_data_read
1126 :
1127 : type(drv_input_data_t), intent(inout) :: indata
1128 : character(len=32), intent(in) :: name
1129 : integer, intent(in) :: idx
1130 : type(physics_state),target, intent(inout) :: state(begchunk:endchunk)
1131 : type(physics_buffer_desc), pointer, dimension(:,:) :: pbuf2d
1132 : integer, intent(in) :: recno
1133 :
1134 0 : type(physics_buffer_desc), pointer, dimension(:) :: phys_buffer_chunk
1135 : character(len=32) :: radname
1136 : integer :: c, ncol
1137 0 : real(r8) :: mass(pcols,pver,begchunk:endchunk)
1138 0 : real(r8), pointer :: mmr_ptr(:,:)
1139 :
1140 0 : radname = 'rad_'//trim(name)
1141 0 : mass = drv_input_data_read( indata, radname, 'lev', pver, recno )
1142 :
1143 0 : do c = begchunk,endchunk
1144 0 : ncol = state(c)%ncol
1145 0 : phys_buffer_chunk => pbuf_get_chunk(pbuf2d, c)
1146 0 : call rad_cnst_get_aer_mmr(0, idx, state(c), phys_buffer_chunk, mmr_ptr)
1147 0 : mmr_ptr(:ncol,:) = mass(:ncol,:,c)
1148 : enddo
1149 :
1150 0 : end subroutine read_rad_aer_data
1151 :
1152 : !=================================================================================
1153 : !=================================================================================
1154 0 : subroutine read_rad_mam_data(indata, name, mode_idx, spec_idx, state, pbuf2d, recno )
1155 0 : use rad_constituents, only: rad_cnst_get_aer_mmr
1156 : use drv_input_data, only: drv_input_data_read
1157 :
1158 : type(drv_input_data_t), intent(inout) :: indata
1159 : character(len=32), intent(in) :: name
1160 : integer, intent(in) :: mode_idx
1161 : integer, intent(in) :: spec_idx
1162 : type(physics_state),target, intent(inout) :: state(begchunk:endchunk)
1163 : type(physics_buffer_desc), pointer, dimension(:,:) :: pbuf2d
1164 : integer, intent(in) :: recno
1165 :
1166 0 : type(physics_buffer_desc), pointer, dimension(:) :: phys_buffer_chunk
1167 : character(len=32) :: radname
1168 : integer :: c, ncol
1169 0 : real(r8) :: mass(pcols,pver,begchunk:endchunk)
1170 0 : real(r8), pointer :: mmr_ptr(:,:)
1171 :
1172 0 : radname = 'rad_'//trim(name)
1173 0 : mass = drv_input_data_read( indata, radname, 'lev', pver, recno )
1174 :
1175 0 : do c = begchunk,endchunk
1176 0 : ncol = state(c)%ncol
1177 0 : phys_buffer_chunk => pbuf_get_chunk(pbuf2d, c)
1178 0 : call rad_cnst_get_aer_mmr(0, mode_idx, spec_idx, 'a', state(c), phys_buffer_chunk, mmr_ptr)
1179 0 : mmr_ptr(:ncol,:) = mass(:ncol,:,c)
1180 : enddo
1181 :
1182 0 : end subroutine read_rad_mam_data
1183 :
1184 : end module radiation_data
|