Line data Source code
1 : module atm_stream_ndep
2 :
3 : !-----------------------------------------------------------------------
4 : ! Contains methods for reading in nitrogen deposition data file
5 : ! Also includes functions for dynamic ndep file handling and
6 : ! interpolation.
7 : !-----------------------------------------------------------------------
8 : !
9 : use ESMF , only : ESMF_Clock, ESMF_Mesh
10 : use ESMF , only : ESMF_SUCCESS, ESMF_LOGERR_PASSTHRU, ESMF_END_ABORT
11 : use ESMF , only : ESMF_Finalize, ESMF_LogFoundError
12 : use nuopc_shr_methods , only : chkerr
13 : use dshr_strdata_mod , only : shr_strdata_type
14 : use shr_kind_mod , only : r8 => shr_kind_r8, CL => shr_kind_cl, CS => shr_kind_cs
15 : use shr_log_mod , only : errMsg => shr_log_errMsg
16 : use spmd_utils , only : mpicom, masterproc, iam
17 : use spmd_utils , only : mpi_character, mpi_integer
18 : use cam_logfile , only : iulog
19 : use cam_abortutils , only : endrun
20 :
21 : implicit none
22 : private
23 :
24 : public :: stream_ndep_readnl ! read runtime options
25 : public :: stream_ndep_init ! position datasets for dynamic ndep
26 : public :: stream_ndep_interp ! interpolates between two years of ndep file data
27 :
28 : private :: stream_ndep_check_units ! Check the units and make sure they can be used
29 :
30 : ! The ndep stream is not needed for aquaplanet or simple model configurations. It
31 : ! is disabled by setting the namelist variable stream_ndep_data_filename to 'UNSET' or empty string.
32 : logical, public, protected :: ndep_stream_active = .false.
33 :
34 : type(shr_strdata_type) :: sdat_ndep ! input data stream
35 : logical, public :: stream_ndep_is_initialized = .false.
36 : character(len=CS) :: stream_varlist_ndep(2)
37 : type(ESMF_Clock) :: model_clock
38 :
39 : character(len=*), parameter :: sourcefile = __FILE__
40 :
41 : character(len=CL) :: stream_ndep_data_filename
42 : character(len=CL) :: stream_ndep_mesh_filename
43 : integer :: stream_ndep_year_first ! first year in stream to use
44 : integer :: stream_ndep_year_last ! last year in stream to use
45 : integer :: stream_ndep_year_align ! align stream_year_firstndep with
46 :
47 : !==============================================================================
48 : contains
49 : !==============================================================================
50 :
51 18432 : subroutine stream_ndep_readnl(nlfile)
52 :
53 : ! Uses:
54 : use shr_nl_mod, only: shr_nl_find_group_name
55 :
56 : ! input/output variables
57 : character(len=*), intent(in) :: nlfile
58 :
59 : ! local variables
60 : integer :: nu_nml ! unit for namelist file
61 : integer :: nml_error ! namelist i/o error flag
62 : integer :: ierr
63 : character(*), parameter :: subName = "('stream_ndep_readnl')"
64 : !-----------------------------------------------------------------------
65 :
66 : namelist /ndep_stream_nl/ &
67 : stream_ndep_data_filename, &
68 : stream_ndep_mesh_filename, &
69 : stream_ndep_year_first, &
70 : stream_ndep_year_last, &
71 : stream_ndep_year_align
72 :
73 : ! Default values for namelist
74 1536 : stream_ndep_data_filename = ' '
75 1536 : stream_ndep_mesh_filename = ' '
76 1536 : stream_ndep_year_first = 1 ! first year in stream to use
77 1536 : stream_ndep_year_last = 1 ! last year in stream to use
78 1536 : stream_ndep_year_align = 1 ! align stream_ndep_year_first with this model year
79 :
80 : ! For now variable list in stream data file is hard-wired
81 4608 : stream_varlist_ndep = (/'NDEP_NHx_month', 'NDEP_NOy_month'/)
82 :
83 : ! Read ndep_stream namelist
84 1536 : if (masterproc) then
85 2 : open( newunit=nu_nml, file=trim(nlfile), status='old', iostat=nml_error )
86 2 : if (nml_error /= 0) then
87 0 : call endrun(subName//': ERROR opening '//trim(nlfile)//errMsg(sourcefile, __LINE__))
88 : end if
89 2 : call shr_nl_find_group_name(nu_nml, 'ndep_stream_nl', status=nml_error)
90 2 : if (nml_error == 0) then
91 2 : read(nu_nml, nml=ndep_stream_nl, iostat=nml_error)
92 2 : if (nml_error /= 0) then
93 0 : call endrun(' ERROR reading ndep_stream_nl namelist'//errMsg(sourcefile, __LINE__))
94 : end if
95 : end if
96 2 : close(nu_nml)
97 : endif
98 1536 : call mpi_bcast(stream_ndep_mesh_filename, len(stream_ndep_mesh_filename), mpi_character, 0, mpicom, ierr)
99 1536 : if (ierr /= 0) call endrun(trim(subname)//": FATAL: mpi_bcast: stream_ndep_mesh_filename")
100 1536 : call mpi_bcast(stream_ndep_data_filename, len(stream_ndep_data_filename), mpi_character, 0, mpicom, ierr)
101 1536 : if (ierr /= 0) call endrun(trim(subname)//": FATAL: mpi_bcast: stream_ndep_data_filename")
102 1536 : call mpi_bcast(stream_ndep_year_first, 1, mpi_integer, 0, mpicom, ierr)
103 1536 : if (ierr /= 0) call endrun(trim(subname)//": FATAL: mpi_bcast: stream_ndep_year_first")
104 1536 : call mpi_bcast(stream_ndep_year_last, 1, mpi_integer, 0, mpicom, ierr)
105 1536 : if (ierr /= 0) call endrun(trim(subname)//": FATAL: mpi_bcast: stream_ndep_year_last")
106 1536 : call mpi_bcast(stream_ndep_year_align, 1, mpi_integer, 0, mpicom, ierr)
107 1536 : if (ierr /= 0) call endrun(trim(subname)//": FATAL: mpi_bcast: stream_ndep_year_align")
108 :
109 1536 : ndep_stream_active = len_trim(stream_ndep_data_filename)>0 .and. stream_ndep_data_filename/='UNSET'
110 :
111 : ! Check whether the stream is being used.
112 1536 : if (.not.ndep_stream_active) then
113 0 : if (masterproc) then
114 0 : write(iulog,'(a)') ' '
115 0 : write(iulog,'(a)') 'NDEP STREAM IS NOT USED.'
116 0 : write(iulog,'(a)') ' '
117 : endif
118 0 : return
119 : endif
120 :
121 1536 : if (masterproc) then
122 2 : write(iulog,'(a)' ) ' '
123 2 : write(iulog,'(a,i8)') 'stream ndep settings:'
124 2 : write(iulog,'(a,a)' ) ' stream_ndep_data_filename = ',trim(stream_ndep_data_filename)
125 2 : write(iulog,'(a,a)' ) ' stream_ndep_mesh_filename = ',trim(stream_ndep_mesh_filename)
126 2 : write(iulog,'(a,a,a)') ' stream_varlist_ndep = ',trim(stream_varlist_ndep(1)), trim(stream_varlist_ndep(2))
127 2 : write(iulog,'(a,i8)') ' stream_ndep_year_first = ',stream_ndep_year_first
128 2 : write(iulog,'(a,i8)') ' stream_ndep_year_last = ',stream_ndep_year_last
129 2 : write(iulog,'(a,i8)') ' stream_ndep_year_align = ',stream_ndep_year_align
130 2 : write(iulog,'(a)' ) ' '
131 : endif
132 :
133 : end subroutine stream_ndep_readnl
134 :
135 1536 : subroutine stream_ndep_init(model_mesh, model_clock, rc)
136 : use dshr_strdata_mod, only: shr_strdata_init_from_inline
137 :
138 : ! input/output variables
139 : type(ESMF_CLock), intent(in) :: model_clock
140 : type(ESMF_Mesh) , intent(in) :: model_mesh
141 : integer , intent(out) :: rc
142 :
143 : ! local variables
144 : character(*), parameter :: subName = "('stream_ndep_init')"
145 :
146 1536 : rc = ESMF_SUCCESS
147 1536 : if (.not.ndep_stream_active) then
148 : return
149 : end if
150 : !
151 : ! Initialize data stream information.
152 : ! Read in units
153 1536 : call stream_ndep_check_units(stream_ndep_data_filename)
154 :
155 : ! Initialize the cdeps data type sdat_ndep
156 : call shr_strdata_init_from_inline(sdat_ndep, &
157 : my_task = iam, &
158 : logunit = iulog, &
159 : compname = 'ATM', &
160 : model_clock = model_clock, &
161 : model_mesh = model_mesh, &
162 : stream_meshfile = trim(stream_ndep_mesh_filename), &
163 : stream_filenames = (/trim(stream_ndep_data_filename)/), &
164 : stream_yearFirst = stream_ndep_year_first, &
165 : stream_yearLast = stream_ndep_year_last, &
166 : stream_yearAlign = stream_ndep_year_align, &
167 : stream_fldlistFile = stream_varlist_ndep, &
168 : stream_fldListModel = stream_varlist_ndep, &
169 : stream_lev_dimname = 'null', &
170 : stream_mapalgo = 'bilinear', &
171 : stream_offset = 0, &
172 : stream_taxmode = 'cycle', &
173 : stream_dtlimit = 1.0e30_r8, &
174 : stream_tintalgo = 'linear', &
175 : stream_name = 'Nitrogen deposition data ', &
176 3072 : rc = rc)
177 1536 : if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, line=__LINE__, file=__FILE__)) then
178 0 : call ESMF_Finalize(endflag=ESMF_END_ABORT)
179 : end if
180 :
181 1536 : end subroutine stream_ndep_init
182 :
183 : !================================================================
184 1536 : subroutine stream_ndep_check_units( stream_fldFileName_ndep)
185 :
186 : !--------------------------------------------------------
187 : ! Check that units are correct on the file and if need any conversion
188 : !--------------------------------------------------------
189 :
190 1536 : use cam_pio_utils , only : cam_pio_createfile, cam_pio_openfile, cam_pio_closefile, pio_subsystem
191 : use pio , only : file_desc_t, io_desc_t, var_desc_t, pio_double, pio_def_dim
192 : use pio , only : pio_bcast_error, pio_seterrorhandling, pio_inq_varid, pio_get_att
193 : use pio , only : PIO_NOERR, PIO_NOWRITE
194 :
195 : ! Arguments
196 : character(len=*), intent(in) :: stream_fldFileName_ndep ! ndep filename
197 : !
198 : ! Local variables
199 : type(file_desc_t) :: File ! NetCDF filehandle for ndep file
200 : type(var_desc_t) :: vardesc ! variable descriptor
201 : integer :: ierr ! error status
202 : integer :: err_handling ! temporary
203 : character(len=CS) :: ndepunits! ndep units
204 : !-----------------------------------------------------------------------
205 :
206 1536 : call cam_pio_openfile( File, trim(stream_fldFileName_ndep), PIO_NOWRITE)
207 1536 : call pio_seterrorhandling(File, PIO_BCAST_ERROR, err_handling)
208 1536 : ierr = pio_inq_varid(File, stream_varlist_ndep(1), vardesc)
209 1536 : if (ierr /= PIO_NOERR) then
210 : call endrun(' ERROR finding variable: '//trim(stream_varlist_ndep(1))//" in file: "// &
211 0 : trim(stream_fldFileName_ndep)//errMsg(sourcefile, __LINE__))
212 : else
213 1536 : ierr = PIO_get_att(File, vardesc, "units", ndepunits)
214 : end if
215 1536 : call pio_seterrorhandling(File, err_handling)
216 1536 : call cam_pio_closefile(File)
217 :
218 : ! Now check to make sure they are correct
219 1536 : if (.not. trim(ndepunits) == "g(N)/m2/s" )then
220 : call endrun(' ERROR in units for nitrogen deposition equal to: '//trim(ndepunits)//" not units expected"// &
221 0 : errMsg(sourcefile, __LINE__))
222 : end if
223 :
224 3072 : end subroutine stream_ndep_check_units
225 :
226 : !================================================================
227 16896 : subroutine stream_ndep_interp(cam_out, rc)
228 :
229 1536 : use dshr_methods_mod , only : dshr_fldbun_getfldptr
230 : use dshr_strdata_mod , only : shr_strdata_advance
231 : use camsrfexch , only : cam_out_t
232 : use ppgrid , only : begchunk, endchunk
233 : use time_manager , only : get_curr_date
234 : use phys_grid , only : get_ncols_p
235 :
236 : ! input/output variables
237 : type(cam_out_t) , intent(inout) :: cam_out(begchunk:endchunk)
238 : integer , intent(out) :: rc
239 :
240 : ! local variables
241 : integer :: i,c,g
242 : integer :: year ! year (0, ...) for nstep+1
243 : integer :: mon ! month (1, ..., 12) for nstep+1
244 : integer :: day ! day of month (1, ..., 31) for nstep+1
245 : integer :: sec ! seconds into current date for nstep+1
246 : integer :: mcdate ! Current model date (yyyymmdd)
247 16896 : real(r8), pointer :: dataptr1d_nhx(:)
248 16896 : real(r8), pointer :: dataptr1d_noy(:)
249 :
250 : ! NDEP read from forcing is expected to be in units of gN/m2/sec - but the mediator
251 : ! expects units of kgN/m2/sec
252 : real(r8), parameter :: scale_ndep = .001_r8
253 :
254 : !-----------------------------------------------------------------------
255 :
256 : ! Advance sdat stream
257 16896 : call get_curr_date(year, mon, day, sec)
258 16896 : mcdate = year*10000 + mon*100 + day
259 16896 : call shr_strdata_advance(sdat_ndep, ymd=mcdate, tod=sec, logunit=iulog, istr='ndepdyn', rc=rc)
260 16896 : if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, line=__LINE__, file=__FILE__)) then
261 0 : call ESMF_Finalize(endflag=ESMF_END_ABORT)
262 : end if
263 :
264 : ! Get pointer for stream data that is time and spatially interpolated to model time and grid
265 16896 : call dshr_fldbun_getFldPtr(sdat_ndep%pstrm(1)%fldbun_model, stream_varlist_ndep(1), fldptr1=dataptr1d_nhx, rc=rc)
266 16896 : if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, line=__LINE__, file=__FILE__)) then
267 0 : call ESMF_Finalize(endflag=ESMF_END_ABORT)
268 : end if
269 16896 : call dshr_fldbun_getFldPtr(sdat_ndep%pstrm(1)%fldbun_model, stream_varlist_ndep(2), fldptr1=dataptr1d_noy, rc=rc)
270 16896 : if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, line=__LINE__, file=__FILE__)) then
271 0 : call ESMF_Finalize(endflag=ESMF_END_ABORT)
272 : end if
273 :
274 16896 : g = 1
275 85008 : do c = begchunk,endchunk
276 1154208 : do i = 1,get_ncols_p(c)
277 1069200 : cam_out(c)%nhx_nitrogen_flx(i) = dataptr1d_nhx(g) * scale_ndep
278 1069200 : cam_out(c)%noy_nitrogen_flx(i) = dataptr1d_noy(g) * scale_ndep
279 1137312 : g = g + 1
280 : end do
281 : end do
282 :
283 33792 : end subroutine stream_ndep_interp
284 :
285 : end module atm_stream_ndep
|