Line data Source code
1 :
2 :
3 : module mo_sulf
4 : !---------------------------------------------------------------
5 : ! ... Annual cycle for sulfur
6 : !---------------------------------------------------------------
7 :
8 : use shr_kind_mod, only : r8 => shr_kind_r8
9 :
10 : use cam_abortutils, only : endrun
11 : use cam_logfile, only : iulog
12 : use tracer_data, only : trfld,trfile
13 : use physics_types, only : physics_state
14 : use ppgrid, only : begchunk, endchunk
15 : use physics_buffer, only : physics_buffer_desc
16 : use ppgrid, only : pcols, pver
17 :
18 : use spmd_utils, only : masterproc
19 :
20 : implicit none
21 :
22 : private
23 : public :: sulf_inti, set_sulf_time, sulf_interp, sulf_readnl
24 :
25 : save
26 :
27 : type(trfld), pointer :: fields(:) => null()
28 : type(trfile) :: file
29 :
30 : logical :: read_sulf = .false.
31 :
32 : character(len=16) :: fld_name = 'SULFATE'
33 : character(len=256) :: filename = 'NONE'
34 : character(len=256) :: filelist = ' '
35 : character(len=256) :: datapath = ' '
36 : character(len=32) :: datatype = 'CYCLICAL'
37 : logical :: rmv_file = .false.
38 : integer :: cycle_yr = 0
39 : integer :: fixed_ymd = 0
40 : integer :: fixed_tod = 0
41 :
42 : logical :: has_sulf_file = .false.
43 :
44 : contains
45 :
46 : !-------------------------------------------------------------------
47 : !-------------------------------------------------------------------
48 1490712 : subroutine sulf_readnl(nlfile)
49 :
50 : use namelist_utils, only: find_group_name
51 : use units, only: getunit, freeunit
52 : use mpishorthand
53 :
54 : character(len=*), intent(in) :: nlfile ! filepath for file containing namelist input
55 :
56 : ! Local variables
57 : integer :: unitn, ierr
58 : character(len=*), parameter :: subname = 'sulf_readnl'
59 :
60 : character(len=16) :: sulf_name
61 : character(len=256) :: sulf_file
62 : character(len=256) :: sulf_filelist
63 : character(len=256) :: sulf_datapath
64 : character(len=32) :: sulf_type
65 : logical :: sulf_rmfile
66 : integer :: sulf_cycle_yr
67 : integer :: sulf_fixed_ymd
68 : integer :: sulf_fixed_tod
69 :
70 : namelist /sulf_nl/ &
71 : sulf_name, &
72 : sulf_file, &
73 : sulf_filelist, &
74 : sulf_datapath, &
75 : sulf_type, &
76 : sulf_rmfile, &
77 : sulf_cycle_yr, &
78 : sulf_fixed_ymd, &
79 : sulf_fixed_tod
80 : !-----------------------------------------------------------------------------
81 :
82 : ! Initialize namelist variables from local module variables.
83 1536 : sulf_name = fld_name
84 1536 : sulf_file = filename
85 1536 : sulf_filelist = filelist
86 1536 : sulf_datapath = datapath
87 1536 : sulf_type = datatype
88 1536 : sulf_rmfile = rmv_file
89 1536 : sulf_cycle_yr = cycle_yr
90 1536 : sulf_fixed_ymd= fixed_ymd
91 1536 : sulf_fixed_tod= fixed_tod
92 :
93 : ! Read namelist
94 1536 : if (masterproc) then
95 2 : unitn = getunit()
96 2 : open( unitn, file=trim(nlfile), status='old' )
97 2 : call find_group_name(unitn, 'sulf_nl', status=ierr)
98 2 : if (ierr == 0) then
99 0 : read(unitn, sulf_nl, iostat=ierr)
100 0 : if (ierr /= 0) then
101 0 : call endrun(subname // ':: ERROR reading namelist')
102 : end if
103 : end if
104 2 : close(unitn)
105 2 : call freeunit(unitn)
106 : end if
107 :
108 : #ifdef SPMD
109 : ! Broadcast namelist variables
110 1536 : call mpibcast(sulf_name, len(sulf_name), mpichar, 0, mpicom)
111 1536 : call mpibcast(sulf_file, len(sulf_file), mpichar, 0, mpicom)
112 1536 : call mpibcast(sulf_filelist, len(sulf_filelist), mpichar, 0, mpicom)
113 1536 : call mpibcast(sulf_datapath, len(sulf_datapath), mpichar, 0, mpicom)
114 1536 : call mpibcast(sulf_type, len(sulf_type), mpichar, 0, mpicom)
115 1536 : call mpibcast(sulf_rmfile, 1, mpilog, 0, mpicom)
116 1536 : call mpibcast(sulf_cycle_yr, 1, mpiint, 0, mpicom)
117 1536 : call mpibcast(sulf_fixed_ymd,1, mpiint, 0, mpicom)
118 1536 : call mpibcast(sulf_fixed_tod,1, mpiint, 0, mpicom)
119 : #endif
120 :
121 : ! Update module variables with user settings.
122 1536 : fld_name = sulf_name
123 1536 : filename = sulf_file
124 1536 : filelist = sulf_filelist
125 1536 : datapath = sulf_datapath
126 1536 : datatype = sulf_type
127 1536 : rmv_file = sulf_rmfile
128 1536 : cycle_yr = sulf_cycle_yr
129 1536 : fixed_ymd = sulf_fixed_ymd
130 1536 : fixed_tod = sulf_fixed_tod
131 :
132 : ! Turn on prescribed volcanics if user has specified an input dataset.
133 1536 : if (len_trim(filename) > 0 .and. filename.ne.'NONE') has_sulf_file = .true.
134 :
135 1536 : end subroutine sulf_readnl
136 :
137 : !-----------------------------------------------------------------------
138 : !-----------------------------------------------------------------------
139 1536 : subroutine sulf_inti()
140 : !-----------------------------------------------------------------------
141 : ! ... Open netCDF file containing annual sulfur data. Initialize
142 : ! arrays with the data to be interpolated to the current time.
143 : !
144 : ! It is assumed that the time coordinate is increasing
145 : ! and represents calendar days; range = [1.,366.).
146 : !-----------------------------------------------------------------------
147 : use spmd_utils, only : masterproc
148 : use mo_chem_utls, only : get_spc_ndx, get_rxt_ndx
149 : use interpolate_data, only : lininterp_init, lininterp, lininterp_finish, interp_type
150 : use tracer_data, only : trcdata_init
151 : use cam_history, only : addfld
152 :
153 : implicit none
154 :
155 : !-----------------------------------------------------------------------
156 : ! ... Dummy args
157 : !-----------------------------------------------------------------------
158 :
159 : !-----------------------------------------------------------------------
160 : ! ... Local variables
161 : !-----------------------------------------------------------------------
162 : integer :: ndxs(5), so4_ndx
163 :
164 : character(len=8), parameter :: fld_names(1) = (/'SULFATE '/)
165 :
166 1536 : ndxs(1) = get_rxt_ndx( 'usr_N2O5_aer' )
167 1536 : ndxs(2) = get_rxt_ndx( 'usr_NO3_aer' )
168 1536 : ndxs(3) = get_rxt_ndx( 'usr_NO2_aer' )
169 1536 : ndxs(4) = get_rxt_ndx( 'usr_HO2_aer' )
170 1536 : ndxs(5) = get_rxt_ndx( 'het1' )
171 1536 : so4_ndx = get_spc_ndx( 'SO4' )
172 1536 : if ( so4_ndx < 1 ) so4_ndx = get_spc_ndx( 'so4_a1' )
173 :
174 9216 : read_sulf = any( ndxs > 0) .and. (so4_ndx < 0) .and. has_sulf_file
175 :
176 1536 : if ( .not. read_sulf ) return
177 :
178 0 : allocate(file%in_pbuf(size(fld_names)))
179 0 : file%in_pbuf(:) = .false.
180 : call trcdata_init( (/ fld_name /), filename, filelist, datapath, fields, file, &
181 0 : rmv_file, cycle_yr, fixed_ymd, fixed_tod, datatype)
182 :
183 0 : call addfld('SULFATE', (/ 'lev' /), 'I','VMR', 'sulfate data' )
184 :
185 1536 : end subroutine sulf_inti
186 :
187 741888 : subroutine set_sulf_time( pbuf2d, state )
188 : !--------------------------------------------------------------------
189 : ! ... Check and set time interpolation indicies
190 : !--------------------------------------------------------------------
191 1536 : use tracer_data, only : advance_trcdata
192 :
193 : implicit none
194 :
195 : !--------------------------------------------------------------------
196 : ! ... Dummy args
197 : !--------------------------------------------------------------------
198 : type(physics_buffer_desc), pointer :: pbuf2d(:,:)
199 : type(physics_state), intent(in):: state(begchunk:endchunk)
200 :
201 370944 : if ( .not. read_sulf ) return
202 :
203 0 : call advance_trcdata( fields, file, state, pbuf2d )
204 :
205 370944 : end subroutine set_sulf_time
206 :
207 1489176 : subroutine sulf_interp( ncol, lchnk, ccm_sulf )
208 : !-----------------------------------------------------------------------
209 : ! ... Time interpolate sulfatei to current time
210 : !-----------------------------------------------------------------------
211 370944 : use cam_history, only : outfld
212 :
213 : implicit none
214 :
215 : !-----------------------------------------------------------------------
216 : ! ... Dummy arguments
217 : !-----------------------------------------------------------------------
218 : integer, intent(in) :: ncol ! columns in chunk
219 : integer, intent(in) :: lchnk ! chunk number
220 : real(r8), intent(out) :: ccm_sulf(:,:) ! output sulfate
221 :
222 : !-----------------------------------------------------------------------
223 : ! ... Local variables
224 : !-----------------------------------------------------------------------
225 :
226 2314006344 : ccm_sulf(:,:) = 0._r8
227 :
228 1489176 : if ( .not. read_sulf ) return
229 :
230 0 : ccm_sulf(:ncol,:) = fields(1)%data(:ncol,:,lchnk)
231 :
232 0 : call outfld( 'SULFATE', ccm_sulf(:ncol,:), ncol, lchnk )
233 :
234 1489176 : end subroutine sulf_interp
235 :
236 : end module mo_sulf
|