FMS  2025.03
Flexible Modeling System
data_override.inc
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 
20 ! This file contains the body of the data_override_r4 and data_override_r8
21 ! modules. These modules are not intended to be used directly - they should be
22 ! used through the data_override_mod API. See data_override.F90 for details.
23 
24 use platform_mod, only: r4_kind, r8_kind, fms_path_len, fms_file_len
25 use yaml_parser_mod
26 use constants_mod, only: deg_to_rad
27 use mpp_mod, only : mpp_error, fatal, warning, note, stdout, stdlog, mpp_max
28 use mpp_mod, only : input_nml_file
29 use horiz_interp_mod, only : horiz_interp_init, horiz_interp_new, horiz_interp_type, &
31 use time_interp_external2_mod, only: time_interp_external_init, &
36  set_override_region, &
38  no_region, inside_region, outside_region, &
39  get_external_fileobj
40 use fms_mod, only: write_version_number, lowercase, check_nml_error
41 use axis_utils2_mod, only : nearest_index, axis_edges
42 use mpp_domains_mod, only : domain2d, mpp_get_compute_domain, null_domain2d,operator(.NE.),operator(.EQ.)
43 use mpp_domains_mod, only : mpp_get_global_domain, mpp_get_data_domain
44 use mpp_domains_mod, only : domainug, mpp_pass_sg_to_ug, mpp_get_ug_sg_domain, null_domainug
45 use time_manager_mod, only: time_type, OPERATOR(>), OPERATOR(<)
46 use fms2_io_mod, only : fmsnetcdffile_t, open_file, close_file, &
47  read_data, fms2_io_init, variable_exists, &
48  get_mosaic_tile_file, file_exists, get_instance_filename
49 use get_grid_version_mod, only: get_grid_version_1, get_grid_version_2
50 use fms_string_utils_mod, only: string
51 
52 implicit none
53 private
54 
55 ! Include variable "version" to be written to log file.
56 #include<file_version.h>
57 
58 !> Private type for holding field and grid information from a data table
59 !> @ingroup data_override_mod
60 type data_type
61  character(len=3) :: gridname
62  character(len=128) :: fieldname_code !< fieldname used in user's code (model)
63  character(len=128) :: fieldname_file !< fieldname used in the netcdf data file
64  character(len=FMS_PATH_LEN) :: file_name !< name of netCDF data file
65  character(len=128) :: interpol_method !< interpolation method (default "bilinear")
66  logical :: ext_weights
67  character(len=128) :: ext_weights_file_name
68  character(len=128) :: ext_weights_source
69  real(FMS_DATA_OVERRIDE_KIND_) :: factor !< For unit conversion, default=1, see OVERVIEW above
70  real(FMS_DATA_OVERRIDE_KIND_) :: lon_start, lon_end, lat_start, lat_end
71  integer :: region_type
72  logical :: multifile = .false.
73  character(len=FMS_PATH_LEN) :: prev_file_name !< name of netCDF data file for previous segment
74  character(len=FMS_PATH_LEN) :: next_file_name !< name of netCDF data file for next segment
75  type(time_type), dimension(:), allocatable :: time_records
76  type(time_type), dimension(:), allocatable :: time_prev_records
77  type(time_type), dimension(:), allocatable :: time_next_records
78 end type data_type
79 
80 !> Private type for holding various data fields for performing data overrides
81 !> @ingroup data_override_mod
82 type override_type
83  character(len=3) :: gridname
84  character(len=128) :: fieldname
85  integer :: t_index !< index for time interp
86  integer :: pt_index !< previous index for time interp
87  integer :: nt_index !< next index for time interp
88  type(horiz_interp_type), allocatable :: horz_interp(:) !< index for horizontal spatial interp
89  integer :: dims(4) !< dimensions(x,y,z,t) of the field in filename
90  integer :: comp_domain(4) !< istart,iend,jstart,jend for compute domain
91  integer :: numthreads
92  real(FMS_DATA_OVERRIDE_KIND_), allocatable :: lon_in(:)
93  real(FMS_DATA_OVERRIDE_KIND_), allocatable :: lat_in(:)
94  logical, allocatable :: need_compute(:)
95  integer :: numwindows
96  integer :: window_size(2)
97  integer :: is_src, ie_src, js_src, je_src
98 end type override_type
99 
100 !> Private type for holding horiz_interp_type for a weight file
101 !! This is needed so that if variables use the same weight file,
102 !! then we won't have to read the weight file again
103 !> @ingroup data_override_mod
104 type fmsexternalweights_type
105  character(len=:), allocatable :: weight_filename !< Name of the weight file
106  type(horiz_interp_type) :: horiz_interp !< Horiz interp type read in from the weight file
107 end type fmsexternalweights_type
108 
109 integer, parameter :: lkind = fms_data_override_kind_
110 integer, parameter :: max_table=100, max_array=100
111 
112 integer :: table_size !< actual size of data table
113 integer :: nweight_files !< Number of weight files that have been used
114 type(fmsExternalWeights_type), allocatable, target :: external_weights(:) !< External weights types
115 logical :: module_is_initialized = .false.
116 
117 type(domain2D) :: ocn_domain,atm_domain,lnd_domain, ice_domain
118 type(domainUG) :: lnd_domainUG
119 
120 real(FMS_DATA_OVERRIDE_KIND_), dimension(:,:), target, allocatable :: lon_local_ocn, lat_local_ocn
121 real(FMS_DATA_OVERRIDE_KIND_), dimension(:,:), target, allocatable :: lon_local_atm, lat_local_atm
122 real(FMS_DATA_OVERRIDE_KIND_), dimension(:,:), target, allocatable :: lon_local_ice, lat_local_ice
123 real(FMS_DATA_OVERRIDE_KIND_), dimension(:,:), target, allocatable :: lon_local_lnd, lat_local_lnd
124 real(FMS_DATA_OVERRIDE_KIND_) :: min_glo_lon_ocn, max_glo_lon_ocn
125 real(FMS_DATA_OVERRIDE_KIND_) :: min_glo_lon_atm, max_glo_lon_atm
126 real(FMS_DATA_OVERRIDE_KIND_) :: min_glo_lon_lnd, max_glo_lon_lnd
127 real(FMS_DATA_OVERRIDE_KIND_) :: min_glo_lon_ice, max_glo_lon_ice
128 integer :: num_fields = 0 !< number of fields in override_array already processed
129 
130 type(data_type), dimension(:), allocatable :: data_table !< user-provided data table
131 
132 type(data_type) :: default_table
133 type(override_type), dimension(max_array) :: override_array !< to store processed fields
134 type(override_type) :: default_array
135 logical :: debug_data_override
136 logical :: grid_center_bug = .false.
137 logical :: reproduce_null_char_bug = .false. !> Flag indicating
138  !! to reproduce the mpp_io bug where lat/lon_bnd were
139  !! not read correctly if null characters are present in
140  !! the netcdf file
141 logical :: use_center_grid_points=.false. !< Flag indicating
142  !! whether or not to use the centroid values of the
143  !! supergrid from the grid file as opposed to calculating it
144  !! by taking the average of the four corner points.
145  !! This is only relevant to OCN and ICE grids.
146 logical :: use_data_table_yaml = .false.
147 
148 namelist /data_override_nml/ debug_data_override, grid_center_bug, reproduce_null_char_bug, use_data_table_yaml, &
149  use_center_grid_points
150 
151 public :: data_override_init_impl_, data_override_unset_atm_, data_override_unset_ocn_, &
152  & data_override_unset_lnd_, data_override_unset_ice_, data_override_0d_, &
153  & data_override_2d_, data_override_3d_, data_override_ug_1d_, &
154  & data_override_ug_2d_
155 
156 contains
157 
158 !> @brief Assign default values for default_table, get domain of component models,
159 !! get global grids of component models.
160 !! Users should call data_override_init before calling data_override
161 !!
162 !! This subroutine should be called by data_override_init.
163 !!
164 !! Data_table is initialized here with default values. Users should provide "real" values
165 !! that will override the default values. Real values can be specified in either data_table
166 !! or data_table.yaml. Each line of data_table contains one data_entry. Items of data_entry
167 !! are comma-separated.
168 subroutine data_override_init_impl_(Atm_domain_in, Ocean_domain_in, Ice_domain_in, Land_domain_in, Land_domainUG_in)
169  type (domain2d), intent(in), optional :: Atm_domain_in !> Atmosphere domain
170  type (domain2d), intent(in), optional :: Ocean_domain_in !> Ocean domain
171  type (domain2d), intent(in), optional :: Ice_domain_in !> Ice domain
172  type (domain2d), intent(in), optional :: Land_domain_in !> Land domain
173  type(domainUG) , intent(in), optional :: Land_domainUG_in !> Land domain, unstructured grid
174 
175  character(len=18), parameter :: grid_file = 'INPUT/grid_spec.nc'
176  integer :: i, iunit, io_status, ierr, use_get_grid_version
177  logical :: atm_on, ocn_on, lnd_on, ice_on, lndUG_on
178  logical :: file_open
179  type(FmsNetcdfFile_t) :: fileobj
180 
181  debug_data_override = .false.
182 
183  read (input_nml_file, data_override_nml, iostat=io_status)
184  ierr = check_nml_error(io_status, 'data_override_nml')
185  iunit = stdlog()
186  write(iunit, data_override_nml)
187 
188 ! grid_center_bug is no longer supported.
189 if (grid_center_bug) then
190  call mpp_error(fatal, "data_override_init: You have overridden the default value of " // &
191  "grid_center_bug and set it to .true. in data_override_nml. This was a temporary workaround " // &
192  "that is no longer supported. Please remove this namelist variable.")
193 endif
194 
195 if (use_data_table_yaml) then
196  call mpp_error(note, "data_override_init: You are using the yaml version of the data table.")
197 else
198  call mpp_error(note, "data_override_init: You are using the legacy version of the data table.")
199 end if
200 
201  atm_on = PRESENT(atm_domain_in)
202  ocn_on = PRESENT(ocean_domain_in)
203  lnd_on = PRESENT(land_domain_in)
204  ice_on = PRESENT(ice_domain_in)
205  lndug_on = PRESENT(land_domainug_in)
206  if(.not. module_is_initialized) then
207  atm_domain = null_domain2d
208  ocn_domain = null_domain2d
209  lnd_domain = null_domain2d
210  ice_domain = null_domain2d
211  lnd_domainug = null_domainug
212  end if
213  if (atm_on) atm_domain = atm_domain_in
214  if (ocn_on) ocn_domain = ocean_domain_in
215  if (lnd_on) lnd_domain = land_domain_in
216  if (ice_on) ice_domain = ice_domain_in
217  if (lndug_on) lnd_domainug = land_domainug_in
218 
219  if(.not. module_is_initialized) then
220  call horiz_interp_init
221  call write_version_number("DATA_OVERRIDE_MOD", version)
222 
223 ! Initialize user-provided data table
224  default_table%gridname = 'non'
225  default_table%fieldname_code = 'none'
226  default_table%fieldname_file = 'none'
227  default_table%file_name = 'none'
228  default_table%factor = 1._lkind
229  default_table%interpol_method = 'bilinear'
230  default_table%multifile = .false.
231  default_table%prev_file_name = ''
232  default_table%next_file_name = ''
233 
234 #ifdef use_yaml
235  if (use_data_table_yaml) then
236  if (file_exists("data_table")) &
237  call mpp_error(fatal, "You cannot have the legacy data_table if use_data_table_yaml=.true.")
238  call read_table_yaml(data_table)
239  allocate(external_weights(table_size))
240  nweight_files = 0
241  else
242  if (file_exists("data_table.yaml"))&
243  call mpp_error(fatal, "You cannot have the yaml data_table if use_data_table_yaml=.false.")
244  allocate(data_table(max_table))
245  do i = 1, max_table
246  data_table(i) = default_table
247  enddo
248  call read_table(data_table)
249  end if
250 #else
251  if (file_exists("data_table.yaml"))&
252  call mpp_error(fatal, "You cannot have the yaml data_table if use_data_table_yaml=.false.")
253 
254  if (use_data_table_yaml) then
255  call mpp_error(fatal, "You cannot have use_data_table_yaml=.true. without compiling with -Duse_yaml")
256  else
257 
258  allocate(data_table(max_table))
259  do i = 1, max_table
260  data_table(i) = default_table
261  enddo
262  call read_table(data_table)
263  end if
264 #endif
265 
266 ! Initialize override array
267  default_array%gridname = 'NONE'
268  default_array%fieldname = 'NONE'
269  default_array%t_index = -1
270  default_array%dims = -1
271  default_array%comp_domain = -1
272  do i = 1, max_array
273  override_array(i) = default_array
274  enddo
276  end if
277 
278  module_is_initialized = .true.
279 
280  if ( .NOT. (atm_on .or. ocn_on .or. lnd_on .or. ice_on .or. lndug_on)) return
281  if (table_size .eq. 0) then
282  call mpp_error(note, "data_table is empty, not doing any data_overrides")
283  return
284  endif
285  call fms2_io_init
286 
287 ! Test if grid_file is already opened
288  inquire (file=trim(grid_file), opened=file_open)
289  if(file_open) call mpp_error(fatal, trim(grid_file)//' already opened')
290 
291  if(.not. open_file(fileobj, grid_file, 'read' )) then
292  call mpp_error(fatal, 'data_override_mod: Error in opening file '//trim(grid_file))
293  endif
294 
295  if(variable_exists(fileobj, "x_T" ) .OR. variable_exists(fileobj, "geolon_t" ) ) then
296  use_get_grid_version = 1
297  call close_file(fileobj)
298  else if(variable_exists(fileobj, "ocn_mosaic_file" ) .OR. variable_exists(fileobj, "gridfiles" ) ) then
299  use_get_grid_version = 2
300  if(variable_exists(fileobj, "gridfiles" ) ) then
301  if(count_ne_1((ocn_on .OR. ice_on), lnd_on, atm_on)) call mpp_error(fatal, 'data_override_mod: the grid file ' //&
302  'is a solo mosaic, one and only one of atm_on, lnd_on or ice_on/ocn_on should be true')
303  end if
304  else
305  call mpp_error(fatal, 'data_override_mod: none of x_T, geolon_t, ocn_mosaic_file or gridfiles exist in '// &
306  & trim(grid_file))
307  endif
308 
309  if(use_get_grid_version .EQ. 1) then
310  if (atm_on .and. .not. allocated(lon_local_atm) ) then
311  call get_grid_version_1(grid_file, 'atm', atm_domain, lon_local_atm, lat_local_atm, &
312  min_glo_lon_atm, max_glo_lon_atm )
313  endif
314  if (ocn_on .and. .not. allocated(lon_local_ocn) ) then
315  call get_grid_version_1(grid_file, 'ocn', ocn_domain, lon_local_ocn, lat_local_ocn, &
316  min_glo_lon_ocn, max_glo_lon_ocn )
317  endif
318 
319  if (lnd_on .and. .not. allocated(lon_local_lnd) ) then
320  call get_grid_version_1(grid_file, 'lnd', lnd_domain, lon_local_lnd, lat_local_lnd, &
321  min_glo_lon_lnd, max_glo_lon_lnd )
322  endif
323 
324  if (ice_on .and. .not. allocated(lon_local_ice) ) then
325  call get_grid_version_1(grid_file, 'ice', ice_domain, lon_local_ice, lat_local_ice, &
326  min_glo_lon_ice, max_glo_lon_ice )
327  endif
328  else
329  if (atm_on .and. .not. allocated(lon_local_atm) ) then
330  call get_grid_version_2(fileobj, 'atm', atm_domain, lon_local_atm, lat_local_atm, &
331  min_glo_lon_atm, max_glo_lon_atm )
332  endif
333 
334  if (ocn_on .and. .not. allocated(lon_local_ocn) ) then
335  call get_grid_version_2(fileobj, 'ocn', ocn_domain, lon_local_ocn, lat_local_ocn, &
336  min_glo_lon_ocn, max_glo_lon_ocn, use_center_grid_points)
337  endif
338 
339  if (lnd_on .and. .not. allocated(lon_local_lnd) ) then
340  call get_grid_version_2(fileobj, 'lnd', lnd_domain, lon_local_lnd, lat_local_lnd, &
341  min_glo_lon_lnd, max_glo_lon_lnd )
342  endif
343 
344  if (ice_on .and. .not. allocated(lon_local_ice) ) then
345  call get_grid_version_2(fileobj, 'ocn', ice_domain, lon_local_ice, lat_local_ice, &
346  min_glo_lon_ice, max_glo_lon_ice, use_center_grid_points)
347  endif
348  end if
349  if(use_get_grid_version .EQ. 2) then
350  call close_file(fileobj)
351  end if
352 end subroutine data_override_init_impl_
353 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
354 !> @brief Implementation of the following truth table:
355 !!
356 !! Arg 1 Arg 2 Arg 3 | Result !!
357 !!
358 !! .true. .true. .true. | .true. !!
359 !! .true. .true. .false. | .true. !!
360 !! .true. .false. .true. | .true. !!
361 !! .true. .false. .false. | .false. !!
362 !! .false. .true. .true. | .true. !!
363 !! .false. .true. .false. | .false. !!
364 !! .false. .false. .true. | .false. !!
365 !! .false. .false. .false. | .true. !!
366 function count_ne_1(in_1, in_2, in_3)
367  logical, intent(in) :: in_1 !< Argument 1
368  logical, intent(in) :: in_2 !< Argument 2
369  logical, intent(in) :: in_3 !< Argument 3
370  logical :: count_ne_1
371 
372  count_ne_1 = .not.(in_1.NEQV.in_2.NEQV.in_3) .OR. (in_1.AND.in_2.AND.in_3)
373 end function count_ne_1
374 
375 subroutine read_table(data_table)
376  type(data_type), dimension(max_table), intent(inout) :: data_table
377 
378  integer :: ntable
379  integer :: ntable_lima
380  integer :: ntable_new
381 
382  integer :: iunit
383  integer :: io_status
384  integer :: index_1col, index_2col
385  character(len=256) :: record
386  type(data_type) :: data_entry
387 
388  logical :: ongrid
389  logical :: table_exists !< Flag indicating existence of data_table
390  character(len=128) :: region, region_type
391 
392  integer :: sunit
393 
394 ! Read coupler_table
395  inquire(file='data_table', exist=table_exists)
396  if (.not. table_exists) then
397  call mpp_error(note, 'data_override_mod: File data_table does not exist.')
398  table_size = 0
399  return
400  end if
401 
402  open(newunit=iunit, file='data_table', action='READ', iostat=io_status)
403  if(io_status/=0) call mpp_error(fatal, 'data_override_mod: Error in opening file data_table.')
404 
405  ntable = 0
406  ntable_lima = 0
407  ntable_new = 0
408 
409  do while (ntable <= max_table)
410  read(iunit,'(a)',end=100) record
411  if (record(1:1) == '#') cycle
412  if (record(1:10) == ' ') cycle
413  ntable=ntable+1
414  if(index(lowercase(record), "inside_region") .ne. 0 .or. index(lowercase(record), "outside_region") .ne. 0) then
415  if(index(lowercase(record), ".false.") .ne. 0 .or. index(lowercase(record), ".true.") .ne. 0 ) then
416  ntable_lima = ntable_lima + 1
417  read(record,*,err=99) data_entry%gridname, data_entry%fieldname_code, data_entry%fieldname_file, &
418  data_entry%file_name, ongrid, data_entry%factor, region, region_type
419  if(ongrid) then
420  data_entry%interpol_method = 'none'
421  else
422  data_entry%interpol_method = 'bilinear'
423  endif
424  else
425  ntable_new=ntable_new+1
426  read(record,*,err=99) data_entry%gridname, data_entry%fieldname_code, data_entry%fieldname_file, &
427  data_entry%file_name, data_entry%interpol_method, data_entry%factor, region, &
428  & region_type
429 
430  if (index(data_entry%file_name, ":") .ne. 0) then
431  data_entry%multifile = .true.
432  index_1col = index(data_entry%file_name, ":")
433  index_2col = index(data_entry%file_name(index_1col+1:), ":")
434  if (index_2col .eq. 0) call mpp_error(fatal, "data_override: when bridging over forcing files, " &
435  // "central forcing files must be preceded AND followed by the column (:) separator")
436  data_entry%prev_file_name = data_entry%file_name(1:index_1col-1)
437  data_entry%next_file_name = data_entry%file_name(index_1col+index_2col+1:)
438  ! once previous/next files are filled in, overwrite current
439  data_entry%file_name = data_entry%file_name(index_1col+1:index_1col+index_2col-1)
440  else
441  data_entry%multifile = .false.
442  data_entry%prev_file_name = ""
443  data_entry%next_file_name = ""
444  endif
445  if (data_entry%interpol_method == 'default') then
446  data_entry%interpol_method = default_table%interpol_method
447  endif
448  if (.not.(data_entry%interpol_method == 'default' .or. &
449  data_entry%interpol_method == 'bicubic' .or. &
450  data_entry%interpol_method == 'bilinear' .or. &
451  data_entry%interpol_method == 'none')) then
452  sunit = stdout()
453  write(sunit,*)" gridname is ", trim(data_entry%gridname)
454  write(sunit,*)" fieldname_code is ", trim(data_entry%fieldname_code)
455  write(sunit,*)" fieldname_file is ", trim(data_entry%fieldname_file)
456  write(sunit,*)" file_name is ", trim(data_entry%file_name)
457  write(sunit,*)" factor is ", data_entry%factor
458  write(sunit,*)" interpol_method is ", trim(data_entry%interpol_method)
459  call mpp_error(fatal, 'data_override_mod: invalid last entry in data_override_table, ' &
460  //'its value should be "default", "bicubic", "bilinear" or "none" ')
461  endif
462  endif
463  if( trim(region_type) == "inside_region" ) then
464  data_entry%region_type = inside_region
465  else if( trim(region_type) == "outside_region" ) then
466  data_entry%region_type = outside_region
467  else
468  call mpp_error(fatal, 'data_override_mod: region type should be inside_region or outside_region')
469  endif
470  if (data_entry%file_name == "") call mpp_error(fatal, &
471  "data_override: filename not given in data_table when region_type is not NO_REGION")
472  if(data_entry%fieldname_file == "") call mpp_error(fatal, &
473  "data_override: fieldname_file must be specified in data_table when region_type is not NO_REGION")
474  if( trim(data_entry%interpol_method) == 'none') call mpp_error(fatal, &
475  "data_override(data_override_init): ongrid must be false when region_type is not NO_REGION")
476  read(region,*) data_entry%lon_start, data_entry%lon_end, data_entry%lat_start, data_entry%lat_end
477  !--- make sure data_entry%lon_end > data_entry%lon_start and data_entry%lat_end > data_entry%lat_start
478  if(data_entry%lon_end .LE. data_entry%lon_start) call mpp_error(fatal, &
479  "data_override: lon_end should be greater than lon_start")
480  if(data_entry%lat_end .LE. data_entry%lat_start) call mpp_error(fatal, &
481  "data_override: lat_end should be greater than lat_start")
482  ! old format
483  else if (index(lowercase(record), ".false.") .ne. 0 .or. index(lowercase(record), ".true.") .ne. 0 ) then
484  ntable_lima = ntable_lima + 1
485  read(record,*,err=99) data_entry%gridname, data_entry%fieldname_code, data_entry%fieldname_file, &
486  data_entry%file_name, ongrid, data_entry%factor
487  if (index(data_entry%file_name, ":") .ne. 0) then
488  data_entry%multifile = .true.
489  index_1col = index(data_entry%file_name, ":")
490  index_2col = index(data_entry%file_name(index_1col+1:), ":")
491  if (index_2col .eq. 0) call mpp_error(fatal, "data_override: when bridging over forcing files, " &
492  // "central forcing files must be preceded AND followed by the column (:) separator")
493  data_entry%prev_file_name = data_entry%file_name(1:index_1col-1)
494  data_entry%next_file_name = data_entry%file_name(index_1col+index_2col+1:)
495  ! once previous/next files are filled in, overwrite current
496  data_entry%file_name = data_entry%file_name(index_1col+1:index_1col+index_2col-1)
497  else
498  data_entry%multifile = .false.
499  data_entry%prev_file_name = ""
500  data_entry%next_file_name = ""
501  endif
502  if(ongrid) then
503  data_entry%interpol_method = 'none'
504  else
505  data_entry%interpol_method = 'bilinear'
506  endif
507  data_entry%lon_start = 0.0_lkind
508  data_entry%lon_end = -1.0_lkind
509  data_entry%lat_start = 0.0_lkind
510  data_entry%lat_end = -1.0_lkind
511  data_entry%region_type = no_region
512  else ! new format
513  ntable_new=ntable_new+1
514  read(record,*,err=99) data_entry%gridname, data_entry%fieldname_code, data_entry%fieldname_file, &
515  data_entry%file_name, data_entry%interpol_method, data_entry%factor
516  if (index(data_entry%file_name, ":") .ne. 0) then
517  index_1col = index(data_entry%file_name, ":")
518  index_2col = index(data_entry%file_name(index_1col+1:), ":")
519  data_entry%multifile = .true.
520  if (index_2col .eq. 0) call mpp_error(fatal, "data_override: when bridging over forcing files, " &
521  // "central forcing files must be preceded AND followed by the column (:) separator")
522  data_entry%prev_file_name = data_entry%file_name(1:index_1col-1)
523  data_entry%next_file_name = data_entry%file_name(index_1col+index_2col+1:)
524  ! once previous/next files are filled in, overwrite current
525  data_entry%file_name = data_entry%file_name(index_1col+1:index_1col+index_2col-1)
526  else
527  data_entry%multifile = .false.
528  data_entry%prev_file_name = ""
529  data_entry%next_file_name = ""
530  endif
531  if (data_entry%interpol_method == 'default') then
532  data_entry%interpol_method = default_table%interpol_method
533  endif
534  if (.not.(data_entry%interpol_method == 'default' .or. &
535  data_entry%interpol_method == 'bicubic' .or. &
536  data_entry%interpol_method == 'bilinear' .or. &
537  data_entry%interpol_method == 'none')) then
538  sunit = stdout()
539  write(sunit,*)" gridname is ", trim(data_entry%gridname)
540  write(sunit,*)" fieldname_code is ", trim(data_entry%fieldname_code)
541  write(sunit,*)" fieldname_file is ", trim(data_entry%fieldname_file)
542  write(sunit,*)" file_name is ", trim(data_entry%file_name)
543  write(sunit,*)" factor is ", data_entry%factor
544  write(sunit,*)" interpol_method is ", trim(data_entry%interpol_method)
545  call mpp_error(fatal, 'data_override_mod: invalid last entry in data_override_table, ' &
546  //'its value should be "default", "bicubic", "bilinear" or "none" ')
547  endif
548  data_entry%lon_start = 0.0_lkind
549  data_entry%lon_end = -1.0_lkind
550  data_entry%lat_start = 0.0_lkind
551  data_entry%lat_end = -1.0_lkind
552  data_entry%region_type = no_region
553  endif
554  data_entry%ext_weights = .false.
555  data_table(ntable) = data_entry
556  enddo
557  call mpp_error(fatal,'too many enries in data_table')
558 99 call mpp_error(fatal,'error in data_table format')
559 100 continue
560  table_size = ntable
561  if(ntable_new*ntable_lima /= 0) call mpp_error(fatal, &
562  'data_override_mod: New and old formats together in same data_table not supported')
563  close(iunit, iostat=io_status)
564  if(io_status/=0) call mpp_error(fatal, 'data_override_mod: Error in closing file data_table')
565 end subroutine read_table
566 
567 #ifdef use_yaml
568 !> @brief Read and parse the data_table.yaml
569 subroutine read_table_yaml(data_table)
570  type(data_type), dimension(:), allocatable, intent(out) :: data_table !< Contents of the data_table.yaml
571 
572  integer, allocatable :: entry_id(:)
573  integer :: sub_block_id(1), sub2_block_id(1)
574  integer :: nentries, mentries
575  integer :: i
576  character(len=50) :: buffer
577  character(len=FMS_FILE_LEN) :: filename !< Name of the expected data_table.yaml
578  integer :: file_id
579 
580  ! If doing and ensemble or nest run add the filename appendix (ens_XX or nest_XX) to the filename
581  call get_instance_filename("data_table.yaml", filename)
582  if (index(trim(filename), "ens_") .ne. 0) then
583  if (file_exists(filename) .and. file_exists("data_table.yaml")) &
584  call mpp_error(fatal, "Both data_table.yaml and "//trim(filename)//" exists, pick one!")
585 
586  !< If the end_* file does not exist, revert back to "data_table.yaml"
587  !! where every ensemble is using the same yaml
588  if (.not. file_exists(filename)) filename = "data_table.yaml"
589  endif
590 
591  file_id = open_and_parse_file(trim(filename))
592 
593  if (file_id==999) then
594  nentries = 0
595  else
596  nentries = get_num_blocks(file_id, "data_table")
597  allocate(data_table(nentries))
598  allocate(entry_id(nentries))
599  call get_block_ids(file_id, "data_table", entry_id)
600 
601  do i = 1, nentries
602  call get_value_from_key(file_id, entry_id(i), "factor", data_table(i)%factor)
603  call get_value_from_key(file_id, entry_id(i), "grid_name", data_table(i)%gridname)
604  call check_for_valid_gridname(data_table(i)%gridname)
605  call get_value_from_key(file_id, entry_id(i), "fieldname_in_model", data_table(i)%fieldname_code)
606 
607  mentries = get_num_blocks(file_id, "override_file", parent_block_id=entry_id(i))
608  data_table(i)%file_name = ""
609  data_table(i)%fieldname_file = ""
610  data_table(i)%interpol_method = "none"
611  data_table(i)%multifile = .false.
612  data_table(i)%ext_weights = .false.
613  data_table(i)%region_type = no_region
614  data_table(i)%prev_file_name = ""
615  data_table(i)%next_file_name = ""
616  data_table(i)%ext_weights_file_name = ""
617  data_table(i)%ext_weights_source = ""
618 
619  ! If there is no override_file block, then not overriding from file, so move on to the next entry
620  if (mentries .eq. 0) cycle
621 
622  if(mentries.gt.1) call mpp_error(fatal, "Too many override_file blocks in data table. "//&
623  "Check your data_table.yaml entry for field:"//trim(data_table(i)%gridname)//":"//&
624  trim(data_table(i)%fieldname_code))
625  call get_block_ids(file_id, "override_file", sub_block_id, parent_block_id=entry_id(i))
626 
627  call get_value_from_key(file_id, sub_block_id(1), "file_name", data_table(i)%file_name)
628  call get_value_from_key(file_id, sub_block_id(1), "fieldname_in_file", data_table(i)%fieldname_file)
629  call get_value_from_key(file_id, sub_block_id(1), "interp_method", data_table(i)%interpol_method)
630  call check_interpol_method(data_table(i)%interpol_method, data_table(i)%file_name, &
631  & data_table(i)%fieldname_file)
632 
633  mentries = get_num_blocks(file_id, "multi_file", parent_block_id=sub_block_id(1))
634  if(mentries.gt.1) call mpp_error(fatal, "Too many multi_file blocks in tata table. "//&
635  "Check your data_table.yaml entry for field:"//trim(data_table(i)%gridname)//":"//&
636  trim(data_table(i)%fieldname_code))
637 
638  if(mentries.gt.0) data_table(i)%multifile = .true.
639 
640  if (data_table(i)%multifile) then
641  call get_block_ids(file_id, "multi_file", sub2_block_id, parent_block_id=sub_block_id(1))
642  call get_value_from_key(file_id, sub2_block_id(1), "prev_file_name", data_table(i)%prev_file_name)
643  call get_value_from_key(file_id, sub2_block_id(1), "next_file_name", data_table(i)%next_file_name)
644  if (trim(data_table(i)%prev_file_name) .eq. "" .and. trim(data_table(i)%next_file_name) .eq. "") &
645  call mpp_error(fatal, "The prev_file_name and next_file_name must be present if is_multi_file. "//&
646  "Check your data_table.yaml entry for field:"//trim(data_table(i)%gridname)//":"//&
647  trim(data_table(i)%fieldname_code))
648  endif
649 
650  mentries = get_num_blocks(file_id, "external_weights", parent_block_id=sub_block_id(1))
651  if(mentries.gt.1) call mpp_error(fatal, "Too many external_weight blocks in data table. "//&
652  "Check your data_table.yaml entry for field:"//trim(data_table(i)%gridname)//":"//&
653  trim(data_table(i)%fieldname_code))
654 
655  if(mentries.gt.0) data_table(i)%ext_weights = .true.
656 
657  if (data_table(i)%ext_weights) then
658  call get_block_ids(file_id, "external_weights", sub2_block_id, parent_block_id=sub_block_id(1))
659  call get_value_from_key(file_id, sub2_block_id(1), "file_name", data_table(i)%ext_weights_file_name)
660  call get_value_from_key(file_id, sub2_block_id(1), "source", data_table(i)%ext_weights_source)
661  if (trim(data_table(i)%ext_weights_file_name) .eq. "" .and. trim(data_table(i)%ext_weights_source) .eq. "") &
662  call mpp_error(fatal, "The file_name and source must be present when using external weights"//&
663  "Check your data_table.yaml entry for field:"//trim(data_table(i)%gridname)//":"//&
664  trim(data_table(i)%fieldname_code))
665  endif
666 
667  mentries = get_num_blocks(file_id, "subregion", parent_block_id=entry_id(i))
668  if(mentries.gt.1) call mpp_error(fatal, "Too many subregion blocks in data table. "//&
669  "Check your data_table.yaml entry for field:"//trim(data_table(i)%gridname)//":"//&
670  trim(data_table(i)%fieldname_code))
671 
672  buffer = ""
673  if(mentries.gt.0) then
674  call get_block_ids(file_id, "subregion", sub_block_id, parent_block_id=entry_id(i))
675  call get_value_from_key(file_id, sub_block_id(1), "type", buffer)
676  endif
677 
678  call check_and_set_region_type(buffer, data_table(i)%region_type)
679  if (data_table(i)%region_type .ne. no_region) then
680  call get_value_from_key(file_id, sub_block_id(1), "lon_start", data_table(i)%lon_start)
681  call get_value_from_key(file_id, sub_block_id(1), "lon_end", data_table(i)%lon_end)
682  call get_value_from_key(file_id, sub_block_id(1), "lat_start", data_table(i)%lat_start)
683  call get_value_from_key(file_id, sub_block_id(1), "lat_end", data_table(i)%lat_end)
684  call check_valid_lat_lon(data_table(i)%lon_start, data_table(i)%lon_end, &
685  data_table(i)%lat_start, data_table(i)%lat_end)
686  endif
687  end do
688 
689  end if
690  table_size = nentries !< Because one variable is not enough
691 end subroutine read_table_yaml
692 
693 !> @brief Check if a grid name is valid, crashes if it is not
694 subroutine check_for_valid_gridname(gridname)
695  character(len=*), intent(in) :: gridname !< Gridname
696 
697  select case(trim(gridname))
698  case ("OCN", "ATM", "LND", "ICE")
699  case default
700  call mpp_error(fatal, trim(gridname)//" is not a valid gridname. "//&
701  "The acceptable values are OCN ATM LND and ICE. Check your data_table.yaml")
702  end select
703 end subroutine check_for_valid_gridname
704 
705 !> @brief Check if the interpol method is correct, crashes if it is not
706 subroutine check_interpol_method(interp_method, filename, fieldname)
707  character(len=*), intent(in) :: interp_method !< The interpo_method
708  character(len=*), intent(in) :: filename !< The filename
709  character(len=*), intent(in) :: fieldname !< The fieldname in the file
710 
711  select case(trim(interp_method))
712  case ("bicubic", "bilinear")
713  if (trim(filename) .eq. "" .or. trim(fieldname) .eq. "") call mpp_error(fatal, &
714  "The file_name and the fieldname_file must be set if using the bicubic or bilinear interpolation method."//&
715  " Check your data_table.yaml")
716  case ("none")
717  if (trim(filename) .ne. "" ) then
718  if (trim(fieldname) .eq. "") call mpp_error(fatal, &
719  "If the interpol_method is none and file_name is specified (ongrid case), "//&
720  "you must also specify the fieldname_file")
721  endif
722  case default
723  call mpp_error(fatal, trim(interp_method)//" is not a valid interp method. "//&
724  "The acceptable values are bilinear and bicubic")
725  end select
726 end subroutine check_interpol_method
727 
728 !> @brief Check if a region_type is valid, crashes if it is not. Otherwise it sets the
729 !! correct integer parameter.
730 subroutine check_and_set_region_type(region_type_str, region_type_int)
731  character(len=*), intent(in) :: region_type_str !< The region type as defined in the data.yaml
732  integer, intent(out) :: region_type_int !< The region type as an integer parameter
733 
734  select case(trim(region_type_str))
735  case ("inside_region")
736  region_type_int = inside_region
737  case ("outside_region")
738  region_type_int = outside_region
739  case ("")
740  region_type_int = no_region
741  case default
742  call mpp_error(fatal, trim(region_type_str)//" is not a valid region type. "//&
743  "The acceptable values are inside_region and outside_regioon. Check your data_table.yaml")
744  end select
745 end subroutine check_and_set_region_type
746 
747 !> @brief Check if a region lon_start, lon_end, lat_start and lat_end is valid.
748 !! Crashes if it is not.
749 subroutine check_valid_lat_lon(lon_start, lon_end, lat_start, lat_end)
750  real(FMS_DATA_OVERRIDE_KIND_), intent(in) :: lon_start !< Starting longitude of the data_override region
751  real(FMS_DATA_OVERRIDE_KIND_), intent(in) :: lon_end !< Ending longitude of the data_override region
752  real(FMS_DATA_OVERRIDE_KIND_), intent(in) :: lat_start !< Starting lattiude of the data_override region
753  real(FMS_DATA_OVERRIDE_KIND_), intent(in) :: lat_end !< Ending lattiude of the data_override region
754 
755  if (lon_start > lon_end) call mpp_error(fatal, &
756  "lon_start:"//string(lon_start)//" is greater than lon_end"//string(lon_end)//&
757  ". Check your data_table.yaml.")
758 
759  if (lat_start > lat_end) call mpp_error(fatal, &
760  "lat_start:"//string(lat_start)//" is greater than lat_end:"//string(lat_end)//&
761  ". Check your data_table.yaml.")
762 end subroutine check_valid_lat_lon
763 #endif
764 
765 subroutine data_override_unset_atm_
766  atm_domain = null_domain2d
767  if (allocated(lon_local_atm)) deallocate(lon_local_atm)
768  if (allocated(lat_local_atm)) deallocate(lat_local_atm)
769 end subroutine
770 
771 subroutine data_override_unset_ocn_
772  ocn_domain = null_domain2d
773  if (allocated(lon_local_ocn)) deallocate(lon_local_ocn)
774  if (allocated(lat_local_ocn)) deallocate(lat_local_ocn)
775 end subroutine
776 
777 subroutine data_override_unset_lnd_
778  lnd_domain = null_domain2d
779  if (allocated(lon_local_lnd)) deallocate(lon_local_lnd)
780  if (allocated(lat_local_lnd)) deallocate(lat_local_lnd)
781 end subroutine
782 
783 subroutine data_override_unset_ice_
784  ice_domain = null_domain2d
785  if (allocated(lon_local_ice)) deallocate(lon_local_ice)
786  if (allocated(lat_local_ice)) deallocate(lat_local_ice)
787 end subroutine
788 
789 !> @brief Given a gridname, this routine returns the working domain associated with this gridname
790 subroutine get_domain(gridname, domain, comp_domain)
791  character(len=3), intent(in) :: gridname
792  type(domain2D), intent(inout) :: domain
793  integer, intent(out), optional :: comp_domain(4) !< istart,iend,jstart,jend for compute domain
794 
795  domain = null_domain2d
796  select case (gridname)
797  case('OCN')
798  domain = ocn_domain
799  case('ATM')
800  domain = atm_domain
801  case('LND')
802  domain = lnd_domain
803  case('ICE')
804  domain = ice_domain
805  case default
806  call mpp_error(fatal,'error in data_override get_domain')
807  end select
808  if(domain .EQ. null_domain2d) call mpp_error(fatal,'data_override: failure in get_domain')
809  if(present(comp_domain)) &
810  call mpp_get_compute_domain(domain,comp_domain(1),comp_domain(2),comp_domain(3),comp_domain(4))
811 end subroutine get_domain
812 
813 !> @brief Given a gridname, this routine returns the working domain associated with this gridname
814 subroutine get_domainug(gridname, UGdomain, comp_domain)
815  character(len=3), intent(in) :: gridname
816  type(domainUG), intent(inout) :: UGdomain
817  integer, intent(out), optional :: comp_domain(4) !< istart,iend,jstart,jend for compute domain
818  type(domain2D), pointer :: SGdomain => null()
819 
820  ugdomain = null_domainug
821  select case (gridname)
822  case('LND')
823  ugdomain = lnd_domainug
824  case default
825  call mpp_error(fatal,'error in data_override get_domain')
826  end select
827 ! if(UGdomain .EQ. NULL_DOMAINUG) call mpp_error(FATAL,'data_override: failure in get_domain')
828  if(present(comp_domain)) &
829  call mpp_get_ug_sg_domain(ugdomain,sgdomain)
830  call mpp_get_compute_domain(sgdomain,comp_domain(1),comp_domain(2),comp_domain(3),comp_domain(4))
831 end subroutine get_domainug
832 !===============================================================================================
833 
834 !> @brief Routine to perform data override for scalar fields
835 subroutine data_override_0d_(gridname,fieldname_code,data_out,time,override,data_index)
836  character(len=3), intent(in) :: gridname !< model grid ID (ocn,ice,atm,lnd)
837  character(len=*), intent(in) :: fieldname_code !< field name as used in the model (may be
838  !! different from the name in NetCDF data file)
839  logical, intent(out), optional :: override !< true if the field has been overridden succesfully
840  type(time_type), intent(in) :: time !< (target) model time
841  real(FMS_DATA_OVERRIDE_KIND_), intent(out) :: data_out !< output data array returned by this call
842  integer, intent(in), optional :: data_index
843 
844  type(time_type) :: first_record !< first record of "current" file
845  type(time_type) :: last_record !< last record of "current" file
846  character(len=FMS_PATH_LEN) :: filename !< file containing source data
847  character(len=FMS_PATH_LEN) :: prevfilename !< file containing previous source data, when using multiple files
848  character(len=FMS_PATH_LEN) :: nextfilename !< file containing next source data, when using multiple files
849  character(len=128) :: fieldname !< fieldname used in the data file
850  integer :: index1 !< field index in data_table
851  integer :: dims(4)
852  integer :: prev_dims(4) !< dimensions of previous source data, when using multiple files
853  integer :: next_dims(4) !< dimensions of next source data, when using multiple files
854  integer :: id_time !< index for time interp in override array
855  integer :: id_time_prev=-1 !< time index for previous file, when using multiple files
856  integer :: id_time_next=-1 !< time index for next file, when using multiple files
857  integer :: curr_position !< position of the field currently processed in override_array
858  integer :: i
859  real(FMS_DATA_OVERRIDE_KIND_) :: factor
860  logical :: multifile !< use multiple consecutive files for override
861 
862  if(.not.module_is_initialized) &
863  call mpp_error(fatal,'Error: need to call data_override_init first')
864 
865 !1 Look for the data file in data_table
866  if(PRESENT(override)) override = .false.
867  if (present(data_index)) then
868  index1 = data_index
869  else
870  index1 = -1
871  do i = 1, table_size
872  if( trim(gridname) /= trim(data_table(i)%gridname)) cycle
873  if( trim(fieldname_code) /= trim(data_table(i)%fieldname_code)) cycle
874  index1 = i ! field found
875  exit
876  enddo
877  if(index1 .eq. -1) then
878  if(debug_data_override) &
879  call mpp_error(warning,'this field is NOT found in data_table: '//trim(fieldname_code))
880  return ! NO override was performed
881  endif
882  endif
883 
884  fieldname = data_table(index1)%fieldname_file ! fieldname in netCDF data file
885  factor = data_table(index1)%factor
886  multifile = data_table(index1)%multifile
887 
888  if(fieldname == "") then
889  data_out = factor
890  if(PRESENT(override)) override = .true.
891  return
892  else
893  filename = data_table(index1)%file_name
894  if (filename == "") call mpp_error(fatal,'data_override: filename not given in data_table')
895  if (multifile) prevfilename = data_table(index1)%prev_file_name
896  if (multifile) nextfilename = data_table(index1)%next_file_name
897  endif
898 
899 !3 Check if fieldname has been previously processed
900 !$OMP SINGLE
901  curr_position = -1
902  if(num_fields > 0 ) then
903  do i = 1, num_fields
904  if(trim(override_array(i)%gridname) /= trim(gridname)) cycle
905  if(trim(override_array(i)%fieldname) /= trim(fieldname_code)) cycle
906  curr_position = i
907  exit
908  enddo
909  endif
910 
911  if(curr_position < 0) then ! the field has not been processed previously
912  num_fields = num_fields + 1
913  curr_position = num_fields
914  ! record fieldname, gridname in override_array
915  override_array(curr_position)%fieldname = fieldname_code
916  override_array(curr_position)%gridname = gridname
917  id_time = init_external_field(filename,fieldname,verbose=debug_data_override)
918  if(id_time<0) call mpp_error(fatal,'data_override:field not found in init_external_field 1')
919  override_array(curr_position)%t_index = id_time
920  else !curr_position >0
921  !9 Get id_time previously stored in override_array
922  id_time = override_array(curr_position)%t_index
923  endif !if curr_position < 0
924 
925 
926  ! if using consecutive files for data_override, get time axis for previous and next files
927  ! and check spatial dims for consistency
928  if_multi1: if (multifile) then
929  id_time_prev = -1
930  if_prev1: if (trim(prevfilename) /= '') then
931  id_time_prev = init_external_field(prevfilename,fieldname,verbose=debug_data_override)
932  dims = get_external_field_size(id_time)
933  prev_dims = get_external_field_size(id_time_prev)
934  ! check consistency of spatial dims
935  if ((prev_dims(1) .ne. dims(1)) .or. (prev_dims(2) .ne. dims(2)) .or. &
936  (prev_dims(3) .ne. dims(3))) then
937  call mpp_error(fatal, 'data_override: dimensions mismatch between consecutive forcing files')
938  endif
939  allocate(data_table(index1)%time_prev_records(prev_dims(4)))
940  call get_time_axis(id_time_prev,data_table(index1)%time_prev_records)
941  endif if_prev1
942  id_time_next = -1
943  if_next1: if (trim(nextfilename) /= '') then
944  id_time_next = init_external_field(nextfilename,fieldname,verbose=debug_data_override)
945  dims = get_external_field_size(id_time)
946  next_dims = get_external_field_size(id_time_next)
947  ! check consistency of spatial dims
948  if ((next_dims(1) .ne. dims(1)) .or. (next_dims(2) .ne. dims(2)) .or. &
949  (next_dims(3) .ne. dims(3))) then
950  call mpp_error(fatal, 'data_override: dimensions mismatch between consecutive forcing files')
951  endif
952  allocate(data_table(index1)%time_next_records(next_dims(4)))
953  call get_time_axis(id_time_next,data_table(index1)%time_next_records)
954  endif if_next1
955  endif if_multi1
956 
957 
958  !10 do time interp to get data in compute_domain
959 
960  ! if using consecutive files, allow to perform time interpolation between the last record of previous
961  ! file and first record of current file OR between the last record of current file and first record of
962  ! next file hence "bridging" over files.
963  if_multi2: if (multifile) then
964  dims = get_external_field_size(id_time)
965  if (.not. allocated(data_table(index1)%time_records)) allocate(data_table(index1)%time_records(dims(4)))
966  call get_time_axis(id_time,data_table(index1)%time_records)
967 
968  first_record = data_table(index1)%time_records(1)
969  last_record = data_table(index1)%time_records(dims(4))
970 
971  if_time2: if (time<first_record) then
972  if (id_time_prev<0) call mpp_error(fatal,'data_override:previous file needed with multifile')
973  prev_dims = get_external_field_size(id_time_prev)
974  if (time<data_table(index1)%time_prev_records(prev_dims(4))) call mpp_error(fatal, &
975  'data_override: time_interp_external_bridge should only be called to bridge with previous file')
976  call time_interp_external_bridge(id_time_prev, id_time,time,data_out,verbose=debug_data_override)
977  elseif (time>last_record) then
978  if (id_time_next<0) call mpp_error(fatal,'data_override:next file needed with multifile')
979  if (time>data_table(index1)%time_next_records(1)) call mpp_error(fatal, &
980  'data_override: time_interp_external_bridge should only be called to bridge with next file')
981  call time_interp_external_bridge(id_time, id_time_next,time,data_out,verbose=debug_data_override)
982  else ! first_record < time < last_record, do not use bridge
983  call time_interp_external(id_time,time,data_out,verbose=debug_data_override)
984  endif if_time2
985  else ! standard behavior
986  call time_interp_external(id_time,time,data_out,verbose=debug_data_override)
987  endif if_multi2
988 
989 
990  data_out = data_out*factor
991 !$OMP END SINGLE
992 
993  if(PRESENT(override)) override = .true.
994 
995 end subroutine data_override_0d_
996 
997 !> @brief This routine performs data override for 2D fields.
998 subroutine data_override_2d_(gridname,fieldname,data_2D,time,override, is_in, ie_in, js_in, je_in)
999  character(len=3), intent(in) :: gridname !< model grid ID
1000  character(len=*), intent(in) :: fieldname !< field to override
1001  logical, intent(out), optional :: override !< true if the field has been overridden succesfully
1002  type(time_type), intent(in) :: time !< model time
1003  real(FMS_DATA_OVERRIDE_KIND_), dimension(:,:), intent(inout) :: data_2D !< data returned by this call
1004  integer, optional, intent(in) :: is_in, ie_in, js_in, je_in
1005  real(FMS_DATA_OVERRIDE_KIND_), dimension(:,:,:), allocatable :: data_3D
1006  integer :: index1
1007  integer :: i
1008 
1009 !1 Look for the data file in data_table
1010  if(PRESENT(override)) override = .false.
1011  index1 = -1
1012  do i = 1, table_size
1013  if( trim(gridname) /= trim(data_table(i)%gridname)) cycle
1014  if( trim(fieldname) /= trim(data_table(i)%fieldname_code)) cycle
1015  index1 = i ! field found
1016  exit
1017  enddo
1018  if(index1 .eq. -1) return ! NO override was performed
1019 
1020  allocate(data_3d(size(data_2d,1),size(data_2d,2),1))
1021  data_3d(:,:,1) = data_2d
1022  call data_override_3d_(gridname,fieldname,data_3d,time,override,data_index=index1,&
1023  is_in=is_in,ie_in=ie_in,js_in=js_in,je_in=je_in)
1024 
1025  data_2d(:,:) = data_3d(:,:,1)
1026  deallocate(data_3d)
1027 end subroutine data_override_2d_
1028 
1029 !> @brief This routine performs data override for 3D fields
1030 subroutine data_override_3d_(gridname,fieldname_code,return_data,time,override,data_index, is_in, ie_in, js_in, je_in)
1031  character(len=3), intent(in) :: gridname !< model grid ID
1032  character(len=*), intent(in) :: fieldname_code !< field name as used in the model
1033  logical, optional, intent(out) :: override !< true if the field has been overridden succesfully
1034  type(time_type), intent(in) :: time !< (target) model time
1035  integer, optional, intent(in) :: data_index
1036  real(FMS_DATA_OVERRIDE_KIND_), dimension(:,:,:), intent(inout) :: return_data !< data returned by this call
1037  integer, optional, intent(in) :: is_in, ie_in, js_in, je_in
1038  logical, dimension(:,:,:), allocatable :: mask_out
1039 
1040  character(len=FMS_PATH_LEN) :: filename !< file containing source data
1041  character(len=FMS_PATH_LEN) :: filename2 !< file containing source data
1042  character(len=FMS_PATH_LEN) :: prevfilename !< file containing source data for previous file
1043  character(len=FMS_PATH_LEN) :: prevfilename2 !< file containing source data for previous file
1044  character(len=FMS_PATH_LEN) :: nextfilename !< file containing source data for next file
1045  character(len=FMS_PATH_LEN) :: nextfilename2 !< file containing source data for next file
1046  character(len=128) :: fieldname !< fieldname used in the data file
1047  integer :: i,j
1048  integer :: dims(4)
1049  integer :: prev_dims(4) !< dimensions of previous source data, when using multiple files
1050  integer :: next_dims(4) !< dimensions of next source data, when using multiple files
1051  integer :: index1 !< field index in data_table
1052  integer :: id_time !< index for time interp in override array
1053  integer :: id_time_prev=-1 !< time index for previous file, when using multiple files
1054  integer :: id_time_next=-1 !< time index for next file, when using multiple files
1055  integer :: axis_sizes(4)
1056  character(len=32) :: axis_names(4)
1057  real(FMS_DATA_OVERRIDE_KIND_), dimension(:,:), pointer :: lon_local =>null() !< of output (target) grid cells
1058  real(FMS_DATA_OVERRIDE_KIND_), dimension(:,:), pointer :: lat_local =>null() !< of output (target) grid cells
1059  real(FMS_DATA_OVERRIDE_KIND_), dimension(:), allocatable :: lon_tmp, lat_tmp
1060 
1061  logical :: data_file_is_2D = .false. !< data in netCDF file is 2D
1062  logical :: ongrid, use_comp_domain
1063  type(domain2D) :: domain
1064  type(time_type) :: first_record !< first record of "current" file
1065  type(time_type) :: last_record !< last record of "current" file
1066  integer :: curr_position !< position of the field currently processed in override_array
1067  real(FMS_DATA_OVERRIDE_KIND_) :: factor
1068  integer, dimension(4) :: comp_domain = 0 !< istart,iend,jstart,jend for compute domain
1069  integer :: nxd, nyd, nxc, nyc, nwindows
1070  integer :: nwindows_x, ipos, jpos, window_size(2)
1071  integer :: istart, iend, jstart, jend
1072  integer :: isw, iew, jsw, jew
1073  integer :: omp_get_num_threads, window_id
1074  logical :: need_compute
1075  real(FMS_DATA_OVERRIDE_KIND_) :: lat_min, lat_max
1076  integer :: is_src, ie_src, js_src, je_src
1077  logical :: exists
1078  logical :: multifile !< use multiple consecutive files for override
1079  type(FmsNetcdfFile_t) :: fileobj
1080  integer :: startingi !< Starting x index for the compute domain relative to the input buffer
1081  integer :: endingi !< Ending x index for the compute domain relative to the input buffer
1082  integer :: startingj !< Starting y index for the compute domain relative to the input buffer
1083  integer :: endingj !< Ending y index for the compute domain relative to the input buffer
1084  integer :: nhalox !< Number of halos in the x direction
1085  integer :: nhaloy !< Number of halos in the y direction
1086  logical :: found_weight_file !< .True. if the weight file has already been read
1087  integer :: nglat !< Number of latitudes in the global domain
1088  integer :: nglon !< Number of longitudes in the global domain
1089 
1090  use_comp_domain = .false.
1091  if(.not.module_is_initialized) &
1092  call mpp_error(fatal,'Error: need to call data_override_init first')
1093 
1094 !1 Look for the data file in data_table
1095  if(PRESENT(override)) override = .false.
1096  if (present(data_index)) then
1097  index1 = data_index
1098  else
1099  index1 = -1
1100  do i = 1, table_size
1101  if( trim(gridname) /= trim(data_table(i)%gridname)) cycle
1102  if( trim(fieldname_code) /= trim(data_table(i)%fieldname_code)) cycle
1103  index1 = i ! field found
1104  exit
1105  enddo
1106  if(index1 .eq. -1) then
1107  if(debug_data_override) &
1108  call mpp_error(warning,'this field is NOT found in data_table: '//trim(fieldname_code))
1109  return ! NO override was performed
1110  endif
1111  endif
1112 
1113  fieldname = data_table(index1)%fieldname_file ! fieldname in netCDF data file
1114  factor = data_table(index1)%factor
1115  multifile = data_table(index1)%multifile
1116 
1117  if(fieldname == "") then
1118  return_data = factor
1119  if(PRESENT(override)) override = .true.
1120  return
1121  else
1122  filename = data_table(index1)%file_name
1123  if (filename == "") call mpp_error(fatal,'data_override: filename not given in data_table')
1124  if (multifile) prevfilename = data_table(index1)%prev_file_name
1125  if (multifile) nextfilename = data_table(index1)%next_file_name
1126  endif
1127 
1128  ongrid = (data_table(index1)%interpol_method == 'none')
1129 
1130 !3 Check if fieldname has been previously processed
1131 !$OMP CRITICAL
1132  curr_position = -1
1133  if(num_fields > 0 ) then
1134  do i = 1, num_fields
1135  if(trim(override_array(i)%gridname) /= trim(gridname)) cycle
1136  if(trim(override_array(i)%fieldname) /= trim(fieldname_code)) cycle
1137  curr_position = i
1138  exit
1139  enddo
1140  endif
1141 
1142  if(curr_position < 0) then ! the field has not been processed previously
1143  num_fields = num_fields + 1
1144  curr_position = num_fields
1145 
1146 ! Get working domain from model's gridname
1147  call get_domain(gridname,domain,comp_domain)
1148  call mpp_get_data_domain(domain, xsize=nxd, ysize=nyd)
1149  nxc = comp_domain(2)-comp_domain(1) + 1
1150  nyc = comp_domain(4)-comp_domain(3) + 1
1151 
1152 ! record fieldname, gridname in override_array
1153  override_array(curr_position)%fieldname = fieldname_code
1154  override_array(curr_position)%gridname = gridname
1155  override_array(curr_position)%comp_domain = comp_domain
1156 ! get number of threads
1157  override_array(curr_position)%numthreads = 1
1158 #if defined(_OPENMP)
1159  override_array(curr_position)%numthreads = omp_get_num_threads()
1160 #endif
1161 !--- data_override may be called from physics windows. The following are possible situations
1162 !--- 1. size(return_data,1) == nxd and size(return_data,2) == nyd
1163 !--- (on return_data domain and there is only one window).
1164 !--- 2. nxc is divisible by size(return_data,1), nyc is divisible by size(return_data,2),
1165 !--- nwindow = (nxc/size(return_data(1))*(nyc/size(return_data,2)),
1166 !--- also we require nwindows is divisible by nthreads.
1167 !--- The another restrition is that size(return_data,1) == ie_in - is_in + 1,
1168 !--- size(return_data,2) == je_in - js_in + 1
1169  nwindows = 1
1170  if( nxd == size(return_data,1) .AND. nyd == size(return_data,2) ) then !
1171  use_comp_domain = .false.
1172  else if ( mod(nxc, size(return_data,1)) ==0 .AND. mod(nyc, size(return_data,2)) ==0 ) then
1173  use_comp_domain = .true.
1174  nwindows = (nxc/size(return_data,1))*(nyc/size(return_data,2))
1175  else
1176  call mpp_error(fatal, &
1177  & "data_override: data is not on data domain and compute domain is not divisible by size(data)")
1178  endif
1179  override_array(curr_position)%window_size(1) = size(return_data,1)
1180  override_array(curr_position)%window_size(2) = size(return_data,2)
1181 
1182  window_size = override_array(curr_position)%window_size
1183  override_array(curr_position)%numwindows = nwindows
1184  if( mod(nwindows, override_array(curr_position)%numthreads) .NE. 0 ) then
1185  call mpp_error(fatal, "data_override: nwindow is not divisible by nthreads")
1186  endif
1187  allocate(override_array(curr_position)%need_compute(nwindows))
1188  override_array(curr_position)%need_compute = .true.
1189 
1190 !4 get index for time interp
1191  if(ongrid) then
1192  if( data_table(index1)%region_type .NE. no_region ) then
1193  call mpp_error(fatal,.NE.'data_override: ongrid must be false when region_type NO_REGION')
1194  endif
1195 
1196 ! Allow on-grid data_overrides on cubed sphere grid
1197  inquire(file=trim(filename),exist=exists)
1198  if (.not. exists) then
1199  call get_mosaic_tile_file(filename,filename2,.false.,domain)
1200  filename = filename2
1201  endif
1202 
1203  ! if using consecutive files for data_override, get file names
1204  if_multi3: if (multifile) then
1205  if_prev3: if (trim(prevfilename) /= '') then
1206  inquire(file=trim(prevfilename),exist=exists)
1207  if (.not. exists) then
1208  call get_mosaic_tile_file(prevfilename,prevfilename2,.false.,domain)
1209  prevfilename = prevfilename2
1210  endif
1211  endif if_prev3
1212  if_next3: if (trim(nextfilename) /= '') then
1213  inquire(file=trim(nextfilename),exist=exists)
1214  if (.not. exists) then
1215  call get_mosaic_tile_file(nextfilename,nextfilename2,.false.,domain)
1216  nextfilename = nextfilename2
1217  endif
1218  endif if_next3
1219  endif if_multi3
1220 
1221  !--- we always only pass data on compute domain
1222  id_time = init_external_field(filename,fieldname,domain=domain,verbose=debug_data_override, &
1223  use_comp_domain=use_comp_domain, nwindows=nwindows, ongrid=ongrid)
1224 
1225  ! if using consecutive files for data_override, get time axis for previous and next files
1226  ! and check spatial dims for consistency
1227  if_multi4: if (multifile) then
1228  id_time_prev = -1
1229  if_prev4:if (trim(prevfilename) /= '') then
1230  id_time_prev = init_external_field(prevfilename,fieldname,domain=domain, &
1231  verbose=debug_data_override,use_comp_domain=use_comp_domain, &
1232  nwindows = nwindows, ongrid=ongrid)
1233  dims = get_external_field_size(id_time)
1234  prev_dims = get_external_field_size(id_time_prev)
1235  ! check consistency of spatial dims
1236  if ((prev_dims(1) .ne. dims(1)) .or. (prev_dims(2) .ne. dims(2)) .or. &
1237  (prev_dims(3) .ne. dims(3))) then
1238  call mpp_error(fatal, 'data_override: dimensions mismatch between consecutive forcing files')
1239  endif
1240  allocate(data_table(index1)%time_prev_records(prev_dims(4)))
1241  call get_time_axis(id_time_prev,data_table(index1)%time_prev_records)
1242  endif if_prev4
1243  id_time_next = -1
1244  if_next4: if (trim(nextfilename) /= '') then
1245  id_time_next = init_external_field(nextfilename,fieldname,domain=domain, &
1246  verbose=debug_data_override,use_comp_domain=use_comp_domain, &
1247  nwindows = nwindows, ongrid=ongrid)
1248  dims = get_external_field_size(id_time)
1249  next_dims = get_external_field_size(id_time_next)
1250  ! check consistency of spatial dims
1251  if ((next_dims(1) .ne. dims(1)) .or. (next_dims(2) .ne. dims(2)) .or. &
1252  (next_dims(3) .ne. dims(3))) then
1253  call mpp_error(fatal, 'data_override: dimensions mismatch between consecutive forcing files')
1254  endif
1255  allocate(data_table(index1)%time_next_records(next_dims(4)))
1256  call get_time_axis(id_time_next,data_table(index1)%time_next_records)
1257  endif if_next4
1258  endif if_multi4
1259 
1260  dims = get_external_field_size(id_time)
1261  override_array(curr_position)%dims = dims
1262  if(id_time<0) call mpp_error(fatal,'data_override:field not found in init_external_field 1')
1263  override_array(curr_position)%t_index = id_time
1264  override_array(curr_position)%pt_index = id_time_prev
1265  override_array(curr_position)%nt_index = id_time_next
1266  else !ongrid=false
1267  id_time = init_external_field(filename,fieldname,domain=domain, axis_names=axis_names,&
1268  axis_sizes=axis_sizes, verbose=debug_data_override,override=.true.,use_comp_domain=use_comp_domain, &
1269  nwindows = nwindows)
1270 
1271  ! if using consecutive files for data_override, get time axis for previous and next files
1272  ! and check spatial dims for consistency
1273  if_multi5: if (multifile) then
1274  id_time_prev = -1
1275  if_prev5: if (trim(prevfilename) /= '') then
1276  id_time_prev = init_external_field(prevfilename,fieldname,domain=domain, axis_names=axis_names,&
1277  axis_sizes=axis_sizes, verbose=debug_data_override,override=.true.,use_comp_domain=use_comp_domain, &
1278  nwindows = nwindows)
1279  prev_dims = get_external_field_size(id_time_prev)
1280  allocate(data_table(index1)%time_prev_records(prev_dims(4)))
1281  call get_time_axis(id_time_prev,data_table(index1)%time_prev_records)
1282  endif if_prev5
1283  id_time_next = -1
1284  if_next5: if (trim(nextfilename) /= '') then
1285  id_time_next = init_external_field(nextfilename,fieldname,domain=domain, axis_names=axis_names,&
1286  axis_sizes=axis_sizes, verbose=debug_data_override,override=.true.,use_comp_domain=use_comp_domain, &
1287  nwindows = nwindows)
1288  next_dims = get_external_field_size(id_time_next)
1289  allocate(data_table(index1)%time_next_records(next_dims(4)))
1290  call get_time_axis(id_time_next,data_table(index1)%time_next_records)
1291  endif if_next5
1292  endif if_multi5
1293 
1294  dims = get_external_field_size(id_time)
1295  override_array(curr_position)%dims = dims
1296  if(id_time<0) call mpp_error(fatal,'data_override:field not found in init_external_field 2')
1297  override_array(curr_position)%t_index = id_time
1298  override_array(curr_position)%pt_index = id_time_prev
1299  override_array(curr_position)%nt_index = id_time_next
1300 
1301  ! get lon and lat of the input (source) grid, assuming that axis%data contains
1302  ! lat and lon of the input grid (in degrees)
1303 
1304  allocate(override_array(curr_position)%horz_interp(nwindows))
1305  allocate(override_array(curr_position)%lon_in(axis_sizes(1)+1))
1306  allocate(override_array(curr_position)%lat_in(axis_sizes(2)+1))
1307  if(get_external_fileobj(filename, fileobj)) then
1308  call axis_edges(fileobj, axis_names(1), override_array(curr_position)%lon_in, &
1309  reproduce_null_char_bug_flag=reproduce_null_char_bug)
1310  call axis_edges(fileobj, axis_names(2), override_array(curr_position)%lat_in, &
1311  reproduce_null_char_bug_flag=reproduce_null_char_bug)
1312  else
1313  call mpp_error(fatal,'data_override: file '//trim(filename)//' is not opened in time_interp_external')
1314  end if
1315 ! convert lon_in and lat_in from deg to radian
1316  override_array(curr_position)%lon_in = override_array(curr_position)%lon_in * real(deg_to_rad, lkind)
1317  override_array(curr_position)%lat_in = override_array(curr_position)%lat_in * real(deg_to_rad, lkind)
1318 
1319  !--- find the region of the source grid that cover the local model grid.
1320  !--- currently we only find the index range for j-direction because
1321  !--- of the cyclic condition in i-direction. The purpose of this is to
1322  !--- decrease the memory usage and increase the IO performance.
1323  select case(gridname)
1324  case('OCN')
1325  lon_local => lon_local_ocn; lat_local => lat_local_ocn
1326  case('ICE')
1327  lon_local => lon_local_ice; lat_local => lat_local_ice
1328  case('ATM')
1329  lon_local => lon_local_atm; lat_local => lat_local_atm
1330  case('LND')
1331  lon_local => lon_local_lnd; lat_local => lat_local_lnd
1332  case default
1333  call mpp_error(fatal,'error: gridname not recognized in data_override')
1334  end select
1335 
1336  lat_min = minval(lat_local)
1337  lat_max = maxval(lat_local)
1338  is_src = 1
1339  ie_src = axis_sizes(1)
1340  js_src = 1
1341  je_src = axis_sizes(2)
1342  ! do j = 1, axis_sizes(2)+1
1343  ! if( override_array(curr_position)%lat_in(j) > lat_min ) exit
1344  ! js_src = j
1345  ! enddo
1346  ! do j = 1, axis_sizes(2)+1
1347  ! je_src = j
1348  ! if( override_array(curr_position)%lat_in(j) > lat_max ) exit
1349  ! enddo
1350 
1351  !--- bicubic interpolation need one extra point in each direction. Also add
1352  !--- one more point for because lat_in is in the corner but the interpolation
1353  !--- use center points.
1354  select case (data_table(index1)%interpol_method)
1355  case ('bilinear')
1356  js_src = max(1, js_src-1)
1357  je_src = min(axis_sizes(2), je_src+1)
1358  case ('bicubic')
1359  js_src = max(1, js_src-2)
1360  je_src = min(axis_sizes(2), je_src+2)
1361  end select
1362  override_array(curr_position)%is_src = is_src
1363  override_array(curr_position)%ie_src = ie_src
1364  override_array(curr_position)%js_src = js_src
1365  override_array(curr_position)%je_src = je_src
1366  call reset_src_data_region(id_time, is_src, ie_src, js_src, je_src)
1367  if (multifile) then
1368  if (trim(prevfilename) /= '') then
1369  call reset_src_data_region(id_time_prev, is_src, ie_src, js_src, je_src)
1370  endif
1371  if (trim(nextfilename) /= '') then
1372  call reset_src_data_region(id_time_next, is_src, ie_src, js_src, je_src)
1373  endif
1374  endif
1375 
1376 ! Find the index of lon_start, lon_end, lat_start and lat_end in the input grid (nearest points)
1377  if( data_table(index1)%region_type .NE. no_region ) then
1378  allocate( lon_tmp(axis_sizes(1)), lat_tmp(axis_sizes(2)) )
1379  call read_data(fileobj, axis_names(1), lon_tmp)
1380  call read_data(fileobj, axis_names(2), lat_tmp)
1381  ! limit lon_start, lon_end are inside lon_in
1382  ! lat_start, lat_end are inside lat_in
1383  if(data_table(index1)%lon_start < lon_tmp(1) .OR. data_table(index1)%lon_start .GT. lon_tmp(axis_sizes(1)))&
1384  call mpp_error(fatal, "data_override: lon_start is outside lon_T")
1385  if( data_table(index1)%lon_end < lon_tmp(1) .OR. data_table(index1)%lon_end .GT. lon_tmp(axis_sizes(1))) &
1386  call mpp_error(fatal, "data_override: lon_end is outside lon_T")
1387  if(data_table(index1)%lat_start < lat_tmp(1) .OR. data_table(index1)%lat_start .GT. lat_tmp(axis_sizes(2)))&
1388  call mpp_error(fatal, "data_override: lat_start is outside lat_T")
1389  if( data_table(index1)%lat_end < lat_tmp(1) .OR. data_table(index1)%lat_end .GT. lat_tmp(axis_sizes(2))) &
1390  call mpp_error(fatal, "data_override: lat_end is outside lat_T")
1391  istart = nearest_index(data_table(index1)%lon_start, lon_tmp)
1392  iend = nearest_index(data_table(index1)%lon_end, lon_tmp)
1393  jstart = nearest_index(data_table(index1)%lat_start, lat_tmp)
1394  jend = nearest_index(data_table(index1)%lat_end, lat_tmp)
1395  ! adjust the index according to is_src and js_src
1396  istart = istart - is_src + 1
1397  iend = iend - is_src + 1
1398  jstart = jstart - js_src + 1
1399  jend = jend - js_src + 1
1400  call set_override_region(id_time, data_table(index1)%region_type, istart, iend, jstart, jend)
1401  if (multifile) then
1402  if (trim(prevfilename) /= '') then
1403  call set_override_region(id_time_prev, data_table(index1)%region_type, istart, iend, jstart, jend)
1404  endif
1405  if (trim(nextfilename) /= '') then
1406  call set_override_region(id_time_next, data_table(index1)%region_type, istart, iend, jstart, jend)
1407  endif
1408  endif
1409  deallocate(lon_tmp, lat_tmp)
1410  endif
1411 
1412  endif !ongrid
1413  else !curr_position >0
1414  dims = override_array(curr_position)%dims
1415  comp_domain = override_array(curr_position)%comp_domain
1416  nxc = comp_domain(2)-comp_domain(1) + 1
1417  nyc = comp_domain(4)-comp_domain(3) + 1
1418  is_src = override_array(curr_position)%is_src
1419  ie_src = override_array(curr_position)%ie_src
1420  js_src = override_array(curr_position)%js_src
1421  je_src = override_array(curr_position)%je_src
1422  window_size = override_array(curr_position)%window_size
1423  !---make sure data size match window_size
1424  if( window_size(1) .NE. size(return_data,1) .OR. window_size(2) .NE. size(return_data,2) ) then
1425  call mpp_error(fatal, "data_override: window_size does not match size(data)")
1426  endif
1427 !9 Get id_time previously stored in override_array
1428  id_time = override_array(curr_position)%t_index
1429  id_time_prev = override_array(curr_position)%pt_index
1430  id_time_next = override_array(curr_position)%nt_index
1431  endif
1432 !$OMP END CRITICAL
1433 
1434  if( override_array(curr_position)%numwindows > 1 ) then
1435  if( .NOT. PRESENT(is_in) .OR. .NOT. PRESENT(is_in) .OR. .NOT. PRESENT(is_in) .OR. .NOT. PRESENT(is_in) ) then
1436  call mpp_error(fatal, "data_override: is_in, ie_in, js_in, je_in must be present when nwindows > 1")
1437  endif
1438  endif
1439 
1440  isw = comp_domain(1)
1441  iew = comp_domain(2)
1442  jsw = comp_domain(3)
1443  jew = comp_domain(4)
1444  window_id = 1
1445  if( override_array(curr_position)%numwindows > 1 ) then
1446  nxc = comp_domain(2) - comp_domain(1) + 1
1447  nwindows_x = nxc/window_size(1)
1448  ipos = (is_in-1)/window_size(1) + 1
1449  jpos = (js_in-1)/window_size(2)
1450 
1451  window_id = jpos*nwindows_x + ipos
1452  isw = isw + is_in - 1
1453  iew = isw + ie_in - is_in
1454  jsw = jsw + js_in - 1
1455  jew = jsw + je_in - js_in
1456  endif
1457 
1458  if( ongrid ) then
1459  need_compute = .false.
1460  else
1461  !--- find the index for windows.
1462  need_compute=override_array(curr_position)%need_compute(window_id)
1463  endif
1464 
1465  !--- call horiz_interp_new is not initialized
1466 
1467  if( need_compute ) then
1468  select case(gridname)
1469  case('OCN')
1470  lon_local => lon_local_ocn; lat_local => lat_local_ocn
1471  case('ICE')
1472  lon_local => lon_local_ice; lat_local => lat_local_ice
1473  case('ATM')
1474  lon_local => lon_local_atm; lat_local => lat_local_atm
1475  case('LND')
1476  lon_local => lon_local_lnd; lat_local => lat_local_lnd
1477  case default
1478  call mpp_error(fatal,'error: gridname not recognized in data_override')
1479  end select
1480 
1481  if (data_table(index1)%ext_weights) then
1482  found_weight_file = .false.
1483  do i = 1, nweight_files
1484  if (external_weights(i)%weight_filename .eq. trim(data_table(index1)%ext_weights_file_name)) then
1485  override_array(curr_position)%horz_interp(window_id) = external_weights(i)%horiz_interp
1486  found_weight_file = .true.
1487  exit
1488  endif
1489  enddo
1490 
1491  if (.not. found_weight_file) then
1492  nweight_files = nweight_files + 1
1493  external_weights(nweight_files)%weight_filename = trim(data_table(index1)%ext_weights_file_name)
1494 
1495  call mpp_get_global_domain(domain, xsize=nglon, ysize=nglat)
1496  call horiz_interp_read_weights(external_weights(nweight_files)%horiz_interp, &
1497  external_weights(nweight_files)%weight_filename, &
1498  lon_local(isw:iew,jsw:jew), lat_local(isw:iew,jsw:jew), &
1499  override_array(curr_position)%lon_in(is_src:ie_src+1), &
1500  override_array(curr_position)%lat_in(js_src:je_src+1), &
1501  data_table(index1)%ext_weights_source, &
1502  data_table(index1)%interpol_method, isw, iew, jsw, jew, nglon, nglat)
1503 
1504  override_array(curr_position)%horz_interp(window_id) = external_weights(nweight_files)%horiz_interp
1505  endif
1506  else
1507  select case (data_table(index1)%interpol_method)
1508  case ('bilinear')
1509  call horiz_interp_new (override_array(curr_position)%horz_interp(window_id), &
1510  override_array(curr_position)%lon_in(is_src:ie_src+1), &
1511  override_array(curr_position)%lat_in(js_src:je_src+1), &
1512  lon_local(isw:iew,jsw:jew), lat_local(isw:iew,jsw:jew), interp_method="bilinear")
1513  case ('bicubic')
1514  call horiz_interp_new (override_array(curr_position)%horz_interp(window_id), &
1515  override_array(curr_position)%lon_in(is_src:ie_src+1), &
1516  override_array(curr_position)%lat_in(js_src:je_src+1), &
1517  lon_local(isw:iew,jsw:jew), lat_local(isw:iew,jsw:jew), interp_method="bicubic")
1518  end select
1519  endif
1520  override_array(curr_position)%need_compute(window_id) = .false.
1521  endif
1522 
1523  ! Determine if data in netCDF file is 2D or not
1524  data_file_is_2d = .false.
1525  if((dims(3) == 1) .and. (size(return_data,3)>1)) data_file_is_2d = .true.
1526 
1527  if(dims(3) .NE. 1 .and. (size(return_data,3) .NE. dims(3))) &
1528  call mpp_error(fatal, .NE..NE."data_override: dims(3) 1 and size(return_data,3) dims(3)")
1529 
1530 
1531  dims = get_external_field_size(id_time)
1532  if (.not. allocated(data_table(index1)%time_records)) allocate(data_table(index1)%time_records(dims(4)))
1533  call get_time_axis(id_time,data_table(index1)%time_records)
1534 
1535  first_record = data_table(index1)%time_records(1)
1536  last_record = data_table(index1)%time_records(dims(4))
1537 
1538  if(ongrid) then
1539  if (.not. use_comp_domain) then
1540  !< Determine the size of the halox and the part of `data` that is in the compute domain
1541  nhalox = (size(return_data,1) - nxc)/2
1542  nhaloy = (size(return_data,2) - nyc)/2
1543  startingi = lbound(return_data,1) + nhalox
1544  startingj = lbound(return_data,2) + nhaloy
1545  endingi = ubound(return_data,1) - nhalox
1546  endingj = ubound(return_data,2) - nhaloy
1547  end if
1548 
1549 !10 do time interp to get data in compute_domain
1550  if(data_file_is_2d) then
1551  if (use_comp_domain) then
1552 
1553  ! if using consecutive files, allow to perform time interpolation between the last record of previous
1554  ! file and first record of current file OR between the last record of current file and first record of
1555  ! next file hence "bridging" over files.
1556  if_multi6: if (multifile) then
1557  if_time6: if (time<first_record) then
1558  ! previous file must be init and time must be between last record of previous file and
1559  ! first record of current file
1560  if (id_time_prev<0) call mpp_error(fatal,'data_override:previous file needed with multifile')
1561  prev_dims = get_external_field_size(id_time_prev)
1562  if (time<data_table(index1)%time_prev_records(prev_dims(4))) call mpp_error(fatal, &
1563  'data_override: time_interp_external_bridge should only be called to bridge with previous file')
1564  ! bridge with previous file
1565  call time_interp_external_bridge(id_time_prev,id_time,time,return_data(:,:,1),verbose=debug_data_override, &
1566  is_in=is_in,ie_in=ie_in,js_in=js_in,je_in=je_in,window_id=window_id)
1567  elseif (time>last_record) then
1568  ! next file must be init and time must be between last record of current file and
1569  ! first record of next file
1570  if (id_time_next<0) call mpp_error(fatal,'data_override:next file needed with multifile')
1571  if (time>data_table(index1)%time_next_records(1)) call mpp_error(fatal, &
1572  'data_override: time_interp_external_bridge should only be called to bridge with next file')
1573  ! bridge with next file
1574  call time_interp_external_bridge(id_time,id_time_next,time,return_data(:,:,1),verbose=debug_data_override, &
1575  is_in=is_in,ie_in=ie_in,js_in=js_in,je_in=je_in,window_id=window_id)
1576  else ! first_record <= time <= last_record, do not use bridge
1577  call time_interp_external(id_time,time,return_data(:,:,1),verbose=debug_data_override, &
1578  is_in=is_in,ie_in=ie_in,js_in=js_in,je_in=je_in,window_id=window_id)
1579  endif if_time6
1580  else ! standard behavior
1581  call time_interp_external(id_time,time,return_data(:,:,1),verbose=debug_data_override, &
1582  is_in=is_in,ie_in=ie_in,js_in=js_in,je_in=je_in,window_id=window_id)
1583  endif if_multi6
1584 
1585  else
1586  !> If this in an ongrid case and you are not in the compute domain, send in `data` to be the correct
1587  !! size
1588 
1589  ! if using consecutive files, allow to perform time interpolation between the last record of previous
1590  ! file and first record of current file OR between the last record of current file and first record of
1591  ! next file hence "bridging" over files.
1592  if_multi7: if (multifile) then
1593  if_time7: if (time<first_record) then
1594  ! previous file must be init and time must be between last record of previous file and
1595  ! first record of current file
1596  if (id_time_prev<0) call mpp_error(fatal,'data_override:previous file needed with multifile')
1597  prev_dims = get_external_field_size(id_time_prev)
1598  if (time<data_table(index1)%time_prev_records(prev_dims(4))) call mpp_error(fatal, &
1599  'data_override: time_interp_external_bridge should only be called to bridge with previous file')
1600  ! bridge with previous file
1601  call time_interp_external_bridge(id_time_prev,id_time,time,&
1602  return_data(startingi:endingi,startingj:endingj,1), &
1603  verbose=debug_data_override,is_in=is_in,ie_in=ie_in, &
1604  js_in=js_in,je_in=je_in,window_id=window_id)
1605  elseif (time>last_record) then
1606  ! next file must be init and time must be between last record of current file and
1607  ! first record of next file
1608  if (id_time_next<0) call mpp_error(fatal,'data_override:next file needed with multifile')
1609  if (time>data_table(index1)%time_next_records(1)) call mpp_error(fatal, &
1610  'data_override: time_interp_external_bridge should only be called to bridge with next file')
1611  ! bridge with next file
1612  call time_interp_external_bridge(id_time,id_time_next,time,&
1613  return_data(startingi:endingi,startingj:endingj,1), &
1614  verbose=debug_data_override,is_in=is_in,ie_in=ie_in, &
1615  js_in=js_in,je_in=je_in,window_id=window_id)
1616  else ! first_record <= time <= last_record, do not use bridge
1617  call time_interp_external(id_time,time,return_data(startingi:endingi,startingj:endingj,1), &
1618  verbose=debug_data_override,is_in=is_in,ie_in=ie_in, &
1619  js_in=js_in,je_in=je_in,window_id=window_id)
1620  endif if_time7
1621  else ! standard behavior
1622  call time_interp_external(id_time,time,return_data(startingi:endingi,startingj:endingj,1), &
1623  verbose=debug_data_override,is_in=is_in,ie_in=ie_in, &
1624  js_in=js_in,je_in=je_in,window_id=window_id)
1625  endif if_multi7
1626 
1627  end if
1628  return_data(:,:,1) = return_data(:,:,1)*factor
1629  do i = 2, size(return_data,3)
1630  return_data(:,:,i) = return_data(:,:,1)
1631  end do
1632  else
1633  if (use_comp_domain) then
1634 
1635  ! if using consecutive files, allow to perform time interpolation between the last record of previous
1636  ! file and first record of current file OR between the last record of current file and first record of
1637  ! next file hence "bridging" over files.
1638  if_multi8: if (multifile) then
1639  if_time8: if (time<first_record) then
1640  if (id_time_prev<0) call mpp_error(fatal,'data_override:previous file needed with multifile')
1641  prev_dims = get_external_field_size(id_time_prev)
1642  if (time<data_table(index1)%time_prev_records(prev_dims(4))) call mpp_error(fatal, &
1643  'data_override: time_interp_external_bridge should only be called to bridge with previous file')
1644  call time_interp_external_bridge(id_time_prev,id_time,time,return_data,verbose=debug_data_override, &
1645  is_in=is_in,ie_in=ie_in,js_in=js_in,je_in=je_in,window_id=window_id)
1646  elseif (time>last_record) then
1647  if (id_time_next<0) call mpp_error(fatal,'data_override:next file needed with multifile')
1648  if (time>data_table(index1)%time_next_records(1)) call mpp_error(fatal, &
1649  'data_override: time_interp_external_bridge should only be called to bridge with next file')
1650  call time_interp_external_bridge(id_time,id_time_next,time,return_data,verbose=debug_data_override, &
1651  is_in=is_in,ie_in=ie_in,js_in=js_in,je_in=je_in,window_id=window_id)
1652  else ! first_record <= time <= last_record, do not use bridge
1653  call time_interp_external(id_time,time,return_data,verbose=debug_data_override, &
1654  is_in=is_in,ie_in=ie_in,js_in=js_in,je_in=je_in,window_id=window_id)
1655  endif if_time8
1656  else ! standard behavior
1657  call time_interp_external(id_time,time,return_data,verbose=debug_data_override, &
1658  is_in=is_in,ie_in=ie_in,js_in=js_in,je_in=je_in,window_id=window_id)
1659  endif if_multi8
1660 
1661  else
1662  !> If this in an ongrid case and you are not in the compute domain, send in `data` to be the correct
1663  !! size
1664 
1665  ! if using consecutive files, allow to perform time interpolation between the last record of previous
1666  ! file and first record of current file OR between the last record of current file and first record of
1667  ! next file hence "bridging" over files.
1668  if_multi9: if (multifile) then
1669  if_time9: if (time<first_record) then
1670  if (id_time_prev<0) call mpp_error(fatal,'data_override:previous file needed with multifile')
1671  prev_dims = get_external_field_size(id_time_prev)
1672  if (time<data_table(index1)%time_prev_records(prev_dims(4))) call mpp_error(fatal, &
1673  'data_override: time_interp_external_bridge should only be called to bridge with previous file')
1674  call time_interp_external_bridge(id_time_prev,id_time,time,&
1675  return_data(startingi:endingi,startingj:endingj,:), &
1676  verbose=debug_data_override,is_in=is_in,ie_in=ie_in, &
1677  js_in=js_in,je_in=je_in,window_id=window_id)
1678  elseif (time>last_record) then
1679  if (id_time_next<0) call mpp_error(fatal,'data_override:next file needed with multifile')
1680  if (time>data_table(index1)%time_next_records(1)) call mpp_error(fatal, &
1681  'data_override: time_interp_external_bridge should only be called to bridge with next file')
1682  call time_interp_external_bridge(id_time,id_time_next,time,&
1683  return_data(startingi:endingi,startingj:endingj,:), &
1684  verbose=debug_data_override,is_in=is_in,ie_in=ie_in, &
1685  js_in=js_in,je_in=je_in,window_id=window_id)
1686  else ! first_record <= time <= last_record, do not use bridge
1687  call time_interp_external(id_time,time,return_data(startingi:endingi,startingj:endingj,:), &
1688  verbose=debug_data_override,is_in=is_in,ie_in=ie_in, &
1689  js_in=js_in,je_in=je_in,window_id=window_id)
1690  endif if_time9
1691  else ! standard behavior
1692  call time_interp_external(id_time,time,return_data(startingi:endingi,startingj:endingj,:), &
1693  verbose=debug_data_override,is_in=is_in,ie_in=ie_in, &
1694  js_in=js_in,je_in=je_in,window_id=window_id)
1695  endif if_multi9
1696 
1697  end if
1698  return_data = return_data*factor
1699  endif
1700  else ! off grid case
1701 ! do time interp to get global data
1702  if(data_file_is_2d) then
1703  if( data_table(index1)%region_type == no_region ) then
1704 
1705  ! if using consecutive files, allow to perform time interpolation between the last record of previous
1706  ! file and first record of current file OR between the last record of current file and first record of
1707  ! next file hence "bridging" over files.
1708  if_multi10: if (multifile) then
1709  if_time10: if (time<first_record) then
1710  if (id_time_prev<0) call mpp_error(fatal,'data_override:previous file needed with multifile')
1711  prev_dims = get_external_field_size(id_time_prev)
1712  if (time<data_table(index1)%time_prev_records(prev_dims(4))) call mpp_error(fatal, &
1713  'data_override: time_interp_external_bridge should only be called to bridge with previous file')
1714  call time_interp_external_bridge(id_time_prev,id_time,time,return_data(:,:,1), &
1715  verbose=debug_data_override, &
1716  horz_interp=override_array(curr_position)%horz_interp(window_id), &
1717  is_in=is_in,ie_in=ie_in,js_in=js_in,je_in=je_in,window_id=window_id)
1718  elseif (time>last_record) then
1719  if (id_time_next<0) call mpp_error(fatal,'data_override:next file needed with multifile')
1720  if (time>data_table(index1)%time_next_records(1)) call mpp_error(fatal, &
1721  'data_override: time_interp_external_bridge should only be called to bridge with next file')
1722  call time_interp_external_bridge(id_time,id_time_next,time,return_data(:,:,1), &
1723  verbose=debug_data_override, &
1724  horz_interp=override_array(curr_position)%horz_interp(window_id), &
1725  is_in=is_in,ie_in=ie_in,js_in=js_in,je_in=je_in,window_id=window_id)
1726  else ! first_record <= time <= last_record, do not use bridge
1727  call time_interp_external(id_time,time,return_data(:,:,1),verbose=debug_data_override, &
1728  horz_interp=override_array(curr_position)%horz_interp(window_id), &
1729  is_in=is_in,ie_in=ie_in,js_in=js_in,je_in=je_in,window_id=window_id)
1730  endif if_time10
1731  else ! standard behavior
1732  call time_interp_external(id_time,time,return_data(:,:,1),verbose=debug_data_override, &
1733  horz_interp=override_array(curr_position)%horz_interp(window_id), &
1734  is_in=is_in,ie_in=ie_in,js_in=js_in,je_in=je_in,window_id=window_id)
1735  endif if_multi10
1736 
1737  return_data(:,:,1) = return_data(:,:,1)*factor
1738  do i = 2, size(return_data,3)
1739  return_data(:,:,i) = return_data(:,:,1)
1740  enddo
1741  else
1742  allocate(mask_out(size(return_data,1), size(return_data,2),1))
1743  mask_out = .false.
1744 
1745  ! if using consecutive files, allow to perform time interpolation between the last record of previous
1746  ! file and first record of current file OR between the last record of current file and first record of
1747  ! next file hence "bridging" over files.
1748  if_multi11: if (multifile) then
1749  if_time11: if (time<first_record) then
1750  if (id_time_prev<0) call mpp_error(fatal,'data_override:previous file needed with multifile')
1751  prev_dims = get_external_field_size(id_time_prev)
1752  if (time<data_table(index1)%time_prev_records(prev_dims(4))) call mpp_error(fatal, &
1753  'data_override: time_interp_external_bridge should only be called to bridge with previous file')
1754  call time_interp_external_bridge(id_time_prev,id_time,time,return_data(:,:,1), &
1755  verbose=debug_data_override, &
1756  horz_interp=override_array(curr_position)%horz_interp(window_id), &
1757  mask_out =mask_out(:,:,1), &
1758  is_in=is_in,ie_in=ie_in,js_in=js_in,je_in=je_in,window_id=window_id)
1759  elseif (time>last_record) then
1760  if (id_time_next<0) call mpp_error(fatal,'data_override:next file needed with multifile')
1761  if (time>data_table(index1)%time_next_records(1)) call mpp_error(fatal, &
1762  'data_override: time_interp_external_bridge should only be called to bridge with next file')
1763  call time_interp_external_bridge(id_time,id_time_next,time,return_data(:,:,1), &
1764  verbose=debug_data_override, &
1765  horz_interp=override_array(curr_position)%horz_interp(window_id), &
1766  mask_out =mask_out(:,:,1), &
1767  is_in=is_in,ie_in=ie_in,js_in=js_in,je_in=je_in,window_id=window_id)
1768  else ! first_record <= time <= last_record, do not use bridge
1769  call time_interp_external(id_time,time,return_data(:,:,1),verbose=debug_data_override, &
1770  horz_interp=override_array(curr_position)%horz_interp(window_id), &
1771  mask_out =mask_out(:,:,1), &
1772  is_in=is_in,ie_in=ie_in,js_in=js_in,je_in=je_in,window_id=window_id)
1773  endif if_time11
1774  else ! standard behavior
1775  call time_interp_external(id_time,time,return_data(:,:,1),verbose=debug_data_override, &
1776  horz_interp=override_array(curr_position)%horz_interp(window_id), &
1777  mask_out =mask_out(:,:,1), &
1778  is_in=is_in,ie_in=ie_in,js_in=js_in,je_in=je_in,window_id=window_id)
1779  endif if_multi11
1780 
1781  where(mask_out(:,:,1))
1782  return_data(:,:,1) = return_data(:,:,1)*factor
1783  end where
1784  do i = 2, size(return_data,3)
1785  where(mask_out(:,:,1))
1786  return_data(:,:,i) = return_data(:,:,1)
1787  end where
1788  enddo
1789  deallocate(mask_out)
1790  endif
1791  else
1792  if( data_table(index1)%region_type == no_region ) then
1793 
1794  ! if using consecutive files, allow to perform time interpolation between the last record of previous
1795  ! file and first record of current file OR between the last record of current file and first record of
1796  ! next file hence "bridging" over files.
1797  if_multi12: if (multifile) then
1798  if_time12: if (time<first_record) then
1799  if (id_time_prev<0) call mpp_error(fatal,'data_override:previous file needed with multifile')
1800  prev_dims = get_external_field_size(id_time_prev)
1801  if (time<data_table(index1)%time_prev_records(prev_dims(4))) call mpp_error(fatal, &
1802  'data_override: time_interp_external_bridge should only be called to bridge with previous file')
1803  call time_interp_external_bridge(id_time_prev,id_time,time,return_data,verbose=debug_data_override, &
1804  horz_interp=override_array(curr_position)%horz_interp(window_id), &
1805  is_in=is_in,ie_in=ie_in,js_in=js_in,je_in=je_in,window_id=window_id)
1806  elseif (time>last_record) then
1807  if (id_time_next<0) call mpp_error(fatal,'data_override:next file needed with multifile')
1808  if (time>data_table(index1)%time_next_records(1)) call mpp_error(fatal, &
1809  'data_override: time_interp_external_bridge should only be called to bridge with next file')
1810  call time_interp_external_bridge(id_time,id_time_next,time,return_data,verbose=debug_data_override, &
1811  horz_interp=override_array(curr_position)%horz_interp(window_id), &
1812  is_in=is_in,ie_in=ie_in,js_in=js_in,je_in=je_in,window_id=window_id)
1813  else ! first_record <= time <= last_record, do not use bridge
1814  call time_interp_external(id_time,time,return_data,verbose=debug_data_override, &
1815  horz_interp=override_array(curr_position)%horz_interp(window_id), &
1816  is_in=is_in,ie_in=ie_in,js_in=js_in,je_in=je_in,window_id=window_id)
1817  endif if_time12
1818  else ! standard behavior
1819  call time_interp_external(id_time,time,return_data,verbose=debug_data_override, &
1820  horz_interp=override_array(curr_position)%horz_interp(window_id), &
1821  is_in=is_in,ie_in=ie_in,js_in=js_in,je_in=je_in,window_id=window_id)
1822  endif if_multi12
1823 
1824  return_data = return_data*factor
1825  else
1826  allocate(mask_out(size(return_data,1), size(return_data,2), size(return_data,3)) )
1827  mask_out = .false.
1828 
1829  ! if using consecutive files, allow to perform time interpolation between the last record of previous
1830  ! file and first record of current file OR between the last record of current file and first record of
1831  ! next file hence "bridging" over files.
1832  if_multi13: if (multifile) then
1833  if_time13: if (time<first_record) then
1834  if (id_time_prev<0) call mpp_error(fatal,'data_override:previous file needed with multifile')
1835  prev_dims = get_external_field_size(id_time_prev)
1836  if (time<data_table(index1)%time_prev_records(prev_dims(4))) call mpp_error(fatal, &
1837  'data_override: time_interp_external_bridge should only be called to bridge with previous file')
1838  call time_interp_external_bridge(id_time_prev,id_time,time,return_data,verbose=debug_data_override, &
1839  horz_interp=override_array(curr_position)%horz_interp(window_id), &
1840  mask_out =mask_out, &
1841  is_in=is_in,ie_in=ie_in,js_in=js_in,je_in=je_in,window_id=window_id)
1842  elseif (time>last_record) then
1843  if (id_time_next<0) call mpp_error(fatal,'data_override:next file needed with multifile')
1844  if (time>data_table(index1)%time_next_records(1)) call mpp_error(fatal, &
1845  'data_override: time_interp_external_bridge should only be called to bridge with next file')
1846  call time_interp_external_bridge(id_time,id_time_next,time,return_data,verbose=debug_data_override, &
1847  horz_interp=override_array(curr_position)%horz_interp(window_id), &
1848  mask_out =mask_out, &
1849  is_in=is_in,ie_in=ie_in,js_in=js_in,je_in=je_in,window_id=window_id)
1850  else ! first_record <= time <= last_record, do not use bridge
1851  call time_interp_external(id_time,time,return_data,verbose=debug_data_override, &
1852  horz_interp=override_array(curr_position)%horz_interp(window_id), &
1853  mask_out =mask_out, &
1854  is_in=is_in,ie_in=ie_in,js_in=js_in,je_in=je_in,window_id=window_id)
1855  endif if_time13
1856  else ! standard behavior
1857  call time_interp_external(id_time,time,return_data,verbose=debug_data_override, &
1858  horz_interp=override_array(curr_position)%horz_interp(window_id), &
1859  mask_out =mask_out, &
1860  is_in=is_in,ie_in=ie_in,js_in=js_in,je_in=je_in,window_id=window_id)
1861  endif if_multi13
1862 
1863  where(mask_out)
1864  return_data = return_data*factor
1865  end where
1866  deallocate(mask_out)
1867  endif
1868  endif
1869 
1870  endif
1871 
1872  if(PRESENT(override)) override = .true.
1873 end subroutine data_override_3d_
1874 
1875 !> @brief Data override for 1D unstructured grids
1876 subroutine data_override_ug_1d_(gridname,fieldname,return_data,time,override)
1877  character(len=3), intent(in) :: gridname !< model grid ID
1878  character(len=*), intent(in) :: fieldname !< field to override
1879  real(FMS_DATA_OVERRIDE_KIND_), dimension(:), intent(inout) :: return_data !< data returned by this call
1880  type(time_type), intent(in) :: time !< model time
1881  logical, intent(out), optional :: override !< true if the field has been overridden succesfully
1882  !local vars
1883  real(FMS_DATA_OVERRIDE_KIND_), dimension(:,:), allocatable :: data_SG
1884  type(domainUG) :: UG_domain
1885  integer :: index1
1886  integer :: i
1887  integer, dimension(4) :: comp_domain = 0 !< istart,iend,jstart,jend for compute domain
1888 
1889  !1 Look for the data file in data_table
1890  if(PRESENT(override)) override = .false.
1891  index1 = -1
1892  do i = 1, table_size
1893  if( trim(gridname) /= trim(data_table(i)%gridname)) cycle
1894  if( trim(fieldname) /= trim(data_table(i)%fieldname_code)) cycle
1895  index1 = i ! field found
1896  exit
1897  enddo
1898  if(index1 .eq. -1) return ! NO override was performed
1899 
1900  call get_domainug(gridname,ug_domain,comp_domain)
1901  allocate(data_sg(comp_domain(1):comp_domain(2),comp_domain(3):comp_domain(4)))
1902 
1903  call data_override_2d_(gridname,fieldname,data_sg,time,override)
1904 
1905  call mpp_pass_sg_to_ug(ug_domain, data_sg(:,:), return_data(:))
1906 
1907  deallocate(data_sg)
1908 end subroutine data_override_ug_1d_
1909 
1910 !> @brief Data override for 2D unstructured grids
1911 subroutine data_override_ug_2d_(gridname,fieldname,return_data,time,override)
1912  character(len=3), intent(in) :: gridname !< model grid ID
1913  character(len=*), intent(in) :: fieldname !< field to override
1914  real(FMS_DATA_OVERRIDE_KIND_), dimension(:,:), intent(inout) :: return_data !< data returned by this call
1915  type(time_type), intent(in) :: time !< model time
1916  logical, intent(out), optional :: override !< true if the field has been overridden succesfully
1917  !local vars
1918  real(FMS_DATA_OVERRIDE_KIND_), dimension(:,:,:), allocatable :: data_SG
1919  real(FMS_DATA_OVERRIDE_KIND_), dimension(:,:), allocatable :: data_UG
1920  type(domainUG) :: UG_domain
1921  integer :: index1
1922  integer :: i, nlevel, nlevel_max
1923  integer, dimension(4) :: comp_domain = 0 !< istart,iend,jstart,jend for compute domain
1924 
1925 !1 Look for the data file in data_table
1926  if(PRESENT(override)) override = .false.
1927  index1 = -1
1928  do i = 1, table_size
1929  if( trim(gridname) /= trim(data_table(i)%gridname)) cycle
1930  if( trim(fieldname) /= trim(data_table(i)%fieldname_code)) cycle
1931  index1 = i ! field found
1932  exit
1933  enddo
1934  if(index1 .eq. -1) return ! NO override was performed
1935 
1936  nlevel = size(return_data,2)
1937  nlevel_max = nlevel
1938  call mpp_max(nlevel_max)
1939 
1940  call get_domainug(gridname,ug_domain,comp_domain)
1941  allocate(data_sg(comp_domain(1):comp_domain(2),comp_domain(3):comp_domain(4),nlevel_max))
1942  allocate(data_ug(size(return_data,1), nlevel_max))
1943  data_sg = 0._lkind
1944  call data_override_3d_(gridname,fieldname,data_sg,time,override)
1945 
1946  call mpp_pass_sg_to_ug(ug_domain, data_sg(:,:,:), data_ug(:,:))
1947  return_data(:,1:nlevel) = data_ug(:,1:nlevel)
1948 
1949  deallocate(data_sg, data_ug)
1950 end subroutine data_override_ug_2d_
Perform 1D interpolation between grids.
Definition: axis_utils2.F90:60
subroutine, public fms2_io_init()
Reads the fms2_io_nml.
Definition: fms2_io.F90:396
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
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:524
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:702
character(:) function, allocatable, public string(v, fmt)
Converts a number or a Boolean value to a string.
subroutine, public horiz_interp_init
Initialize module and writes version number to logfile.out.
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....
Passes data from a structured grid to an unstructured grid Example usage:
The domain2D type contains all the necessary information to define the global, compute and data domai...
Domain information for managing data on unstructured grids.
integer function stdout()
This function returns the current standard fortran unit numbers for output.
Definition: mpp_util.inc:43
integer function stdlog()
This function returns the current standard fortran unit numbers for log messages. Log messages,...
Definition: mpp_util.inc:59
Error handler.
Definition: mpp.F90:382
Reduction operations. Find the max of scalar a from the PEs in pelist result is also automatically br...
Definition: mpp.F90:538
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.
subroutine, public time_interp_external_init()
Initialize the time_interp_external_mod module.
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, public get_time_axis(index, time)
Provide data from external file interpolated to current model time. Data may be local to current proc...
Type to represent amounts of time. Implemented as seconds and days to allow for larger intervals.
integer function, public open_and_parse_file(filename)
Opens and parses a yaml file.
subroutine, public get_block_ids(file_id, block_name, block_ids, parent_block_id)
Gets the ids of the blocks with block_name in the yaml file If parent_block_id is present,...
integer function, public get_num_blocks(file_id, block_name, parent_block_id)
Determines the number of blocks with block_name in the yaml file If parent_block_id is present,...
Dermine the value of a key from a keyname.
Definition: yaml_parser.F90:60
Allocates space and initializes a derived-type variable that contains pre-computed interpolation indi...
Subroutines for reading in weight files and using that to fill in the horiz_interp type instead calcu...