FV3 Bundle
data_override.F90
Go to the documentation of this file.
1 !***********************************************************************
2 !* GNU Lesser General Public License
3 !*
4 !* This file is part of the GFDL Flexible Modeling System (FMS).
5 !*
6 !* FMS is free software: you can redistribute it and/or modify it under
7 !* the terms of the GNU Lesser General Public License as published by
8 !* the Free Software Foundation, either version 3 of the License, or (at
9 !* your option) any later version.
10 !*
11 !* FMS is distributed in the hope that it will be useful, but WITHOUT
12 !* ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
13 !* FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
14 !* for more details.
15 !*
16 !* You should have received a copy of the GNU Lesser General Public
17 !* License along with FMS. If not, see <http://www.gnu.org/licenses/>.
18 !***********************************************************************
19 
20 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
21 !! !!
22 !! GNU General Public License !!
23 !! !!
24 !! This file is part of the Flexible Modeling System (FMS). !!
25 !! !!
26 !! FMS is free software; you can redistribute it and/or modify !!
27 !! it and are expected to follow the terms of the GNU General Public !!
28 !! License as published by the Free Software Foundation. !!
29 !! !!
30 !! FMS is distributed in the hope that it will be useful, !!
31 !! but WITHOUT ANY WARRANTY; without even the implied warranty of !!
32 !! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the !!
33 !! GNU General Public License for more details. !!
34 !! !!
35 !! You should have received a copy of the GNU General Public License !!
36 !! along with FMS; if not, write to: !!
37 !! Free Software Foundation, Inc. !!
38 !! 59 Temple Place, Suite 330 !!
39 !! Boston, MA 02111-1307 USA !!
40 !! or see: !!
41 !! http://www.gnu.org/licenses/gpl.txt !!
42 !! !!
43 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
45 !
46 ! <CONTACT EMAIL="Zhi.Liang@noaa.gov">
47 ! Z. Liang
48 ! </CONTACT>
49 !
50 ! <CONTACT EMAIL="Matthew.Harrison@noaa.gov">
51 ! M.J. Harrison
52 ! </CONTACT>
53 !
54 ! <CONTACT EMAIL="Michael.Winton@noaa.gov">
55 ! M. Winton
56 ! </CONTACT>
57 
58 !<OVERVIEW>
59 ! Given a gridname, fieldname and model time this routine will get data in a file whose
60 ! path is described in a user-provided data_table, do spatial and temporal interpolation if
61 ! necessary to convert data to model's grid and time.
62 !
63 ! Before using data_override a data_table must be created with the following entries:
64 ! gridname, fieldname_code, fieldname_file, file_name, ongrid, factor.
65 !
66 ! More explainations about data_table entries can be found in the source code (defining data_type)
67 !
68 ! If user wants to override fieldname_code with a const, set fieldname_file in data_table = ""
69 ! and factor = const
70 !
71 ! If user wants to override fieldname_code with data from a file, set fieldname_file = name in
72 ! the netCDF data file, factor then will be for unit conversion (=1 if no conversion required)
73 !
74 ! A field can be overriden globally (by default) or users can specify one or two regions in which
75 ! data_override will take place, field values outside the region will not be affected.
76 !</OVERVIEW>
77 #include <fms_platform.h>
78 use platform_mod, only: r8_kind
79 use constants_mod, only: pi
80 use mpp_io_mod, only: axistype,mpp_close,mpp_open,mpp_get_axis_data,mpp_rdonly,mpp_ascii
81 use mpp_mod, only : mpp_error,fatal,warning,mpp_pe,stdout,stdlog,mpp_root_pe, note, mpp_min, mpp_max, mpp_chksum
82 use mpp_mod, only : input_nml_file
83 use horiz_interp_mod, only : horiz_interp_init, horiz_interp_new, horiz_interp_type, &
84  assignment(=), horiz_interp_del
90 use fms_mod, only: write_version_number, field_exist, lowercase, file_exist, open_namelist_file, check_nml_error, close_file
92 use mpp_domains_mod, only : domain2d, mpp_get_compute_domain, null_domain2d,operator(.NE.),operator(.EQ.)
96 use mpp_domains_mod, only : domainug, mpp_pass_sg_to_ug, mpp_get_ug_sg_domain, null_domainug
97 
98 use time_manager_mod, only: time_type
99 
100 implicit none
101 private
102 
103 ! Include variable "version" to be written to log file.
104 #include<file_version.h>
105 
107  character(len=3) :: gridname
108  character(len=128) :: fieldname_code !fieldname used in user's code (model)
109  character(len=128) :: fieldname_file ! fieldname used in the netcdf data file
110  character(len=512) :: file_name ! name of netCDF data file
111  character(len=128) :: interpol_method ! interpolation method (default "bilinear")
112  real :: factor ! For unit conversion, default=1, see OVERVIEW above
113  real :: lon_start, lon_end, lat_start, lat_end
114  integer :: region_type
115 end type data_type
116 
117 
119  character(len=3) :: gridname
120  character(len=128) :: fieldname
121  integer :: t_index !index for time interp
122  type(horiz_interp_type), pointer :: horz_interp(:) =>null() ! index for horizontal spatial interp
123  integer :: dims(4) ! dimensions(x,y,z,t) of the field in filename
124  integer :: comp_domain(4) ! istart,iend,jstart,jend for compute domain
125  integer :: numthreads
126  real, _allocatable :: lon_in(:) _null
127  real, _allocatable :: lat_in(:) _null
128  logical, _allocatable :: need_compute(:) _null
129  integer :: numwindows
130  integer :: window_size(2)
131  integer :: is_src, ie_src, js_src, je_src
132 end type override_type
133 
134  integer, parameter :: max_table=100, max_array=100
135  integer :: table_size ! actual size of data table
136  integer, parameter :: annual=1, monthly=2, daily=3, hourly=4, undef=-1
137  real, parameter :: tpi=2*pi
139  logical :: module_is_initialized = .false.
140 
142 type(domainug),save :: lnd_domainug
143 
144 real, dimension(:,:), target, allocatable :: lon_local_ocn, lat_local_ocn
145 real, dimension(:,:), target, allocatable :: lon_local_atm, lat_local_atm
146 real, dimension(:,:), target, allocatable :: lon_local_ice, lat_local_ice
147 real, dimension(:,:), target, allocatable :: lon_local_lnd, lat_local_lnd
152 integer:: num_fields = 0 ! number of fields in override_array already processed
153 type(data_type), dimension(max_table) :: data_table ! user-provided data table
155 type(override_type), dimension(max_array), save :: override_array ! to store processed fields
157 logical :: atm_on, ocn_on, lnd_on, ice_on
158 logical :: lndug_on
160 logical :: grid_center_bug = .false.
161 
162 namelist /data_override_nml/ debug_data_override, grid_center_bug
163 
164 interface data_override
165  module procedure data_override_0d
166  module procedure data_override_2d
167  module procedure data_override_3d
168 end interface
169 
171  module procedure data_override_ug_1d
172  module procedure data_override_ug_2d
173 end interface
174 
176 public :: data_override_ug
177 
178 contains
179 !===============================================================================================
180 ! <SUBROUTINE NAME="data_override_init">
181 ! <DESCRIPTION>
182 ! Assign default values for default_table, get domain of component models,
183 ! get global grids of component models.
184 ! Users should call data_override_init before calling data_override
185 ! </DESCRIPTION>
186 ! <TEMPLATE>
187 ! call data_override_init
188 ! </TEMPLATE>
189 subroutine data_override_init(Atm_domain_in, Ocean_domain_in, Ice_domain_in, Land_domain_in, Land_domainUG_in)
190  type(domain2d), intent(in), optional :: atm_domain_in
191  type(domain2d), intent(in), optional :: ocean_domain_in, ice_domain_in
192  type(domain2d), intent(in), optional :: land_domain_in
193  type(domainug) , intent(in), optional :: land_domainug_in
194 
195 ! <NOTE>
196 ! This subroutine should be called in coupler_init after
197 ! (ocean/atmos/land/ice)_model_init have been called.
198 !
199 ! data_override_init can be called more than once, in one call the user can pass
200 ! up to 4 domains of component models, at least one domain must be present in
201 ! any call
202 !
203 ! Data_table is initialized here with default values. Users should provide "real" values
204 ! that will override the default values. Real values can be given using data_table, each
205 ! line of data_table contains one data_entry. Items of data_entry are comma separated.
206 !
207 ! </NOTE>
208  character(len=128) :: grid_file = 'INPUT/grid_spec.nc'
209  integer :: is,ie,js,je,count
210  integer :: i, iunit, ntable, ntable_lima, ntable_new, unit,io_status, ierr
211  character(len=256) :: record
212  logical :: file_open
213  logical :: ongrid
214  character(len=128) :: region, region_type
215 
216 
217  type(data_type) :: data_entry
218 
219  debug_data_override = .false.
220 
221 #ifdef INTERNAL_FILE_NML
222  read (input_nml_file, data_override_nml, iostat=io_status)
223  ierr = check_nml_error(io_status, 'data_override_nml')
224 #else
225  iunit = open_namelist_file()
226  ierr=1; do while (ierr /= 0)
227  read (iunit, nml=data_override_nml, iostat=io_status, end=10)
228  ierr = check_nml_error(io_status, 'data_override_nml')
229  enddo
230 10 call close_file (iunit)
231 #endif
232  unit = stdlog()
233  write(unit, data_override_nml)
234 
235 ! if(module_is_initialized) return
236 
237  atm_on = PRESENT(atm_domain_in)
238  ocn_on = PRESENT(ocean_domain_in)
239  lnd_on = PRESENT(land_domain_in)
240  ice_on = PRESENT(ice_domain_in)
241  lndug_on = PRESENT(land_domainug_in)
242  if(.not. module_is_initialized) then
243  atm_domain = null_domain2d
244  ocn_domain = null_domain2d
245  lnd_domain = null_domain2d
246  ice_domain = null_domain2d
248  end if
249  if (atm_on) atm_domain = atm_domain_in
250  if (ocn_on) ocn_domain = ocean_domain_in
251  if (lnd_on) lnd_domain = land_domain_in
252  if (ice_on) ice_domain = ice_domain_in
253  if (lndug_on) lnd_domainug = land_domainug_in
254 
255  if(.not. module_is_initialized) then
256  call horiz_interp_init
257  radian_to_deg = 180./pi
258  deg_to_radian = pi/180.
259 
260  call write_version_number("DATA_OVERRIDE_MOD", version)
261 
262 ! Initialize user-provided data table
263  default_table%gridname = 'none'
264  default_table%fieldname_code = 'none'
265  default_table%fieldname_file = 'none'
266  default_table%file_name = 'none'
267  default_table%factor = 1.
268  default_table%interpol_method = 'bilinear'
269  do i = 1,max_table
271  enddo
272 
273 ! Read coupler_table
274  call mpp_open(iunit, 'data_table', action=mpp_rdonly)
275  ntable = 0
276  ntable_lima = 0
277  ntable_new = 0
278 
279  do while (ntable <= max_table)
280  read(iunit,'(a)',end=100) record
281  if (record(1:1) == '#') cycle
282  if (record(1:10) == ' ') cycle
283  ntable=ntable+1
284  if (index(lowercase(record), "inside_region") .ne. 0 .or. index(lowercase(record), "outside_region") .ne. 0) then
285  if(index(lowercase(record), ".false.") .ne. 0 .or. index(lowercase(record), ".true.") .ne. 0 ) then
286  ntable_lima = ntable_lima + 1
287  read(record,*,err=99) data_entry%gridname, data_entry%fieldname_code, data_entry%fieldname_file, &
288  data_entry%file_name, ongrid, data_entry%factor, region, region_type
289  if(ongrid) then
290  data_entry%interpol_method = 'none'
291  else
292  data_entry%interpol_method = 'bilinear'
293  endif
294  else
295  ntable_new=ntable_new+1
296  read(record,*,err=99) data_entry%gridname, data_entry%fieldname_code, data_entry%fieldname_file, &
297  data_entry%file_name, data_entry%interpol_method, data_entry%factor, region, region_type
298  if (data_entry%interpol_method == 'default') then
299  data_entry%interpol_method = default_table%interpol_method
300  endif
301  if (.not.(data_entry%interpol_method == 'default' .or. &
302  data_entry%interpol_method == 'bicubic' .or. &
303  data_entry%interpol_method == 'bilinear' .or. &
304  data_entry%interpol_method == 'none')) then
305  unit = stdout()
306  write(unit,*)" gridname is ", trim(data_entry%gridname)
307  write(unit,*)" fieldname_code is ", trim(data_entry%fieldname_code)
308  write(unit,*)" fieldname_file is ", trim(data_entry%fieldname_file)
309  write(unit,*)" file_name is ", trim(data_entry%file_name)
310  write(unit,*)" factor is ", data_entry%factor
311  write(unit,*)" interpol_method is ", trim(data_entry%interpol_method)
312  call mpp_error(fatal, 'data_override_mod: invalid last entry in data_override_table, ' &
313  //'its value should be "default", "bicubic", "bilinear" or "none" ')
314  endif
315  endif
316  if( trim(region_type) == "inside_region" ) then
317  data_entry%region_type = inside_region
318  else if( trim(region_type) == "outside_region" ) then
319  data_entry%region_type = outside_region
320  else
321  call mpp_error(fatal, 'data_override_mod: region type should be inside_region or outside_region')
322  endif
323  if (data_entry%file_name == "") call mpp_error(fatal, &
324  "data_override: filename not given in data_table when region_type is not NO_REGION")
325  if(data_entry%fieldname_file == "") call mpp_error(fatal, &
326  "data_override: fieldname_file must be specified in data_table when region_type is not NO_REGION")
327  if( trim(data_entry%interpol_method) == 'none') call mpp_error(fatal, &
328  "data_override(data_override_init): ongrid must be false when region_type is not NO_REGION")
329  read(region,*) data_entry%lon_start, data_entry%lon_end, data_entry%lat_start, data_entry%lat_end
330  !--- make sure data_entry%lon_end > data_entry%lon_start and data_entry%lat_end > data_entry%lat_start
331  if(data_entry%lon_end .LE. data_entry%lon_start) call mpp_error(fatal, &
332  "data_override: lon_end should be greater than lon_start")
333  if(data_entry%lat_end .LE. data_entry%lat_start) call mpp_error(fatal, &
334  "data_override: lat_end should be greater than lat_start")
335  else if (index(lowercase(record), ".false.") .ne. 0 .or. index(lowercase(record), ".true.") .ne. 0 ) then ! old format
336  ntable_lima = ntable_lima + 1
337  read(record,*,err=99) data_entry%gridname, data_entry%fieldname_code, data_entry%fieldname_file, &
338  data_entry%file_name, ongrid, data_entry%factor
339  if(ongrid) then
340  data_entry%interpol_method = 'none'
341  else
342  data_entry%interpol_method = 'bilinear'
343  endif
344  data_entry%lon_start = 0.0
345  data_entry%lon_end = -1.0
346  data_entry%lat_start = 0.0
347  data_entry%lat_end = -1.0
348  data_entry%region_type = no_region
349  else ! new format
350  ntable_new=ntable_new+1
351  read(record,*,err=99) data_entry%gridname, data_entry%fieldname_code, data_entry%fieldname_file, &
352  data_entry%file_name, data_entry%interpol_method, data_entry%factor
353  if (data_entry%interpol_method == 'default') then
354  data_entry%interpol_method = default_table%interpol_method
355  endif
356  if (.not.(data_entry%interpol_method == 'default' .or. &
357  data_entry%interpol_method == 'bicubic' .or. &
358  data_entry%interpol_method == 'bilinear' .or. &
359  data_entry%interpol_method == 'none')) then
360  unit = stdout()
361  write(unit,*)" gridname is ", trim(data_entry%gridname)
362  write(unit,*)" fieldname_code is ", trim(data_entry%fieldname_code)
363  write(unit,*)" fieldname_file is ", trim(data_entry%fieldname_file)
364  write(unit,*)" file_name is ", trim(data_entry%file_name)
365  write(unit,*)" factor is ", data_entry%factor
366  write(unit,*)" interpol_method is ", trim(data_entry%interpol_method)
367  call mpp_error(fatal, 'data_override_mod: invalid last entry in data_override_table, ' &
368  //'its value should be "default", "bicubic", "bilinear" or "none" ')
369  endif
370  data_entry%lon_start = 0.0
371  data_entry%lon_end = -1.0
372  data_entry%lat_start = 0.0
373  data_entry%lat_end = -1.0
374  data_entry%region_type = no_region
375  endif
376  data_table(ntable) = data_entry
377  enddo
378  call mpp_error(fatal,'too many enries in data_table')
379 99 call mpp_error(fatal,'error in data_table format')
380 100 continue
381  table_size = ntable
382  if(ntable_new*ntable_lima /= 0) call mpp_error(fatal, &
383  'data_override_mod: New and old formats together in same data_table not supported')
384  call mpp_close(iunit)
385 ! Initialize override array
386  default_array%gridname = 'NONE'
387  default_array%fieldname = 'NONE'
388  default_array%t_index = -1
389  default_array%dims = -1
390  default_array%comp_domain = -1
391  do i = 1, max_array
393  enddo
395  endif
396 
397  module_is_initialized = .true.
398 
399  if ( .NOT. (atm_on .or. ocn_on .or. lnd_on .or. ice_on .or. lndug_on)) return
400  call fms_io_init
401 
402 ! Test if grid_file is already opened
403  inquire (file=trim(grid_file), opened=file_open)
404  if(file_open) call mpp_error(fatal, trim(grid_file)//' already opened')
405 
406  if(field_exist(grid_file, "x_T" ) .OR. field_exist(grid_file, "geolon_t" ) ) then
407  if (atm_on .and. .not. allocated(lon_local_atm) ) then
408  call mpp_get_compute_domain( atm_domain,is,ie,js,je)
409  allocate(lon_local_atm(is:ie,js:je), lat_local_atm(is:ie,js:je))
410  call get_grid_version_1(grid_file, 'atm', atm_domain, is, ie, js, je, lon_local_atm, lat_local_atm, &
412  endif
413  if (ocn_on .and. .not. allocated(lon_local_ocn) ) then
414  call mpp_get_compute_domain( ocn_domain,is,ie,js,je)
415  allocate(lon_local_ocn(is:ie,js:je), lat_local_ocn(is:ie,js:je))
416  call get_grid_version_1(grid_file, 'ocn', ocn_domain, is, ie, js, je, lon_local_ocn, lat_local_ocn, &
418  endif
419 
420  if (lnd_on .and. .not. allocated(lon_local_lnd) ) then
421  call mpp_get_compute_domain( lnd_domain,is,ie,js,je)
422  allocate(lon_local_lnd(is:ie,js:je), lat_local_lnd(is:ie,js:je))
423  call get_grid_version_1(grid_file, 'lnd', lnd_domain, is, ie, js, je, lon_local_lnd, lat_local_lnd, &
425  endif
426 
427  if (ice_on .and. .not. allocated(lon_local_ice) ) then
428  call mpp_get_compute_domain( ice_domain,is,ie,js,je)
429  allocate(lon_local_ice(is:ie,js:je), lat_local_ice(is:ie,js:je))
430  call get_grid_version_1(grid_file, 'ice', ice_domain, is, ie, js, je, lon_local_ice, lat_local_ice, &
432  endif
433  else if(field_exist(grid_file, "ocn_mosaic_file" ) .OR. field_exist(grid_file, "gridfiles" ) ) then
434  if(field_exist(grid_file, "gridfiles" ) ) then
435  count = 0
436  if (atm_on) count = count + 1
437  if (lnd_on) count = count + 1
438  if ( ocn_on .OR. ice_on ) count = count + 1
439  if(count .NE. 1) call mpp_error(fatal, 'data_override_mod: the grid file is a solo mosaic, ' // &
440  'one and only one of atm_on, lnd_on or ice_on/ocn_on should be true')
441  endif
442  if (atm_on .and. .not. allocated(lon_local_atm) ) then
443  call mpp_get_compute_domain(atm_domain,is,ie,js,je)
444  allocate(lon_local_atm(is:ie,js:je), lat_local_atm(is:ie,js:je))
445  call get_grid_version_2(grid_file, 'atm', atm_domain, is, ie, js, je, lon_local_atm, lat_local_atm, &
447  endif
448 
449  if (ocn_on .and. .not. allocated(lon_local_ocn) ) then
450  call mpp_get_compute_domain( ocn_domain,is,ie,js,je)
451  allocate(lon_local_ocn(is:ie,js:je), lat_local_ocn(is:ie,js:je))
452  call get_grid_version_2(grid_file, 'ocn', ocn_domain, is, ie, js, je, lon_local_ocn, lat_local_ocn, &
454  endif
455 
456  if (lnd_on .and. .not. allocated(lon_local_lnd) ) then
457  call mpp_get_compute_domain( lnd_domain,is,ie,js,je)
458  allocate(lon_local_lnd(is:ie,js:je), lat_local_lnd(is:ie,js:je))
459  call get_grid_version_2(grid_file, 'lnd', lnd_domain, is, ie, js, je, lon_local_lnd, lat_local_lnd, &
461  endif
462 
463  if (ice_on .and. .not. allocated(lon_local_ice) ) then
464  call mpp_get_compute_domain( ice_domain,is,ie,js,je)
465  allocate(lon_local_ice(is:ie,js:je), lat_local_ice(is:ie,js:je))
466  call get_grid_version_2(grid_file, 'ocn', ice_domain, is, ie, js, je, lon_local_ice, lat_local_ice, &
468  endif
469  else
470  call mpp_error(fatal, 'data_override_mod: none of x_T, geolon_t, ocn_mosaic_file or gridfiles exist in '//trim(grid_file))
471  end if
472 
473 end subroutine data_override_init
474 ! </SUBROUTINE>
475 !===============================================================================================
476 
477 !===============================================================================================
478 ! <SUBROUTINE NAME="data_override_unset_domain">
479 ! <DESCRIPTION>
480 ! Unset domains that had previously been set for use by data_override.
481 ! </DESCRIPTION>
482 ! <TEMPLATE>
483 ! call data_override_unset_domain
484 ! </TEMPLATE>
485 subroutine data_override_unset_domains(unset_Atm, unset_Ocean, &
486  unset_Ice, unset_Land, must_be_set)
487  logical, intent(in), optional :: unset_atm, unset_ocean, unset_ice, unset_land
488  logical, intent(in), optional :: must_be_set
489 
490 ! <NOTE>
491 ! This subroutine deallocates any data override domains that have been set.
492 ! </NOTE>
493  logical :: fail_if_not_set
494 
495  fail_if_not_set = .true. ; if (present(must_be_set)) fail_if_not_set = must_be_set
496 
497  if (.not.module_is_initialized) call mpp_error(fatal, &
498  "data_override_unset_domains called with an unititialized data_override module.")
499 
500  if (PRESENT(unset_atm)) then ; if (unset_atm) then
501  if (fail_if_not_set .and. .not.atm_on) call mpp_error(fatal, &
502  "data_override_unset_domains attempted to work on an Atm_domain that had not been set.")
503  atm_domain = null_domain2d
504  atm_on = .false.
505  if (allocated(lon_local_atm)) deallocate(lon_local_atm)
506  if (allocated(lat_local_atm)) deallocate(lat_local_atm)
507  endif ; endif
508  if (PRESENT(unset_ocean)) then ; if (unset_ocean) then
509  if (fail_if_not_set .and. .not.ocn_on) call mpp_error(fatal, &
510  "data_override_unset_domains attempted to work on an Ocn_domain that had not been set.")
511  ocn_domain = null_domain2d
512  ocn_on = .false.
513  if (allocated(lon_local_ocn)) deallocate(lon_local_ocn)
514  if (allocated(lat_local_ocn)) deallocate(lat_local_ocn)
515  endif ; endif
516  if (PRESENT(unset_land)) then ; if (unset_land) then
517  if (fail_if_not_set .and. .not.lnd_on) call mpp_error(fatal, &
518  "data_override_unset_domains attempted to work on a Land_domain that had not been set.")
519  lnd_domain = null_domain2d
520  lnd_on = .false.
521  if (allocated(lon_local_lnd)) deallocate(lon_local_lnd)
522  if (allocated(lat_local_lnd)) deallocate(lat_local_lnd)
523  endif ; endif
524  if (PRESENT(unset_ice)) then ; if (unset_ice) then
525  if (fail_if_not_set .and. .not.ice_on) call mpp_error(fatal, &
526  "data_override_unset_domains attempted to work on an Ice_domain that had not been set.")
527  ice_domain = null_domain2d
528  ice_on = .false.
529  if (allocated(lon_local_ice)) deallocate(lon_local_ice)
530  if (allocated(lat_local_ice)) deallocate(lat_local_ice)
531  endif ; endif
532 
533 end subroutine data_override_unset_domains
534 ! </SUBROUTINE>
535 !===============================================================================================
536 
537 
538 subroutine check_grid_sizes(domain_name, Domain, nlon, nlat)
539 character(len=12), intent(in) :: domain_name
540 type(domain2d), intent(in) :: domain
541 integer, intent(in) :: nlon, nlat
542 
543 character(len=184) :: error_message
544 integer :: xsize, ysize
545 
546 call mpp_get_global_domain(domain, xsize=xsize, ysize=ysize)
547 if(nlon .NE. xsize .OR. nlat .NE. ysize) then
548  error_message = 'Error in data_override_init. Size of grid as specified by '// &
549  ' does not conform to that specified by grid_spec.nc.'// &
550  ' From : by From grid_spec.nc: by '
551  error_message( 59: 70) = domain_name
552  error_message(130:141) = domain_name
553  write(error_message(143:146),'(i4)') xsize
554  write(error_message(150:153),'(i4)') ysize
555  write(error_message(174:177),'(i4)') nlon
556  write(error_message(181:184),'(i4)') nlat
557  call mpp_error(fatal,error_message)
558 endif
559 
560 end subroutine check_grid_sizes
561 !===============================================================================================
562 subroutine get_domain(gridname, domain, comp_domain)
563 ! Given a gridname, this routine returns the working domain associated with this gridname
564  character(len=3), intent(in) :: gridname
565  type(domain2D), intent(inout) :: domain
566  integer, intent(out), optional :: comp_domain(4) ! istart,iend,jstart,jend for compute domain
567 
568  domain = null_domain2d
569  select case (gridname)
570  case('OCN')
571  domain = ocn_domain
572  case('ATM')
573  domain = atm_domain
574  case('LND')
575  domain = lnd_domain
576  case('ICE')
577  domain = ice_domain
578  case default
579  call mpp_error(fatal,'error in data_override get_domain')
580  end select
581  if(domain .EQ. null_domain2d) call mpp_error(fatal,'data_override: failure in get_domain')
582  if(present(comp_domain)) &
583  call mpp_get_compute_domain(domain,comp_domain(1),comp_domain(2),comp_domain(3),comp_domain(4))
584 end subroutine get_domain
585 
586 subroutine get_domainug(gridname, UGdomain, comp_domain)
587 ! Given a gridname, this routine returns the working domain associated with this gridname
588  character(len=3), intent(in) :: gridname
589  type(domainUG), intent(inout) :: UGdomain
590  integer, intent(out), optional :: comp_domain(4) ! istart,iend,jstart,jend for compute domain
591  type(domain2D), pointer :: SGdomain => null()
592 
593  ugdomain = null_domainug
594  select case (gridname)
595  case('LND')
596  ugdomain = lnd_domainug
597  case default
598  call mpp_error(fatal,'error in data_override get_domain')
599  end select
600 ! if(UGdomain .EQ. NULL_DOMAINUG) call mpp_error(FATAL,'data_override: failure in get_domain')
601  if(present(comp_domain)) &
602  call mpp_get_ug_sg_domain(ugdomain,sgdomain)
603  call mpp_get_compute_domain(sgdomain,comp_domain(1),comp_domain(2),comp_domain(3),comp_domain(4))
604 end subroutine get_domainug
605 !===============================================================================================
606 
607 ! <SUBROUTINE NAME="data_override_2d">
608 ! <DESCRIPTION>
609 ! This routine performs data override for 2D fields; for usage, see data_override_3d.
610 ! </DESCRIPTION>
611 subroutine data_override_2d(gridname,fieldname,data_2D,time,override, is_in, ie_in, js_in, je_in)
612  character(len=3), intent(in) :: gridname ! model grid ID
613  character(len=*), intent(in) :: fieldname ! field to override
614  logical, intent(out), optional :: override ! true if the field has been overriden succesfully
615  type(time_type), intent(in) :: time ! model time
616  real, dimension(:,:), intent(inout) :: data_2D !data returned by this call
617  integer, optional, intent(in) :: is_in, ie_in, js_in, je_in
618 ! real, dimension(size(data_2D,1),size(data_2D,2),1) :: data_3D
619  real, dimension(:,:,:), allocatable :: data_3D
620  integer :: index1
621  integer :: i
622 
623 !1 Look for the data file in data_table
624  if(PRESENT(override)) override = .false.
625  index1 = -1
626  do i = 1, table_size
627  if( trim(gridname) /= trim(data_table(i)%gridname)) cycle
628  if( trim(fieldname) /= trim(data_table(i)%fieldname_code)) cycle
629  index1 = i ! field found
630  exit
631  enddo
632  if(index1 .eq. -1) return ! NO override was performed
633 
634  allocate(data_3d(size(data_2d,1),size(data_2d,2),1))
635  data_3d(:,:,1) = data_2d
636  call data_override_3d(gridname,fieldname,data_3d,time,override,data_index=index1,&
637  is_in=is_in,ie_in=ie_in,js_in=js_in,je_in=je_in)
638 
639  data_2d(:,:) = data_3d(:,:,1)
640  deallocate(data_3d)
641 end subroutine data_override_2d
642 ! </SUBROUTINE>
643 !===============================================================================================
644 
645 ! <SUBROUTINE NAME="data_override_3d">
646 ! <DESCRIPTION>
647 ! This routine performs data override for 3D fields
648 ! <TEMPLATE>
649 ! call data_override(gridname,fieldname,data,time,override)
650 ! </TEMPLATE>
651 ! </DESCRIPTION>
652 
653 ! <IN NAME="gridname" TYPE="character" DIM="(*)">
654 ! Grid name (Ocean, Ice, Atmosphere, Land)
655 ! </IN>
656 ! <IN NAME="fieldname_code" TYPE="character" DIM="(*)">
657 ! Field name as used in the code (may be different from the name in NetCDF data file)
658 ! </IN>
659 ! <OUT NAME="data" TYPE="real" DIM="(:,:,:)">
660 ! array containing output data
661 ! </OUT>
662 ! <IN NAME="time" TYPE="time_type">
663 ! model time
664 ! </IN>
665 ! <OUT NAME="override" TYPE="logical">
666 ! TRUE if the field is overriden, FALSE otherwise
667 ! </OUT>
668 ! <IN NAME="data_index" TYPE="integer">
669 ! </IN>
670 subroutine data_override_3d(gridname,fieldname_code,data,time,override,data_index, is_in, ie_in, js_in, je_in)
671  character(len=3), intent(in) :: gridname ! model grid ID
672  character(len=*), intent(in) :: fieldname_code ! field name as used in the model
673  logical, optional, intent(out) :: override ! true if the field has been overriden succesfully
674  type(time_type), intent(in) :: time !(target) model time
675  integer, optional, intent(in) :: data_index
676  real, dimension(:,:,:), intent(inout) :: data !data returned by this call
677  integer, optional, intent(in) :: is_in, ie_in, js_in, je_in
678  logical, dimension(:,:,:), allocatable :: mask_out
679 
680  character(len=512) :: filename, filename2 !file containing source data
681  character(len=128) :: fieldname ! fieldname used in the data file
682  integer :: i,j
683  integer :: dims(4)
684  integer :: index1 ! field index in data_table
685  integer :: id_time !index for time interp in override array
686  integer :: axis_sizes(4)
687  real, dimension(:,:), pointer :: lon_local =>null(), &
688  lat_local =>null() !of output (target) grid cells
689  real, dimension(:), allocatable :: lon_tmp, lat_tmp
690 
691  type(axistype) :: axis_centers(4), axis_bounds(4)
692  logical :: data_file_is_2D = .false. !data in netCDF file is 2D
693  logical :: ongrid, use_comp_domain
694  type(domain2D) :: domain
695  integer :: curr_position ! position of the field currently processed in override_array
696  real :: factor
697  integer, dimension(4) :: comp_domain = 0 ! istart,iend,jstart,jend for compute domain
698  integer :: nxd, nyd, nxc, nyc, nwindows
699  integer :: nwindows_x, ipos, jpos, window_size(2)
700  integer :: istart, iend, jstart, jend
701  integer :: isw, iew, jsw, jew, n
702  integer :: omp_get_num_threads, omp_get_thread_num, thread_id, window_id
703  logical :: need_compute
704  real :: lat_min, lat_max
705  integer :: is_src, ie_src, js_src, je_src
706  logical :: exists
707 
708  use_comp_domain = .false.
709  if(.not.module_is_initialized) &
710  call mpp_error(fatal,'Error: need to call data_override_init first')
711 
712 !1 Look for the data file in data_table
713  if(PRESENT(override)) override = .false.
714  if (present(data_index)) then
715  index1 = data_index
716  else
717  index1 = -1
718  do i = 1, table_size
719  if( trim(gridname) /= trim(data_table(i)%gridname)) cycle
720  if( trim(fieldname_code) /= trim(data_table(i)%fieldname_code)) cycle
721  index1 = i ! field found
722  exit
723  enddo
724  if(index1 .eq. -1) then
725  if(mpp_pe() == mpp_root_pe() .and. debug_data_override) &
726  call mpp_error(warning,'this field is NOT found in data_table: '//trim(fieldname_code))
727  return ! NO override was performed
728  endif
729  endif
730 
731  fieldname = data_table(index1)%fieldname_file ! fieldname in netCDF data file
732  factor = data_table(index1)%factor
733 
734  if(fieldname == "") then
735  data = factor
736  if(PRESENT(override)) override = .true.
737  return
738  else
739  filename = data_table(index1)%file_name
740  if (filename == "") call mpp_error(fatal,'data_override: filename not given in data_table')
741  endif
742 
743  ongrid = (data_table(index1)%interpol_method == 'none')
744 
745 !3 Check if fieldname has been previously processed
746 !$OMP CRITICAL
747  curr_position = -1
748  if(num_fields > 0 ) then
749  do i = 1, num_fields
750  if(trim(override_array(i)%gridname) /= trim(gridname)) cycle
751  if(trim(override_array(i)%fieldname) /= trim(fieldname_code)) cycle
752  curr_position = i
753  exit
754  enddo
755  endif
756 
757  if(curr_position < 0) then ! the field has not been processed previously
758  num_fields = num_fields + 1
759  curr_position = num_fields
760 
761 ! Get working domain from model's gridname
762  call get_domain(gridname,domain,comp_domain)
763  call mpp_get_data_domain(domain, xsize=nxd, ysize=nyd)
764  nxc = comp_domain(2)-comp_domain(1) + 1
765  nyc = comp_domain(4)-comp_domain(3) + 1
766 
767 ! record fieldname, gridname in override_array
768  override_array(curr_position)%fieldname = fieldname_code
769  override_array(curr_position)%gridname = gridname
770  override_array(curr_position)%comp_domain = comp_domain
771 ! get number of threads
772  override_array(curr_position)%numthreads = 1
773 #if defined(_OPENMP)
774  override_array(curr_position)%numthreads = omp_get_num_threads()
775 #endif
776 !--- data_override may be called from physics windows. The following are possible situations
777 !--- 1. size(data,1) == nxd and size(data,2) == nyd ( on data domain and there is only one window).
778 !--- 2. nxc is divisible by size(data,1), nyc is divisible by size(data,2),
779 !--- nwindow = (nxc/size(data(1))*(nyc/size(data,2)), also we require nwindows is divisible by nthreads.
780 !--- The another restrition is that size(data,1) == ie_in - is_in + 1,
781 !--- size(data,2) == je_in - js_in + 1
782  nwindows = 1
783  if( nxd == size(data,1) .AND. nyd == size(data,2) ) then !
784  use_comp_domain = .false.
785  else if ( mod(nxc, size(data,1)) ==0 .AND. mod(nyc, size(data,2)) ==0 ) then
786  use_comp_domain = .true.
787  nwindows = (nxc/size(data,1))*(nyc/size(data,2))
788  else
789  call mpp_error(fatal, "data_override: data is not on data domain and compute domain is not divisible by size(data)")
790  endif
791  override_array(curr_position)%window_size(1) = size(data,1)
792  override_array(curr_position)%window_size(2) = size(data,2)
793 
794  window_size = override_array(curr_position)%window_size
795  override_array(curr_position)%numwindows = nwindows
796  if( mod(nwindows, override_array(curr_position)%numthreads) .NE. 0 ) then
797  call mpp_error(fatal, "data_override: nwindow is not divisible by nthreads")
798  endif
799  allocate(override_array(curr_position)%need_compute(nwindows))
800  override_array(curr_position)%need_compute = .true.
801 
802 !4 get index for time interp
803  if(ongrid) then
804  if( data_table(index1)%region_type .NE. no_region ) then
805  call mpp_error(fatal,.NE.'data_override: ongrid must be false when region_type NO_REGION')
806  endif
807 
808 ! Allow on-grid data_overrides on cubed sphere grid
809  inquire(file=trim(filename),exist=exists)
810  if (.not. exists) then
811  call get_mosaic_tile_file(filename,filename2,.false.,domain)
812  filename = filename2
813  endif
814 
815  !--- we always only pass data on compute domain
816  id_time = init_external_field(filename,fieldname,domain=domain,verbose=.false., &
817  use_comp_domain=use_comp_domain, nwindows=nwindows)
818  dims = get_external_field_size(id_time)
819  override_array(curr_position)%dims = dims
820  if(id_time<0) call mpp_error(fatal,'data_override:field not found in init_external_field 1')
821  override_array(curr_position)%t_index = id_time
822  else !ongrid=false
823  id_time = init_external_field(filename,fieldname,domain=domain, axis_centers=axis_centers,&
824  axis_sizes=axis_sizes, verbose=.false.,override=.true.,use_comp_domain=use_comp_domain, &
825  nwindows = nwindows)
826  dims = get_external_field_size(id_time)
827  override_array(curr_position)%dims = dims
828  if(id_time<0) call mpp_error(fatal,'data_override:field not found in init_external_field 2')
829  override_array(curr_position)%t_index = id_time
830 
831  ! get lon and lat of the input (source) grid, assuming that axis%data contains
832  ! lat and lon of the input grid (in degrees)
833  call get_axis_bounds(axis_centers(1),axis_bounds(1), axis_centers)
834  call get_axis_bounds(axis_centers(2),axis_bounds(2), axis_centers)
835 
836  allocate(override_array(curr_position)%horz_interp(nwindows))
837  allocate(override_array(curr_position)%lon_in(axis_sizes(1)+1))
838  allocate(override_array(curr_position)%lat_in(axis_sizes(2)+1))
839  call mpp_get_axis_data(axis_bounds(1),override_array(curr_position)%lon_in)
840  call mpp_get_axis_data(axis_bounds(2),override_array(curr_position)%lat_in)
841 
842 ! convert lon_in and lat_in from deg to radian
843  override_array(curr_position)%lon_in = override_array(curr_position)%lon_in * deg_to_radian
844  override_array(curr_position)%lat_in = override_array(curr_position)%lat_in * deg_to_radian
845 
846  !--- find the region of the source grid that cover the local model grid.
847  !--- currently we only find the index range for j-direction because
848  !--- of the cyclic condition in i-direction. The purpose of this is to
849  !--- decrease the memory usage and increase the IO performance.
850  select case(gridname)
851  case('OCN')
852  lon_local => lon_local_ocn; lat_local => lat_local_ocn
853  case('ICE')
854  lon_local => lon_local_ice; lat_local => lat_local_ice
855  case('ATM')
856  lon_local => lon_local_atm; lat_local => lat_local_atm
857  case('LND')
858  lon_local => lon_local_lnd; lat_local => lat_local_lnd
859  case default
860  call mpp_error(fatal,'error: gridname not recognized in data_override')
861  end select
862 
863  lat_min = minval(lat_local)
864  lat_max = maxval(lat_local)
865  is_src = 1
866  ie_src = axis_sizes(1)
867  js_src = 1
868  je_src = axis_sizes(2)
869  do j = 1, axis_sizes(2)+1
870  if( override_array(curr_position)%lat_in(j) > lat_min ) exit
871  js_src = j
872  enddo
873  do j = 1, axis_sizes(2)+1
874  je_src = j
875  if( override_array(curr_position)%lat_in(j) > lat_max ) exit
876  enddo
877 
878  !--- bicubic interpolation need one extra point in each direction. Also add
879  !--- one more point for because lat_in is in the corner but the interpolation
880  !--- use center points.
881  select case (data_table(index1)%interpol_method)
882  case ('bilinear')
883  js_src = max(1, js_src-1)
884  je_src = min(axis_sizes(2), je_src+1)
885  case ('bicubic')
886  js_src = max(1, js_src-2)
887  je_src = min(axis_sizes(2), je_src+2)
888  end select
889  override_array(curr_position)%is_src = is_src
890  override_array(curr_position)%ie_src = ie_src
891  override_array(curr_position)%js_src = js_src
892  override_array(curr_position)%je_src = je_src
893  call reset_src_data_region(id_time, is_src, ie_src, js_src, je_src)
894 
895 ! Find the index of lon_start, lon_end, lat_start and lat_end in the input grid (nearest points)
896  if( data_table(index1)%region_type .NE. no_region ) then
897  allocate( lon_tmp(axis_sizes(1)), lat_tmp(axis_sizes(2)) )
898  call mpp_get_axis_data(axis_centers(1), lon_tmp)
899  call mpp_get_axis_data(axis_centers(2), lat_tmp)
900  ! limit lon_start, lon_end are inside lon_in
901  ! lat_start, lat_end are inside lat_in
902  if( data_table(index1)%lon_start < lon_tmp(1) .OR. data_table(index1)%lon_start .GT. lon_tmp(axis_sizes(1))) &
903  call mpp_error(fatal, "data_override: lon_start is outside lon_T")
904  if( data_table(index1)%lon_end < lon_tmp(1) .OR. data_table(index1)%lon_end .GT. lon_tmp(axis_sizes(1))) &
905  call mpp_error(fatal, "data_override: lon_end is outside lon_T")
906  if( data_table(index1)%lat_start < lat_tmp(1) .OR. data_table(index1)%lat_start .GT. lat_tmp(axis_sizes(2))) &
907  call mpp_error(fatal, "data_override: lat_start is outside lat_T")
908  if( data_table(index1)%lat_end < lat_tmp(1) .OR. data_table(index1)%lat_end .GT. lat_tmp(axis_sizes(2))) &
909  call mpp_error(fatal, "data_override: lat_end is outside lat_T")
910  istart = nearest_index(data_table(index1)%lon_start, lon_tmp)
911  iend = nearest_index(data_table(index1)%lon_end, lon_tmp)
912  jstart = nearest_index(data_table(index1)%lat_start, lat_tmp)
913  jend = nearest_index(data_table(index1)%lat_end, lat_tmp)
914  ! adjust the index according to is_src and js_src
915  istart = istart - is_src + 1
916  iend = iend - is_src + 1
917  jstart = jstart - js_src + 1
918  jend = jend - js_src + 1
919  call set_override_region(id_time, data_table(index1)%region_type, istart, iend, jstart, jend)
920  deallocate(lon_tmp, lat_tmp)
921  endif
922 
923  endif
924  else !curr_position >0
925  dims = override_array(curr_position)%dims
926  comp_domain = override_array(curr_position)%comp_domain
927  is_src = override_array(curr_position)%is_src
928  ie_src = override_array(curr_position)%ie_src
929  js_src = override_array(curr_position)%js_src
930  je_src = override_array(curr_position)%je_src
931  window_size = override_array(curr_position)%window_size
932  !---make sure data size match window_size
933  if( window_size(1) .NE. size(data,1) .OR. window_size(2) .NE. size(data,2) ) then
934  call mpp_error(fatal, "data_override: window_size does not match size(data)")
935  endif
936 !9 Get id_time previously stored in override_array
937  id_time = override_array(curr_position)%t_index
938  endif
939 !$OMP END CRITICAL
940 
941  if( override_array(curr_position)%numwindows > 1 ) then
942  if( .NOT. PRESENT(is_in) .OR. .NOT. PRESENT(is_in) .OR. .NOT. PRESENT(is_in) .OR. .NOT. PRESENT(is_in) ) then
943  call mpp_error(fatal, "data_override: is_in, ie_in, js_in, je_in must be present when nwindows > 1")
944  endif
945  endif
946 
947  isw = comp_domain(1)
948  iew = comp_domain(2)
949  jsw = comp_domain(3)
950  jew = comp_domain(4)
951  window_id = 1
952  if( override_array(curr_position)%numwindows > 1 ) then
953  nxc = comp_domain(2) - comp_domain(1) + 1
954  nwindows_x = nxc/window_size(1)
955  ipos = (is_in-1)/window_size(1) + 1
956  jpos = (js_in-1)/window_size(2)
957 
958  window_id = jpos*nwindows_x + ipos
959  isw = isw + is_in - 1
960  iew = isw + ie_in - is_in
961  jsw = jsw + js_in - 1
962  jew = jsw + je_in - js_in
963  endif
964 
965  if( ongrid ) then
966  need_compute = .false.
967  else
968  !--- find the index for windows.
969  need_compute=override_array(curr_position)%need_compute(window_id)
970  endif
971 
972  !--- call horiz_interp_new is not initialized
973 
974  if( need_compute ) then
975  select case(gridname)
976  case('OCN')
977  lon_local => lon_local_ocn; lat_local => lat_local_ocn
978  case('ICE')
979  lon_local => lon_local_ice; lat_local => lat_local_ice
980  case('ATM')
981  lon_local => lon_local_atm; lat_local => lat_local_atm
982  case('LND')
983  lon_local => lon_local_lnd; lat_local => lat_local_lnd
984  case default
985  call mpp_error(fatal,'error: gridname not recognized in data_override')
986  end select
987 
988  select case (data_table(index1)%interpol_method)
989  case ('bilinear')
990  call horiz_interp_new (override_array(curr_position)%horz_interp(window_id), &
991  override_array(curr_position)%lon_in(is_src:ie_src+1), &
992  override_array(curr_position)%lat_in(js_src:je_src+1), &
993  lon_local(isw:iew,jsw:jew), lat_local(isw:iew,jsw:jew), interp_method="bilinear")
994  case ('bicubic')
995  call horiz_interp_new (override_array(curr_position)%horz_interp(window_id), &
996  override_array(curr_position)%lon_in(is_src:ie_src+1), &
997  override_array(curr_position)%lat_in(js_src:je_src+1), &
998  lon_local(isw:iew,jsw:jew), lat_local(isw:iew,jsw:jew), interp_method="bicubic")
999  end select
1000  override_array(curr_position)%need_compute(window_id) = .false.
1001  endif
1002 
1003  ! Determine if data in netCDF file is 2D or not
1004  data_file_is_2d = .false.
1005  if((dims(3) == 1) .and. (size(data,3)>1)) data_file_is_2d = .true.
1006 
1007  if(dims(3) .NE. 1 .and. (size(data,3) .NE. dims(3))) &
1008  call mpp_error(fatal, .NE..NE."data_override: dims(3) 1 and size(data,3) dims(3)")
1009 
1010  if(ongrid) then
1011 !10 do time interp to get data in compute_domain
1012  if(data_file_is_2d) then
1013  call time_interp_external(id_time,time,data(:,:,1),verbose=.false., &
1014  is_in=is_in,ie_in=ie_in,js_in=js_in,je_in=je_in,window_id=window_id)
1015  data(:,:,1) = data(:,:,1)*factor
1016  do i = 2, size(data,3)
1017  data(:,:,i) = data(:,:,1)
1018  enddo
1019  else
1020  call time_interp_external(id_time,time,data,verbose=.false., &
1021  is_in=is_in,ie_in=ie_in,js_in=js_in,je_in=je_in,window_id=window_id)
1022  data = data*factor
1023  endif
1024  else ! off grid case
1025 ! do time interp to get global data
1026  if(data_file_is_2d) then
1027  if( data_table(index1)%region_type == no_region ) then
1028  call time_interp_external(id_time,time,data(:,:,1),verbose=.false., &
1029  horz_interp=override_array(curr_position)%horz_interp(window_id), &
1030  is_in=is_in,ie_in=ie_in,js_in=js_in,je_in=je_in,window_id=window_id)
1031  data(:,:,1) = data(:,:,1)*factor
1032  do i = 2, size(data,3)
1033  data(:,:,i) = data(:,:,1)
1034  enddo
1035  else
1036  allocate(mask_out(size(data,1), size(data,2),1))
1037  mask_out = .false.
1038  call time_interp_external(id_time,time,data(:,:,1),verbose=.false., &
1039  horz_interp=override_array(curr_position)%horz_interp(window_id), &
1040  mask_out =mask_out(:,:,1), &
1041  is_in=is_in,ie_in=ie_in,js_in=js_in,je_in=je_in,window_id=window_id)
1042  where(mask_out(:,:,1))
1043  data(:,:,1) = data(:,:,1)*factor
1044  end where
1045  do i = 2, size(data,3)
1046  where(mask_out(:,:,1))
1047  data(:,:,i) = data(:,:,1)
1048  end where
1049  enddo
1050  deallocate(mask_out)
1051  endif
1052  else
1053  if( data_table(index1)%region_type == no_region ) then
1054  call time_interp_external(id_time,time,data,verbose=.false., &
1055  horz_interp=override_array(curr_position)%horz_interp(window_id), &
1056  is_in=is_in,ie_in=ie_in,js_in=js_in,je_in=je_in,window_id=window_id)
1057  data = data*factor
1058  else
1059  allocate(mask_out(size(data,1), size(data,2), size(data,3)) )
1060  mask_out = .false.
1061  call time_interp_external(id_time,time,data,verbose=.false., &
1062  horz_interp=override_array(curr_position)%horz_interp(window_id), &
1063  mask_out =mask_out, &
1064  is_in=is_in,ie_in=ie_in,js_in=js_in,je_in=je_in,window_id=window_id)
1065 
1066  where(mask_out)
1067  data = data*factor
1068  end where
1069  deallocate(mask_out)
1070  endif
1071  endif
1072 
1073  endif
1074 
1075  if(PRESENT(override)) override = .true.
1076 
1077 end subroutine data_override_3d
1078 ! </SUBROUTINE>
1079 
1080 ! <SUBROUTINE NAME="data_override_0d">
1081 ! <DESCRIPTION>
1082 ! This routine performs data override for scalar fields
1083 ! <TEMPLATE>
1084 ! call data_override(fieldname,data,time,override)
1085 ! </TEMPLATE>
1086 ! </DESCRIPTION>
1087 ! <IN NAME="gridname" TYPE="character" DIM="(*)">
1088 ! Grid name (Ocean, Ice, Atmosphere, Land)
1089 ! </IN>
1090 ! <IN NAME="fieldname_code" TYPE="character" DIM="(*)">
1091 ! Field name as used in the code (may be different from the name in NetCDF data file)
1092 ! </IN>
1093 ! <OUT NAME="data" TYPE="real" DIM="(:,:,:)">
1094 ! array containing output data
1095 ! </OUT>
1096 ! <IN NAME="time" TYPE="time_type">
1097 ! model time
1098 ! </IN>
1099 ! <OUT NAME="override" TYPE="logical">
1100 ! TRUE if the field is overriden, FALSE otherwise
1101 ! </OUT>
1102 ! <IN NAME="data_index" TYPE="integer">
1103 ! </IN>
1104 subroutine data_override_0d(gridname,fieldname_code,data,time,override,data_index)
1105  character(len=3), intent(in) :: gridname ! model grid ID
1106  character(len=*), intent(in) :: fieldname_code ! field name as used in the model
1107  logical, intent(out), optional :: override ! true if the field has been overriden succesfully
1108  type(time_type), intent(in) :: time !(target) model time
1109  real, intent(out) :: data !data returned by this call
1110  integer, intent(in), optional :: data_index
1111 
1112  character(len=512) :: filename !file containing source data
1113  character(len=128) :: fieldname ! fieldname used in the data file
1114  integer :: index1 ! field index in data_table
1115  integer :: id_time !index for time interp in override array
1116  integer :: curr_position ! position of the field currently processed in override_array
1117  integer :: i
1118  real :: factor
1119 
1120  if(.not.module_is_initialized) &
1121  call mpp_error(fatal,'Error: need to call data_override_init first')
1122 
1123 !1 Look for the data file in data_table
1124  if(PRESENT(override)) override = .false.
1125  if (present(data_index)) then
1126  index1 = data_index
1127  else
1128  index1 = -1
1129  do i = 1, table_size
1130  if( trim(gridname) /= trim(data_table(i)%gridname)) cycle
1131  if( trim(fieldname_code) /= trim(data_table(i)%fieldname_code)) cycle
1132  index1 = i ! field found
1133  exit
1134  enddo
1135  if(index1 .eq. -1) then
1136  if(mpp_pe() == mpp_root_pe() .and. debug_data_override) &
1137  call mpp_error(warning,'this field is NOT found in data_table: '//trim(fieldname_code))
1138  return ! NO override was performed
1139  endif
1140  endif
1141 
1142  fieldname = data_table(index1)%fieldname_file ! fieldname in netCDF data file
1143  factor = data_table(index1)%factor
1144 
1145  if(fieldname == "") then
1146  data = factor
1147  if(PRESENT(override)) override = .true.
1148  return
1149  else
1150  filename = data_table(index1)%file_name
1151  if (filename == "") call mpp_error(fatal,'data_override: filename not given in data_table')
1152  endif
1153 
1154 !3 Check if fieldname has been previously processed
1155 !$OMP SINGLE
1156  curr_position = -1
1157  if(num_fields > 0 ) then
1158  do i = 1, num_fields
1159  if(trim(override_array(i)%gridname) /= trim(gridname)) cycle
1160  if(trim(override_array(i)%fieldname) /= trim(fieldname_code)) cycle
1161  curr_position = i
1162  exit
1163  enddo
1164  endif
1165 
1166  if(curr_position < 0) then ! the field has not been processed previously
1167  num_fields = num_fields + 1
1168  curr_position = num_fields
1169  ! record fieldname, gridname in override_array
1170  override_array(curr_position)%fieldname = fieldname_code
1171  override_array(curr_position)%gridname = gridname
1172  id_time = init_external_field(filename,fieldname,verbose=.false.)
1173  if(id_time<0) call mpp_error(fatal,'data_override:field not found in init_external_field 1')
1174  override_array(curr_position)%t_index = id_time
1175  else !curr_position >0
1176  !9 Get id_time previously stored in override_array
1177  id_time = override_array(curr_position)%t_index
1178  endif !if curr_position < 0
1179 
1180  !10 do time interp to get data in compute_domain
1181  call time_interp_external(id_time, time, data, verbose=.false.)
1182  data = data*factor
1183 !$OMP END SINGLE
1184 
1185  if(PRESENT(override)) override = .true.
1186 
1187 end subroutine data_override_0d
1188 ! </SUBROUTINE>
1189 
1190 subroutine data_override_ug_1d(gridname,fieldname,data,time,override)
1191  character(len=3), intent(in) :: gridname ! model grid ID
1192  character(len=*), intent(in) :: fieldname ! field to override
1193  real, dimension(:), intent(inout) :: data !data returned by this call
1194  type(time_type), intent(in) :: time ! model time
1195  logical, intent(out), optional :: override ! true if the field has been overriden succesfully
1196  !local vars
1197  real, dimension(:,:), allocatable :: data_SG
1198  type(domainUG) :: UG_domain
1199  integer :: index1
1200  integer :: i
1201  integer, dimension(4) :: comp_domain = 0 ! istart,iend,jstart,jend for compute domain
1202 
1203  !1 Look for the data file in data_table
1204  if(PRESENT(override)) override = .false.
1205  index1 = -1
1206  do i = 1, table_size
1207  if( trim(gridname) /= trim(data_table(i)%gridname)) cycle
1208  if( trim(fieldname) /= trim(data_table(i)%fieldname_code)) cycle
1209  index1 = i ! field found
1210  exit
1211  enddo
1212  if(index1 .eq. -1) return ! NO override was performed
1213 
1214  call get_domainug(gridname,ug_domain,comp_domain)
1215  allocate(data_sg(comp_domain(1):comp_domain(2),comp_domain(3):comp_domain(4)))
1216 
1217  call data_override_2d(gridname,fieldname,data_sg,time,override)
1218 
1219  call mpp_pass_sg_to_ug(ug_domain, data_sg(:,:), data(:))
1220 
1221  deallocate(data_sg)
1222 
1223 end subroutine data_override_ug_1d
1224 
1225 subroutine data_override_ug_2d(gridname,fieldname,data,time,override)
1226  character(len=3), intent(in) :: gridname ! model grid ID
1227  character(len=*), intent(in) :: fieldname ! field to override
1228  real, dimension(:,:), intent(inout) :: data !data returned by this call
1229  type(time_type), intent(in) :: time ! model time
1230  logical, intent(out), optional :: override ! true if the field has been overriden succesfully
1231  !local vars
1232  real, dimension(:,:,:), allocatable :: data_SG
1233  real, dimension(:,:), allocatable :: data_UG
1234  type(domainUG) :: UG_domain
1235  integer :: index1
1236  integer :: i, nlevel, nlevel_max
1237  integer, dimension(4) :: comp_domain = 0 ! istart,iend,jstart,jend for compute domain
1238 
1239 !1 Look for the data file in data_table
1240  if(PRESENT(override)) override = .false.
1241  index1 = -1
1242  do i = 1, table_size
1243  if( trim(gridname) /= trim(data_table(i)%gridname)) cycle
1244  if( trim(fieldname) /= trim(data_table(i)%fieldname_code)) cycle
1245  index1 = i ! field found
1246  exit
1247  enddo
1248  if(index1 .eq. -1) return ! NO override was performed
1249 
1250  nlevel = size(data,2)
1251  nlevel_max = nlevel
1252  call mpp_max(nlevel_max)
1253 
1254  call get_domainug(gridname,ug_domain,comp_domain)
1255  allocate(data_sg(comp_domain(1):comp_domain(2),comp_domain(3):comp_domain(4),nlevel_max))
1256  allocate(data_ug(size(data,1), nlevel_max))
1257  data_sg = 0.0
1258  call data_override_3d(gridname,fieldname,data_sg,time,override)
1259 
1260  call mpp_pass_sg_to_ug(ug_domain, data_sg(:,:,:), data_ug(:,:))
1261  data(:,1:nlevel) = data_ug(:,1:nlevel)
1262 
1263  deallocate(data_sg, data_ug)
1264 
1265 end subroutine data_override_ug_2d
1266 
1267 !===============================================================================================
1268 
1269 ! Get lon and lat of three model (target) grids from grid_spec.nc
1270 subroutine get_grid_version_1(grid_file, mod_name, domain, isc, iec, jsc, jec, lon, lat, min_lon, max_lon)
1271  character(len=*), intent(in) :: grid_file
1272  character(len=*), intent(in) :: mod_name
1273  type(domain2d), intent(in) :: domain
1274  integer, intent(in) :: isc, iec, jsc, jec
1275  real, dimension(isc:,jsc:), intent(out) :: lon, lat
1276  real, intent(out) :: min_lon, max_lon
1277 
1278  integer :: i, j, siz(4)
1279  integer :: nlon, nlat ! size of global lon and lat
1280  real, dimension(:,:,:), allocatable :: lon_vert, lat_vert !of OCN grid vertices
1281  real, dimension(:), allocatable :: glon, glat ! lon and lat of 1-D grid of atm/lnd
1282  logical :: is_new_grid
1283  integer :: is, ie, js, je
1284  integer :: isd, ied, jsd, jed
1285  integer :: isg, ieg, jsg, jeg
1286  type(domain2d) :: domain2
1287  character(len=3) :: xname, yname
1288 
1289  call mpp_get_data_domain(domain, isd, ied, jsd, jed)
1290  call mpp_get_global_domain(domain, isg, ieg, jsg, jeg)
1291 
1292  select case(mod_name)
1293  case('ocn', 'ice')
1294  is_new_grid = .false.
1295  if(field_exist(grid_file, 'x_T')) then
1296  is_new_grid = .true.
1297  else if(field_exist(grid_file, 'geolon_t')) then
1298  is_new_grid = .false.
1299  else
1300  call mpp_error(fatal,'data_override: both x_T and geolon_t is not in the grid file '//trim(grid_file) )
1301  endif
1302 
1303  if(is_new_grid) then
1304  call field_size(grid_file, 'x_T', siz)
1305  nlon = siz(1); nlat = siz(2)
1306  call check_grid_sizes(trim(mod_name)//'_domain ', domain, nlon, nlat)
1307  allocate(lon_vert(isc:iec,jsc:jec,4), lat_vert(isc:iec,jsc:jec,4) )
1308  call read_data(trim(grid_file), 'x_vert_T', lon_vert, domain)
1309  call read_data(trim(grid_file), 'y_vert_T', lat_vert, domain)
1310 
1311 !2 Global lon and lat of ocean grid cell centers are determined from adjacent vertices
1312  lon(:,:) = (lon_vert(:,:,1) + lon_vert(:,:,2) + lon_vert(:,:,3) + lon_vert(:,:,4))*0.25
1313  lat(:,:) = (lat_vert(:,:,1) + lat_vert(:,:,2) + lat_vert(:,:,3) + lat_vert(:,:,4))*0.25
1314  else
1315  if(grid_center_bug) call mpp_error(note, &
1316  'data_override: grid_center_bug is set to true, the grid center location may be incorrect')
1317  call field_size(grid_file, 'geolon_vert_t', siz)
1318  nlon = siz(1) - 1; nlat = siz(2) - 1;
1319  call check_grid_sizes(trim(mod_name)//'_domain ', domain, nlon, nlat)
1320  call mpp_copy_domain(domain, domain2)
1321  call mpp_set_compute_domain(domain2, isc, iec+1, jsc, jec+1, iec-isc+2, jec-jsc+2 )
1322  call mpp_set_data_domain (domain2, isd, ied+1, jsd, jed+1, ied-isd+2, jed-jsd+2 )
1323  call mpp_set_global_domain (domain2, isg, ieg+1, jsg, jeg+1, ieg-isg+2, jeg-jsg+2 )
1324  allocate(lon_vert(isc:iec+1,jsc:jec+1,1))
1325  allocate(lat_vert(isc:iec+1,jsc:jec+1,1))
1326  call read_data(trim(grid_file), 'geolon_vert_t', lon_vert, domain2)
1327  call read_data(trim(grid_file), 'geolat_vert_t', lat_vert, domain2)
1328 
1329  if(grid_center_bug) then
1330  do j = jsc, jec
1331  do i = isc, iec
1332  lon(i,j) = (lon_vert(i,j,1) + lon_vert(i+1,j,1))/2.
1333  lat(i,j) = (lat_vert(i,j,1) + lat_vert(i,j+1,1))/2.
1334  enddo
1335  enddo
1336  else
1337  do j = jsc, jec
1338  do i = isc, iec
1339  lon(i,j) = (lon_vert(i,j,1) + lon_vert(i+1,j,1) + &
1340  lon_vert(i+1,j+1,1) + lon_vert(i,j+1,1))*0.25
1341  lat(i,j) = (lat_vert(i,j,1) + lat_vert(i+1,j,1) + &
1342  lat_vert(i+1,j+1,1) + lat_vert(i,j+1,1))*0.25
1343  enddo
1344  enddo
1345  end if
1346  call mpp_deallocate_domain(domain2)
1347  endif
1348  deallocate(lon_vert)
1349  deallocate(lat_vert)
1350  case('atm', 'lnd')
1351  if(trim(mod_name) == 'atm') then
1352  xname = 'xta'; yname = 'yta'
1353  else
1354  xname = 'xtl'; yname = 'ytl'
1355  endif
1356  call field_size(grid_file, xname, siz)
1357  nlon = siz(1); allocate(glon(nlon))
1358  call read_data(grid_file, xname, glon, no_domain = .true.)
1359 
1360  call field_size(grid_file, yname, siz)
1361  nlat = siz(1); allocate(glat(nlat))
1362  call read_data(grid_file, yname, glat, no_domain = .true.)
1363  call check_grid_sizes(trim(mod_name)//'_domain ', domain, nlon, nlat)
1364 
1365  is = isc - isg + 1; ie = iec - isg + 1
1366  js = jsc - jsg + 1; je = jec - jsg + 1
1367  do j = js, jec
1368  do i = is, ie
1369  lon(i,j) = glon(i)
1370  lat(i,j) = glat(j)
1371  enddo
1372  enddo
1373  deallocate(glon)
1374  deallocate(glat)
1375  case default
1376  call mpp_error(fatal, "data_override_mod: mod_name should be 'atm', 'ocn', 'ice' or 'lnd' ")
1377  end select
1378 
1379  ! convert from degree to radian
1380  lon = lon * deg_to_radian
1381  lat = lat* deg_to_radian
1382  min_lon = minval(lon)
1383  max_lon = maxval(lon)
1384  call mpp_min(min_lon)
1385  call mpp_max(max_lon)
1386 
1387 
1388 end subroutine get_grid_version_1
1389 
1390 ! Get global lon and lat of three model (target) grids from mosaic.nc
1391 ! z1l: currently we assume the refinement ratio is 2 and there is one tile on each pe.
1392 subroutine get_grid_version_2(mosaic_file, mod_name, domain, isc, iec, jsc, jec, lon, lat, min_lon, max_lon)
1393  character(len=*), intent(in) :: mosaic_file
1394  character(len=*), intent(in) :: mod_name
1395  type(domain2d), intent(in) :: domain
1396  integer, intent(in) :: isc, iec, jsc, jec
1397  real, dimension(isc:,jsc:), intent(out) :: lon, lat
1398  real, intent(out) :: min_lon, max_lon
1399 
1400  integer :: i, j, siz(4)
1401  integer :: nlon, nlat ! size of global grid
1402  integer :: nlon_super, nlat_super ! size of global supergrid.
1403  integer :: isd, ied, jsd, jed
1404  integer :: isg, ieg, jsg, jeg
1405  integer :: isc2, iec2, jsc2, jec2
1406  character(len=256) :: solo_mosaic_file, grid_file
1407  real, allocatable :: tmpx(:,:), tmpy(:,:)
1408  type(domain2d) :: domain2
1409 
1410  if(trim(mod_name) .NE. 'atm' .AND. trim(mod_name) .NE. 'ocn' .AND. &
1411  trim(mod_name) .NE. 'ice' .AND. trim(mod_name) .NE. 'lnd' ) call mpp_error(fatal, &
1412  "data_override_mod: mod_name should be 'atm', 'ocn', 'ice' or 'lnd' ")
1413 
1414  call mpp_get_data_domain(domain, isd, ied, jsd, jed)
1415  call mpp_get_global_domain(domain, isg, ieg, jsg, jeg)
1416 
1417  ! get the grid file to read
1418  if(field_exist(mosaic_file, trim(mod_name)//'_mosaic_file' )) then
1419  call read_data(mosaic_file, trim(mod_name)//'_mosaic_file', solo_mosaic_file)
1420  solo_mosaic_file = 'INPUT/'//trim(solo_mosaic_file)
1421  else
1422  solo_mosaic_file = mosaic_file
1423  end if
1424  call get_mosaic_tile_grid(grid_file, solo_mosaic_file, domain)
1425 
1426  call field_size(grid_file, 'area', siz)
1427  nlon_super = siz(1); nlat_super = siz(2)
1428  if( mod(nlon_super,2) .NE. 0) call mpp_error(fatal, &
1429  'data_override_mod: '//trim(mod_name)//' supergrid longitude size can not be divided by 2')
1430  if( mod(nlat_super,2) .NE. 0) call mpp_error(fatal, &
1431  'data_override_mod: '//trim(mod_name)//' supergrid latitude size can not be divided by 2')
1432  nlon = nlon_super/2;
1433  nlat = nlat_super/2;
1434  call check_grid_sizes(trim(mod_name)//'_domain ', domain, nlon, nlat)
1435 
1436  !--- setup the domain for super grid.
1437  call mpp_copy_domain(domain, domain2)
1438  call mpp_set_compute_domain(domain2, 2*isc-1, 2*iec+1, 2*jsc-1, 2*jec+1, 2*iec-2*isc+3, 2*jec-2*jsc+3 )
1439  call mpp_set_data_domain (domain2, 2*isd-1, 2*ied+1, 2*jsd-1, 2*jed+1, 2*ied-2*isd+3, 2*jed-2*jsd+3 )
1440  call mpp_set_global_domain (domain2, 2*isg-1, 2*ieg+1, 2*jsg-1, 2*jeg+1, 2*ieg-2*isg+3, 2*jeg-2*jsg+3 )
1441 
1442  call mpp_get_compute_domain(domain2, isc2, iec2, jsc2, jec2)
1443  if(isc2 .NE. 2*isc-1 .OR. iec2 .NE. 2*iec+1 .OR. jsc2 .NE. 2*jsc-1 .OR. jec2 .NE. 2*jec+1) then
1444  call mpp_error(fatal, 'data_override_mod: '//trim(mod_name)//' supergrid domain is not set properly')
1445  endif
1446 
1447  allocate(tmpx(isc2:iec2, jsc2:jec2), tmpy(isc2:iec2, jsc2:jec2) )
1448  call read_data( grid_file, 'x', tmpx, domain2)
1449  call read_data( grid_file, 'y', tmpy, domain2)
1450  ! copy data onto model grid
1451  if(trim(mod_name) == 'ocn' .OR. trim(mod_name) == 'ice') then
1452  do j = jsc, jec
1453  do i = isc, iec
1454  lon(i,j) = (tmpx(i*2-1,j*2-1)+tmpx(i*2+1,j*2-1)+tmpx(i*2+1,j*2+1)+tmpx(i*2-1,j*2+1))*0.25
1455  lat(i,j) = (tmpy(i*2-1,j*2-1)+tmpy(i*2+1,j*2-1)+tmpy(i*2+1,j*2+1)+tmpy(i*2-1,j*2+1))*0.25
1456  end do
1457  end do
1458  else
1459  do j = jsc, jec
1460  do i = isc, iec
1461  lon(i,j) = tmpx(i*2,j*2)
1462  lat(i,j) = tmpy(i*2,j*2)
1463  end do
1464  end do
1465  endif
1466 
1467  ! convert to radian
1468  lon = lon * deg_to_radian
1469  lat = lat * deg_to_radian
1470 
1471  deallocate(tmpx, tmpy)
1472  min_lon = minval(lon)
1473  max_lon = maxval(lon)
1474  call mpp_min(min_lon)
1475  call mpp_max(max_lon)
1476 
1477  call mpp_deallocate_domain(domain2)
1478 
1479 end subroutine get_grid_version_2
1480 
1481 !===============================================================================================
1482 end module data_override_mod
1483 
1484 #ifdef test_data_override
1485 
1486  program test
1487 
1488  ! Input data and path_names file for this program is in:
1489  ! /archive/pjp/unit_tests/test_data_override/lima/exp1
1490 
1491  use mpp_mod, only: input_nml_file, stdout, mpp_chksum
1493  use fms_mod, only: fms_init, fms_end, mpp_npes, file_exist, open_namelist_file, check_nml_error, close_file
1494  use fms_mod, only: error_mesg, fatal, file_exist, field_exist, field_size
1495  use fms_io_mod, only: read_data, fms_io_exit
1496  use constants_mod, only: constants_init, pi
1499  use diag_manager_mod, only: send_data, diag_axis_init
1501  use mpp_mod, only : fatal, warning, mpp_debug, note, mpp_clock_sync,mpp_clock_detailed
1502  use mpp_mod, only : mpp_pe, mpp_npes, mpp_node, mpp_root_pe, mpp_error, mpp_set_warn_level
1503  use mpp_mod, only : mpp_declare_pelist, mpp_set_current_pelist, mpp_sync, mpp_sync_self
1504  use mpp_mod, only : mpp_clock_begin, mpp_clock_end, mpp_clock_id
1505  use mpp_mod, only : mpp_init, mpp_exit, mpp_chksum, stdout, stderr
1506  use mpp_mod, only : input_nml_file
1507  use mpp_mod, only : mpp_get_current_pelist, mpp_broadcast
1508  use mpp_domains_mod, only : global_data_domain, bitwise_exact_sum, bgrid_ne, cgrid_ne, dgrid_ne
1509  use mpp_domains_mod, only : fold_south_edge, fold_north_edge, fold_west_edge, fold_east_edge
1510  use mpp_domains_mod, only : mpp_domain_time, cyclic_global_domain, nupdate,eupdate, xupdate, yupdate, scalar_pair
1511  use mpp_domains_mod, only : domain1d, domain2d, domaincommunicator2d, bitwise_efp_sum
1512  use mpp_domains_mod, only : mpp_get_compute_domain, mpp_get_data_domain, mpp_domains_set_stack_size
1514  use mpp_domains_mod, only : mpp_domains_init, mpp_domains_exit, mpp_broadcast_domain
1517  use mpp_domains_mod, only : mpp_get_neighbor_pe, mpp_define_mosaic, mpp_nullify_domain_list
1518  use mpp_domains_mod, only : north, north_east, east, south_east, corner, center
1519  use mpp_domains_mod, only : south, south_west, west, north_west, mpp_define_mosaic_pelist
1520  use mpp_domains_mod, only : mpp_get_global_domain, zero, ninety, minus_ninety
1522  use mpp_domains_mod, only : mpp_define_nest_domains, nest_domain_type
1523  use mpp_domains_mod, only : mpp_get_c2f_index, mpp_update_nest_fine
1524  use mpp_domains_mod, only : mpp_get_f2c_index, mpp_update_nest_coarse
1525  use mpp_domains_mod, only : mpp_get_domain_shift, edgeupdate, mpp_deallocate_domain
1527  use mpp_domains_mod, only : mpp_do_group_update, mpp_clear_group_update
1529  use mpp_domains_mod, only : wupdate, supdate, mpp_get_compute_domains
1530  use mpp_domains_mod, only : domainug, mpp_define_unstruct_domain, mpp_get_ug_domain_tile_id
1531  use mpp_domains_mod, only : mpp_get_ug_compute_domain, mpp_pass_sg_to_ug, mpp_pass_ug_to_sg
1532  use mpp_domains_mod, only : mpp_get_ug_global_domain, mpp_global_field_ug
1534 
1535  implicit none
1536 
1537  integer :: stdoutunit
1538  integer :: num_threads = 1
1539  integer :: omp_get_thread_num
1540  integer :: isw, iew, jsw, jew
1541  integer, allocatable :: is_win(:), js_win(:)
1542  integer :: nx_dom, ny_dom, nx_win, ny_win
1543  type(domain2d) :: domain
1544  integer :: nlon, nlat, siz(4)
1545  real, allocatable, dimension(:) :: x, y
1546  real, allocatable, dimension(:,:) :: lon, lat
1547  real, allocatable, dimension(:,:) :: sst, ice
1548  integer :: id_x, id_y, id_lon, id_lat, id_sst, id_ice
1549  integer :: i, j, is, ie, js, je, unit, io, ierr, n
1550  real :: rad_to_deg
1551  character(len=36) :: message
1552  type(time_type) :: time
1553  logical :: used
1554  logical, allocatable :: ov_sst(:), ov_ice(:)
1555  integer, dimension(2) :: layout = (/0,0/)
1556  character(len=256) :: solo_mosaic_file, tile_file
1557  character(len=128) :: grid_file = "INPUT/grid_spec.nc"
1558  integer :: window(2) = (/1,1/)
1559  integer :: get_cpu_affinity, base_cpu
1560  integer :: nthreads=1
1561  integer :: nwindows
1562  integer :: nx_cubic=0, ny_cubic=0, nx_latlon=0, ny_latlon=0
1563 
1564  namelist / test_data_override_nml / layout, window, nthreads, nx_cubic, ny_cubic, nx_latlon, ny_latlon
1565 
1566  call fms_init
1567  call constants_init
1569  call diag_manager_init
1570 
1571  rad_to_deg = 180./pi
1572 
1573 #ifdef INTERNAL_FILE_NML
1574  read (input_nml_file, test_data_override_nml, iostat=io)
1575  ierr = check_nml_error(io, 'test_data_override_nml')
1576 #else
1577  if (file_exist('input.nml')) then
1578  unit = open_namelist_file( )
1579  ierr=1
1580  do while (ierr /= 0)
1581  read(unit, nml=test_data_override_nml, iostat=io, end=10)
1582  ierr = check_nml_error(io, 'test_data_override_nml')
1583  enddo
1584 10 call close_file (unit)
1585  endif
1586 #endif
1587 
1588  if(field_exist(grid_file, "x_T" ) ) then
1589  call field_size(grid_file, 'x_T', siz)
1590  nlon = siz(1)
1591  nlat = siz(2)
1592  else if(field_exist(grid_file, "geolon_t" ) ) then
1593  call field_size(grid_file, 'geolon_t', siz)
1594  nlon = siz(1)
1595  nlat = siz(2)
1596  else if (field_exist(grid_file, "ocn_mosaic_file" )) then
1597  call read_data(grid_file, 'ocn_mosaic_file', solo_mosaic_file)
1598  solo_mosaic_file = 'INPUT/'//trim(solo_mosaic_file)
1599  call field_size(solo_mosaic_file, 'gridfiles', siz)
1600  if( siz(2) .NE. 1) call error_mesg('test_data_override', 'only support single tile mosaic, contact developer', fatal)
1601  call read_data(solo_mosaic_file, 'gridfiles', tile_file)
1602  tile_file = 'INPUT/'//trim(tile_file)
1603  call field_size(tile_file, 'area', siz)
1604  if(mod(siz(1),2) .NE. 0 .OR. mod(siz(2),2) .NE. 0 ) call error_mesg('test_data_override', &
1605  "test_data_override: supergrid size can not be divided by 2", fatal)
1606  nlon = siz(1)/2
1607  nlat = siz(2)/2
1608  else
1609  call error_mesg('test_data_override', 'x_T, geolon_t and ocn_mosaic_file does not exist', fatal)
1610  end if
1611 
1612  if(layout(1)*layout(2) .NE. mpp_npes() ) then
1613  call mpp_define_layout( (/1,nlon,1,nlat/), mpp_npes(), layout )
1614  end if
1615 
1616 
1617  call mpp_define_domains( (/1,nlon,1,nlat/), layout, domain, name='test_data_override')
1618  call data_override_init(ice_domain_in=domain, ocean_domain_in=domain)
1619  call data_override_init(ice_domain_in=domain, ocean_domain_in=domain)
1620  call mpp_get_compute_domain(domain, is, ie, js, je)
1621  call get_grid
1622 
1623  allocate(x(nlon), y(nlat))
1624 
1625  do i=1,nlon
1626  x(i) = i
1627  enddo
1628  do j=1,nlat
1629  y(j) = j
1630  enddo
1631 
1632  time = set_date(2,1,1,0,0,0)
1633 
1634  allocate(sst(is:ie,js:je), ice(is:ie,js:je))
1635  sst = 0
1636  ice = 0
1637 
1638  id_x = diag_axis_init('x', x, 'point_E', 'x', long_name='point_E', domain2=domain)
1639  id_y = diag_axis_init('y', y, 'point_N', 'y', long_name='point_N', domain2=domain)
1640 
1641  time = time + set_time(0,1)
1642 
1643  id_lon = register_static_field('test_data_override_mod', 'lon', (/id_x,id_y/), 'longitude', 'Degrees')
1644  id_lat = register_static_field('test_data_override_mod', 'lat', (/id_x,id_y/), 'longitude', 'Degrees')
1645  id_sst = register_diag_field('test_data_override_mod', 'sst', (/id_x,id_y/), time, 'SST', 'K')
1646  id_ice = register_diag_field('test_data_override_mod', 'ice', (/id_x,id_y/), time, 'ICE', ' ')
1647  used = send_data(id_lon, lon, time)
1648  used = send_data(id_lat, lat, time)
1649 
1650  sst = 0.
1651  ice = 0.
1652 
1653 ! get number of threads
1654 
1655 nx_dom = ie - is + 1
1656 ny_dom = je - js + 1
1657 if( mod( nx_dom, window(1) ) .NE. 0 ) call error_mesg('test_data_override', &
1658  "nx_dom is not divisible by window(1)", fatal)
1659 if( mod( ny_dom, window(2) ) .NE. 0 ) call error_mesg('test_data_override', &
1660  "ny_dom is not divisible by window(2)", fatal)
1661 
1662 nwindows = window(1)*window(2)
1663 !$ call omp_set_num_threads(nthreads)
1664 !$ base_cpu = get_cpu_affinity()
1665 !$OMP PARALLEL
1666 !$ call set_cpu_affinity( base_cpu + omp_get_thread_num() )
1667 !$OMP END PARALLEL
1668 
1669 nx_win = nx_dom/window(1)
1670 ny_win = ny_dom/window(2)
1671 allocate(is_win(nwindows), js_win(nwindows))
1672 
1673 i = 1
1674 do jsw = js,je,ny_win
1675  do isw = is,ie,nx_win
1676  is_win(i) = isw
1677  js_win(i) = jsw
1678  i = i + 1
1679  enddo
1680 enddo
1681 
1682 allocate(ov_sst(nwindows), ov_ice(nwindows))
1683 !$OMP parallel do schedule(static) default(shared) private(isw, iew, jsw, jew)
1684 do n = 1, nwindows
1685  isw = is_win(n)
1686  iew = isw + nx_win - 1
1687  jsw = js_win(n)
1688  jew = jsw + ny_win - 1
1689  call data_override('OCN','sst_obs',sst(isw:iew,jsw:jew),time,override=ov_sst(n), &
1690  is_in=isw-is+1, ie_in=iew-is+1, js_in=jsw-js+1, je_in=jew-js+1)
1691  call data_override('ICE', 'sic_obs', ice(isw:iew,jsw:jew), time, override=ov_ice(n), &
1692  is_in=isw-is+1, ie_in=iew-is+1, js_in=jsw-js+1, je_in=jew-js+1)
1693 enddo
1694 
1695  if(any(.NOT. ov_sst) .or. any(.not.ov_ice)) then
1696  if(any(.NOT. ov_sst)) then
1697  message = 'override failed for ice'
1698  else if(any(.NOT. ov_ice)) then
1699  message = 'override failed for sst'
1700  else
1701  message = 'override failed for both sst and ice'
1702  endif
1703  call error_mesg('test_data_override', trim(message), fatal)
1704  endif
1705 
1706  stdoutunit = stdout()
1707  write(stdoutunit,*)"===>NOTE from test_data_override: sst chksum = ", mpp_chksum(sst)
1708  write(stdoutunit,*)"===>NOTE from test_data_override: ice chksum = ", mpp_chksum(ice)
1709 
1710 
1711  if(id_sst > 0) used = send_data(id_sst, sst, time)
1712  if(id_ice > 0) used = send_data(id_ice, ice, time)
1713 
1714  if(nx_cubic > 0 .and. ny_cubic > 0) then
1715  call test_unstruct_grid( 'Cubic-Grid', time )
1716  endif
1717  if(nx_latlon > 0 .and. ny_latlon > 0) then
1718  call test_unstruct_grid( 'Latlon-Grid', time )
1719  endif
1720 
1721 
1722 
1723 
1724 !-------------------------------------------------------------------------------------------------------
1725 ! What follows is a test of calendar conversion
1726 
1727 !Time = set_date(1980,2,27,0,0,0)
1728 !call print_time(Time)
1729 !call data_override('OCN','sst_obs',sst,Time)
1730 !if(id_sst > 0) used = send_data(id_sst, sst, Time)
1731 
1732 !Time = set_date(1980,2,28,0,0,0)
1733 !call print_time(Time)
1734 !call data_override('OCN','sst_obs',sst,Time)
1735 !if(id_sst > 0) used = send_data(id_sst, sst, Time)
1736 
1737 !Time = set_date(1980,2,29,0,0,0)
1738 !call print_time(Time)
1739 !call data_override('OCN','sst_obs',sst,Time)
1740 !if(id_sst > 0) used = send_data(id_sst, sst, Time)
1741 
1742 !Time = set_date(1980,3,1,0,0,0)
1743 !call print_time(Time)
1744 !call data_override('OCN','sst_obs',sst,Time)
1745 !if(id_sst > 0) used = send_data(id_sst, sst, Time)
1746 
1747 !Time = set_date(1980,3,2,0,0,0)
1748 !call print_time(Time)
1749 !call data_override('OCN','sst_obs',sst,Time)
1750 !if(id_sst > 0) used = send_data(id_sst, sst, Time)
1751 !-------------------------------------------------------------------------------------------------------
1752 
1753  call diag_manager_end(time)
1754  call fms_io_exit
1755  call fms_end
1756 
1757 contains
1758 
1759 !=================================================================================================================================
1760  subroutine get_grid
1761  real, allocatable, dimension(:,:,:) :: lon_vert_glo, lat_vert_glo
1762  real, allocatable, dimension(:,:) :: lon_global, lat_global
1763  integer, dimension(4) :: siz
1764  character(len=128) :: message
1765 
1766 
1767  if(field_exist(grid_file, 'x_T')) then
1768  call field_size(grid_file, 'x_T', siz)
1769  if(siz(1) /= nlon .or. siz(2) /= nlat) then
1770  write(message,'(a,2i4)') 'x_T is wrong shape. shape(x_T)=',siz(1:2)
1771  call error_mesg('test_data_override', trim(message), fatal)
1772  endif
1773  allocate(lon_vert_glo(nlon,nlat,4), lat_vert_glo(nlon,nlat,4) )
1774  allocate(lon_global(nlon,nlat ), lat_global(nlon,nlat ) )
1775  call read_data(trim(grid_file), 'x_vert_T', lon_vert_glo, no_domain=.true.)
1776  call read_data(trim(grid_file), 'y_vert_T', lat_vert_glo, no_domain=.true.)
1777  lon_global(:,:) = (lon_vert_glo(:,:,1) + lon_vert_glo(:,:,2) + lon_vert_glo(:,:,3) + lon_vert_glo(:,:,4))*0.25
1778  lat_global(:,:) = (lat_vert_glo(:,:,1) + lat_vert_glo(:,:,2) + lat_vert_glo(:,:,3) + lat_vert_glo(:,:,4))*0.25
1779  else if(field_exist(grid_file, "geolon_t" ) ) then
1780  call field_size(grid_file, 'geolon_vert_t', siz)
1781  if(siz(1) /= nlon+1 .or. siz(2) /= nlat+1) then
1782  write(message,'(a,2i4)') 'geolon_vert_t is wrong shape. shape(geolon_vert_t)=',siz(1:2)
1783  call error_mesg('test_data_override', trim(message), fatal)
1784  endif
1785  allocate(lon_vert_glo(nlon+1,nlat+1,1), lat_vert_glo(nlon+1,nlat+1,1))
1786  allocate(lon_global(nlon, nlat ), lat_global(nlon, nlat ))
1787  call read_data(trim(grid_file), 'geolon_vert_t', lon_vert_glo, no_domain=.true.)
1788  call read_data(trim(grid_file), 'geolat_vert_t', lat_vert_glo, no_domain=.true.)
1789 
1790  do i = 1, nlon
1791  do j = 1, nlat
1792  lon_global(i,j) = (lon_vert_glo(i,j,1) + lon_vert_glo(i+1,j,1) + &
1793  lon_vert_glo(i+1,j+1,1) + lon_vert_glo(i,j+1,1))*0.25
1794  lat_global(i,j) = (lat_vert_glo(i,j,1) + lat_vert_glo(i+1,j,1) + &
1795  lat_vert_glo(i+1,j+1,1) + lat_vert_glo(i,j+1,1))*0.25
1796  enddo
1797  enddo
1798  else if( field_exist(grid_file, "ocn_mosaic_file") ) then ! reading from mosaic file
1799  call field_size(tile_file, 'area', siz)
1800  if(siz(1) /= nlon*2 .or. siz(2) /= nlat*2) then
1801  write(message,'(a,2i4)') 'area is wrong shape. shape(area)=',siz(1:2)
1802  call error_mesg('test_data_override', trim(message), fatal)
1803  endif
1804  allocate(lon_vert_glo(siz(1)+1,siz(2)+1,1), lat_vert_glo(siz(1)+1,siz(2)+1,1))
1805  allocate(lon_global(nlon, nlat ), lat_global(nlon, nlat ))
1806  call read_data( tile_file, 'x', lon_vert_glo, no_domain=.true.)
1807  call read_data( tile_file, 'y', lat_vert_glo, no_domain=.true.)
1808  do j = 1, nlat
1809  do i = 1, nlon
1810  lon_global(i,j) = lon_vert_glo(i*2,j*2,1)
1811  lat_global(i,j) = lat_vert_glo(i*2,j*2,1)
1812  end do
1813  end do
1814  end if
1815 
1816  allocate(lon(is:ie,js:je), lat(is:ie,js:je))
1817  lon = lon_global(is:ie,js:je)
1818  lat = lat_global(is:ie,js:je)
1819 
1820  deallocate(lon_vert_glo)
1821  deallocate(lat_vert_glo)
1822  deallocate(lon_global)
1823  deallocate(lat_global)
1824 
1825  end subroutine get_grid
1826 
1827  subroutine test_unstruct_grid( type, Time )
1829  character(len=*), intent(in) :: type
1830  type(time_type), intent(in) :: Time !(target) model time
1831 
1832  integer :: pe, npes
1833  integer :: nx, ny, nz=40, stackmax=4000000
1834  integer :: unit=7
1835  integer :: stdunit = 6
1836  logical :: debug=.false., opened
1837 
1838  integer :: mpes = 0
1839  integer :: whalo = 2, ehalo = 2, shalo = 2, nhalo = 2
1840  character(len=32) :: warn_level = "fatal"
1841  integer :: layout_cubic(2) = (/0,0/)
1842  integer :: layout_tripolar(2) = (/0,0/)
1843  integer :: layout_ensemble(2) = (/0,0/)
1844  logical :: do_sleep = .false.
1845  integer :: num_iter = 1
1846  integer :: num_fields = 4
1847 
1848  type(domain2D) :: SG_domain
1849  type(domainUG) :: UG_domain
1850  integer :: num_contact, ntiles, npes_per_tile
1851  integer :: isc, iec, jsc, jec, isd, ied, jsd, jed
1852  integer :: ism, iem, jsm, jem, lsg, leg
1853 
1854  integer, allocatable, dimension(:) :: pe_start, pe_end, npts_tile, grid_index, ntiles_grid
1855  integer, allocatable, dimension(:,:) :: layout2D, global_indices
1856  real, allocatable, dimension(:,:) :: x1, x2, g1, g2
1857  real, allocatable, dimension(:,:,:) :: a1, a2, gdata
1858  real, allocatable, dimension(:,:) :: rmask
1859  real, allocatable, dimension(:) :: frac_crit
1860  logical, allocatable, dimension(:,:,:) :: lmask,msk
1861  integer, allocatable, dimension(:) :: isl, iel, jsl, jel
1862  character(len=3) :: text
1863  integer :: tile
1864  integer :: ntotal_land, istart, iend, pos
1865  integer :: outunit, errunit, k, l
1866 
1867  call mpp_memuse_begin()
1868  npes = mpp_npes()
1869 
1870  outunit = stdout()
1871  errunit = stderr()
1872  !--- check the type
1873  select case(type)
1874  case ( 'Cubic-Grid' )
1875  if( nx_cubic == 0 ) then
1876  call mpp_error(note,'test_unstruct_update: for Cubic_grid mosaic, nx_cubic is zero, '//&
1877  'No test is done for Cubic-Grid mosaic. ' )
1878  return
1879  endif
1880  if( nx_cubic .NE. ny_cubic ) then
1881  call mpp_error(note,'test_unstruct_update: for Cubic_grid mosaic, nx_cubic does not equal ny_cubic, '//&
1882  'No test is done for Cubic-Grid mosaic. ' )
1883  return
1884  endif
1885  nx = nx_cubic
1886  ny = ny_cubic
1887  ntiles = 6
1888  num_contact = 12
1889  if( mod(npes, ntiles) == 0 ) then
1890  npes_per_tile = npes/ntiles
1891  write(outunit,*)'NOTE from test_unstruct_update ==> For Mosaic "', trim(type), &
1892  '", each tile will be distributed over ', npes_per_tile, ' processors.'
1893  else
1894  call mpp_error(note,'test_unstruct_update: npes should be multiple of ntiles No test is done for '//trim(type))
1895  return
1896  endif
1897  if(layout_cubic(1)*layout_cubic(2) == npes_per_tile) then
1898  layout = layout_cubic
1899  else
1900  call mpp_define_layout( (/1,nx,1,ny/), npes_per_tile, layout )
1901  endif
1902  allocate(frac_crit(ntiles))
1903  frac_crit(1) = 0.3; frac_crit(2) = 0.1; frac_crit(3) = 0.6
1904  frac_crit(4) = 0.2; frac_crit(5) = 0.4; frac_crit(6) = 0.5
1905  allocate(layout2d(2,ntiles), global_indices(4,ntiles), pe_start(ntiles), pe_end(ntiles) )
1906  do n = 1, ntiles
1907  pe_start(n) = (n-1)*npes_per_tile
1908  pe_end(n) = n*npes_per_tile-1
1909  end do
1910 
1911  do n = 1, ntiles
1912  global_indices(:,n) = (/1,nx,1,ny/)
1913  layout2d(:,n) = layout
1914  end do
1915 
1916  call define_cubic_mosaic(type, SG_domain, (/nx,nx,nx,nx,nx,nx/), (/ny,ny,ny,ny,ny,ny/), &
1917  global_indices, layout2D, pe_start, pe_end )
1918  case ( 'Latlon-Grid' )
1919  if(nx_latlon == 0 .OR. ny_latlon == 0 ) then
1920  call mpp_error(note,'test_unstruct_update: for latlon mosaic, nx_latlon and ny_latlon are zero, '//&
1921  'No test is done for Lalton-Grid mosaic. ' )
1922  return
1923  endif
1924  nx = nx_latlon
1925  ny = ny_latlon
1926  ntiles = 1
1927  npes_per_tile = npes
1928  allocate(frac_crit(ntiles))
1929  frac_crit(1) = 0.3
1930  call mpp_define_layout((/1,nx,1,ny/), npes, layout)
1931  call mpp_define_domains((/1,nx,1,ny/), layout, sg_domain, xflags = cyclic_global_domain)
1932  case default
1933  call mpp_error(fatal, 'test_group_update: no such test: '//type)
1934  end select
1935 
1936  !--- setup data
1937  call mpp_get_compute_domain( sg_domain, isc, iec, jsc, jec )
1938  call mpp_get_data_domain ( sg_domain, isd, ied, jsd, jed )
1939 
1940  allocate(lmask(nx,ny,ntiles))
1941  allocate(npts_tile(ntiles))
1942  lmask = .false.
1943  if(mpp_pe() == mpp_root_pe() ) then
1944  allocate(rmask(nx,ny))
1945  !--- construct gmask.
1946  do n = 1, ntiles
1947  call random_number(rmask)
1948  do j = 1, ny
1949  do i = 1, nx
1950  if(rmask(i,j) > frac_crit(n)) then
1951  lmask(i,j,n) = .true.
1952  endif
1953  enddo
1954  enddo
1955  npts_tile(n) = count(lmask(:,:,n))
1956  enddo
1957  ntotal_land = sum(npts_tile)
1958  allocate(grid_index(ntotal_land))
1959  l = 0
1960  allocate(isl(0:mpp_npes()-1), iel(0:mpp_npes()-1))
1961  allocate(jsl(0:mpp_npes()-1), jel(0:mpp_npes()-1))
1962  call mpp_get_compute_domains(sg_domain,xbegin=isl,xend=iel,ybegin=jsl,yend=jel)
1963 
1964  do n = 1, ntiles
1965  do j = 1, ny
1966  do i = 1, nx
1967  if(lmask(i,j,n)) then
1968  l = l + 1
1969  grid_index(l) = (j-1)*nx+i
1970  endif
1971  enddo
1972  enddo
1973  enddo
1974  deallocate(rmask, isl, iel, jsl, jel)
1975  endif
1976  call mpp_broadcast(npts_tile, ntiles, mpp_root_pe())
1977  if(mpp_pe() .NE. mpp_root_pe()) then
1978  ntotal_land = sum(npts_tile)
1979  allocate(grid_index(ntotal_land))
1980  endif
1981  call mpp_broadcast(grid_index, ntotal_land, mpp_root_pe())
1982 
1983  allocate(ntiles_grid(ntotal_land))
1984  ntiles_grid = 1
1985  !--- define the unstructured grid domain
1986  call mpp_define_unstruct_domain(ug_domain, sg_domain, npts_tile, ntiles_grid, mpp_npes(), 1, grid_index, name="LAND unstruct")
1987  call mpp_get_ug_compute_domain(ug_domain, istart, iend)
1988 
1989  !--- figure out lmask according to grid_index
1990  pos = 0
1991  do n = 1, ntiles
1992  do l = 1, npts_tile(n)
1993  pos = pos + 1
1994  j = (grid_index(pos)-1)/nx + 1
1995  i = mod((grid_index(pos)-1),nx) + 1
1996  lmask(i,j,n) = .true.
1997  enddo
1998  enddo
1999 
2000  !--- test the 2-D data is on computing domain
2001  allocate( a1(isc:iec, jsc:jec,1), a2(isc:iec,jsc:jec,1 ) )
2002  allocate(msk(isc:iec, jsc:jec,1)); msk = .false.
2003 
2004  tile = mpp_pe()/npes_per_tile + 1
2005  do j = jsc, jec
2006  do i = isc, iec
2007  msk(i,j,1) = lmask(i,j,tile)
2008  enddo
2009  enddo
2010  !First override the test SG data from file/field
2011  call data_override_init(land_domain_in=sg_domain)
2012  call data_override('LND','sst_obs',a1(:,:,1),time)
2013 
2014  !Create the test UG data
2015  a2 = -9999
2016  do j = jsc, jec
2017  do i = isc, iec
2018  if(.NOT. msk(i,j,1)) a2(i,j,1)=a1(i,j,1)
2019  enddo
2020  enddo
2021 
2022  allocate(x1(istart:iend,1), x2(istart:iend,1))
2023  x1 = -99999
2024  x2 = -999999
2025  !--- fill the value of x2
2026 
2027  !Now override the test UG data from the same file/field
2028  call data_override_init(land_domainug_in=ug_domain)
2029  call data_override_ug('LND','sst_obs',x2(:,1),time)
2030 
2031  !Ensure you get the same UG data from the SG data
2032  call mpp_pass_sg_to_ug(ug_domain, a1(:,:,1), x1(:,1))
2033  call compare_checksums_2d(x1, x2, type//' SG2UG 2-D compute domain')
2034 
2035  !Ensure you get the same SG data from the UG data if you transform back
2036  call mpp_pass_ug_to_sg(ug_domain, x1(:,1), a2(:,:,1))
2037  call compare_checksums(a1(:,:,1:1),a2(:,:,1:1),type//' UG2SG 2-D compute domain')
2038 
2039  deallocate(a1,a2,x1,x2)
2040 
2041  !--- test the 3-D data is on computing domain
2042  allocate( a1(isc:iec, jsc:jec,nz), a2(isc:iec,jsc:jec,nz ) )
2043 
2044 ! tile = mpp_pe()/npes_per_tile + 1
2045 ! do k = 1, nz
2046 ! do j = jsc, jec
2047 ! do i = isc, iec
2048 ! a1(i,j,k) = gdata(i,j,tile)
2049 ! if(a1(i,j,k) .NE. -999) a1(i,j,k) = a1(i,j,k) + k*1.e-6
2050 ! enddo
2051 ! enddo
2052 ! enddo
2053 
2054  !First override the test SG data from file/field
2055  call data_override_init(land_domain_in=sg_domain)
2056  call data_override('LND','sst_obs',a1,time)
2057 
2058  a2 = -999
2059  !For this test on non-Land points a2 must match a1
2060  do k = 1, nz
2061  do j = jsc, jec
2062  do i = isc, iec
2063  if(.NOT. msk(i,j,1)) a2(i,j,k)=a1(i,j,k)
2064  enddo
2065  enddo
2066  enddo
2067 
2068  allocate(x1(istart:iend,nz), x2(istart:iend,nz))
2069  x1 = -999
2070  x2 = -999
2071  !--- fill the value of x2
2072  ! tile = mpp_get_UG_domain_tile_id(UG_domain)
2073  ! pos = 0
2074  ! do n = 1, tile-1
2075  ! pos = pos + npts_tile(n)
2076  ! enddo
2077  ! do l = istart, iend
2078  ! i = mod((grid_index(pos+l)-1), nx) + 1
2079  ! j = (grid_index(pos+l)-1)/nx + 1
2080  ! do k = 1, nz
2081  ! x2(l,k) = gdata(i,j,tile) + k*1.e-6
2082  ! enddo
2083  ! enddo
2084 
2085  !Now override the test UG data from the same file/field
2086  call data_override_init(land_domainug_in=ug_domain)
2087  call data_override_ug('LND','sst_obs',x2,time)
2088 
2089  !Ensure you get the same UG data from the SG data
2090  call mpp_pass_sg_to_ug(ug_domain, a1, x1)
2091  call compare_checksums_2d(x1, x2, type//' SG2UG 3-D compute domain')
2092  !Ensure you get the same SG data from the UG data if you transform back
2093  call mpp_pass_ug_to_sg(ug_domain, x1, a2)
2094  call compare_checksums(a1,a2,type//' UG2SG 3-D compute domain')
2095  deallocate(a1,a2,x1,x2)
2096 
2097 
2098  end subroutine test_unstruct_grid
2099 
2100  subroutine compare_checksums( a, b, string )
2101  real, intent(in), dimension(:,:,:) :: a, b
2102  character(len=*), intent(in) :: string
2103  integer(LONG_KIND) :: sum1, sum2
2104  integer :: i, j, k,pe
2105 
2106  ! z1l can not call mpp_sync here since there might be different number of tiles on each pe.
2107  ! mpp_sync()
2108  call mpp_sync_self()
2109  pe = mpp_pe()
2110 
2111  if(size(a,1) .ne. size(b,1) .or. size(a,2) .ne. size(b,2) .or. size(a,3) .ne. size(b,3) ) &
2112  call mpp_error(fatal,'compare_chksum: size of a and b does not match')
2113 
2114 !print*, "pe,a(1,1,1),b(1,1,1)=", pe,a(1,1,1),b(1,1,1)
2115 
2116  do k = 1, size(a,3)
2117  do j = 1, size(a,2)
2118  do i = 1, size(a,1)
2119  if(a(i,j,k) .ne. b(i,j,k)) then
2120  print*, "pe,i,j,k", pe,i,j,k
2121  print*, "a =", a(i,j,k)
2122  print*, "b =", b(i,j,k)
2123 ! write(stdunit,'(a,i3,a,i3,a,i3,a,i3,a,f20.9,a,f20.9)')" at pe ", mpp_pe(), &
2124 ! ", at point (",i,", ", j, ", ", k, "), a = ", a(i,j,k), ", b = ", b(i,j,k)
2125  call mpp_error(fatal, trim(string)//': point by point comparison are not OK.')
2126  endif
2127  enddo
2128  enddo
2129  enddo
2130 
2131  sum1 = mpp_chksum( a, (/pe/) )
2132  sum2 = mpp_chksum( b, (/pe/) )
2133 
2134  if( sum1.EQ.sum2 )then
2135  if( pe.EQ.mpp_root_pe() )call mpp_error( note, trim(string)//': OK.' )
2136  !--- in some case, even though checksum agree, the two arrays
2137  ! actually are different, like comparing (1.1,-1.2) with (-1.1,1.2)
2138  !--- hence we need to check the value point by point.
2139  else
2140  call mpp_error( fatal, trim(string)//': chksums are not OK.' )
2141  end if
2142  end subroutine compare_checksums
2143 
2144  !###########################################################################
2145  subroutine compare_checksums_2d( a, b, string )
2146  real, intent(in), dimension(:,:) :: a, b
2147  character(len=*), intent(in) :: string
2148  integer(LONG_KIND) :: sum1, sum2
2149  integer :: i, j,pe
2150 
2151  ! z1l can not call mpp_sync here since there might be different number of tiles on each pe.
2152  ! mpp_sync()
2153  call mpp_sync_self()
2154  pe = mpp_pe()
2155 
2156  if(size(a,1) .ne. size(b,1) .or. size(a,2) .ne. size(b,2) ) &
2157  call mpp_error(fatal,'compare_chksum_2D: size of a and b does not match')
2158 
2159 !print*, "a(1,1),b(1,1)=", a(1,1),b(1,1)
2160 
2161  do j = 1, size(a,2)
2162  do i = 1, size(a,1)
2163  if(a(i,j) .ne. b(i,j)) then
2164  print*, "i,j= ", i,j
2165  print*, "a =", a(i,j)
2166  print*, "b =", b(i,j)
2167 ! write(stdunit,'(a,i3,a,i3,a,i3,a,f20.9,a,f20.9)')"at pe ", mpp_pe(), &
2168 ! ", at point (",i,", ", j, "),a=", a(i,j), ",b=", b(i,j)
2169  call mpp_error(fatal, trim(string)//': point by point comparison are not OK.')
2170  endif
2171  enddo
2172  enddo
2173 
2174  sum1 = mpp_chksum( a, (/pe/) )
2175  sum2 = mpp_chksum( b, (/pe/) )
2176 
2177  if( sum1.EQ.sum2 )then
2178  if( pe.EQ.mpp_root_pe() )call mpp_error( note, trim(string)//': OK.' )
2179  !--- in some case, even though checksum agree, the two arrays
2180  ! actually are different, like comparing (1.1,-1.2) with (-1.1,1.2)
2181  !--- hence we need to check the value point by point.
2182  else
2183  call mpp_error( fatal, trim(string)//': chksums are not OK.' )
2184  end if
2185  end subroutine compare_checksums_2d
2186  subroutine define_cubic_mosaic(type, domain, ni, nj, global_indices, layout, pe_start, pe_end)
2187  character(len=*), intent(in) :: type
2188  type(domain2d), intent(inout) :: domain
2189  integer, intent(in) :: global_indices(:,:), layout(:,:)
2190  integer, intent(in) :: ni(:), nj(:)
2191  integer, intent(in) :: pe_start(:), pe_end(:)
2192  integer, dimension(12) :: istart1, iend1, jstart1, jend1, tile1
2193  integer, dimension(12) :: istart2, iend2, jstart2, jend2, tile2
2194  integer :: ntiles, num_contact, msize(2)
2195  integer :: whalo = 2, ehalo = 2, shalo = 2, nhalo = 2
2196 
2197 
2198  ntiles = 6
2199  num_contact = 12
2200  if(size(pe_start(:)) .NE. 6 .OR. size(pe_end(:)) .NE. 6 ) call mpp_error(fatal, &
2201  "define_cubic_mosaic: size of pe_start and pe_end should be 6")
2202  if(size(global_indices,1) .NE. 4) call mpp_error(fatal, &
2203  "define_cubic_mosaic: size of first dimension of global_indices should be 4")
2204  if(size(global_indices,2) .NE. 6) call mpp_error(fatal, &
2205  "define_cubic_mosaic: size of second dimension of global_indices should be 6")
2206  if(size(layout,1) .NE. 2) call mpp_error(fatal, &
2207  "define_cubic_mosaic: size of first dimension of layout should be 2")
2208  if(size(layout,2) .NE. 6) call mpp_error(fatal, &
2209  "define_cubic_mosaic: size of second dimension of layout should be 6")
2210  if(size(ni(:)) .NE. 6 .OR. size(nj(:)) .NE. 6) call mpp_error(fatal, &
2211  "define_cubic_mosaic: size of ni and nj should be 6")
2212 
2213  !--- Contact line 1, between tile 1 (EAST) and tile 2 (WEST)
2214  tile1(1) = 1; tile2(1) = 2
2215  istart1(1) = ni(1); iend1(1) = ni(1); jstart1(1) = 1; jend1(1) = nj(1)
2216  istart2(1) = 1; iend2(1) = 1; jstart2(1) = 1; jend2(1) = nj(2)
2217  !--- Contact line 2, between tile 1 (NORTH) and tile 3 (WEST)
2218  tile1(2) = 1; tile2(2) = 3
2219  istart1(2) = 1; iend1(2) = ni(1); jstart1(2) = nj(1); jend1(2) = nj(1)
2220  istart2(2) = 1; iend2(2) = 1; jstart2(2) = nj(3); jend2(2) = 1
2221  !--- Contact line 3, between tile 1 (WEST) and tile 5 (NORTH)
2222  tile1(3) = 1; tile2(3) = 5
2223  istart1(3) = 1; iend1(3) = 1; jstart1(3) = 1; jend1(3) = nj(1)
2224  istart2(3) = ni(5); iend2(3) = 1; jstart2(3) = nj(5); jend2(3) = nj(5)
2225  !--- Contact line 4, between tile 1 (SOUTH) and tile 6 (NORTH)
2226  tile1(4) = 1; tile2(4) = 6
2227  istart1(4) = 1; iend1(4) = ni(1); jstart1(4) = 1; jend1(4) = 1
2228  istart2(4) = 1; iend2(4) = ni(6); jstart2(4) = nj(6); jend2(4) = nj(6)
2229  !--- Contact line 5, between tile 2 (NORTH) and tile 3 (SOUTH)
2230  tile1(5) = 2; tile2(5) = 3
2231  istart1(5) = 1; iend1(5) = ni(2); jstart1(5) = nj(2); jend1(5) = nj(2)
2232  istart2(5) = 1; iend2(5) = ni(3); jstart2(5) = 1; jend2(5) = 1
2233  !--- Contact line 6, between tile 2 (EAST) and tile 4 (SOUTH)
2234  tile1(6) = 2; tile2(6) = 4
2235  istart1(6) = ni(2); iend1(6) = ni(2); jstart1(6) = 1; jend1(6) = nj(2)
2236  istart2(6) = ni(4); iend2(6) = 1; jstart2(6) = 1; jend2(6) = 1
2237  !--- Contact line 7, between tile 2 (SOUTH) and tile 6 (EAST)
2238  tile1(7) = 2; tile2(7) = 6
2239  istart1(7) = 1; iend1(7) = ni(2); jstart1(7) = 1; jend1(7) = 1
2240  istart2(7) = ni(6); iend2(7) = ni(6); jstart2(7) = nj(6); jend2(7) = 1
2241  !--- Contact line 8, between tile 3 (EAST) and tile 4 (WEST)
2242  tile1(8) = 3; tile2(8) = 4
2243  istart1(8) = ni(3); iend1(8) = ni(3); jstart1(8) = 1; jend1(8) = nj(3)
2244  istart2(8) = 1; iend2(8) = 1; jstart2(8) = 1; jend2(8) = nj(4)
2245  !--- Contact line 9, between tile 3 (NORTH) and tile 5 (WEST)
2246  tile1(9) = 3; tile2(9) = 5
2247  istart1(9) = 1; iend1(9) = ni(3); jstart1(9) = nj(3); jend1(9) = nj(3)
2248  istart2(9) = 1; iend2(9) = 1; jstart2(9) = nj(5); jend2(9) = 1
2249  !--- Contact line 10, between tile 4 (NORTH) and tile 5 (SOUTH)
2250  tile1(10) = 4; tile2(10) = 5
2251  istart1(10) = 1; iend1(10) = ni(4); jstart1(10) = nj(4); jend1(10) = nj(4)
2252  istart2(10) = 1; iend2(10) = ni(5); jstart2(10) = 1; jend2(10) = 1
2253  !--- Contact line 11, between tile 4 (EAST) and tile 6 (SOUTH)
2254  tile1(11) = 4; tile2(11) = 6
2255  istart1(11) = ni(4); iend1(11) = ni(4); jstart1(11) = 1; jend1(11) = nj(4)
2256  istart2(11) = ni(6); iend2(11) = 1; jstart2(11) = 1; jend2(11) = 1
2257  !--- Contact line 12, between tile 5 (EAST) and tile 6 (WEST)
2258  tile1(12) = 5; tile2(12) = 6
2259  istart1(12) = ni(5); iend1(12) = ni(5); jstart1(12) = 1; jend1(12) = nj(5)
2260  istart2(12) = 1; iend2(12) = 1; jstart2(12) = 1; jend2(12) = nj(6)
2261  msize(1) = maxval(ni(:)/layout(1,:)) + whalo + ehalo + 1 ! make sure memory domain size is no smaller than
2262  msize(2) = maxval(nj(:)/layout(2,:)) + shalo + nhalo + 1 ! data domain size
2263  call mpp_define_mosaic(global_indices, layout, domain, ntiles, num_contact, tile1, tile2, &
2264  istart1, iend1, jstart1, jend1, istart2, iend2, jstart2, jend2, &
2265  pe_start, pe_end, symmetry = .true., whalo=whalo, ehalo=ehalo, &
2266  shalo=shalo, nhalo=nhalo, name = trim(type), memory_size = msize )
2267 
2268  return
2269 
2270  end subroutine define_cubic_mosaic
2271 
2272 !=================================================================================================================================
2273  end program test
2274 #endif
Definition: fms.F90:20
integer, parameter undef
type(domain2d), save ocn_domain
integer, parameter annual
real, dimension(:,:), allocatable, target lon_local_lnd
integer, parameter max_array
integer, parameter, public noleap
subroutine data_override_3d(gridname, fieldname_code, data, time, override, data_index, is_in, ie_in, js_in, je_in)
subroutine define_cubic_mosaic(type, domain, ni, nj, global_indices, layout, pe_start, pe_end)
subroutine, public mpp_memuse_begin
subroutine, public data_override_init(Atm_domain_in, Ocean_domain_in, Ice_domain_in, Land_domain_in, Land_domainUG_in)
real, dimension(:,:), allocatable, target lat_local_ice
subroutine, public time_interp_external_init()
integer, parameter daily
subroutine, public horiz_interp_del(Interp)
real, dimension(:,:), allocatable, target lon_local_ice
type(override_type), dimension(max_array), save override_array
integer, parameter, public inside_region
type(domain2d), save lnd_domain
subroutine compare_checksums(a, b, string)
integer function, public register_static_field(module_name, field_name, axes, long_name, units, missing_value, range, mask_variant, standard_name, DYNAMIC, do_not_log, interp_method, tile_count, area, volume, realm)
integer function, public nearest_index(value, array)
Definition: axis_utils.F90:443
integer, parameter max_table
subroutine data_override_ug_1d(gridname, fieldname, data, time, override)
subroutine check_grid_sizes(domain_name, Domain, nlon, nlat)
subroutine, public diag_manager_end(time)
int get_cpu_affinity(void)
Definition: affinity.c:44
type(data_type) default_table
Definition: mpp.F90:39
integer function, public check_nml_error(IOSTAT, NML_NAME)
Definition: fms.F90:658
real, dimension(:,:), allocatable, target lat_local_ocn
subroutine data_override_0d(gridname, fieldname_code, data, time, override, data_index)
type(override_type), save default_array
type(domainug), save lnd_domainug
subroutine get_grid_version_2(mosaic_file, mod_name, domain, isc, iec, jsc, jec, lon, lat, min_lon, max_lon)
character(len=input_str_length), dimension(:), allocatable, target, public input_nml_file
Definition: mpp.F90:1378
subroutine get_grid_version_1(grid_file, mod_name, domain, isc, iec, jsc, jec, lon, lat, min_lon, max_lon)
subroutine, public set_calendar_type(type, err_msg)
integer, parameter hourly
real, parameter, public pi
Ratio of circle circumference to diameter [N/A].
Definition: constants.F90:74
int field_exist(const char *file, const char *name)
Definition: read_mosaic.c:69
subroutine, public fms_init(localcomm)
Definition: fms.F90:353
integer, parameter monthly
subroutine test_unstruct_grid(type, Time)
type(domain2d), save, public null_domain2d
real, dimension(:,:), allocatable, target lon_local_atm
integer function, public init_external_field(file, fieldname, format, threading, domain, desired_units, verbose, axis_centers, axis_sizes, override, correct_leap_year_inconsistency, permit_calendar_conversion, use_comp_domain, ierr, nwindows, ignore_axis_atts)
subroutine, public diag_manager_init(diag_model_subset, time_init, err_msg)
subroutine, public fms_io_init()
Definition: fms_io.F90:638
subroutine, public field_size(filename, fieldname, siz, field_found, domain, no_domain)
Definition: fms_io.F90:4941
integer, parameter, public julian
subroutine data_override_2d(gridname, fieldname, data_2D, time, override, is_in, ie_in, js_in, je_in)
subroutine get_grid
integer function, dimension(4), public get_external_field_size(index)
subroutine, public mpp_memuse_end(text, unit)
real, dimension(:,:), allocatable, target lat_local_lnd
subroutine, public horiz_interp_init
integer, parameter, public no_region
real, dimension(:,:), allocatable, target lat_local_atm
type(domainug), save, public null_domainug
subroutine, public reset_src_data_region(index, is, ie, js, je)
integer, parameter r8_kind
Definition: platform.F90:24
subroutine, public fms_end()
Definition: fms.F90:476
subroutine, public fms_io_exit()
Definition: fms_io.F90:750
subroutine, public set_override_region(index, region_type, is_region, ie_region, js_region, je_region)
subroutine, public get_axis_bounds(axis, axis_bound, axes, bnd_name, err_msg)
Definition: axis_utils.F90:149
#define max(a, b)
Definition: mosaic_util.h:33
subroutine, public data_override_unset_domains(unset_Atm, unset_Ocean, unset_Ice, unset_Land, must_be_set)
subroutine, public get_mosaic_tile_grid(grid_file, mosaic_file, domain, tile_count)
Definition: fms_io.F90:7850
subroutine compare_checksums_2d(a, b, string)
real, parameter tpi
subroutine get_domainug(gridname, UGdomain, comp_domain)
type(data_type), dimension(max_table) data_table
type(domain2d), save atm_domain
real, dimension(:,:), allocatable, target lon_local_ocn
#define min(a, b)
Definition: mosaic_util.h:32
subroutine get_domain(gridname, domain, comp_domain)
subroutine, public error_mesg(routine, message, level)
Definition: fms.F90:529
subroutine, public print_time(Time, str, unit)
type(domain2d), save ice_domain
subroutine, public constants_init
dummy routine.
Definition: constants.F90:141
subroutine data_override_ug_2d(gridname, fieldname, data, time, override)
logical module_is_initialized
integer, parameter, public outside_region