Line data Source code
1 : !-------------------------------------------------------------------
2 : ! manages reading and interpolation of prescribed volcanic aerosol
3 : ! Created by: Francis Vitt
4 : !-------------------------------------------------------------------
5 : module prescribed_volcaero
6 :
7 : use shr_kind_mod, only : r8 => shr_kind_r8
8 : use cam_abortutils, only : endrun
9 : use spmd_utils, only : masterproc
10 : use tracer_data, only : trfld, trfile
11 : use cam_logfile, only : iulog
12 :
13 : implicit none
14 : private
15 : save
16 :
17 : type(trfld), pointer :: fields(:)
18 : type(trfile) :: file
19 :
20 : public :: prescribed_volcaero_readnl
21 : public :: prescribed_volcaero_register
22 : public :: prescribed_volcaero_init
23 : public :: prescribed_volcaero_adv
24 : public :: write_prescribed_volcaero_restart
25 : public :: read_prescribed_volcaero_restart
26 : public :: has_prescribed_volcaero
27 : public :: init_prescribed_volcaero_restart
28 :
29 :
30 : logical :: has_prescribed_volcaero = .false.
31 : character(len=8), parameter :: volcaero_name = 'VOLC_MMR'
32 : character(len=13), parameter :: volcrad_name = 'VOLC_RAD_GEOM'
33 : character(len=9), parameter :: volcmass_name = 'VOLC_MASS'
34 : character(len=11), parameter :: volcmass_column_name = 'VOLC_MASS_C'
35 :
36 : ! These variables are settable via the namelist (with longer names)
37 : character(len=16) :: fld_name = 'MMRVOLC'
38 : character(len=256) :: filename = 'NONE'
39 : character(len=256) :: filelist = ''
40 : character(len=256) :: datapath = ''
41 : character(len=32) :: data_type = 'SERIAL'
42 : logical :: rmv_file = .false.
43 : integer :: cycle_yr = 0
44 : integer :: fixed_ymd = 0
45 : integer :: fixed_tod = 0
46 : integer :: radius_ndx
47 :
48 : contains
49 :
50 : !-------------------------------------------------------------------
51 : !-------------------------------------------------------------------
52 2304 : subroutine prescribed_volcaero_readnl(nlfile)
53 :
54 : use namelist_utils, only: find_group_name
55 : use units, only: getunit, freeunit
56 : use mpishorthand
57 :
58 : character(len=*), intent(in) :: nlfile ! filepath for file containing namelist input
59 :
60 : ! Local variables
61 : integer :: unitn, ierr
62 : character(len=*), parameter :: subname = 'prescribed_volcaero_readnl'
63 :
64 : character(len=16) :: prescribed_volcaero_name
65 : character(len=256) :: prescribed_volcaero_file
66 : character(len=256) :: prescribed_volcaero_filelist
67 : character(len=256) :: prescribed_volcaero_datapath
68 : character(len=32) :: prescribed_volcaero_type
69 : logical :: prescribed_volcaero_rmfile
70 : integer :: prescribed_volcaero_cycle_yr
71 : integer :: prescribed_volcaero_fixed_ymd
72 : integer :: prescribed_volcaero_fixed_tod
73 :
74 : namelist /prescribed_volcaero_nl/ &
75 : prescribed_volcaero_name, &
76 : prescribed_volcaero_file, &
77 : prescribed_volcaero_filelist, &
78 : prescribed_volcaero_datapath, &
79 : prescribed_volcaero_type, &
80 : prescribed_volcaero_rmfile, &
81 : prescribed_volcaero_cycle_yr, &
82 : prescribed_volcaero_fixed_ymd, &
83 : prescribed_volcaero_fixed_tod
84 : !-----------------------------------------------------------------------------
85 :
86 : ! Initialize namelist variables from local module variables.
87 1536 : prescribed_volcaero_name = fld_name
88 1536 : prescribed_volcaero_file = filename
89 1536 : prescribed_volcaero_filelist = filelist
90 1536 : prescribed_volcaero_datapath = datapath
91 1536 : prescribed_volcaero_type = data_type
92 1536 : prescribed_volcaero_rmfile = rmv_file
93 1536 : prescribed_volcaero_cycle_yr = cycle_yr
94 1536 : prescribed_volcaero_fixed_ymd= fixed_ymd
95 1536 : prescribed_volcaero_fixed_tod= fixed_tod
96 :
97 : ! Read namelist
98 1536 : if (masterproc) then
99 2 : unitn = getunit()
100 2 : open( unitn, file=trim(nlfile), status='old' )
101 2 : call find_group_name(unitn, 'prescribed_volcaero_nl', status=ierr)
102 2 : if (ierr == 0) then
103 0 : read(unitn, prescribed_volcaero_nl, iostat=ierr)
104 0 : if (ierr /= 0) then
105 0 : call endrun(subname // ':: ERROR reading namelist')
106 : end if
107 : end if
108 2 : close(unitn)
109 2 : call freeunit(unitn)
110 : end if
111 :
112 : #ifdef SPMD
113 : ! Broadcast namelist variables
114 1536 : call mpibcast(prescribed_volcaero_name, len(prescribed_volcaero_name), mpichar, 0, mpicom)
115 1536 : call mpibcast(prescribed_volcaero_file, len(prescribed_volcaero_file), mpichar, 0, mpicom)
116 1536 : call mpibcast(prescribed_volcaero_filelist, len(prescribed_volcaero_filelist), mpichar, 0, mpicom)
117 1536 : call mpibcast(prescribed_volcaero_datapath, len(prescribed_volcaero_datapath), mpichar, 0, mpicom)
118 1536 : call mpibcast(prescribed_volcaero_type, len(prescribed_volcaero_type), mpichar, 0, mpicom)
119 1536 : call mpibcast(prescribed_volcaero_rmfile, 1, mpilog, 0, mpicom)
120 1536 : call mpibcast(prescribed_volcaero_cycle_yr, 1, mpiint, 0, mpicom)
121 1536 : call mpibcast(prescribed_volcaero_fixed_ymd,1, mpiint, 0, mpicom)
122 1536 : call mpibcast(prescribed_volcaero_fixed_tod,1, mpiint, 0, mpicom)
123 : #endif
124 :
125 : ! Update module variables with user settings.
126 1536 : fld_name = prescribed_volcaero_name
127 1536 : filename = prescribed_volcaero_file
128 1536 : filelist = prescribed_volcaero_filelist
129 1536 : datapath = prescribed_volcaero_datapath
130 1536 : data_type = prescribed_volcaero_type
131 1536 : rmv_file = prescribed_volcaero_rmfile
132 1536 : cycle_yr = prescribed_volcaero_cycle_yr
133 1536 : fixed_ymd = prescribed_volcaero_fixed_ymd
134 1536 : fixed_tod = prescribed_volcaero_fixed_tod
135 :
136 : ! Turn on prescribed volcanics if user has specified an input dataset.
137 1536 : if (len_trim(filename) > 0 .and. filename.ne.'NONE') has_prescribed_volcaero = .true.
138 :
139 1536 : end subroutine prescribed_volcaero_readnl
140 :
141 : !-------------------------------------------------------------------
142 : !-------------------------------------------------------------------
143 1536 : subroutine prescribed_volcaero_register()
144 : use ppgrid, only: pver,pcols
145 : use physics_buffer, only : pbuf_add_field, dtype_r8
146 :
147 : integer :: idx
148 :
149 1536 : if (has_prescribed_volcaero) then
150 0 : call pbuf_add_field(volcaero_name,'physpkg',dtype_r8,(/pcols,pver/),idx)
151 0 : call pbuf_add_field(volcrad_name, 'physpkg',dtype_r8,(/pcols,pver/),idx)
152 :
153 : endif
154 :
155 1536 : endsubroutine prescribed_volcaero_register
156 :
157 : !-------------------------------------------------------------------
158 : !-------------------------------------------------------------------
159 1536 : subroutine prescribed_volcaero_init()
160 :
161 1536 : use tracer_data, only : trcdata_init
162 : use cam_history, only : addfld, horiz_only
163 : use physics_buffer, only : pbuf_get_index
164 :
165 : implicit none
166 :
167 : integer :: ndx, istat
168 : integer :: errcode
169 : character(len=32) :: specifier(1)
170 :
171 1536 : if ( has_prescribed_volcaero ) then
172 0 : if ( masterproc ) then
173 0 : write(iulog,*) 'volcanic aerosol is prescribed in :'//trim(filename)
174 : endif
175 : else
176 : return
177 : endif
178 :
179 0 : specifier(1) = trim(volcaero_name)//':'//trim(fld_name)
180 :
181 :
182 0 : allocate(file%in_pbuf(size(specifier)))
183 0 : file%in_pbuf(:) = .true.
184 : call trcdata_init( specifier, filename, filelist, datapath, fields, file, &
185 0 : rmv_file, cycle_yr, fixed_ymd, fixed_tod, data_type)
186 :
187 :
188 0 : call addfld(volcaero_name, (/ 'lev' /), 'I','kg/kg', 'prescribed volcanic aerosol dry mass mixing ratio' )
189 0 : call addfld(volcrad_name, (/ 'lev' /), 'I','m', 'volcanic aerosol geometric-mean radius' )
190 0 : call addfld(volcmass_name, (/ 'lev' /), 'I','kg/m^2', 'volcanic aerosol vertical mass path in layer' )
191 0 : call addfld(volcmass_column_name, horiz_only, 'I','kg/m^2', 'volcanic aerosol column mass' )
192 :
193 0 : radius_ndx = pbuf_get_index(volcrad_name, errcode)
194 :
195 1536 : end subroutine prescribed_volcaero_init
196 :
197 : !-------------------------------------------------------------------
198 : !-------------------------------------------------------------------
199 741888 : subroutine prescribed_volcaero_adv( state, pbuf2d)
200 :
201 1536 : use tracer_data, only : advance_trcdata
202 : use physics_types,only : physics_state
203 : use ppgrid, only : begchunk, endchunk
204 : use ppgrid, only : pcols, pver
205 : use string_utils, only : to_lower, GLC
206 : use cam_history, only : outfld
207 : use physconst, only : mwdry ! molecular weight dry air ~ kg/kmole
208 : use physconst, only : boltz, gravit ! J/K/molecule
209 : use tropopause, only : tropopause_find_cam
210 :
211 : use physics_buffer, only : physics_buffer_desc, pbuf_get_field, pbuf_get_chunk
212 :
213 : implicit none
214 :
215 : type(physics_state), intent(in) :: state(begchunk:endchunk)
216 :
217 : type(physics_buffer_desc), pointer :: pbuf2d(:,:)
218 :
219 370944 : type(physics_buffer_desc), pointer :: pbuf_chnk(:)
220 :
221 : integer :: c,ncol,i,k
222 : real(r8) :: to_mmr(pcols,pver)
223 : real(r8), parameter :: molmass = 47.9981995_r8
224 : real(r8) :: ptrop
225 : real(r8) :: concvolc ! micrograms of wetted aerosol per cubic centimeter
226 : real(r8) :: volcmass(pcols,pver)
227 : real(r8) :: columnmass(pcols)
228 : real(r8) :: mmrvolc
229 : integer :: tropLev(pcols)
230 :
231 : real(r8) :: outdata(pcols,pver)
232 370944 : real(r8), pointer :: data(:,:)
233 370944 : real(r8), pointer :: radius(:,:)
234 :
235 : !WACCM-derived relation between mass concentration and wet aerosol radius in meters
236 : real(r8),parameter :: radius_conversion = 1.9e-4_r8
237 :
238 370944 : if( .not. has_prescribed_volcaero ) return
239 :
240 0 : call advance_trcdata( fields, file, state, pbuf2d )
241 :
242 : ! copy prescribed tracer fields into state svariable with the correct units
243 0 : do c = begchunk,endchunk
244 0 : pbuf_chnk => pbuf_get_chunk(pbuf2d, c)
245 0 : call pbuf_get_field(pbuf_chnk, radius_ndx, radius)
246 0 : radius(:,:) = 0._r8
247 0 : ncol = state(c)%ncol
248 0 : select case ( to_lower(trim(fields(1)%units(:GLC(fields(1)%units)))) )
249 : case ("molec/cm3","/cm3","molecules/cm3","cm^-3","cm**-3")
250 0 : to_mmr(:ncol,:) = (molmass*1.e6_r8*boltz*state(c)%t(:ncol,:))/(mwdry*state(c)%pmiddry(:ncol,:))
251 : case ('kg/kg','mmr','kg kg-1')
252 0 : to_mmr(:ncol,:) = 1._r8
253 : case ('mol/mol','mole/mole','vmr','fraction')
254 0 : to_mmr(:ncol,:) = molmass/mwdry
255 : case default
256 0 : write(iulog,*) 'prescribed_volcaero_adv: units = ',trim(fields(1)%units) ,' are not recognized'
257 0 : call endrun('prescribed_volcaero_adv: units are not recognized')
258 : end select
259 :
260 0 : call pbuf_get_field(pbuf_chnk, fields(1)%pbuf_ndx, data)
261 0 : data(:ncol,:) = to_mmr(:ncol,:) * data(:ncol,:) ! mmr
262 :
263 : !REMOVECAM - no longer need this when CAM is retired and pcols no longer exists
264 0 : tropLev(:) = 0
265 : !REMOVECAM_END
266 0 : call tropopause_find_cam(state(c), tropLev)
267 0 : do i = 1,ncol
268 0 : do k = 1,pver
269 : ! set to zero below tropopause
270 0 : if ( k >= tropLev(i) ) then
271 0 : data(i,k) = 0._r8
272 : endif
273 0 : mmrvolc = data(i,k)
274 0 : if (mmrvolc > 0._r8) then
275 0 : concvolc = (mmrvolc * state(c)%pdel(i,k))/(gravit * state(c)%zm(i,k))
276 0 : radius(i,k) = radius_conversion*(concvolc**(1._r8/3._r8))
277 : endif
278 : enddo
279 : enddo
280 :
281 0 : volcmass(:ncol,:) = data(:ncol,:)*state(c)%pdel(:ncol,:)/gravit
282 0 : columnmass(:ncol) = sum(volcmass(:ncol,:), 2)
283 :
284 0 : call outfld( volcaero_name, data(:,:), pcols, state(c)%lchnk)
285 0 : call outfld( volcrad_name, radius(:,:), pcols, state(c)%lchnk)
286 0 : call outfld( volcmass_name, volcmass(:,:), pcols, state(c)%lchnk)
287 0 : call outfld( volcmass_column_name, columnmass(:), pcols, state(c)%lchnk)
288 :
289 : enddo
290 :
291 370944 : end subroutine prescribed_volcaero_adv
292 :
293 : !-------------------------------------------------------------------
294 1536 : subroutine init_prescribed_volcaero_restart( piofile )
295 370944 : use pio, only : file_desc_t
296 : use tracer_data, only : init_trc_restart
297 : implicit none
298 : type(file_desc_t),intent(inout) :: pioFile ! pio File pointer
299 :
300 1536 : call init_trc_restart( 'prescribed_volcaero', piofile, file )
301 :
302 1536 : end subroutine init_prescribed_volcaero_restart
303 : !-------------------------------------------------------------------
304 1536 : subroutine write_prescribed_volcaero_restart( piofile )
305 1536 : use tracer_data, only : write_trc_restart
306 : use pio, only : file_desc_t
307 : implicit none
308 :
309 : type(file_desc_t) :: piofile
310 :
311 1536 : call write_trc_restart( piofile, file )
312 :
313 1536 : end subroutine write_prescribed_volcaero_restart
314 : !-------------------------------------------------------------------
315 : !-------------------------------------------------------------------
316 768 : subroutine read_prescribed_volcaero_restart( pioFile )
317 1536 : use tracer_data, only : read_trc_restart
318 : use pio, only : file_desc_t
319 : implicit none
320 :
321 : type(file_desc_t) :: piofile
322 :
323 768 : call read_trc_restart( 'prescribed_volcaero', piofile, file )
324 :
325 768 : end subroutine read_prescribed_volcaero_restart
326 :
327 : end module prescribed_volcaero
|