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