Line data Source code
1 : module waccm_forcing
2 : !================================================================================================
3 : !
4 : ! Provides WACCM forcing data for use without interactive chemistry
5 : ! -- for GHG chemistry
6 : !
7 : ! FVITT 21 Mar 2011 -- creation
8 : !
9 : !================================================================================================
10 :
11 : use shr_kind_mod, only: r8 => shr_kind_r8
12 : use cam_abortutils, only: endrun
13 :
14 : use tracer_data, only : trfld, trfile
15 : use ppgrid, only : pcols, pver
16 : use ppgrid, only : begchunk, endchunk
17 : use spmd_utils, only : masterproc
18 :
19 : implicit none
20 : private
21 : save
22 :
23 : ! Public interfaces
24 : public :: waccm_forcing_init
25 : public :: waccm_forcing_adv
26 : public :: get_cnst ! return prescribed constituents for nlte
27 : public :: get_solar ! return prescribed net solar heating rate
28 : public :: waccm_forcing_readnl
29 :
30 : ! Private module data
31 :
32 : type(trfld), pointer :: fields(:) => null()
33 : type(trfile) :: file
34 :
35 : integer, parameter :: N_FLDS = 7
36 : integer, parameter :: N_MMRS = 6
37 :
38 : character(len=256) :: filename = ''
39 : character(len=256) :: filelist = ''
40 : character(len=256) :: datapath = ''
41 : character(len=32) :: datatype = 'CYCLICAL'
42 : logical :: rmv_file = .false.
43 :
44 : integer :: cycle_yr = 0
45 : integer :: fixed_ymd = 0
46 : integer :: fixed_tod = 0
47 :
48 : integer :: o1_ndx=1, o2_ndx=2, o3_ndx=3, no_ndx=4, h_ndx=5, co2_ndx=6, qrs_ndx=7
49 : character(len=16) :: specifier(N_FLDS) = (/ 'O ','O2 ','O3 ', 'NO ', 'H ','CO2 ', 'QRS_TOT' /)
50 : real(r8), parameter :: molmass(N_MMRS) = (/ 15.99940_r8, 31.99880_r8, 47.99820_r8, 30.00614_r8, 1.007400_r8, 44.00980_r8 /)
51 :
52 : !================================================================================================
53 : contains
54 : !================================================================================================
55 :
56 : !-------------------------------------------------------------------
57 : !-------------------------------------------------------------------
58 1536 : subroutine waccm_forcing_readnl(nlfile)
59 :
60 : use namelist_utils, only: find_group_name
61 : use units, only: getunit, freeunit
62 : use mpishorthand
63 :
64 : character(len=*), intent(in) :: nlfile ! filepath for file containing namelist input
65 :
66 : ! Local variables
67 : integer :: unitn, ierr
68 : character(len=*), parameter :: subname = 'prescribed_aero_readnl'
69 :
70 : character(len=16) :: waccm_forcing_specifier(N_FLDS)
71 : character(len=256) :: waccm_forcing_file
72 : character(len=256) :: waccm_forcing_filelist
73 : character(len=256) :: waccm_forcing_datapath
74 : character(len=32) :: waccm_forcing_type
75 : logical :: waccm_forcing_rmfile
76 : integer :: waccm_forcing_cycle_yr
77 : integer :: waccm_forcing_fixed_ymd
78 : integer :: waccm_forcing_fixed_tod
79 :
80 : namelist /waccm_forcing_nl/ &
81 : waccm_forcing_specifier, &
82 : waccm_forcing_file, &
83 : waccm_forcing_filelist, &
84 : waccm_forcing_datapath, &
85 : waccm_forcing_type, &
86 : waccm_forcing_rmfile, &
87 : waccm_forcing_cycle_yr, &
88 : waccm_forcing_fixed_ymd, &
89 : waccm_forcing_fixed_tod
90 : !-----------------------------------------------------------------------------
91 :
92 : ! Initialize namelist variables from local module variables.
93 12288 : waccm_forcing_specifier= specifier
94 1536 : waccm_forcing_file = filename
95 1536 : waccm_forcing_filelist = filelist
96 1536 : waccm_forcing_datapath = datapath
97 1536 : waccm_forcing_type = datatype
98 1536 : waccm_forcing_rmfile = rmv_file
99 1536 : waccm_forcing_cycle_yr = cycle_yr
100 1536 : waccm_forcing_fixed_ymd= fixed_ymd
101 1536 : waccm_forcing_fixed_tod= fixed_tod
102 :
103 : ! Read namelist
104 1536 : if (masterproc) then
105 2 : unitn = getunit()
106 2 : open( unitn, file=trim(nlfile), status='old' )
107 2 : call find_group_name(unitn, 'waccm_forcing_nl', status=ierr)
108 2 : if (ierr == 0) then
109 0 : read(unitn, waccm_forcing_nl, iostat=ierr)
110 0 : if (ierr /= 0) then
111 0 : call endrun(subname // ':: ERROR reading namelist')
112 : end if
113 : end if
114 2 : close(unitn)
115 2 : call freeunit(unitn)
116 : end if
117 :
118 : #ifdef SPMD
119 : ! Broadcast namelist variables
120 1536 : call mpibcast(waccm_forcing_specifier,len(waccm_forcing_specifier(1))*N_FLDS, mpichar, 0, mpicom)
121 1536 : call mpibcast(waccm_forcing_file, len(waccm_forcing_file), mpichar, 0, mpicom)
122 1536 : call mpibcast(waccm_forcing_filelist, len(waccm_forcing_filelist), mpichar, 0, mpicom)
123 1536 : call mpibcast(waccm_forcing_datapath, len(waccm_forcing_datapath), mpichar, 0, mpicom)
124 1536 : call mpibcast(waccm_forcing_type, len(waccm_forcing_type), mpichar, 0, mpicom)
125 1536 : call mpibcast(waccm_forcing_rmfile, 1, mpilog, 0, mpicom)
126 1536 : call mpibcast(waccm_forcing_cycle_yr, 1, mpiint, 0, mpicom)
127 1536 : call mpibcast(waccm_forcing_fixed_ymd,1, mpiint, 0, mpicom)
128 1536 : call mpibcast(waccm_forcing_fixed_tod,1, mpiint, 0, mpicom)
129 : #endif
130 :
131 : ! Update module variables with user settings.
132 12288 : specifier = waccm_forcing_specifier
133 1536 : filename = waccm_forcing_file
134 1536 : filelist = waccm_forcing_filelist
135 1536 : datapath = waccm_forcing_datapath
136 1536 : datatype = waccm_forcing_type
137 1536 : rmv_file = waccm_forcing_rmfile
138 1536 : cycle_yr = waccm_forcing_cycle_yr
139 1536 : fixed_ymd = waccm_forcing_fixed_ymd
140 1536 : fixed_tod = waccm_forcing_fixed_tod
141 :
142 1536 : end subroutine waccm_forcing_readnl
143 :
144 : !------------------------------------------------------------------------
145 : !------------------------------------------------------------------------
146 0 : subroutine waccm_forcing_init()
147 :
148 : use tracer_data, only : trcdata_init
149 : use cam_history, only : addfld
150 :
151 : implicit none
152 :
153 : integer :: i
154 :
155 :
156 0 : allocate(file%in_pbuf(size(specifier)))
157 0 : file%in_pbuf(:) = .false.
158 : call trcdata_init( specifier, filename, filelist, datapath, fields, file, &
159 0 : rmv_file, cycle_yr, fixed_ymd, fixed_tod, datatype)
160 :
161 0 : do i = 1,N_FLDS
162 0 : call addfld( 'WFRC_'//trim(fields(i)%fldnam), (/ 'lev' /), 'I', fields(i)%units, 'for waccm forcing' )
163 : enddo
164 :
165 0 : return
166 0 : end subroutine waccm_forcing_init
167 :
168 : !=======================================================================
169 :
170 : !------------------------------------------------------------------------
171 : !------------------------------------------------------------------------
172 0 : subroutine waccm_forcing_adv (state, pbuf2d)
173 :
174 0 : use tracer_data, only : advance_trcdata
175 : use physics_types,only : physics_state
176 : use string_utils, only : to_lower, GLC
177 : use cam_history, only : outfld
178 : use physconst, only : mwdry ! molecular weight dry air ~ kg/kmole
179 : use physconst, only : boltz ! J/K/molecule
180 : use physics_buffer, only : physics_buffer_desc
181 :
182 : implicit none
183 :
184 : type(physics_state), intent(in):: state(begchunk:endchunk)
185 : type(physics_buffer_desc), pointer :: pbuf2d(:,:)
186 :
187 :
188 : integer :: c,ncol,i
189 : real(r8) :: to_mmr(pcols,pver)
190 :
191 0 : call advance_trcdata( fields, file, state, pbuf2d )
192 :
193 : ! set the tracer fields with the correct units
194 0 : do i = 1,N_FLDS
195 :
196 0 : do c = begchunk,endchunk
197 0 : ncol = state(c)%ncol
198 :
199 0 : if ( i<=N_MMRS ) then
200 :
201 0 : select case ( to_lower(trim(fields(i)%units(:GLC(fields(i)%units)))) )
202 : case ("molec/cm3","/cm3","molecules/cm3","cm^-3","cm**-3")
203 0 : to_mmr(:ncol,:) = (molmass(i)*1.e6_r8*boltz*state(c)%t(:ncol,:))/(mwdry*state(c)%pmiddry(:ncol,:))
204 : case ('kg/kg','mmr')
205 0 : to_mmr(:ncol,:) = 1._r8
206 : case ('mol/mol','mole/mole','vmr','fraction')
207 0 : to_mmr(:ncol,:) = molmass(i)/mwdry
208 : case default
209 0 : print*, 'waccm_forcing_adv: units = ',trim(fields(i)%units) ,' are not recognized'
210 0 : call endrun('waccm_forcing_adv: units are not recognized')
211 : end select
212 :
213 0 : call outfld( 'WFRC_'//trim(fields(i)%fldnam), fields(i)%data(:ncol,:,c), ncol, state(c)%lchnk )
214 :
215 0 : fields(i)%data(:ncol,:,c) = to_mmr(:ncol,:) * fields(i)%data(:ncol,:,c)
216 : else
217 0 : call outfld( 'WFRC_'//trim(fields(i)%fldnam), fields(i)%data(:ncol,:,c), ncol, state(c)%lchnk )
218 : endif
219 :
220 : enddo
221 : enddo
222 :
223 0 : return
224 0 : end subroutine waccm_forcing_adv
225 :
226 :
227 : !================================================================================================
228 :
229 0 : subroutine get_cnst (lchnk, co2, o1, o2, no, h, o3)
230 : !
231 : ! Get mass mixing ratios specified from input dataset for used in Fomichev routines
232 : !-------------------------------------------------------------------------
233 :
234 : ! Arguments
235 : integer, intent(in) :: lchnk ! chunk identifier
236 :
237 : real(r8), optional, pointer, dimension(:,:) :: co2
238 : real(r8), optional, pointer, dimension(:,:) :: o1
239 : real(r8), optional, pointer, dimension(:,:) :: o2
240 : real(r8), optional, pointer, dimension(:,:) :: no
241 : real(r8), optional, pointer, dimension(:,:) :: h
242 : real(r8), optional, pointer, dimension(:,:) :: o3
243 :
244 : !------------------------------------------------------------------------
245 :
246 0 : if (present(co2)) then
247 0 : co2 => fields(co2_ndx)%data(:,:,lchnk)
248 : endif
249 0 : if (present(o1)) then
250 0 : o1 => fields(o1_ndx )%data(:,:,lchnk)
251 : endif
252 0 : if (present(o2)) then
253 0 : o2 => fields(o2_ndx )%data(:,:,lchnk)
254 : endif
255 0 : if (present(no)) then
256 0 : no => fields(no_ndx )%data(:,:,lchnk)
257 : endif
258 0 : if (present(h)) then
259 0 : h => fields(h_ndx )%data(:,:,lchnk)
260 : endif
261 0 : if (present(o3)) then
262 0 : o3 => fields(o3_ndx)%data(:,:,lchnk)
263 : endif
264 :
265 0 : end subroutine get_cnst
266 :
267 : !================================================================================================
268 :
269 0 : subroutine get_solar (ncol, lchnk, qrs_mlt)
270 : !
271 : ! Get M/LT solar heating rates specified from input dataset
272 : !-------------------------------------------------------------------------
273 :
274 : ! Arguments
275 : integer, intent(in) :: ncol ! no. of columns in chunk
276 : integer, intent(in) :: lchnk ! chunk identifier
277 : real(r8), intent(out) :: qrs_mlt(:,:) ! M/LT solar heating rates
278 :
279 : ! Local workspace
280 : integer :: k
281 :
282 : !------------------------------------------------------------------------
283 :
284 0 : qrs_mlt(:ncol,:) = fields(qrs_ndx)%data(:ncol,:,lchnk)
285 :
286 0 : end subroutine get_solar
287 :
288 : end module waccm_forcing
|