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