FMS 2025.01-dev
Flexible Modeling System
Loading...
Searching...
No Matches
time_interp_external2.F90
1!***********************************************************************
2!* GNU Lesser General Public License
3!*
4!* This file is part of the GFDL Flexible Modeling System (FMS).
5!*
6!* FMS is free software: you can redistribute it and/or modify it under
7!* the terms of the GNU Lesser General Public License as published by
8!* the Free Software Foundation, either version 3 of the License, or (at
9!* your option) any later version.
10!*
11!* FMS is distributed in the hope that it will be useful, but WITHOUT
12!* ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
13!* FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
14!* for more details.
15!*
16!* You should have received a copy of the GNU Lesser General Public
17!* License along with FMS. If not, see <http://www.gnu.org/licenses/>.
18!***********************************************************************
19!> @defgroup time_interp_external2_mod time_interp_external2_mod
20!> @ingroup time_interp
21!> @brief Perform I/O and time interpolation of external fields (contained in a file), using
22!! fms2_io.
23!!
24!> @author M.J. Harrison
25!!
26!> Perform I/O and time interpolation for external fields.
27!! Uses udunits library to calculate calendar dates and
28!! convert units. Allows for reading data decomposed across
29!! model horizontal grid using optional domain2d argument
30!!
31!! data are defined over data domain for domain2d data
32!! (halo values are NOT updated by this module)
33
34!> @addtogroup time_interp_external2_mod
35!> @{
36module time_interp_external2_mod
37
38!<NAMELIST NAME="time_interp_external_nml">
39! <DATA NAME="num_io_buffers" TYPE="integer">
40! size of record dimension for internal buffer. This is useful for tuning i/o performance
41! particularly for large datasets (e.g. daily flux fields)
42! </DATA>
43!</NAMELIST>
44
45 use platform_mod, only : r8_kind, r4_kind
46 use fms_mod, only : write_version_number
47 use mpp_mod, only : mpp_error,fatal,warning,mpp_pe, stdout, stdlog, note
48 use mpp_mod, only : input_nml_file, mpp_npes, mpp_root_pe, mpp_broadcast, mpp_get_current_pelist
49 use time_manager_mod, only : time_type, get_date, set_date, operator ( >= ) , operator ( + ) , days_in_month, &
50 operator( - ), operator ( / ), operator ( // ) , days_in_year, increment_time, &
51 set_time, get_time, operator( > ), get_calendar_type, no_calendar
52 use get_cal_time_mod, only : get_cal_time
53 use mpp_domains_mod, only : domain2d, mpp_get_compute_domain, mpp_get_data_domain, &
54 mpp_get_global_domain, null_domain2d
55 use time_interp_mod, only : time_interp, time_interp_init
56 use axis_utils2_mod, only : get_axis_cart, get_axis_modulo, get_axis_modulo_times
57 use fms_mod, only : lowercase, check_nml_error
58 use platform_mod, only: r8_kind, fms_path_len, fms_file_len
59 use horiz_interp_mod, only : horiz_interp, horiz_interp_type
60 use fms2_io_mod, only : valid_t, fmsnetcdfdomainfile_t, open_file, get_unlimited_dimension_name, &
61 variable_att_exists, fmsnetcdffile_t, &
62 variable_exists, get_valid, get_variable_num_dimensions, read_data, &
63 is_valid, close_file, get_dimension_size, get_variable_dimension_names, &
64 get_variable_size, get_time_calendar, get_variable_missing, get_variable_units
65
66 implicit none
67 private
68
69! Include variable "version" to be written to log file.
70#include<file_version.h>
71
72 integer, parameter, public :: NO_REGION=0, inside_region=1, outside_region=2
73 integer, parameter, private :: modulo_year= 0001
74 integer, parameter, private :: LINEAR_TIME_INTERP = 1 ! not used currently
75 integer, parameter, public :: SUCCESS = 0, err_field_not_found = 1
76 real(r8_kind), parameter, private :: DEFAULT_MISSING_VALUE = -1e20_r8_kind
77 integer, private :: max_fields = 100, max_files= 40
78 integer, private :: num_fields = 0, num_files=0
79 ! denotes time intervals in file (interpreted from metadata)
80 integer, private :: num_io_buffers = 2 ! set -1 to read all records from disk into memory
81 logical, private :: module_initialized = .false.
82 logical, private :: debug_this_module = .false.
83
86 public set_override_region, reset_src_data_region
87 public get_external_fileobj
89
90 private find_buf_index,&
91 set_time_modulo
92 !> @}
93
94 !> Represents external fields
95 !> @ingroup time_interp_external2_mod
96 type, private :: ext_fieldtype
97 type(fmsnetcdffile_t), pointer :: fileobj=>null() !< keep unit open when not reading all records
98 character(len=128) :: name, units
99 integer :: siz(4), ndim
100 character(len=32) :: axisname(4)
101 type(domain2d) :: domain
102 type(time_type), dimension(:), pointer :: time =>null() !< midpoint of time interval
103 type(time_type), dimension(:), pointer :: start_time =>null(), end_time =>null()
104 type(time_type), dimension(:), pointer :: period =>null()
105 logical :: modulo_time !< denote climatological time axis
106 real(r8_kind), dimension(:,:,:,:), pointer :: domain_data =>null() !< defined over data domain or global domain
107 logical, dimension(:,:,:,:), pointer :: mask =>null() !< defined over data domain or global domain
108 integer, dimension(:), pointer :: ibuf =>null() !< record numbers associated with buffers
109 real(r8_kind), dimension(:,:,:,:), pointer :: src_data =>null() !< input data buffer
110 type(valid_t) :: valid ! data validator
111 integer :: nbuf
112 logical :: domain_present
113 real(r8_kind) :: slope, intercept
114 integer :: isc,iec,jsc,jec
115 type(time_type) :: modulo_time_beg, modulo_time_end
116 logical :: have_modulo_times, correct_leap_year_inconsistency
117 integer :: region_type
118 integer :: is_region, ie_region, js_region, je_region
119 integer :: is_src, ie_src, js_src, je_src
120 integer :: tdim
121 integer :: numwindows
122 logical, dimension(:,:), pointer :: need_compute=>null()
123 real(r8_kind) :: missing ! missing value
124 end type ext_fieldtype
125
126 !> Holds filename and file object
127 !> @ingroup time_interp_external2_mod
128 type, private :: filetype
129 character(len=FMS_FILE_LEN) :: filename = ''
130 type(FmsNetcdfFile_t), pointer :: fileobj => null()
131 end type filetype
132
133 !> Provide data from external file interpolated to current model time.
134 !! Data may be local to current processor or global, depending on
135 !! "init_external_field" flags. Uses @ref fms2_io for I/O.
136 !!
137 !! @param index index of external field from previous call to init_external_field
138 !! @param time target time for data
139 !! @param [inout] data global or local data array
140 !! @param interp time_interp_external defined interpolation method (optional). Currently
141 !! this module only supports LINEAR_TIME_INTERP.
142 !! @param verbose verbose flag for debugging (optional).
143 !!
144 !> @ingroup time_interp_external2_mod
146 module procedure time_interp_external_0d_r4
147 module procedure time_interp_external_2d_r4
148 module procedure time_interp_external_3d_r4
149 module procedure time_interp_external_0d_r8
150 module procedure time_interp_external_2d_r8
151 module procedure time_interp_external_3d_r8
152 end interface
153
155 module procedure time_interp_external_bridge_0d_r4
156 module procedure time_interp_external_bridge_2d_r4
157 module procedure time_interp_external_bridge_3d_r4
158 module procedure time_interp_external_bridge_0d_r8
159 module procedure time_interp_external_bridge_2d_r8
160 module procedure time_interp_external_bridge_3d_r8
161 end interface
162
163 !> @addtogroup time_interp_external2_mod
164 !> @{
165
166 integer :: outunit
167
168 type(ext_fieldtype), save, private, pointer :: loaded_fields(:) => null()
169 type(filetype), save, private, pointer :: opened_files(:) => null()
170!Balaji: really should use field%missing
171 real(r8_kind), private, parameter :: time_interp_missing=-1e99_r8_kind
172 contains
173
174! <SUBROUTINE NAME="time_interp_external_init">
175!
176! <DESCRIPTION>
177! Initialize the time_interp_external module
178! </DESCRIPTION>
179!
180 !> @brief Initialize the @ref time_interp_external_mod module
182
183 integer :: io_status, logunit, ierr
184
185 namelist /time_interp_external_nml/ num_io_buffers, debug_this_module, &
186 max_fields, max_files
187
188 ! open and read namelist
189
190 if(module_initialized) return
191
192 logunit = stdlog()
193 outunit = stdout()
194 call write_version_number("TIME_INTERP_EXTERNAL_MOD", version)
195
196 read (input_nml_file, time_interp_external_nml, iostat=io_status)
197 ierr = check_nml_error(io_status, 'time_interp_external_nml')
198
199 write(logunit,time_interp_external_nml)
200 call realloc_fields(max_fields)
201 call realloc_files(max_files)
202
203 module_initialized = .true.
204
205 call time_interp_init()
206
207 return
208
209 end subroutine time_interp_external_init
210
211
212!<FUNCTION NAME="init_external_field" TYPE="integer">
213!
214!<DESCRIPTION>
215! initialize an external field. Buffer "num_io_buffers" (default=2) in memory to reduce memory allocations.
216! distributed reads are supported using the optional "domain" flag.
217! Units conversion via the optional "desired_units" flag using udunits_mod.
218!
219! Return integer id of field for future calls to time_interp_external.
220!
221! </DESCRIPTION>
222!
223!<IN NAME="file" TYPE="character(len=*)">
224! filename
225!</IN>
226!<IN NAME="fieldname" TYPE="character(len=*)">
227! fieldname (in file)
228!</IN>
229!<IN NAME="format" TYPE="integer">
230! mpp_io flag for format of file (optional). Currently only "MPP_NETCDF" supported
231!</IN>
232!<IN NAME="threading" TYPE="integer">
233! mpp_io flag for threading (optional). "MPP_SINGLE" means root pe reads global field and distributes to other PEs
234! "MPP_MULTI" means all PEs read data
235!</IN>
236!<IN NAME="domain" TYPE="mpp_domains_mod:domain2d">
237! domain flag (optional)
238!</IN>
239!<IN NAME="desired_units" TYPE="character(len=*)">
240! Target units for data (optional), e.g. convert from deg_K to deg_C.
241! Failure to convert using udunits will result in failure of this module.
242!</IN>
243!<IN NAME="verbose" TYPE="logical">
244! verbose flag for debugging (optional).
245!</IN>
246!<INOUT NAME="axis_centers" TYPE="axistype" DIM="(4)">
247! MPP_IO axistype array for grid centers ordered X-Y-Z-T (optional).
248!</INOUT>
249!<INOUT NAME="axis_sizes" TYPE="integer" DIM="(4)">
250! array of axis lengths ordered X-Y-Z-T (optional).
251!</INOUT>
252
253
254 !> Initialize an external field. Buffer "num_io_buffers" (default=2) in memory to reduce memory allocations.
255 !! distributed reads are supported using the optional "domain" flag.
256 !! Units conversion via the optional "desired_units" flag using udunits_mod.
257 !!
258 !> @return integer id of field for future calls to time_interp_external.
259 !> @param file filename
260 !> @param fieldname fieldname (in file)
261 !> @param format mpp_io flag for format of file(optional). Currently only "MPP_NETCDF" supported
262 !> @param threading mpp_io flag for threading (optional). "MPP_SINGLE" means root pe reads
263 !! global field and distributes to other PEs. "MPP_MULTI" means all PEs read data
264 !> @param domain domain flag (optional)
265 !> @param desired_units Target units for data (optional), e.g. convert from deg_K to deg_C.
266 !! Failure to convert using udunits will result in failure of this module.
267 !> @param verbose verbose flag for debugging (optional).
268 !> @param [out] axis_names List of axis names (optional).
269 !> @param [inout] axis_sizes array of axis lengths ordered X-Y-Z-T (optional).
270 function init_external_field(file,fieldname,domain,desired_units,&
271 verbose,axis_names, axis_sizes,override,correct_leap_year_inconsistency,&
272 permit_calendar_conversion,use_comp_domain,ierr, nwindows, ignore_axis_atts, ongrid )
273
274 character(len=*), intent(in) :: file,fieldname
275 logical, intent(in), optional :: verbose
276 character(len=*), intent(in), optional :: desired_units
277 type(domain2d), intent(in), optional :: domain
278 integer, intent(inout), optional :: axis_sizes(4)
279 character(len=*), intent(out), optional :: axis_names(4)
280 logical, intent(in), optional :: override, correct_leap_year_inconsistency,&
281 permit_calendar_conversion,use_comp_domain
282 integer, intent(out), optional :: ierr
283 integer, intent(in), optional :: nwindows
284 logical, optional :: ignore_axis_atts
285 logical, optional :: ongrid !< Optional flag indicating if the data is ongrid
286
287 logical :: ongrid_local !< Flag indicating if the data is ongrid
288
289 integer :: init_external_field
290
291 real(r8_kind) :: slope, intercept
292 integer :: ndim,ntime,i,j
293 integer :: iscomp,iecomp,jscomp,jecomp,isglobal,ieglobal,jsglobal,jeglobal
294 integer :: isdata,iedata,jsdata,jedata, dxsize, dysize,dxsize_max,dysize_max
295 logical :: verb, transpose_xy,use_comp_domain1
296 real(r8_kind), dimension(:), allocatable :: tstamp, tstart, tend, tavg
297 character(len=1) :: cart
298 character(len=1), dimension(4) :: cart_dir
299 character(len=128) :: units, fld_units
300 character(len=128) :: msg, calendar_type, timebeg, timeend
301 character(len=128) :: timename, timeunits
302 character(len=128), allocatable :: axisname(:)
303 integer, allocatable :: axislen(:)
304 integer :: siz(4), siz_in(4), gxsize, gysize,gxsize_max, gysize_max
305 type(time_type) :: tdiff
306 integer :: yr, mon, day, hr, minu, sec
307 integer :: len, nfile, nfields_orig, nbuf, nx,ny
308 integer :: numwindows
309 logical :: ignore_axatts
310 logical :: have_modulo_time
311 type(fmsnetcdffile_t), pointer :: fileobj=>null()
312 integer, dimension(:), allocatable :: pes !< List of ranks in the current pelist
313
314 if (.not. module_initialized) call mpp_error(fatal,'Must call time_interp_external_init first')
315 if(present(ierr)) ierr = success
316 ignore_axatts=.false.
317 cart_dir(1)='X';cart_dir(2)='Y';cart_dir(3)='Z';cart_dir(4)='T'
318 if(present(ignore_axis_atts)) ignore_axatts = ignore_axis_atts
319 use_comp_domain1 = .false.
320 if(PRESENT(use_comp_domain)) use_comp_domain1 = use_comp_domain
321 verb=.false.
322 if (PRESENT(verbose)) verb=verbose
323 if (debug_this_module) verb = .true.
324 numwindows = 1
325 if(present(nwindows)) numwindows = nwindows
326
327 units = 'same'
328 if (PRESENT(desired_units)) then
329 units = desired_units
330 call mpp_error(fatal,'==> Unit conversion via time_interp_external &
331 &has been temporarily deprecated. Previous versions of&
332 &this module used udunits_mod to perform unit conversion.&
333 & Udunits_mod is in the process of being replaced since &
334 &there were portability issues associated with this code.&
335 & Please remove the desired_units argument from calls to &
336 &this routine.')
337 endif
338 nfile = 0
339 do i=1,num_files
340 if(trim(opened_files(i)%filename) == trim(file)) then
341 nfile = i
342 exit ! file is already opened
343 endif
344 enddo
345 if(nfile == 0) then
346 num_files = num_files + 1
347 if(num_files > max_files) then ! not enough space in the file table, reallocate it
348 !--- z1l: For the case of multiple thread, realoc_files will cause memory leak.
349 !--- If multiple threads are working on file A. One of the thread finished first and
350 !--- begin to work on file B, the realloc_files will cause problem for
351 !--- other threads are working on the file A.
352 ! call realloc_files(2*size(opened_files))
353 call mpp_error(fatal, "time_interp_external: num_files is greater than max_files, "// &
354 "increase time_interp_external_nml max_files")
355 endif
356 opened_files(num_files)%filename = trim(file)
357 allocate(opened_files(num_files)%fileobj)
358 fileobj => opened_files(num_files)%fileobj
359
360 if(.not. open_file(opened_files(num_files)%fileobj, trim(file), 'read')) &
361 call mpp_error(fatal, 'time_interp_external_mod: Error in opening file '//trim(file))
362 else
363 fileobj => opened_files(nfile)%fileobj
364 endif
365
366 call get_unlimited_dimension_name(fileobj, timename)
367 call get_dimension_size(fileobj, timename, ntime)
368
369 if (ntime < 1) then
370 write(msg,'(a15,a,a58)') 'external field ',trim(fieldname),&
371 ' does not have an associated record dimension (REQUIRED) '
372 call mpp_error(fatal,trim(msg))
373 endif
374
375 !--- get time calendar_type
376 call get_time_calendar(fileobj, timename, calendar_type)
377
378 !--- get time units
379 call get_variable_units(fileobj, timename, timeunits)
380
381 !--- get timebeg and timeend
382 have_modulo_time = get_axis_modulo_times(fileobj, timename, timebeg, timeend)
383
384 allocate(pes(mpp_npes()))
385 call mpp_get_current_pelist(pes)
386 allocate(tstamp(ntime),tstart(ntime),tend(ntime),tavg(ntime))
387
388 !< Only root reads the unlimited dimension and broadcasts it to the other ranks
389 if (mpp_root_pe() .eq. mpp_pe()) call read_data(fileobj, timename, tstamp)
390 call mpp_broadcast(tstamp, size(tstamp), mpp_root_pe(), pelist=pes)
391 deallocate(pes)
392
393 transpose_xy = .false.
394 isdata=1; iedata=1; jsdata=1; jedata=1
395 gxsize=1; gysize=1
396 siz_in = 1
397
398 if (PRESENT(domain)) then
399 call mpp_get_compute_domain(domain,iscomp,iecomp,jscomp,jecomp)
400 nx = iecomp-iscomp+1; ny = jecomp-jscomp+1
401 call mpp_get_data_domain(domain,isdata,iedata,jsdata,jedata,dxsize,dxsize_max,dysize,dysize_max)
402 call mpp_get_global_domain(domain,isglobal,ieglobal,jsglobal,jeglobal,gxsize,gxsize_max,gysize,gysize_max)
403 ongrid_local = .false.
404 if (present(ongrid)) ongrid_local = ongrid
405 !> If this is an ongrid case, set is[e]js[e]data to be equal to the compute domain.
406 !! This is what it is used to allocate space for the data!
407 if (ongrid_local) then
408 isdata=iscomp
409 iedata=iecomp
410 jsdata=jscomp
411 jedata=jecomp
412 endif
413 elseif(use_comp_domain1) then
414 call mpp_error(fatal,"init_external_field:"//&
415 " use_comp_domain=true but domain is not present")
416 endif
417
419 nfields_orig = num_fields
420
421 if (.not. variable_exists(fileobj, fieldname) ) then
422 if (present(ierr)) then
423 ierr = err_field_not_found
424 return
425 else
426 call mpp_error(fatal,'external field "'//trim(fieldname)//'" not found in file "'//trim(file)//'"')
427 endif
428 endif
429
430 tavg = -1.0_r8_kind
431 tstart = tstamp
432 tend = tstamp
433 if(variable_att_exists(fileobj, fieldname, 'time_avg_info')) then
434 if(variable_exists(fileobj, 'average_T1')) call read_data(fileobj, 'average_T1', tstart)
435 if(variable_exists(fileobj, 'average_T2')) call read_data(fileobj, 'average_T2', tend)
436 if(variable_exists(fileobj, 'average_DT')) call read_data(fileobj, 'average_DT', tavg)
437 endif
438
439 num_fields = num_fields + 1
440 if(num_fields > max_fields) then
441 !--- z1l: For the case of multiple thread, realoc_fields will cause memory leak.
442 !--- If multiple threads are working on field A. One of the thread finished first and
443 !--- begin to work on field B, the realloc_files will cause problem for
444 !--- other threads are working on the field A.
445 !call realloc_fields(size(field)*2)
446 call mpp_error(fatal, "time_interp_external: num_fields is greater than max_fields, "// &
447 "increase time_interp_external_nml max_fields")
448 endif
449
450 !--- get field attribute
451 call get_variable_units(fileobj, fieldname, fld_units)
452
453 init_external_field = num_fields
454 loaded_fields(num_fields)%fileobj => fileobj
455 loaded_fields(num_fields)%name = trim(fieldname)
456 loaded_fields(num_fields)%units = trim(fld_units)
457 loaded_fields(num_fields)%isc = 1
458 loaded_fields(num_fields)%iec = 1
459 loaded_fields(num_fields)%jsc = 1
460 loaded_fields(num_fields)%jec = 1
461 loaded_fields(num_fields)%region_type = no_region
462 loaded_fields(num_fields)%is_region = 0
463 loaded_fields(num_fields)%ie_region = -1
464 loaded_fields(num_fields)%js_region = 0
465 loaded_fields(num_fields)%je_region = -1
466 if (PRESENT(domain)) then
467 loaded_fields(num_fields)%domain_present = .true.
468 loaded_fields(num_fields)%domain = domain
469 loaded_fields(num_fields)%isc=iscomp;loaded_fields(num_fields)%iec = iecomp
470 loaded_fields(num_fields)%jsc=jscomp;loaded_fields(num_fields)%jec = jecomp
471 else
472 loaded_fields(num_fields)%domain_present = .false.
473 endif
474
475 loaded_fields(num_fields)%valid = get_valid(fileobj, fieldname)
476 ndim = get_variable_num_dimensions(fileobj, fieldname)
477 if (ndim > 4) call mpp_error(fatal, &
478 'invalid array rank <=4d fields supported')
479
480 loaded_fields(num_fields)%siz = 1
481 loaded_fields(num_fields)%ndim = ndim
482 loaded_fields(num_fields)%tdim = 4
483 !--- get field missing value
484 loaded_fields(num_fields)%missing = get_variable_missing(fileobj, fieldname)
485
486 allocate(axisname(ndim), axislen(ndim))
487
488 call get_variable_dimension_names(fileobj, fieldname, axisname)
489 call get_variable_size(fileobj, fieldname, axislen)
490 do j=1,loaded_fields(num_fields)%ndim
491 call get_axis_cart(fileobj, axisname(j), cart)
492 len = axislen(j)
493 if (cart == 'N' .and. .not. ignore_axatts) then
494 write(msg,'(a,"/",a)') trim(file),trim(fieldname)
495 call mpp_error(fatal,'file/field '//trim(msg)// &
496 ' couldnt recognize axis atts in time_interp_external')
497 else if (cart == 'N' .and. ignore_axatts) then
498 cart = cart_dir(j)
499 endif
500 select case (cart)
501 case ('X')
502 if (j.eq.2) transpose_xy = .true.
503 if (.not.PRESENT(domain) .and. .not.PRESENT(override)) then
504 isdata=1;iedata=len
505 iscomp=1;iecomp=len
506 gxsize = len
507 dxsize = len
508 loaded_fields(num_fields)%isc=iscomp;loaded_fields(num_fields)%iec=iecomp
509 elseif (PRESENT(override)) then
510 gxsize = len
511 if (PRESENT(axis_sizes)) axis_sizes(1) = len
512 endif
513 loaded_fields(num_fields)%axisname(1) = axisname(j)
514 if(use_comp_domain1) then
515 loaded_fields(num_fields)%siz(1) = nx
516 else
517 loaded_fields(num_fields)%siz(1) = dxsize
518 endif
519 if (len /= gxsize) then
520 write(msg,'(a,"/",a)') trim(file),trim(fieldname)
521 call mpp_error(fatal,'time_interp_ext, file/field '//trim(msg)//' x dim doesnt match model')
522 endif
523 case ('Y')
524 loaded_fields(num_fields)%axisname(2) = axisname(j)
525 if (.not.PRESENT(domain) .and. .not.PRESENT(override)) then
526 jsdata=1;jedata=len
527 jscomp=1;jecomp=len
528 gysize = len
529 dysize = len
530 loaded_fields(num_fields)%jsc=jscomp;loaded_fields(num_fields)%jec=jecomp
531 elseif (PRESENT(override)) then
532 gysize = len
533 if (PRESENT(axis_sizes)) axis_sizes(2) = len
534 endif
535 if(use_comp_domain1) then
536 loaded_fields(num_fields)%siz(2) = ny
537 else
538 loaded_fields(num_fields)%siz(2) = dysize
539 endif
540 if (len /= gysize) then
541 write(msg,'(a,"/",a)') trim(file),trim(fieldname)
542 call mpp_error(fatal,'time_interp_ext, file/field '//trim(msg)//' y dim doesnt match model')
543 endif
544 case ('Z')
545 loaded_fields(num_fields)%axisname(3) = axisname(j)
546 loaded_fields(num_fields)%siz(3) = len
547 case ('T')
548 loaded_fields(num_fields)%axisname(4) = axisname(j)
549 loaded_fields(num_fields)%siz(4) = ntime
550 loaded_fields(num_fields)%tdim = j
551 end select
552 enddo
553 siz = loaded_fields(num_fields)%siz
554 if(PRESENT(axis_names)) axis_names = loaded_fields(num_fields)%axisname
555 if (PRESENT(axis_sizes) .and. .not.PRESENT(override)) then
556 axis_sizes = loaded_fields(num_fields)%siz
557 endif
558
559 if (verb) write(outunit,'(a,4i6)') 'field x,y,z,t local size= ',siz
560 if (verb) write(outunit,*) 'field contains data in units = ',trim(loaded_fields(num_fields)%units)
561 if (transpose_xy) call mpp_error(fatal,'axis ordering not supported')
562 if (num_io_buffers .le. 1) call mpp_error(fatal,'time_interp_ext:num_io_buffers should be at least 2')
563 nbuf = min(num_io_buffers,siz(4))
564
565 loaded_fields(num_fields)%numwindows = numwindows
566 allocate(loaded_fields(num_fields)%need_compute(nbuf, numwindows))
567 loaded_fields(num_fields)%need_compute = .true.
568
569 allocate(loaded_fields(num_fields)%domain_data(isdata:iedata,jsdata:jedata,siz(3),nbuf),&
570 loaded_fields(num_fields)%mask(isdata:iedata,jsdata:jedata,siz(3),nbuf) )
571 loaded_fields(num_fields)%mask = .false.
572 loaded_fields(num_fields)%domain_data = 0.0_r8_kind
573 slope=1.0_r8_kind;intercept=0.0_r8_kind
574! if (units /= 'same') call convert_units(trim(field(num_fields)%units),trim(units),slope,intercept)
575! if (verb.and.units /= 'same') then
576! write(outunit,*) 'attempting to convert data to units = ',trim(units)
577! write(outunit,'(a,f8.3,a,f8.3)') 'factor = ',slope,' offset= ',intercept
578! endif
579 loaded_fields(num_fields)%slope = slope
580 loaded_fields(num_fields)%intercept = intercept
581 allocate(loaded_fields(num_fields)%ibuf(nbuf))
582 loaded_fields(num_fields)%ibuf = -1
583 loaded_fields(num_fields)%nbuf = 0 ! initialize buffer number so that first reading fills data(:,:,:,1)
584 if(PRESENT(override)) then
585 loaded_fields(num_fields)%is_src = 1
586 loaded_fields(num_fields)%ie_src = gxsize
587 loaded_fields(num_fields)%js_src = 1
588 loaded_fields(num_fields)%je_src = gysize
589 allocate(loaded_fields(num_fields)%src_data(gxsize,gysize,siz(3),nbuf))
590 else
591 loaded_fields(num_fields)%is_src = isdata
592 loaded_fields(num_fields)%ie_src = iedata
593 loaded_fields(num_fields)%js_src = jsdata
594 loaded_fields(num_fields)%je_src = jedata
595 allocate(loaded_fields(num_fields)%src_data(isdata:iedata,jsdata:jedata,siz(3),nbuf))
596 endif
597
598 allocate(loaded_fields(num_fields)%time(ntime))
599 allocate(loaded_fields(num_fields)%period(ntime))
600 allocate(loaded_fields(num_fields)%start_time(ntime))
601 allocate(loaded_fields(num_fields)%end_time(ntime))
602
603 do j=1,ntime
604 loaded_fields(num_fields)%time(j) = get_cal_time(tstamp(j),trim(timeunits),trim(calendar_type), &
605 & permit_calendar_conversion)
606 loaded_fields(num_fields)%start_time(j) = get_cal_time(tstart(j),trim(timeunits),trim(calendar_type), &
607 & permit_calendar_conversion)
608 loaded_fields(num_fields)%end_time(j) = get_cal_time( tend(j),trim(timeunits),trim(calendar_type), &
609 & permit_calendar_conversion)
610 enddo
611
612 if (loaded_fields(num_fields)%modulo_time) then
613 call set_time_modulo(loaded_fields(num_fields)%Time)
614 call set_time_modulo(loaded_fields(num_fields)%start_time)
615 call set_time_modulo(loaded_fields(num_fields)%end_time)
616 endif
617
618 if(present(correct_leap_year_inconsistency)) then
619 loaded_fields(num_fields)%correct_leap_year_inconsistency = correct_leap_year_inconsistency
620 else
621 loaded_fields(num_fields)%correct_leap_year_inconsistency = .false.
622 endif
623
624 if(have_modulo_time) then
625 if(get_calendar_type() == no_calendar) then
626 loaded_fields(num_fields)%modulo_time_beg = set_time(timebeg)
627 loaded_fields(num_fields)%modulo_time_end = set_time(timeend)
628 else
629 loaded_fields(num_fields)%modulo_time_beg = set_date(timebeg)
630 loaded_fields(num_fields)%modulo_time_end = set_date(timeend)
631 endif
632 loaded_fields(num_fields)%have_modulo_times = .true.
633 else
634 loaded_fields(num_fields)%have_modulo_times = .false.
635 endif
636 if(ntime == 1) then
637 call mpp_error(note, 'time_interp_external_mod: file '//trim(file)//' has only one time level')
638 else
639 do j= 1, ntime
640 loaded_fields(num_fields)%period(j) = loaded_fields(num_fields)%end_time(j) &
641 - loaded_fields(num_fields)%start_time(j)
642 if (loaded_fields(num_fields)%period(j) > set_time(0,0)) then
643 call get_time(loaded_fields(num_fields)%period(j), sec, day)
644 sec = sec/2+mod(day,2)*43200
645 day = day/2
646 loaded_fields(num_fields)%time(j) = loaded_fields(num_fields)%start_time(j)+&
647 set_time(sec,day)
648 else
649 if (j > 1 .and. j < ntime) then
650 tdiff = loaded_fields(num_fields)%time(j+1) - loaded_fields(num_fields)%time(j-1)
651 call get_time(tdiff, sec, day)
652 sec = sec/2+mod(day,2)*43200
653 day = day/2
654 loaded_fields(num_fields)%period(j) = set_time(sec,day)
655 sec = sec/2+mod(day,2)*43200
656 day = day/2
657 loaded_fields(num_fields)%start_time(j) = loaded_fields(num_fields)%time(j) - set_time(sec,day)
658 loaded_fields(num_fields)%end_time(j) = loaded_fields(num_fields)%time(j) + set_time(sec,day)
659 elseif ( j == 1) then
660 tdiff = loaded_fields(num_fields)%time(2) - loaded_fields(num_fields)%time(1)
661 call get_time(tdiff, sec, day)
662 loaded_fields(num_fields)%period(j) = set_time(sec,day)
663 sec = sec/2+mod(day,2)*43200
664 day = day/2
665 loaded_fields(num_fields)%start_time(j) = loaded_fields(num_fields)%time(j) - set_time(sec,day)
666 loaded_fields(num_fields)%end_time(j) = loaded_fields(num_fields)%time(j) + set_time(sec,day)
667 else
668 tdiff = loaded_fields(num_fields)%time(ntime) - loaded_fields(num_fields)%time(ntime-1)
669 call get_time(tdiff, sec, day)
670 loaded_fields(num_fields)%period(j) = set_time(sec,day)
671 sec = sec/2+mod(day,2)*43200
672 day = day/2
673 loaded_fields(num_fields)%start_time(j) = loaded_fields(num_fields)%time(j) - set_time(sec,day)
674 loaded_fields(num_fields)%end_time(j) = loaded_fields(num_fields)%time(j) + set_time(sec,day)
675 endif
676 endif
677 enddo
678 endif
679
680 do j=1,ntime-1
681 if (loaded_fields(num_fields)%time(j) >= loaded_fields(num_fields)%time(j+1)) then
682 write(msg,'(A,i20)') "times not monotonically increasing. Filename: " &
683 //trim(file)//" field: "//trim(fieldname)//" timeslice: ", j
684 call mpp_error(fatal, trim(msg))
685 endif
686 enddo
687
688 loaded_fields(num_fields)%modulo_time = get_axis_modulo(fileobj, timename)
689
690 if (verb) then
691 if (loaded_fields(num_fields)%modulo_time) write(outunit,*) 'data are being treated as modulo in time'
692 do j= 1, ntime
693 write(outunit,*) 'time index, ', j
694 call get_date(loaded_fields(num_fields)%start_time(j),yr,mon,day,hr,minu,sec)
695 write(outunit,'(a,i4,a,i2,a,i2,1x,i2,a,i2,a,i2)') &
696 'start time: yyyy/mm/dd hh:mm:ss= ',yr,'/',mon,'/',day,hr,':',minu,':',sec
697 call get_date(loaded_fields(num_fields)%time(j),yr,mon,day,hr,minu,sec)
698 write(outunit,'(a,i4,a,i2,a,i2,1x,i2,a,i2,a,i2)') &
699 'mid time: yyyy/mm/dd hh:mm:ss= ',yr,'/',mon,'/',day,hr,':',minu,':',sec
700 call get_date(loaded_fields(num_fields)%end_time(j),yr,mon,day,hr,minu,sec)
701 write(outunit,'(a,i4,a,i2,a,i2,1x,i2,a,i2,a,i2)') &
702 'end time: yyyy/mm/dd hh:mm:ss= ',yr,'/',mon,'/',day,hr,':',minu,':',sec
703 enddo
704 end if
705
706 deallocate(axisname, axislen)
707 deallocate(tstamp, tstart, tend, tavg)
708
709 return
710
711 end function init_external_field
712
713
714 function get_external_fileobj(filename, fileobj)
715 character(len=*), intent(in) :: filename
716 type(fmsnetcdffile_t), intent(out) :: fileobj
717 logical :: get_external_fileobj
718 integer :: i
719
720 get_external_fileobj = .false.
721 do i=1,num_files
722 if(trim(opened_files(i)%filename) == trim(filename)) then
723 fileobj = opened_files(i)%fileobj
724 get_external_fileobj = .true.
725 exit ! file is already opened
726 endif
727 enddo
728
729 end function get_external_fileobj
730
731
732 subroutine set_time_modulo(Time)
733
734 type(time_type), intent(inout), dimension(:) :: time
735
736 integer :: ntime, n
737 integer :: yr, mon, dy, hr, minu, sec
738
739 ntime = size(time(:))
740
741 do n = 1, ntime
742 call get_date(time(n), yr, mon, dy, hr, minu, sec)
743 yr = modulo_year
744 time(n) = set_date(yr, mon, dy, hr, minu, sec)
745 enddo
746
747
748 end subroutine set_time_modulo
749
750! ============================================================================
751!> load specified record from file
752!> Always loads in r8, casts down for horiz_interp if interp argument is already allocated for r4's.
753subroutine load_record(field, rec, interp, is_in, ie_in, js_in, je_in, window_id_in)
754 type(ext_fieldtype), intent(inout) :: field
755 integer , intent(in) :: rec ! record number
756 type(horiz_interp_type), intent(in), optional :: interp
757 integer, intent(in), optional :: is_in, ie_in, js_in, je_in
758 integer, intent(in), optional :: window_id_in
759
760 ! ---- local vars
761 integer :: ib ! index in the array of input buffers
762 integer :: isw,iew,jsw,jew ! boundaries of the domain on each window
763 integer :: is_region, ie_region, js_region, je_region, i, j
764 integer :: start(4), nread(4)
765 logical :: need_compute
766 integer :: window_id
767 real(r8_kind) :: mask_in(size(field%src_data,1),size(field%src_data,2),size(field%src_data,3))
768 real(r8_kind), allocatable :: mask_out(:,:,:)
769 real(r4_kind), allocatable :: hi_tmp_data(:,:,:,:) !< used to hold a copy of field%domain_data if using r4_kind
770 real(r4_kind), allocatable :: hi_tmp_msk_out(:,:,:) !< used return the field mask if using r4_kind
771
772 window_id = 1
773 if( PRESENT(window_id_in) ) window_id = window_id_in
774 need_compute = .true.
775
776!$OMP CRITICAL
777 ib = find_buf_index(rec,field%ibuf)
778
779 if(ib>0) then
780 !--- do nothing
781 need_compute = .false.
782 else
783 ! calculate current buffer number in round-robin fasion
784 field%nbuf = field%nbuf + 1
785 if(field%nbuf > size(field%domain_data,4).or.field%nbuf <= 0) field%nbuf = 1
786 ib = field%nbuf
787 field%ibuf(ib) = rec
788 field%need_compute(ib,:) = .true.
789
790 if (debug_this_module) write(outunit,*) 'reading record without domain for field ',trim(field%name)
791 start = 1; nread = 1
792 start(1) = field%is_src; nread(1) = field%ie_src - field%is_src + 1
793 start(2) = field%js_src; nread(2) = field%je_src - field%js_src + 1
794 start(3) = 1; nread(3) = size(field%src_data,3)
795 start(field%tdim) = rec; nread(field%tdim) = 1
796 call read_data(field%fileobj,field%name,field%src_data(:,:,:,ib:ib),corner=start,edge_lengths=nread)
797 endif
798!$OMP END CRITICAL
799 isw=field%isc;iew=field%iec
800 jsw=field%jsc;jew=field%jec
801
802 if( field%numwindows > 1) then
803 if( .NOT. PRESENT(is_in) .OR. .NOT. PRESENT(ie_in) .OR. .NOT. PRESENT(js_in) .OR. .NOT. PRESENT(je_in) ) then
804 call mpp_error(fatal, &
805 & 'time_interp_external(load_record): is_in, ie_in, js_in, je_in must be present when numwindows>1')
806 endif
807 isw = isw + is_in - 1
808 iew = isw + ie_in - is_in
809 jsw = jsw + js_in - 1
810 jew = jsw + je_in - js_in
811 endif
812
813 ! interpolate to target grid
814
815 need_compute = field%need_compute(ib, window_id)
816 if(need_compute) then
817 if(PRESENT(interp)) then
818 is_region = field%is_region; ie_region = field%ie_region
819 js_region = field%js_region; je_region = field%je_region
820 mask_in = 0.0_r8_kind
821 where (is_valid(field%src_data(:,:,:,ib), field%valid)) mask_in = 1.0_r8_kind
822 if ( field%region_type .NE. no_region ) then
823 if( any(mask_in == 0.0_r8_kind) ) then
824 call mpp_error(fatal, "time_interp_external: mask_in should be all 1 when region_type is not NO_REGION")
825 endif
826 if( field%region_type == outside_region) then
827 do j = js_region, je_region
828 do i = is_region, ie_region
829 mask_in(i,j,:) = 0.0_r8_kind
830 enddo
831 enddo
832 else ! field%region_choice == INSIDE_REGION
833 do j = 1, size(mask_in,2)
834 do i = 1, size(mask_in,1)
835 if( j<js_region .OR. j>je_region .OR. i<is_region .OR. i>ie_region ) mask_in(i,j,:) = 0.0_r8_kind
836 enddo
837 enddo
838 endif
839 endif
840 !! added for mixed mode. Data is always read in as r8 (via ext_fieldtype). if existing horiz_interp_type was
841 !! initialized in r4, needs to cast down in order to match up with saved values in horiz_interp_type.
842 !! creates some temporary arrays since intent(out) vars can't get passed in directly
843 if (interp%horizInterpReals4_type%is_allocated) then
844 ! allocate (there may be a better way to do this, had issues with gnu)
845 allocate(hi_tmp_msk_out(isw:iew,jsw:jew, SIZE(field%src_data,3)))
846 allocate(hi_tmp_data(lbound(field%domain_data,1):ubound(field%domain_data,1), &
847 lbound(field%domain_data,2):ubound(field%domain_data,2), &
848 lbound(field%domain_data,3):ubound(field%domain_data,3), &
849 lbound(field%domain_data,4):ubound(field%domain_data,4)))
850 ! copy over to r4
851 hi_tmp_data = real(field%domain_data, r4_kind)
852 ! do interpolation
853 call horiz_interp(interp, real(field%src_data(:,:,:,ib),r4_kind), hi_tmp_data(isw:iew,jsw:jew,:,ib), &
854 mask_in=real(mask_in,r4_kind), mask_out=hi_tmp_msk_out)
855 ! assign any output
856 field%domain_data = real(hi_tmp_data, r8_kind)
857 field%mask(isw:iew,jsw:jew,:,ib) = hi_tmp_msk_out(isw:iew,jsw:jew,:) > 0.0_r4_kind
858
859 if(allocated(hi_tmp_data)) deallocate(hi_tmp_data)
860 if(allocated(hi_tmp_msk_out)) deallocate(hi_tmp_msk_out)
861 else
862 allocate(mask_out(isw:iew,jsw:jew, size(field%src_data,3)))
863 call horiz_interp(interp, field%src_data(:,:,:,ib),field%domain_data(isw:iew,jsw:jew,:,ib), &
864 mask_in=mask_in, &
865 mask_out=mask_out)
866 field%mask(isw:iew,jsw:jew,:,ib) = mask_out(isw:iew,jsw:jew,:) > 0.0_r8_kind
867 deallocate(mask_out)
868 endif
869
870 else
871 if ( field%region_type .NE. no_region ) then
872 call mpp_error(fatal, "time_interp_external: region_type should be NO_REGION when interp is not present")
873 endif
874 field%domain_data(isw:iew,jsw:jew,:,ib) = field%src_data(isw:iew,jsw:jew,:,ib)
875 field%mask(isw:iew,jsw:jew,:,ib) = is_valid(field%domain_data(isw:iew,jsw:jew,:,ib),field%valid)
876 endif
877 ! convert units
878 where(field%mask(isw:iew,jsw:jew,:,ib)) field%domain_data(isw:iew,jsw:jew,:,ib) = &
879 field%domain_data(isw:iew,jsw:jew,:,ib)*field%slope + field%intercept
880 field%need_compute(ib, window_id) = .false.
881
882 endif
883
884end subroutine load_record
885
886!> Given a initialized ext_fieldtype and record number, loads the given index into field%src_data
887subroutine load_record_0d(field, rec)
888 type(ext_fieldtype), intent(inout) :: field
889 integer , intent(in) :: rec ! record number
890 ! ---- local vars
891 integer :: ib ! index in the array of input buffers
892 integer :: start(4), nread(4)
893
894 ib = find_buf_index(rec,field%ibuf)
895
896 if(ib>0) then
897 return
898 else
899 ! calculate current buffer number in round-robin fasion
900 field%nbuf = field%nbuf + 1
901 if(field%nbuf > size(field%domain_data,4).or.field%nbuf <= 0) field%nbuf = 1
902 ib = field%nbuf
903 field%ibuf(ib) = rec
904
905 if (debug_this_module) write(outunit,*) 'reading record without domain for field ',trim(field%name)
906 start = 1; nread = 1
907 start(3) = 1; nread(3) = size(field%src_data,3)
908 start(field%tdim) = rec; nread(field%tdim) = 1
909 call read_data(field%fileobj,field%name,field%src_data(:,:,:,ib),corner=start,edge_lengths=nread)
910 if ( field%region_type .NE. no_region ) then
911 call mpp_error(fatal, "time_interp_external: region_type should be NO_REGION when field is scalar")
912 endif
913 field%domain_data(1,1,:,ib) = field%src_data(1,1,:,ib)
914 field%mask(1,1,:,ib) = is_valid(field%domain_data(1,1,:,ib),field%valid)
915 ! convert units
916 where(field%mask(1,1,:,ib)) field%domain_data(1,1,:,ib) = &
917 field%domain_data(1,1,:,ib)*field%slope + field%intercept
918 endif
919
920end subroutine load_record_0d
921
922!> Reallocates src_data for field from module level loaded_fields array
923subroutine reset_src_data_region(index, is, ie, js, je)
924 integer, intent(in) :: index
925 integer, intent(in) :: is, ie, js, je
926 integer :: nk, nbuf
927
928 if( is == loaded_fields(index)%is_src .AND. ie == loaded_fields(index)%ie_src .AND. &
929 js == loaded_fields(index)%js_src .AND. ie == loaded_fields(index)%je_src ) return
930
931 if( .NOT. ASSOCIATED(loaded_fields(index)%src_data) ) call mpp_error(fatal, &
932 "time_interp_external: field(index)%src_data is not associated")
933 nk = size(loaded_fields(index)%src_data,3)
934 nbuf = size(loaded_fields(index)%src_data,4)
935 deallocate(loaded_fields(index)%src_data)
936 allocate(loaded_fields(index)%src_data(is:ie,js:je,nk,nbuf))
937 loaded_fields(index)%is_src = is
938 loaded_fields(index)%ie_src = ie
939 loaded_fields(index)%js_src = js
940 loaded_fields(index)%je_src = je
941
942
943end subroutine reset_src_data_region
944
945subroutine set_override_region(index, region_type, is_region, ie_region, js_region, je_region)
946 integer, intent(in) :: index, region_type
947 integer, intent(in) :: is_region, ie_region, js_region, je_region
948
949 loaded_fields(index)%region_type = region_type
950 loaded_fields(index)%is_region = is_region
951 loaded_fields(index)%ie_region = ie_region
952 loaded_fields(index)%js_region = js_region
953 loaded_fields(index)%je_region = je_region
954
955 return
956
957end subroutine set_override_region
958
959!> reallocates array of fields, increasing its size
960subroutine realloc_files(n)
961 integer, intent(in) :: n ! new size
962
963 type(filetype), pointer :: ptr(:)
964 integer :: i
965
966 if (associated(opened_files)) then
967 if (n <= size(opened_files)) return ! do nothing, if requested size no more than current
968 endif
969
970 allocate(ptr(n))
971 do i = 1, size(ptr)
972 ptr(i)%filename = ''
973 if(Associated(ptr(i)%fileobj)) then
974 call close_file(ptr(i)%fileobj)
975 deallocate(ptr(i)%fileobj)
976 endif
977 enddo
978
979 if (associated(opened_files))then
980 ptr(1:size(opened_files)) = opened_files(:)
981 deallocate(opened_files)
982 endif
983 opened_files => ptr
984
985end subroutine realloc_files
986
987!> reallocates array of fields,increasing its size
988subroutine realloc_fields(n)
989 integer, intent(in) :: n ! new size
990
991 type(ext_fieldtype), pointer :: ptr(:)
992 integer :: i, ier
993
994 if (associated(loaded_fields)) then
995 if (n <= size(loaded_fields)) return ! do nothing if requested size no more then current
996 endif
997
998 allocate(ptr(n))
999 do i=1,size(ptr)
1000 ptr(i)%fileobj => null()
1001 ptr(i)%name=''
1002 ptr(i)%units=''
1003 ptr(i)%siz=-1
1004 ptr(i)%ndim=-1
1005 ptr(i)%domain = null_domain2d
1006 if (ASSOCIATED(ptr(i)%time)) DEALLOCATE(ptr(i)%time, stat=ier)
1007 if (ASSOCIATED(ptr(i)%start_time)) DEALLOCATE(ptr(i)%start_time, stat=ier)
1008 if (ASSOCIATED(ptr(i)%end_time)) DEALLOCATE(ptr(i)%end_time, stat=ier)
1009 if (ASSOCIATED(ptr(i)%period)) DEALLOCATE(ptr(i)%period, stat=ier)
1010 ptr(i)%modulo_time=.false.
1011 if (ASSOCIATED(ptr(i)%domain_data)) DEALLOCATE(ptr(i)%domain_data, stat=ier)
1012 if (ASSOCIATED(ptr(i)%ibuf)) DEALLOCATE(ptr(i)%ibuf, stat=ier)
1013 if (ASSOCIATED(ptr(i)%src_data)) DEALLOCATE(ptr(i)%src_data, stat=ier)
1014 ptr(i)%nbuf=-1
1015 ptr(i)%domain_present=.false.
1016 ptr(i)%slope=1.0_r8_kind
1017 ptr(i)%intercept=0.0_r8_kind
1018 ptr(i)%isc=-1;ptr(i)%iec=-1
1019 ptr(i)%jsc=-1;ptr(i)%jec=-1
1020 enddo
1021 if (associated(loaded_fields)) then
1022 ptr(1:size(loaded_fields)) = loaded_fields(:)
1023 deallocate(loaded_fields)
1024 endif
1025 loaded_fields=>ptr
1026
1027end subroutine realloc_fields
1028
1029!> simple linear search for given value in given list
1030!! TODO should use better search if this list is bigger
1031function find_buf_index(indx,buf)
1032 integer :: indx
1033 integer, dimension(:) :: buf
1034 integer :: find_buf_index
1035
1036 integer :: nbuf, i
1037
1038 nbuf = size(buf(:))
1039
1040 find_buf_index = -1
1041
1042 do i=1,nbuf
1043 if (buf(i) == indx) then
1044 find_buf_index = i
1045 exit
1046 endif
1047 enddo
1048
1049end function find_buf_index
1050
1051!> Returns size of field after call to init_external_field.
1052!! Ordering is X/Y/Z/T.
1053!! This call only makes sense for non-distributed reads.
1055
1056 integer :: index !< returned from previous call to init_external_field.
1057 integer :: get_external_field_size(4)
1058
1059 if (index .lt. 1 .or. index .gt. num_fields) &
1060 call mpp_error(fatal,'invalid index in call to get_external_field_size')
1061
1062
1063 get_external_field_size(1) = loaded_fields(index)%siz(1)
1064 get_external_field_size(2) = loaded_fields(index)%siz(2)
1065 get_external_field_size(3) = loaded_fields(index)%siz(3)
1066 get_external_field_size(4) = loaded_fields(index)%siz(4)
1067
1068end function get_external_field_size
1069
1070!> return missing value
1072
1073 integer :: index !< returned from previous call to init_external_field.
1074 real(r8_kind) :: get_external_field_missing
1075
1076 if (index .lt. 1 .or. index .gt. num_fields) &
1077 call mpp_error(fatal,'invalid index in call to get_external_field_size')
1078
1079
1080 get_external_field_missing = loaded_fields(index)%missing
1081
1082end function get_external_field_missing
1083
1084subroutine get_time_axis(index, time)
1085 integer , intent(in) :: index !< field id
1086 type(time_type), intent(out) :: time(:) !< array of time values to be filled
1087
1088 integer :: n !< size of the data to be assigned
1089
1090 if (index < 1.or.index > num_fields) &
1091 call mpp_error(fatal,'invalid index in call to get_time_axis')
1092
1093 n = min(size(time),size(loaded_fields(index)%time))
1094
1095 time(1:n) = loaded_fields(index)%time(1:n)
1096end subroutine
1097
1098
1099!> exit time_interp_external_mod. Close all open files and
1100!! release storage
1102
1103 integer :: i
1104 !
1105 ! release storage arrays
1106 !
1107 do i=1,num_fields
1108 deallocate(loaded_fields(i)%time,loaded_fields(i)%start_time,loaded_fields(i)%end_time,&
1109 loaded_fields(i)%period,loaded_fields(i)%domain_data,loaded_fields(i)%mask,loaded_fields(i)%ibuf)
1110 if (ASSOCIATED(loaded_fields(i)%src_data)) deallocate(loaded_fields(i)%src_data)
1111 loaded_fields(i)%domain = null_domain2d
1112 loaded_fields(i)%nbuf = 0
1113 loaded_fields(i)%slope = 0.0_r8_kind
1114 loaded_fields(i)%intercept = 0.0_r8_kind
1115 enddo
1116
1117 !-- close all the files opended
1118 do i = 1, num_files
1119 call close_file(opened_files(i)%fileobj)
1120 deallocate(opened_files(i)%fileobj)
1121 enddo
1122
1123 deallocate(loaded_fields)
1124 deallocate(opened_files)
1125
1126 num_fields = 0
1127
1128 module_initialized = .false.
1129
1130end subroutine time_interp_external_exit
1131
1132#include "time_interp_external2_r4.fh"
1133#include "time_interp_external2_r8.fh"
1134
1135#include "time_interp_external2_bridge_r4.fh"
1136#include "time_interp_external2_bridge_r8.fh"
1137
1138end module time_interp_external2_mod
1139!> @}
1140! close documentation grouping
subroutine, public get_axis_cart(fileobj, axisname, cart)
Returns X,Y,Z or T cartesian attribute.
logical function, public get_axis_modulo_times(fileobj, axisname, tbeg, tend)
logical function, public get_axis_modulo(fileobj, axisname)
Checks if 'modulo' variable exists for a given axis.
Close a netcdf or domain file opened with open_file or open_virtual_file.
Definition fms2_io.F90:166
Opens a given netcdf or domain file.
Definition fms2_io.F90:122
Read data from a defined field in a file.
Definition fms2_io.F90:292
Added for mixed precision support. Updates force time_manager math to be done with kind=8 reals _wrap...
Subroutine for performing the horizontal interpolation between two grids.
These routines retrieve the axis specifications associated with the compute domains....
These routines retrieve the axis specifications associated with the data domains. The domain is a der...
These routines retrieve the axis specifications associated with the global domains....
The domain2D type contains all the necessary information to define the global, compute and data domai...
Perform parallel broadcasts.
Definition mpp.F90:1091
Error handler.
Definition mpp.F90:382
subroutine, public time_interp_external_init()
Initialize the time_interp_external_mod module.
subroutine, public get_time_axis(index, time)
integer function, public init_external_field(file, fieldname, domain, desired_units, verbose, axis_names, axis_sizes, override, correct_leap_year_inconsistency, permit_calendar_conversion, use_comp_domain, ierr, nwindows, ignore_axis_atts, ongrid)
Initialize an external field. Buffer "num_io_buffers" (default=2) in memory to reduce memory allocati...
subroutine load_record_0d(field, rec)
Given a initialized ext_fieldtype and record number, loads the given index into fieldsrc_data.
subroutine, public time_interp_external_exit()
exit time_interp_external_mod. Close all open files and release storage
subroutine load_record(field, rec, interp, is_in, ie_in, js_in, je_in, window_id_in)
load specified record from file Always loads in r8, casts down for horiz_interp if interp argument is...
subroutine realloc_files(n)
reallocates array of fields, increasing its size
real(r8_kind) function, public get_external_field_missing(index)
return missing value
integer function, private find_buf_index(indx, buf)
simple linear search for given value in given list TODO should use better search if this list is bigg...
subroutine realloc_fields(n)
reallocates array of fields,increasing its size
integer function, dimension(4), public get_external_field_size(index)
Returns size of field after call to init_external_field. Ordering is X/Y/Z/T. This call only makes se...
subroutine, public reset_src_data_region(index, is, ie, js, je)
Reallocates src_data for field from module level loaded_fields array.
Provide data from external file interpolated to current model time. Data may be local to current proc...
Returns a weight and dates or indices for interpolating between two dates. The interface fraction_of_...
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...
integer function, public get_calendar_type()
Returns default calendar type for mapping from time to date.
type(time_type) function, public increment_time(time, seconds, days, ticks, err_msg, allow_neg_inc)
Increments a time by seconds and days.
integer function, public days_in_month(time, err_msg)
Given a time, computes the corresponding date given the selected date time mapping algorithm.
subroutine, public get_time(time, seconds, days, ticks, err_msg)
Returns days and seconds ( < 86400 ) corresponding to a time. err_msg should be checked for any error...
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...
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.