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