FMS  2025.04
Flexible Modeling System
interpolator.F90
1 !***********************************************************************
2 !* Apache License 2.0
3 !*
4 !* This file is part of the GFDL Flexible Modeling System (FMS).
5 !*
6 !* Licensed under the Apache License, Version 2.0 (the "License");
7 !* you may not use this file except in compliance with the License.
8 !* You may obtain a copy of the License at
9 !*
10 !* http://www.apache.org/licenses/LICENSE-2.0
11 !*
12 !* FMS is distributed in the hope that it will be useful, but WITHOUT
13 !* WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied;
14 !* without even the implied warranty of MERCHANTABILITY or FITNESS FOR A
15 !* PARTICULAR PURPOSE. See the License for the specific language
16 !* governing permissions and limitations under the License.
17 !***********************************************************************
18 !> @defgroup interpolator_mod interpolator_mod
19 !> @ingroup interpolator
20 !> @brief A module to interpolate climatology data to model the grid.
21 !> @author William Cooke <William.Cooke@noaa.gov>
22 
23 module interpolator_mod
24 
25 use mpp_mod, only : mpp_error, &
26  fatal, &
27  mpp_pe, &
28  mpp_init, &
29  mpp_exit, &
30  mpp_npes, &
31  warning, &
32  note, &
33  input_nml_file
34 use mpp_domains_mod, only : mpp_domains_init, &
38  domain2d, &
41 use diag_manager_mod, only : diag_manager_init, get_base_time, &
43  diag_axis_init
44 use fms_mod, only : lowercase, write_version_number, &
45  fms_init, &
46  mpp_root_pe, stdlog, &
48 use fms2_io_mod, only : fmsnetcdffile_t, fms2_io_file_exist => file_exists, dimension_exists, &
49  open_file, fms2_io_read_data=>read_data, &
50  variable_exists, get_variable_num_dimensions, &
51  get_num_variables, get_dimension_size, &
52  get_variable_units, get_variable_names, &
53  get_time_calendar, close_file, &
54  get_variable_dimension_names, get_variable_sense
55 
56 use horiz_interp_mod, only : horiz_interp_type, &
59  assignment(=), &
60  horiz_interp, &
62 use time_manager_mod, only : time_type, &
63  set_time, &
64  set_date, &
66  days_in_year, &
68  leap_year, &
69  julian, noleap, &
70  get_date, &
71  get_date_julian, set_date_no_leap, &
72  set_date_julian, get_date_no_leap, &
73  print_date, &
74  operator(+), &
75  operator(-), &
76  operator(*), &
77  operator(>), &
78  operator(<), &
79  assignment(=), &
81 use time_interp_mod, only : time_interp, year
82 use constants_mod, only : grav, pi, seconds_per_day
83 use platform_mod, only : r4_kind, r8_kind, r16_kind, fms_path_len, fms_file_len
84 
85 !--------------------------------------------------------------------
86 
87 implicit none
88 private
89 
90 !---------------------------------------------------------------------
91 !------- interfaces --------
92 
93 public interpolator_init, &
94  interpolator, &
101  read_data
102 
103 !> Interpolates a field to a model grid
104 !!
105 !> Example usage:
106 !! ~~~~~~~~~~{.f90}
107 !! call interpolator (sulfate, model_time, p_half, model_data, name, is, js, clim_units)
108 !! call interpolator (o3, model_time, p_half, model_data, "ozone", is, js)
109 !! ~~~~~~~~~~
110 !!
111 !! The first option is used to generate sulfate models.
112 !!
113 !! The sulfate data is set by
114 !! ~~~~~~~~~~{.f90}
115 !! type(interpolate_type), intent(inout) :: sulfate
116 !! ~~~~~~~~~~
117 !! The name of the model is set by
118 !! ~~~~~~~~~~{.f90}
119 !! character(len=*), intent(in) :: name
120 !! ~~~~~~~~~~
121 !! The units used in this model are outputted to
122 !! ~~~~~~~~~~{.f90}
123 !! character(len=*), intent(out), optional :: clim_units
124 !! ~~~~~~~~~~
125 !!
126 !! The second option is generate ozone models.
127 !!
128 !! The ozone data is set by
129 !! ~~~~~~~~~~{.f90}
130 !! type(interpolate_type), intent(inout) :: o3
131 !! ~~~~~~~~~~
132 !!
133 !! Both of these options use the following variables in the model.
134 !!
135 !! The time used in the model is set by
136 !!
137 !! ~~~~~~~~~~{.f90}
138 !! type(time_type), intent(in) :: model_time
139 !! ~~~~~~~~~~
140 !! The model pressure field is set by
141 !! ~~~~~~~~~~{.f90}
142 !! real, intent(in), dimension(:,:,:) :: p_half
143 !! ~~~~~~~~~~
144 !!
145 !! @param [inout] <clim_type> The interpolate type previously defined by a call to interpolator_init
146 !! @param [in] <field_name> The name of a field that you wish to interpolate
147 !! @param [in] <Time> The model time that you wish to interpolate to
148 !! @param [in] <phalf> The half level model pressure field
149 !! @param [in] <is> Index for the physics window
150 !! @param [in] <js> Index for the physics window
151 !! @param [out] <interp_data> The model fields with the interpolated climatology data
152 !! @param [out] <clim_units> The units of field_name
153 !> @ingroup interpolator_mod
154 interface interpolator
155  module procedure interpolator_4d_r4, interpolator_4d_r8
156  module procedure interpolator_3d_r4, interpolator_3d_r8
157  module procedure interpolator_2d_r4, interpolator_2d_r8
158  module procedure interpolator_4d_no_time_axis_r4, interpolator_4d_no_time_axis_r8
159  module procedure interpolator_3d_no_time_axis_r4, interpolator_3d_no_time_axis_r8
160  module procedure interpolator_2d_no_time_axis_r4, interpolator_2d_no_time_axis_r8
161 end interface interpolator
162 
163 !> Private assignment override interface for interpolate type
164 !> @ingroup interpolator_mod
165 interface assignment(=)
166  module procedure interpolate_type_eq
167 end interface
168 
170  module procedure interpolator_init_r4
171  module procedure interpolator_init_r8
172 end interface interpolator_init
173 
175  module procedure get_axis_latlon_data_r4
176  module procedure get_axis_latlon_data_r8
177 end interface get_axis_latlon_data
178 
180  module procedure get_axis_level_data_r4
181  module procedure get_axis_level_data_r8
182 end interface get_axis_level_data
183 
184 interface cell_center2
185  module procedure cell_center2_r4
186  module procedure cell_center2_r8
187 end interface cell_center2
188 
189 interface cart_to_latlon
190  module procedure cart_to_latlon_r4
191  module procedure cart_to_latlon_r8
192 end interface cart_to_latlon
193 
194 interface latlon2xyz
195  module procedure latlon2xyz_r4
196  module procedure latlon2xyz_r8
197 end interface latlon2xyz
198 
199 interface diag_read_data
200  module procedure diag_read_data_r4
201  module procedure diag_read_data_r8
202 end interface diag_read_data
203 
204 interface read_data
205  module procedure read_data_r4
206  module procedure read_data_r8
207 end interface read_data
208 
210  module procedure read_data_no_time_axis_r4
211  module procedure read_data_no_time_axis_r8
212 end interface read_data_no_time_axis
213 
214 interface interp_linear
215  module procedure interp_linear_r4
216  module procedure interp_linear_r8
217 end interface interp_linear
218 
219 !> Private interface for weighted scalar interpolation
220 !!
221 !> Example usage:
222 !! ~~~~~~~~~~{.f90}
223 !! call interp_weighted_scalar (pclim, phalf(ilon,j,:),hinterp_data(ilon,j,:,:),interp_data(ilon,j,:,:))
224 !! call interp_weighted_scalar (pclim, phalf(ilon,j,:),hinterp_data(ilon,j,:),interp_data(ilon,j,:))
225 !! ~~~~~~~~~~
226 !!
227 !! @param [in] <grdin> Input grid
228 !! @param [in] <grdout> Output grid
229 !! @param [in] <datin> Input data
230 !! @param [out] <datout> Output data
231 !> @ingroup interpolator_mod
233  module procedure interp_weighted_scalar_1d_r4, interp_weighted_scalar_1d_r8
234  module procedure interp_weighted_scalar_2d_r4, interp_weighted_scalar_2d_r8
235 end interface interp_weighted_scalar
236 
237 !---------------------------------------------------------------------
238 !----------- version number for this module --------------------------
239 
240 ! Include variable "version" to be written to log file.
241 #include<file_version.h>
242 
243 !> Redundant climatology data between fields
244 !> @ingroup interpolate_type
245 
246 type, private :: interpolate_r4_type
247 logical :: is_allocated = .false.
248 real(r4_kind), allocatable :: lat(:) !< No description
249 real(r4_kind), allocatable :: lon(:) !< No description
250 real(r4_kind), allocatable :: latb(:) !< No description
251 real(r4_kind), allocatable :: lonb(:) !< No description
252 real(r4_kind), allocatable :: levs(:) !< No description
253 real(r4_kind), allocatable :: halflevs(:) !< No description
254 real(r4_kind), allocatable :: data5d(:,:,:,:,:) !< (nlatmod,nlonmod,nlevclim,size(time_init,2),nfields)
255 real(r4_kind), allocatable :: pmon_pyear(:,:,:,:) !< No description
256 real(r4_kind), allocatable :: pmon_nyear(:,:,:,:) !< No description
257 real(r4_kind), allocatable :: nmon_nyear(:,:,:,:) !< No description
258 real(r4_kind), allocatable :: nmon_pyear(:,:,:,:) !< No description
259 real(r4_kind) :: tweight !< No description
260 real(r4_kind) :: tweight1 !< The time weight between the climatology years
261 real(r4_kind) :: tweight2 !< No description
262 real(r4_kind) :: tweight3 !< The time weight between the month
263 end type interpolate_r4_type
264 
265 type, private :: interpolate_r8_type
266 logical :: is_allocated = .false.
267 real(r8_kind), allocatable :: lat(:) !< No description
268 real(r8_kind), allocatable :: lon(:) !< No description
269 real(r8_kind), allocatable :: latb(:) !< No description
270 real(r8_kind), allocatable :: lonb(:) !< No description
271 real(r8_kind), allocatable :: levs(:) !< No description
272 real(r8_kind), allocatable :: halflevs(:) !< No description
273 real(r8_kind), allocatable :: data5d(:,:,:,:,:) !< (nlatmod,nlonmod,nlevclim,size(time_init,2),nfields)
274 real(r8_kind), allocatable :: pmon_pyear(:,:,:,:) !< No description
275 real(r8_kind), allocatable :: pmon_nyear(:,:,:,:) !< No description
276 real(r8_kind), allocatable :: nmon_nyear(:,:,:,:) !< No description
277 real(r8_kind), allocatable :: nmon_pyear(:,:,:,:) !< No description
278 real(r8_kind) :: tweight !< No description
279 real(r8_kind) :: tweight1 !< The time weight between the climatology years
280 real(r8_kind) :: tweight2 !< No description
281 real(r8_kind) :: tweight3 !< The time weight between the month
282 end type interpolate_r8_type
283 
284 type, public :: interpolate_type
285 private
286 !Redundant data between fields
287 !All climatology data
288 type(interpolate_r4_type) :: r4_type
289 type(interpolate_r8_type) :: r8_type
290 type(horiz_interp_type) :: interph !< No description
291 type(time_type), allocatable :: time_slice(:) !< An array of the times within the climatology.
292 type(fmsnetcdffile_t) :: fileobj ! object that stores opened file information
293 character(len=FMS_PATH_LEN) :: file_name !< Climatology filename
294 integer :: time_flag !< Linear or seaonal interpolation?
295 integer :: level_type !< Pressure or Sigma level
296 integer :: is,ie,js,je !< No description
297 integer :: vertical_indices !< direction of vertical
298  !! data axis
299 logical :: climatological_year !< Is data for year = 0000?
300 
301 !Field specific data for nfields
302 character(len=64), allocatable :: field_name(:) !< name of this field
303 logical, allocatable :: has_level(:) !< indicate if the variable has level dimension
304 integer, allocatable :: time_init(:,:) !< second index is the number of time_slices being
305  !! kept. 2 or ntime.
306 integer, allocatable :: mr(:) !< Flag for conversion of climatology to mixing ratio.
307 integer, allocatable :: out_of_bounds(:) !< Flag for when surface pressure is out of bounds.
308 !++lwh
309 integer, allocatable :: vert_interp(:) !< Flag for type of vertical interpolation.
310 !--lwh
311 !integer :: indexm, indexp, climatology
312 integer,dimension(:), allocatable :: indexm !< No description
313 integer,dimension(:), allocatable :: indexp !< No description
314 integer,dimension(:), allocatable :: climatology !< No description
315 
316 type(time_type), allocatable :: clim_times(:,:) !< No description
317 logical :: separate_time_vary_calc !< No description
318 integer :: itaum !< No description
319 integer :: itaup !< No description
320 
321 end type interpolate_type
322 
323 !> @addtogroup interpolator_mod
324 !> @{
325 
326 logical :: module_is_initialized = .false.
327 logical :: clim_diag_initialized = .false.
328 
329 integer :: ndim !< No description
330 integer :: nvar !< No description
331 integer :: ntime !< No description
332 integer :: nlat !< No description
333 integer :: nlatb !< No description
334 integer :: nlon !< No description
335 integer :: nlonb !< No description
336 integer :: nlev !< No description
337 integer :: nlevh !< No description
338 integer :: len, ntime_in, num_fields !< No description
339 
340 ! pletzer real, allocatable :: time_in(:)
341 ! sjs real, allocatable :: climdata(:,:,:), climdata2(:,:,:)
342 
343 character(len=64) :: units !< No description
344 integer :: sense !< No description
345 
346 integer, parameter :: max_diag_fields = 30 !< No description
347 
348 ! flags to indicate direction of vertical axis in data file
349 integer, parameter :: increasing_downward = 1, increasing_upward = -1 !< Flags to indicate direction
350  !! of vertical axis in data file
351 !++lwh
352 ! Flags to indicate whether the time interpolation should be linear or some other scheme for seasonal data.
353 ! NOTIME indicates that data file has no time axis.
354 integer, parameter :: linear = 1, seasonal = 2, bilinear = 3, notime = 4 !< Flags to indicate whether the time
355  !! interpolation should be linear or some
356  !! other scheme for seasonal data.
357  !! NOTIME indicates
358  !! that data file has no time axis.
359 
360 ! Flags to indicate where climatology pressure levels are pressure or sigma levels
361 integer, parameter :: pressure = 1, sigma = 2 !< Flags to indicate where climatology pressure
362  !! levels are pressure or sigma levels
363 
364 ! Flags to indicate whether the climatology units are mixing ratio (kg/kg) or column integral (kg/m2).
365 ! Vertical interpolation scheme requires mixing ratio at this time.
366 integer, parameter :: no_conv = 1, kg_m2 = 2 !< Flags to indicate whether the climatology units
367  !! are mixing ratio (kg/kg) or column integral (kg/m2).
368  !! Vertical interpolation scheme requires mixing ratio at
369  !! this time.
370 
371 ! Flags to indicate what to do when the model surface pressure exceeds the climatology surface pressure level.
372 integer, parameter, public :: constant = 1, zero = 2 !< Flags to indicate what to do when the model surface
373  !! pressure exceeds the climatology surface pressure level.
374 
375 ! Flags to indicate the type of vertical interpolation
376 integer, parameter, public :: interp_weighted_p = 10, interp_linear_p = 20, interp_log_p = 30 !< Flags to indicate
377  !! the type of vertical interpolation
378 !--lwh
379 
380 real(r8_kind), parameter :: tpi = (2.0_r8_kind*pi) ! 4.*acos(0.)
381 real(r8_kind), parameter :: dtr = tpi/360._r8_kind
382 
383 
384 
385 integer :: num_clim_diag = 0 !< No description
386 character(len=64) :: climo_diag_name(max_diag_fields) !< No description
387 integer :: climo_diag_id(max_diag_fields), hinterp_id(max_diag_fields) !< No description
388 real(r8_kind) :: missing_value = -1.e10_r8_kind !< No description
389 ! sjs integer :: itaum, itaup
390 
391 #ifdef ENABLE_QUAD_PRECISION
392 ! Higher precision (kind=16) for grid geometrical factors:
393  integer, parameter:: f_p = r16_kind !< Higher precision (kind=16) for grid geometrical factors
394 #else
395 ! 64-bit precision (kind=8)
396  integer, parameter:: f_p = r8_kind !< 64-bit precision (kind=8)
397 #endif
398 
399 logical :: read_all_on_init = .false. !< No description
400 integer :: verbose = 0 !< No description
401 logical :: conservative_interp = .true. !< No description
402 logical :: retain_cm3_bug = .false. !< No description
403 logical :: use_mpp_io = .false. !< Set to true to use mpp_io, otherwise fms2io is used
404 
405 namelist /interpolator_nml/ &
407 
408 contains
409 
410 !> @brief Assignment overload routine for interpolate_type, to be used
411 !! through the assignment interface
412 subroutine interpolate_type_eq (Out, In)
413 
414 type(interpolate_type), intent(in) :: In
415 type(interpolate_type), intent(inout) :: Out
416 
417  if(in%r4_type%is_allocated) then
418  if (allocated(in%r4_type%lat)) out%r4_type%lat = in%r4_type%lat
419  if (allocated(in%r4_type%lon)) out%r4_type%lon = in%r4_type%lon
420  if (allocated(in%r4_type%latb)) out%r4_type%latb = in%r4_type%latb
421  if (allocated(in%r4_type%lonb)) out%r4_type%lonb = in%r4_type%lonb
422  if (allocated(in%r4_type%levs)) out%r4_type%levs = in%r4_type%levs
423  if (allocated(in%r4_type%halflevs)) out%r4_type%halflevs = in%r4_type%halflevs
424  else if(in%r8_type%is_allocated) then
425  if (allocated(in%r8_type%lat)) out%r8_type%lat = in%r8_type%lat
426  if (allocated(in%r8_type%lon)) out%r8_type%lon = in%r8_type%lon
427  if (allocated(in%r8_type%latb)) out%r8_type%latb = in%r8_type%latb
428  if (allocated(in%r8_type%lonb)) out%r8_type%lonb = in%r8_type%lonb
429  if (allocated(in%r8_type%levs)) out%r8_type%levs = in%r8_type%levs
430  if (allocated(in%r8_type%halflevs)) out%r8_type%halflevs = in%r8_type%halflevs
431  end if
432 
433  out%interph = in%interph
434  if (allocated(in%time_slice)) out%time_slice = in%time_slice
435  out%file_name = in%file_name
436  out%time_flag = in%time_flag
437  out%level_type = in%level_type
438  out%is = in%is
439  out%ie = in%ie
440  out%js = in%js
441  out%je = in%je
442  out%vertical_indices = in%vertical_indices
443  out%climatological_year = in%climatological_year
444  if (allocated(in%has_level )) out%has_level = in%has_level
445  if (allocated(in%field_name )) out%field_name = in%field_name
446  if (allocated(in%time_init )) out%time_init = in%time_init
447  if (allocated(in%mr )) out%mr = in%mr
448  if (allocated(in%out_of_bounds)) out%out_of_bounds = in%out_of_bounds
449  if (allocated(in%vert_interp )) out%vert_interp = in%vert_interp
450  if(in%r4_type%is_allocated) then
451  if (allocated(in%r4_type%data5d )) out%r4_type%data5d = in%r4_type%data5d
452  if (allocated(in%r4_type%pmon_pyear )) out%r4_type%pmon_pyear = in%r4_type%pmon_pyear
453  if (allocated(in%r4_type%pmon_nyear )) out%r4_type%pmon_nyear = in%r4_type%pmon_nyear
454  if (allocated(in%r4_type%nmon_nyear )) out%r4_type%nmon_nyear = in%r4_type%nmon_nyear
455  if (allocated(in%r4_type%nmon_pyear )) out%r4_type%nmon_pyear = in%r4_type%nmon_pyear
456  else if(in%r8_type%is_allocated) then
457  if (allocated(in%r8_type%data5d )) out%r8_type%data5d = in%r8_type%data5d
458  if (allocated(in%r8_type%pmon_pyear )) out%r8_type%pmon_pyear = in%r8_type%pmon_pyear
459  if (allocated(in%r8_type%pmon_nyear )) out%r8_type%pmon_nyear = in%r8_type%pmon_nyear
460  if (allocated(in%r8_type%nmon_nyear )) out%r8_type%nmon_nyear = in%r8_type%nmon_nyear
461  if (allocated(in%r8_type%nmon_pyear )) out%r8_type%nmon_pyear = in%r8_type%nmon_pyear
462  end if
463  if (allocated(in%indexm )) out%indexm = in%indexm
464  if (allocated(in%indexp )) out%indexp = in%indexp
465  if (allocated(in%climatology )) out%climatology = in%climatology
466  if (allocated(in%clim_times )) out%clim_times = in%clim_times
467  out%separate_time_vary_calc = in%separate_time_vary_calc
468  if(in%r4_type%is_allocated) then
469  out%r4_type%tweight = in%r4_type%tweight
470  out%r4_type%tweight1 = in%r4_type%tweight1
471  out%r4_type%tweight2 = in%r4_type%tweight2
472  out%r4_type%tweight3 = in%r4_type%tweight3
473  else if(in%r8_type%is_allocated) then
474  out%r8_type%tweight = in%r8_type%tweight
475  out%r8_type%tweight1 = in%r8_type%tweight1
476  out%r8_type%tweight2 = in%r8_type%tweight2
477  out%r8_type%tweight3 = in%r8_type%tweight3
478  end if
479  out%itaum = in%itaum
480  out%itaup = in%itaup
481 
482  out%r4_type%is_allocated = in%r4_type%is_allocated
483  out%r8_type%is_allocated = out%r8_type%is_allocated
484 
485 end subroutine interpolate_type_eq
486 
487 
488 
489 
490 !#######################################################################
491 !
492 !---------------------------------------------------------------------
493 !> @brief check_climo_units checks the units that the climatology
494 !! data is using. This is needed to allow for conversion of
495 !! datasets to mixing ratios which is what the vertical
496 !! interpolation scheme requires. The default is to assume no
497 !! conversion is needed. If the units are those of a column
498 !! burden (kg/m2) then conversion to mixing ratio is flagged.
499 !!
500 !! @param [in] <units> The units which you will be checking
501 function check_climo_units(units)
502 ! Function to check the units that the climatology data is using.
503 ! This is needed to allow for conversion of datasets to mixing ratios which is what the
504 ! vertical interpolation scheme requires
505 ! The default is to assume no conversion is needed.
506 ! If the units are those of a column burden (kg/m2) then conversion to mixing ratio is flagged.
507 !
508 character(len=*), intent(in) :: units
509 
510 integer :: check_climo_units
511 
512 check_climo_units = no_conv
513 select case(chomp(units))
514  case('kg/m2', 'kg/m^2', 'kg/m**2', 'kg m^-2', 'kg m**-2')
516  case('molecules/cm2/s', 'molecule/cm2/s', 'molec/cm2/s')
518  case('kg/m2/s')
520 end select
521 
522 end function check_climo_units
523 !
524 !#######################################################################
525 !
526 !---------------------------------------------------------------------
527 !> @brief init_clim_diag is a routine to register diagnostic fields
528 !! for the climatology file. This routine calculates the domain
529 !! decompostion of the climatology fields for later export
530 !! through send_data. The ids created here are for column
531 !! burdens that will diagnose the vertical interpolation
532 !! routine.
533 !!
534 !! @param [inout] <clim_type> The interpolate type containing the
535 !! names of the fields in the climatology file
536 !! @param [in] <mod_axes> The axes of the model
537 !! @param [in] <init_time> The model initialization time
538 !!
539 !! @throw FATAL, "init_clim_diag : You must call interpolator_init before calling init_clim_diag"
540 !! @throw FATAL, "init_clim_diag : Trying to set up too many diagnostic fields for the climatology data"
541 subroutine init_clim_diag(clim_type, mod_axes, init_time)
542 !
543 ! Routine to register diagnostic fields for the climatology file.
544 ! This routine calculates the domain decompostion of the climatology fields
545 ! for later export through send_data.
546 ! The ids created here are for column burdens that will diagnose the vertical interpolation routine.
547 ! climo_diag_id : 'module_name = climo' is intended for use with the model vertical resolution.
548 ! hinterp_id : 'module_name = 'hinterp' is intended for use with the climatology vertical resolution.
549 
550 ! INTENT INOUT :
551 ! clim_type : The interpolate type containing the names of the fields in the climatology file.
552 !
553 ! INTENT IN :
554 ! mod_axes : The axes of the model.
555 ! init_time : The model initialization time.
556 !
557 type(interpolate_type), intent(inout) :: clim_type
558 integer , intent(in) :: mod_axes(:)
559 type(time_type) , intent(in) :: init_time
560 
561 integer :: axes(2),nxd,nyd,ndivs,i
562 type(domain2d) :: domain
563 integer :: domain_layout(2), iscomp, iecomp,jscomp,jecomp
564 
565 if(clim_type%r4_type%is_allocated) call init_clim_diag_r4(clim_type, mod_axes, init_time)
566 if(clim_type%r8_type%is_allocated) call init_clim_diag_r8(clim_type, mod_axes, init_time)
567 
568 end subroutine init_clim_diag
569 
570 
571 
572 !
573 !---------------------------------------------------------------------
574 !> @brief obtain_interpolator_time_slices makes sure that the
575 !! appropriate time slices are available for interpolation on
576 !! this time step.
577 !!
578 !! @param [inout] <clim_type> The interpolate type previously defined
579 !! by a call to interpolator_init
580 !! @param [in] <Time> The model time that you wish to interpolate to
581 !!
582 !! @throw FATAL "interpolator_timeslice 1: file="
583 !! @throw FATAL "interpolator_timeslice 2: file="
584 !! @throw FATAL "interpolator_timeslice 3: file="
585 !! @throw FATAL "interpolator_timeslice 4: file="
586 !! @throw FATAL "interpolator_timeslice 5: file="
587 !! @throw FATAL "interpolator_timeslice : No data from the previous
588 !! climatology time but we have the next time. How did
589 !! this happen?"
590 subroutine obtain_interpolator_time_slices (clim_type, Time)
591 
592 ! Makes sure that appropriate time slices are available for interpolation
593 ! on this time step
594 !
595 ! INTENT INOUT
596 ! clim_type : The interpolate type previously defined by a call to interpolator_init
597 !
598 ! INTENT IN
599 ! Time : The model time that you wish to interpolate to.
600 
601 type(interpolate_type), intent(inout) :: clim_type
602 type(time_type) , intent(in) :: time
603 
604 integer :: taum, taup
605 integer :: modyear, modmonth, modday, modhour, modminute, modsecond
606 integer :: climyear, climmonth, climday, climhour, climminute, climsecond
607 integer :: year1, month1, day, hour, minute, second
608 integer :: climatology, m
609 type(time_type) :: t_prev, t_next
610 type(time_type), dimension(2) :: month
611 integer :: indexm, indexp, yearm, yearp
612 integer :: i, n
613 character(len=256) :: err_msg
614 
615 if(clim_type%r4_type%is_allocated) call obtain_interpolator_time_slices_r4(clim_type, time)
616 if(clim_type%r8_type%is_allocated) call obtain_interpolator_time_slices_r8(clim_type, time)
617 
618 
619 end subroutine obtain_interpolator_time_slices
620 
621 
622 !#####################################################################
623 !
624 !---------------------------------------------------------------------
625 !> @brief unset_interpolator_time_flag sets a flag in clim_type to
626 !! false.
627 !!
628 !! @param [inout] <clim_type> The interpolate type containing the names of the fields in the climatology file
629 subroutine unset_interpolator_time_flag (clim_type)
630 
631 type(interpolate_type), intent(inout) :: clim_type
632 
633 
634  clim_type%separate_time_vary_calc = .false.
635 
636 
637 end subroutine unset_interpolator_time_flag
638 
639 
640 !#####################################################################
641 !
642 !---------------------------------------------------------------------
643 !> @brief interpolator_end receives interpolate data as input
644 !! and deallocates its memory.
645 !!
646 !! @param [inout] <clim_type> The interpolate type whose components will be deallocated
647 subroutine interpolator_end(clim_type)
648 ! Subroutine to deallocate the interpolate type clim_type.
649 !
650 ! INTENT INOUT
651 ! clim_type : allocate type whose components will be deallocated.
652 !
653 type(interpolate_type), intent(inout) :: clim_type
654 integer :: logunit
655 
656 logunit=stdlog()
657 if ( mpp_pe() == mpp_root_pe() ) then
658  write (logunit,'(/,(a))') 'Exiting interpolator, have a nice day ...'
659 end if
660 
661 if(clim_type%r4_type%is_allocated) then
662  if (allocated (clim_type%r4_type%lat )) deallocate(clim_type%r4_type%lat)
663  if (allocated (clim_type%r4_type%lon )) deallocate(clim_type%r4_type%lon)
664  if (allocated (clim_type%r4_type%latb )) deallocate(clim_type%r4_type%latb)
665  if (allocated (clim_type%r4_type%lonb )) deallocate(clim_type%r4_type%lonb)
666  if (allocated (clim_type%r4_type%levs )) deallocate(clim_type%r4_type%levs)
667  if (allocated (clim_type%r4_type%halflevs)) deallocate(clim_type%r4_type%halflevs)
668  if (allocated (clim_type%r4_type%data5d )) deallocate(clim_type%r4_type%data5d)
669 else if(clim_type%r8_type%is_allocated) then
670  if (allocated (clim_type%r8_type%lat )) deallocate(clim_type%r8_type%lat)
671  if (allocated (clim_type%r8_type%lon )) deallocate(clim_type%r8_type%lon)
672  if (allocated (clim_type%r8_type%latb )) deallocate(clim_type%r8_type%latb)
673  if (allocated (clim_type%r8_type%lonb )) deallocate(clim_type%r8_type%lonb)
674  if (allocated (clim_type%r8_type%levs )) deallocate(clim_type%r8_type%levs)
675  if (allocated (clim_type%r8_type%halflevs)) deallocate(clim_type%r8_type%halflevs)
676  if (allocated (clim_type%r8_type%data5d)) deallocate(clim_type%r8_type%data5d)
677 end if
678 
679 if (allocated (clim_type%time_slice)) deallocate(clim_type%time_slice)
680 if (allocated (clim_type%field_name)) deallocate(clim_type%field_name)
681 if (allocated (clim_type%time_init )) deallocate(clim_type%time_init)
682 if (allocated (clim_type%has_level)) deallocate(clim_type%has_level)
683 if (allocated (clim_type%mr )) deallocate(clim_type%mr)
684 if (allocated (clim_type%out_of_bounds )) deallocate(clim_type%out_of_bounds)
685 if (allocated (clim_type%vert_interp )) deallocate(clim_type%vert_interp)
686 if (allocated(clim_type%indexm)) deallocate(clim_type%indexm)
687 if (allocated(clim_type%indexp)) deallocate(clim_type%indexp)
688 if (allocated(clim_type%clim_times)) deallocate(clim_type%clim_times)
689 if (allocated(clim_type%climatology)) deallocate(clim_type%climatology)
690 
691 call horiz_interp_del(clim_type%interph)
692 
693 if(clim_type%r4_type%is_allocated) then
694  if (allocated(clim_type%r4_type%pmon_pyear)) then
695  deallocate(clim_type%r4_type%pmon_pyear)
696  deallocate(clim_type%r4_type%pmon_nyear)
697  deallocate(clim_type%r4_type%nmon_nyear)
698  deallocate(clim_type%r4_type%nmon_pyear)
699  end if
700 else if(clim_type%r8_type%is_allocated) then
701  if (allocated(clim_type%r8_type%pmon_pyear)) then
702  deallocate(clim_type%r8_type%pmon_pyear)
703  deallocate(clim_type%r8_type%pmon_nyear)
704  deallocate(clim_type%r8_type%nmon_nyear)
705  deallocate(clim_type%r8_type%nmon_pyear)
706  end if
707 endif
708 
709 clim_type%r4_type%is_allocated=.false.
710 clim_type%r8_type%is_allocated=.false.
711 
712 !! RSH mod
713 if( .not.(clim_type%TIME_FLAG .eq. linear .and. read_all_on_init) &
714  .and. (clim_type%TIME_FLAG.ne.notime) ) then
715 ! read_all_on_init)) .or. clim_type%TIME_FLAG .eq. BILINEAR ) then
716  call close_file(clim_type%fileobj)
717 endif
718 
719 
720 module_is_initialized = .false.
721 
722 end subroutine interpolator_end
723 !
724 !#######################################################################
725 !
726 !++lwh
727 !
728 !---------------------------------------------------------------------
729 !> @brief query_interpolator receives an interpolate type as input
730 !! and returns the number of fields and field names.
731 !!
732 !! @param [in] <clim_type> The interpolate type which contains the data
733 !! @param [out] <nfields> OPTIONAL: No description
734 !! @param [out] <field_names> OPTIONAL: No description
735 subroutine query_interpolator( clim_type, nfields, field_names )
736 !
737 ! Query an interpolate_type variable to find the number of fields and field names.
738 !
739 type(interpolate_type), intent(in) :: clim_type
740 integer, intent(out), optional :: nfields
741 character(len=*), dimension(:), intent(out), optional :: field_names
742 
743 if( present( nfields ) ) nfields = SIZE( clim_type%field_name(:) )
744 if( present( field_names ) ) field_names = clim_type%field_name
745 
746 end subroutine query_interpolator
747 !--lwh
748 !
749 !#######################################################################
750 !
751 !---------------------------------------------------------------------
752 !> @brief chomp receives a string from NetCDF files and removes
753 !! CHAR(0) from the end of this string.
754 !!
755 !! @param [in] <string> The string from the NetCDF file
756 function chomp(string)
757 !
758 ! A function to remove CHAR(0) from the end of strings read from NetCDF files.
759 !
760 character(len=*), intent(in) :: string
761 character(len=64) :: chomp
762 
763 integer :: len
764 
765 len = len_trim(string)
766 if (string(len:len) == char(0)) len = len -1
767 
768 chomp = string(:len)
769 
770 end function chomp
771 !
772 !########################################################################
773 
774 #include "interpolator_r4.fh"
775 #include "interpolator_r8.fh"
776 
777 end module interpolator_mod
778 
779 !> @}
780 ! close documentation grouping
781 !
782 !#######################################################################
subroutine, public diag_manager_init(diag_model_subset, time_init, err_msg)
Initialize Diagnostics Manager.
Register a diagnostic field for a given module.
Send data over to output fields.
Close a netcdf or domain file opened with open_file or open_virtual_file.
Definition: fms2_io.F90:165
Opens a given netcdf or domain file.
Definition: fms2_io.F90:121
integer function, public check_nml_error(IOSTAT, NML_NAME)
Checks the iostat argument that is returned after reading a namelist and determines if the error code...
Definition: fms.F90:523
subroutine, public write_version_number(version, tag, unit)
Prints to the log file (or a specified unit) the version id string and tag name.
Definition: fms.F90:701
subroutine, public fms_init(localcomm, alt_input_nml_path)
Initializes the FMS module and also calls the initialization routines for all modules in the MPP pack...
Definition: fms.F90:286
subroutine, public horiz_interp_del(Interp)
Deallocates memory used by "horiz_interp_type" variables. Must be called before reinitializing with h...
subroutine, public horiz_interp_init
Initialize module and writes version number to logfile.out.
Subroutine for performing the horizontal interpolation between two grids.
integer nlat
No description.
subroutine, public query_interpolator(clim_type, nfields, field_names)
query_interpolator receives an interpolate type as input and returns the number of fields and field n...
subroutine, public init_clim_diag(clim_type, mod_axes, init_time)
init_clim_diag is a routine to register diagnostic fields for the climatology file....
integer, parameter increasing_upward
Flags to indicate direction of vertical axis in data file.
subroutine, public interpolate_type_eq(Out, In)
Assignment overload routine for interpolate_type, to be used through the assignment interface.
character(len=64) units
No description.
integer ndim
No description.
integer, parameter f_p
64-bit precision (kind=8)
integer ntime
No description.
integer, parameter sigma
Flags to indicate where climatology pressure levels are pressure or sigma levels.
integer, parameter, public zero
Flags to indicate what to do when the model surface pressure exceeds the climatology surface pressure...
integer, dimension(max_diag_fields) hinterp_id
No description.
integer nlev
No description.
logical use_mpp_io
Set to true to use mpp_io, otherwise fms2io is used.
integer, parameter max_diag_fields
No description.
character(len=64), dimension(max_diag_fields) climo_diag_name
No description.
integer nlon
No description.
subroutine, public unset_interpolator_time_flag(clim_type)
unset_interpolator_time_flag sets a flag in clim_type to false.
integer, parameter kg_m2
Flags to indicate whether the climatology units are mixing ratio (kg/kg) or column integral (kg/m2)....
subroutine, public interpolator_end(clim_type)
interpolator_end receives interpolate data as input and deallocates its memory.
integer sense
No description.
integer nvar
No description.
integer function check_climo_units(units)
check_climo_units checks the units that the climatology data is using. This is needed to allow for co...
logical retain_cm3_bug
No description.
integer verbose
No description.
character(len=64) function chomp(string)
chomp receives a string from NetCDF files and removes CHAR(0) from the end of this string.
logical read_all_on_init
No description.
integer nlatb
No description.
integer, parameter, public interp_log_p
Flags to indicate the type of vertical interpolation.
integer nlonb
No description.
integer nlevh
No description.
integer num_fields
No description.
integer num_clim_diag
No description.
real(r8_kind) missing_value
No description.
integer, parameter notime
Flags to indicate whether the time interpolation should be linear or some other scheme for seasonal d...
logical conservative_interp
No description.
subroutine, public obtain_interpolator_time_slices(clim_type, Time)
obtain_interpolator_time_slices makes sure that the appropriate time slices are available for interpo...
Private interface for weighted scalar interpolation.
Interpolates a field to a model grid.
subroutine mpp_domains_init(flags)
Initialize domain decomp package.
Set up a domain decomposition.
Retrieve layout associated with a domain decomposition. Given a global 2D domain and the number of di...
These routines retrieve the axis specifications associated with the compute domains....
Fill in a global array from domain-decomposed arrays.
Performs halo updates for a given domain.
The domain2D type contains all the necessary information to define the global, compute and data domai...
subroutine mpp_init(flags, localcomm, test_level, alt_input_nml_path)
Initialize the mpp_mod module. Must be called before any usage.
subroutine mpp_exit()
Finalizes process termination. To be called at the end of a run. Certain mpi implementations(openmpi)...
integer function stdlog()
This function returns the current standard fortran unit numbers for log messages. Log messages,...
Definition: mpp_util.inc:58
integer function mpp_npes()
Returns processor count for current pelist.
Definition: mpp_util.inc:420
integer function mpp_pe()
Returns processor ID.
Definition: mpp_util.inc:406
Error handler.
Definition: mpp.F90:381
Returns a weight and dates or indices for interpolating between two dates. The interface fraction_of_...
logical function, public leap_year(Time, err_msg)
Returns true if the year corresponding to the input time is a leap year (for default calendar)....
type(time_type) function, public decrement_time(Time, seconds, days, ticks, err_msg, allow_neg_inc)
Decrements a time by seconds and days.
subroutine, public get_date(time, year, month, day, hour, minute, second, tick, err_msg)
Gets the date for different calendar types. Given a time_interval, returns the corresponding date und...
integer function, public days_in_year(Time)
Returns the number of days in the calendar year corresponding to the date represented by time for the...
subroutine, public print_date(Time, str, unit)
Prints the time to standard output (or optional unit) as a date.
real(kind=r8_kind) function, public time_type_to_real(time)
Converts time to seconds and returns it as a real number.
integer function, public get_calendar_type()
Returns default calendar type for mapping from time to date.
Given an input date in year, month, days, etc., creates a time_type that represents this time interva...
Given some number of seconds and days, returns the corresponding time_type.
Type to represent amounts of time. Implemented as seconds and days to allow for larger intervals.
Allocates space and initializes a derived-type variable that contains pre-computed interpolation indi...
Redundant climatology data between fields.