FMS  2024.03
Flexible Modeling System
mpp_define_nest_domains.inc
1 ! -*-f90-*-
2 
3 
4 !***********************************************************************
5 !* GNU Lesser General Public License
6 !*
7 !* This file is part of the GFDL Flexible Modeling System (FMS).
8 !*
9 !* FMS is free software: you can redistribute it and/or modify it under
10 !* the terms of the GNU Lesser General Public License as published by
11 !* the Free Software Foundation, either version 3 of the License, or (at
12 !* your option) any later version.
13 !*
14 !* FMS is distributed in the hope that it will be useful, but WITHOUT
15 !* ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
16 !* FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
17 !* for more details.
18 !*
19 !* You should have received a copy of the GNU Lesser General Public
20 !* License along with FMS. If not, see <http://www.gnu.org/licenses/>.
21 !***********************************************************************
22 
23 !> @file
24 !> @brief Routines for use in @ref mpp_domains_mod for routines utilizing
25 !! domains with nested grids
26 
27 !> @addtogroup mpp_domains_mod
28 !> @{
29 
30 !#############################################################################
31 !> @brief Set up a domain to pass data between aligned coarse and fine grid of nested model.
32 !!
33 !> Set up a domain to pass data between aligned coarse and fine grid of a nested
34 !! model. Supports multiple and telescoping nests. A telescoping nest is defined as
35 !! a nest within a nest. Nest domains may span multiple tiles, but cannot contain a
36 !! coarse-grid, cube corner. Concurrent nesting is the only supported mechanism,
37 !! i.e. coarse and fine grid are on individual, non-overlapping, processor lists.
38 !! Coarse and fine grid domain need to be defined before calling mpp_define_nest_domains.
39 !! An mpp_broadcast is needed to broadcast both fine and coarse grid domain onto all processors.\n
40 !!\n
41 !! mpp_update_nest_coarse is used to pass data from fine grid to coarse grid computing domain.
42 !! mpp_update_nest_fine is used to pass data from coarse grid to fine grid halo.
43 !! You may call mpp_get_C2F_index before calling mpp_update_nest_fine to get the index for
44 !! passing data from coarse to fine. You may call mpp_get_F2C_index before calling
45 !! mpp_update_nest_coarse to get the index for passing data from coarse to fine.
46 !!
47 !> @note The following tests for nesting of regular lat-lon grids upon a cubed-sphere
48 !! grid are done in test_mpp_domains:\n
49 !! a) a first-level nest spanning multiple cubed-sphere faces (tiles 1, 2, & 4)\n
50 !! b) a first-level nest wholly contained within tile 3\n
51 !! c) a second-level nest contained within the nest mentioned in a)\n
52 !! Tests are done for data at T, E, C, N-cell center.\n
53 !!
54 !!
55 !! Below is an example to pass data between fine and coarse grid (More details on how to
56 !! use the nesting domain update are available in routine test_update_nest_domain of
57 !! test_fms/mpp/test_mpp_domains.F90.\n
58 !!\n
59 !! @code{.F90}
60 !! if( concurrent ) then
61 !! call mpp_broadcast_domain(domain_fine)
62 !! call mpp_broadcast_domain(domain_coarse)
63 !! endif
64 !!
65 !! call mpp_define_nest_domains(nest_domain,domain,num_nest,nest_level(1:num_nest), &
66 !! tile_fine(1:num_nest), tile_coarse(1:num_nest), &
67 !! istart_coarse(1:num_nest), icount_coarse(1:num_nest), &
68 !! jstart_coarse(1:num_nest), jcount_coarse(1:num_nest), &
69 !! npes_nest_tile, x_refine(1:num_nest), y_refine(1:num_nest), &
70 !! extra_halo=extra_halo, name="nest_domain")
71 !!
72 !! call mpp_get_C2F_index(nest_domain, isw_f, iew_f, jsw_f, jew_f, isw_c, iew_c, jsw_c, jew_c, WEST, level)
73 !! call mpp_get_C2F_index(nest_domain, ise_f, iee_f, jse_f, jee_f, ise_c, iee_c, jse_c, jee_c, EAST, level)
74 !! call mpp_get_C2F_index(nest_domain, iss_f, ies_f, jss_f, jes_f, iss_c, ies_c, jss_c, jes_c, SOUTH, level)
75 !! call mpp_get_C2F_index(nest_domain, isn_f, ien_f, jsn_f, jen_f, isn_c, ien_c, jsn_c, jen_c, NORTH, level)
76 !!
77 !! allocate(wbuffer(isw_c:iew_c, jsw_c:jew_c,nz))
78 !! allocate(ebuffer(ise_c:iee_c, jse_c:jee_c,nz))
79 !! allocate(sbuffer(iss_c:ies_c, jss_c:jes_c,nz))
80 !! allocate(nbuffer(isn_c:ien_c, jsn_c:jen_c,nz))
81 !! call mpp_update_nest_fine(x, nest_domain, wbuffer, sbuffer, ebuffer, nbuffer)
82 !!
83 !! call mpp_get_F2C_index(nest_domain, is_c, ie_c, js_c, je_c, is_f, ie_f, js_f, je_f, nest_level=level)
84 !! allocate(buffer (is_f:ie_f, js_f:je_f,nz))
85 !! call mpp_update_nest_coarse(x, nest_domain, buffer)
86 !
87 !! call mpp_define_nest_domains (nest_domain, domain, num_nest, nest_level, tile_fine, tile_coarse,
88 !! istart_coarse, icount_coarse, jstart_coarse, jcount_coarse,
89 !! npes_nest_tile, x_refine, y_refine, extra_halo, name)
90 !!
91 !! @endcode
92 !!
93 !> @note Currently the contact will be limited to overlap contact.
94 !!
95 subroutine mpp_define_nest_domains(nest_domain, domain, num_nest, nest_level, tile_fine, tile_coarse, &
96  istart_coarse, icount_coarse, jstart_coarse, jcount_coarse, npes_nest_tile, &
97  x_refine, y_refine, extra_halo, name)
98  type(nest_domain_type), intent(inout) :: nest_domain !< holds the information to pass data
99  !! between nest and parent grids.
100  type(domain2d), target, intent(in ) :: domain !< domain for the grid defined in the current pelist
101  integer, intent(in ) :: num_nest !< number of nests
102  integer, intent(in ) :: nest_level(:) !< array containing the nest level for each nest
103  !!(>1 implies a telescoping nest)
104  integer, intent(in ) :: tile_fine(:), tile_coarse(:) !< array containing tile number of the
105  !! nest grid (monotonically increasing starting with 7),
106  !! array containing tile number of the parent grid corresponding
107  !! to the lower left corner of a given nest
108  integer, intent(in ) :: istart_coarse(:), icount_coarse(:), jstart_coarse(:), jcount_coarse(:) !<start
109  !! array containing index in the parent grid of the lower left corner of a given
110  !! nest, count: array containing span of the nest on the parent grid
111  integer, intent(in ) :: npes_nest_tile(:) !< array containing number of pes to allocated
112  !! to each defined tile
113  integer, intent(in ) :: x_refine(:), y_refine(:) !< array containing refinement ratio
114  !! for each nest
115  integer, optional, intent(in ) :: extra_halo !< extra halo for passing data from coarse grid to fine grid.
116  !! default is 0 and currently only support extra_halo = 0.
117  character(len=*), optional, intent(in ) :: name !< name of the nest domain
118 
119  integer :: n, l, m, my_tile_coarse
120  integer :: npes_level, prev_tile_coarse
121  integer :: extra_halo_local, npes_nest_top
122  integer, dimension(:), allocatable :: pes, pe_start_pos, pe_end_pos, pelist_level
123  logical, dimension(:), allocatable :: is_nest_fine, is_nest_coarse
124  integer, dimension(num_nest) :: istart_fine, iend_fine, jstart_fine, jend_fine
125  integer, dimension(num_nest) :: iend_coarse, jend_coarse
126  integer :: nnest, nlevels, ntiles_top, ntiles, pos
127  logical :: is_first
128 
129  if(PRESENT(name)) then
130  if(len_trim(name) > name_length) then
131  call mpp_error(fatal, "mpp_domains_define.inc(mpp_define_nest_domain): "// &
132  "the len_trim of optional argument name ="//trim(name)// &
133  " is greater than NAME_LENGTH, change the argument name or increase NAME_LENGTH")
134  endif
135  nest_domain%name = name
136  endif
137 
138  extra_halo_local = 0
139  if(present(extra_halo)) then
140  if(extra_halo .NE. 0) call mpp_error(fatal, &
141  & "mpp_define_nest_domains.inc: only support extra_halo=0, contact developer")
142  extra_halo_local = extra_halo
143  endif
144 
145  !---make sure dimension size is correct
146  if(size(tile_fine(:)) .NE. num_nest) call mpp_error(fatal, &
147  & .NE."mpp_define_nest_domains.inc: size(tile_fine) num_nest")
148  if(size(tile_coarse(:)) .NE. num_nest) call mpp_error(fatal, &
149  & .NE."mpp_define_nest_domains.inc: size(tile_coarse) num_nest")
150  if(size(istart_coarse(:)) .NE. num_nest) call mpp_error(fatal, &
151  & .NE."mpp_define_nest_domains.inc: size(istart_coarse) num_nest")
152  if(size(icount_coarse(:)) .NE. num_nest) call mpp_error(fatal, &
153  & .NE."mpp_define_nest_domains.inc: size(icount_coarse) num_nest")
154  if(size(jstart_coarse(:)) .NE. num_nest) call mpp_error(fatal, &
155  & .NE."mpp_define_nest_domains.inc: size(jstart_coarse) num_nest")
156  if(size(jcount_coarse(:)) .NE. num_nest) call mpp_error(fatal, &
157  & .NE."mpp_define_nest_domains.inc: size(jcount_coarse) num_nest")
158 
159  do n = 1, num_nest
160  if(istart_coarse(n) < 1) call mpp_error(fatal, "mpp_define_nest_domains.inc: istart_coarse < 1")
161  if(icount_coarse(n) < 1) call mpp_error(fatal, "mpp_define_nest_domains.inc: iend_coarse < 1")
162  if(jstart_coarse(n) < 1) call mpp_error(fatal, "mpp_define_nest_domains.inc: jstart_coarse < 1")
163  if(jcount_coarse(n) < 1) call mpp_error(fatal, "mpp_define_nest_domains.inc: jend_coarse < 1")
164  iend_coarse(n) = istart_coarse(n) + icount_coarse(n) - 1
165  jend_coarse(n) = jstart_coarse(n) + jcount_coarse(n) - 1
166  istart_fine(n) = 1 ; iend_fine(n) = icount_coarse(n)*x_refine(n)
167  jstart_fine(n) = 1 ; jend_fine(n) = jcount_coarse(n)*y_refine(n)
168  end do
169 
170  !--- make sure the nest level is monotonic and no jumping
171  if(nest_level(1) .NE. 1) call mpp_error(fatal, .NE."mpp_define_nest_domains.inc: nest_level(1) 1")
172  do n = 2, num_nest
173  if(nest_level(n) < nest_level(n-1)) call mpp_error(fatal, &
174  & "mpp_define_nest_domains.inc: nest_level is not monotone increasing")
175  if(nest_level(n) > nest_level(n-1)+1) call mpp_error(fatal, &
176  & "mpp_define_nest_domains.inc: nest_level(n) > nest_level(n-1)+1")
177  enddo
178  nlevels = nest_level(num_nest)
179 
180  !---make sure tile_fine and tile_nest are monotone increasing.
181  do n = 2, num_nest
182  if(tile_fine(n) < tile_fine(n-1)) call mpp_error(fatal, &
183  & "mpp_define_nest_domains.inc: tile_fine is not monotone increasing")
184  if(tile_coarse(n) < tile_coarse(n-1)) call mpp_error(fatal, "mpp_define_nest_domains.inc: "// &
185  "tile_coarse is not monotone increasing")
186  enddo
187 
188  allocate( pes(mpp_npes()) )
189  call mpp_get_current_pelist(pes)
190 
191  nest_domain%num_nest = num_nest
192  allocate(nest_domain%tile_fine(num_nest), nest_domain%tile_coarse(num_nest) )
193  allocate(nest_domain%istart_fine(num_nest), nest_domain%iend_fine(num_nest) )
194  allocate(nest_domain%jstart_fine(num_nest), nest_domain%jend_fine(num_nest) )
195  allocate(nest_domain%istart_coarse(num_nest), nest_domain%iend_coarse(num_nest) )
196  allocate(nest_domain%jstart_coarse(num_nest), nest_domain%jend_coarse(num_nest) )
197 
198  !---Added to enable moving nests
199  if (associated(nest_domain%nest_level)) deallocate(nest_domain%nest_level) !< Check if allocated
200  allocate(nest_domain%nest_level(num_nest))
201 
202  nest_domain%tile_fine = tile_fine(1:num_nest)
203  nest_domain%tile_coarse = tile_coarse(1:num_nest)
204  nest_domain%istart_fine = istart_fine(1:num_nest)
205  nest_domain%iend_fine = iend_fine(1:num_nest)
206  nest_domain%jstart_fine = jstart_fine(1:num_nest)
207  nest_domain%jend_fine = jend_fine(1:num_nest)
208  nest_domain%istart_coarse = istart_coarse(1:num_nest)
209  nest_domain%iend_coarse = iend_coarse(1:num_nest)
210  nest_domain%jstart_coarse = jstart_coarse(1:num_nest)
211  nest_domain%jend_coarse = jend_coarse(1:num_nest)
212 
213  !---Added to enable moving nests; need this information when shifting the nest position.
214  nest_domain%nest_level = nest_level(1:num_nest)
215 
216  !--- make sure the tile_id of top level of grid is continuous and starting from 1
217  if(mpp_pe()==mpp_root_pe()) then
218  ntiles_top = domain%ntiles
219  do n = 1, ntiles_top
220  if(domain%tile_id_all(n) .NE. n) call mpp_error(fatal, &
221  "mpp_define_nest_domains.inc: top level grid tile_id should be 1, 2, ..,ntiles")
222  enddo
223  endif
224  call mpp_broadcast(ntiles_top, mpp_root_pe())
225  !--- make sure the nest grid tile_ids are continuous
226  do n = 1, num_nest
227  if(tile_fine(n) .NE. ntiles_top+n) then
228  print*, "tile_fine, ntile_top, n=", tile_fine(n), ntiles_top, n, mpp_pe()
229  call mpp_error(fatal, "mpp_define_nest_domains.inc: tile_id is not continuous")
230  endif
231  enddo
232 
233  allocate(pe_start_pos(ntiles_top+num_nest))
234  allocate(pe_end_pos(ntiles_top+num_nest))
235 
236  do n = 2, ntiles_top
237  if(npes_nest_tile(n) .NE. npes_nest_tile(n-1)) call mpp_error(fatal, &
238  "mpp_define_nest_domains.inc: all the tiles in top grid should use same number of MPI ranks")
239  enddo
240 
241  npes_nest_top = npes_nest_tile(1)*ntiles_top
242 
243  !--- get the pe start and end pos for each tile
244  do n = 1, ntiles_top
245  pe_start_pos(n) = 1
246  pe_end_pos(n) = npes_nest_tile(1)*ntiles_top
247  enddo
248  ntiles = tile_fine(num_nest)
249  if(ntiles .NE. ntiles_top + num_nest) call mpp_error(fatal, "mpp_define_nest_domains.inc: "// &
250  .NE."ntiles ntiles_top + num_nest")
251  do n = 1, num_nest
252  pe_start_pos(ntiles_top+n) = pe_end_pos(ntiles_top+n-1) + 1
253  pe_end_pos(ntiles_top+n) = pe_end_pos(ntiles_top+n-1) + npes_nest_tile(tile_fine(n))
254  enddo
255 
256  nest_domain%num_level = nlevels
257  if (associated(nest_domain%nest)) deallocate(nest_domain%nest) !< Check if allocated
258  allocate(nest_domain%nest(nlevels))
259  allocate(pelist_level(mpp_npes()))
260  allocate(is_nest_fine(nlevels))
261  allocate(is_nest_coarse(nlevels))
262 
263  !--- setup pelist for each level
264  pos = 0
265  is_nest_fine(:) = .false.
266  is_nest_coarse(:) = .false.
267  do l = 1, nlevels
268  npes_level = 0
269  pos = 0
270  is_first = .true.
271  prev_tile_coarse = 0
272  !--- first get coarse processor
273  do n = 1, num_nest
274  if(nest_level(n) == l) then
275  if(mpp_pe() .GE. pes(pe_start_pos(tile_fine(n))) .AND. mpp_pe() .LE. pes(pe_end_pos(tile_fine(n)))) then
276  is_nest_fine(l) = .true.
277  endif
278  if(mpp_pe() .GE. pes(pe_start_pos(tile_coarse(n))) .AND. mpp_pe() .LE. pes(pe_end_pos(tile_coarse(n)))) then
279  is_nest_coarse(l) = .true.
280  endif
281  if(pos==0 .OR. (l .NE. 1 .AND. prev_tile_coarse .NE. tile_coarse(n)) ) then
282  do m = pe_start_pos(tile_coarse(n)), pe_end_pos(tile_coarse(n))
283  pos = pos+1
284  pelist_level(pos) = pes(m)
285  enddo
286  npes_level = npes_level + pe_end_pos(tile_coarse(n)) - pe_start_pos(tile_coarse(n)) + 1
287  endif
288  prev_tile_coarse = tile_coarse(n)
289  endif
290  enddo
291  ! fine processor
292  do n = 1, num_nest
293  if(nest_level(n) == l) then
294  do m = pe_start_pos(tile_fine(n)), pe_end_pos(tile_fine(n))
295  pos = pos+1
296  pelist_level(pos) = pes(m)
297  enddo
298  npes_level = npes_level + pe_end_pos(tile_fine(n)) - pe_start_pos(tile_fine(n)) + 1
299  endif
300  enddo
301 
302  if (associated(nest_domain%nest(l)%pelist)) deallocate(nest_domain%nest(l)%pelist) !< Check if allocated
303  allocate(nest_domain%nest(l)%pelist(npes_level))
304  nest_domain%nest(l)%pelist(:) = pelist_level(1:npes_level)
305 
306  call mpp_declare_pelist(nest_domain%nest(l)%pelist)
307  nest_domain%nest(l)%on_level = any(nest_domain%nest(l)%pelist(:)==mpp_pe())
308  nest_domain%nest(l)%is_fine_pe = is_nest_fine(l)
309  nest_domain%nest(l)%is_coarse_pe = is_nest_coarse(l)
310  if(nest_domain%nest(l)%on_level .neqv. (is_nest_fine(l) .OR. is_nest_coarse(l))) then
311  print*, "on_level=", nest_domain%nest(l)%on_level, is_nest_fine(l), is_nest_coarse(l), mpp_pe(),l
312  call mpp_error(fatal, "mpp_define_nest_domains.inc:on_level does not match is_nest_fine/is_nest_coarse")
313  endif
314  if(is_nest_fine(l) .and. is_nest_coarse(l)) then
315  call mpp_error(fatal, "mpp_define_nest_domains.inc: is_nest_fine and is_nest_coarse can not both be true")
316  endif
317  enddo
318 
319  if(count(is_nest_fine)>1) call mpp_error(fatal, "mpp_define_nest_domains.inc: count(is_nest_fine)>1")
320  if(count(is_nest_coarse)>1) call mpp_error(fatal, "mpp_define_nest_domains.inc: count(is_nest_coarse)>1")
321 
322  do l = 1, nlevels
323  !--- setup for each level
324  if(nest_domain%nest(l)%on_level) then
325  call mpp_set_current_pelist(nest_domain%nest(l)%pelist)
326  nnest = count(nest_level==l)
327  nest_domain%nest(l)%num_nest = nnest
328  allocate(nest_domain%nest(l)%tile_fine(nnest), nest_domain%nest(l)%tile_coarse(nnest) )
329  allocate(nest_domain%nest(l)%istart_fine(nnest), nest_domain%nest(l)%iend_fine(nnest) )
330  allocate(nest_domain%nest(l)%jstart_fine(nnest), nest_domain%nest(l)%jend_fine(nnest) )
331  allocate(nest_domain%nest(l)%istart_coarse(nnest), nest_domain%nest(l)%iend_coarse(nnest) )
332  allocate(nest_domain%nest(l)%jstart_coarse(nnest), nest_domain%nest(l)%jend_coarse(nnest) )
333  my_tile_coarse = 0
334 
335  pos=0
336  do n = 1, num_nest
337  if(nest_level(n) ==l) then
338  pos = pos+1
339  nest_domain%nest(l)%tile_fine(pos) = tile_fine(n)
340  nest_domain%nest(l)%tile_coarse(pos) = tile_coarse(n)
341  nest_domain%nest(l)%istart_fine(pos) = istart_fine(n)
342  nest_domain%nest(l)%iend_fine(pos) = iend_fine(n)
343  nest_domain%nest(l)%jstart_fine(pos) = jstart_fine(n)
344  nest_domain%nest(l)%jend_fine(pos) = jend_fine(n)
345  nest_domain%nest(l)%istart_coarse(pos) = istart_coarse(n)
346  nest_domain%nest(l)%iend_coarse(pos) = iend_coarse(n)
347  nest_domain%nest(l)%jstart_coarse(pos) = jstart_coarse(n)
348  nest_domain%nest(l)%jend_coarse(pos) = jend_coarse(n)
349  if(l==1) then
350  my_tile_coarse = 1
351  else if( (mpp_pe() .GE. pes(pe_start_pos(tile_fine(n))) .AND. &
352  & mpp_pe() .LE. pes(pe_end_pos(tile_fine(n)))) .OR. &
353  & (mpp_pe() .GE. pes(pe_start_pos(tile_coarse(n))) .AND. &
354  & mpp_pe() .LE. pes(pe_end_pos(tile_coarse(n)))) ) then
355  my_tile_coarse = tile_coarse(n)
356  endif
357  endif
358  enddo
359  if(my_tile_coarse == 0) call mpp_error(fatal, "mpp_define_nest_domains.inc: my_tile_coarse == 0")
360 
361  if(pos .NE. nest_domain%nest(l)%num_nest) &
362  call mpp_error(fatal, .NE."mpp_define_nest_domains.inc:pos nest_domain%nest(l)%num_nest")
363 
364  if(is_nest_fine(l)) then
365  nest_domain%nest(l)%domain_fine=>domain
366  allocate(nest_domain%nest(l)%domain_coarse)
367  else if(is_nest_coarse(l)) then
368  nest_domain%nest(l)%domain_coarse=>domain
369  allocate(nest_domain%nest(l)%domain_fine)
370  endif
371 !!!! DEBUG CODE ! has problems on coarse domain
372 !!$ print*, 'MPP_BROADCAST_DOMAIN: ', mpp_pe(), l, & !ASSOCIATED(nest_domain%nest(l)%domain_fine), &
373 !!$ nest_domain%nest(l)%domain_fine%tile_id(1), nest_domain%nest(l)%tile_fine, tile_fine
374 !!!! END DEBUG CODE
375  call mpp_broadcast_domain(nest_domain%nest(l)%domain_fine, nest_domain%nest(l)%tile_fine)
376  call mpp_broadcast_domain(nest_domain%nest(l)%domain_coarse, my_tile_coarse)
377  call define_nest_level_type(nest_domain%nest(l), x_refine(l), y_refine(l), extra_halo_local)
378  endif
379  enddo
380 
381 end subroutine mpp_define_nest_domains
382 
383 !> Based on mpp_define_nest_domains, but just resets positioning of nest
384 !> Modifies the parent/coarse start and end indices of the nest location
385 !! Computes new overlaps of nest PEs on parent PEs
386 !! Ramstrom/HRD Moving Nest
387 subroutine mpp_shift_nest_domains(nest_domain, domain, delta_i_coarse, delta_j_coarse, extra_halo)
388  type(nest_domain_type), intent(inout) :: nest_domain !< holds the information to pass data
389  !! between nest and parent grids.
390  type(domain2d), target, intent(in ) :: domain !< domain for the grid defined in the current pelist
391  integer, intent(in ) :: delta_i_coarse(:) !< Array of deltas of coarse grid in y direction
392  integer, intent(in ) :: delta_j_coarse(:) !< Array of deltas of coarse grid in y direction
393  integer, optional, intent(in ) :: extra_halo !< Extra halo size
394 
395  integer :: n, l, my_tile_coarse
396  integer :: num_nest
397  integer :: extra_halo_local
398  integer :: nlevels, pos
399  integer, pointer :: nest_level(:)
400 
401  nest_level => nest_domain%nest_level
402 
403  extra_halo_local = 0
404  if(present(extra_halo)) then
405  if(extra_halo .NE. 0) call mpp_error(fatal, &
406  & "shift mpp_define_nest_domains.inc: only support extra_halo=0, contact developer")
407  extra_halo_local = extra_halo
408  endif
409 
410  num_nest = nest_domain%num_nest
411  nlevels = nest_level(num_nest)
412 
413  !---make sure dimension size is correct
414  if(size(delta_i_coarse(:)) .NE. num_nest) call mpp_error(fatal, &
415  & .NE."shift mpp_define_nest_domains.inc: size(delta_i_coarse) num_nest")
416  if(size(delta_j_coarse(:)) .NE. num_nest) call mpp_error(fatal, &
417  & .NE."shift mpp_define_nest_domains.inc: size(delta_j_coarse) num_nest")
418 
419  ! Step through all nests and apply any modifications
420  ! Should only need to modify {i,j}{start,end}_coarse
421  ! The indices for fine do not change when the nest moves; only if it was changing size (which it cannot)
422  ! The PEs for running each component remain the same
423  ! But the nest may connect to different parent PEs after motion, so that must be recalculated
424  do n = 1, num_nest
425  ! Moving nest code performs validation to ensure nest does not run off edge
426 
427  nest_domain%istart_coarse(n) = nest_domain%istart_coarse(n) + delta_i_coarse(n)
428  nest_domain%iend_coarse(n) = nest_domain%iend_coarse(n) + delta_i_coarse(n)
429 
430  nest_domain%jstart_coarse(n) = nest_domain%jstart_coarse(n) + delta_j_coarse(n)
431  nest_domain%jend_coarse(n) = nest_domain%jend_coarse(n) + delta_j_coarse(n)
432 
433  end do
434 
435  ! Apply the nest motion in the nest_level_type structures
436  do l = 1, nlevels
437  !--- setup for each level
438  if(nest_domain%nest(l)%on_level) then
439  !nnest = count(nest_level==l)
440  my_tile_coarse = 0
441 
442  pos=0
443  do n = 1, num_nest
444  if(nest_level(n) ==l) then
445  pos = pos+1
446  nest_domain%nest(l)%istart_coarse(pos) = nest_domain%istart_coarse(n)
447  nest_domain%nest(l)%iend_coarse(pos) = nest_domain%iend_coarse(n)
448  nest_domain%nest(l)%jstart_coarse(pos) = nest_domain%jstart_coarse(n)
449  nest_domain%nest(l)%jend_coarse(pos) = nest_domain%jend_coarse(n)
450  endif
451  enddo
452 
453  if(pos .NE. nest_domain%nest(l)%num_nest) &
454  call mpp_error(fatal, .NE."shift mpp_define_nest_domains.inc:pos nest_domain%nest(l)%num_nest")
455 
456  ! The nest may connect to different parent PEs after motion, so this must be recalculated
457  call define_nest_level_type(nest_domain%nest(l), nest_domain%nest(l)%x_refine, &
458  & nest_domain%nest(l)%y_refine, extra_halo_local)
459  endif
460  enddo
461 
462 end subroutine mpp_shift_nest_domains
463 
464 
465 subroutine define_nest_level_type(nest_domain, x_refine, y_refine, extra_halo)
466  type(nest_level_type), intent(inout) :: nest_domain !< nest domain to be defined
467  integer, intent(in ) :: extra_halo !< halo value
468  integer, intent(in ) :: x_refine, y_refine !< x and y refinements
469 
470  integer :: n
471  integer :: npes, npes_fine, npes_coarse
472  integer, allocatable :: pes_coarse(:)
473  integer, allocatable :: pes_fine(:)
474  integer, dimension(nest_domain%num_nest) :: my_nest_id
475  integer :: my_num_nest
476 
477  npes = size(nest_domain%pelist(:))
478  npes_coarse = size(nest_domain%domain_coarse%list(:))
479  npes_fine = size(nest_domain%domain_fine%list(:))
480  !--- pes_fine and pes_coarse should be subset of pelist
481  allocate( pes_coarse(npes_coarse) )
482  allocate( pes_fine(npes_fine ) )
483  do n = 1, npes_coarse
484  pes_coarse(n) = nest_domain%domain_coarse%list(n-1)%pe
485  if( .NOT. any(nest_domain%pelist(:) == pes_coarse(n)) ) then
486  call mpp_error(fatal, "mpp_define_nest_domains.inc: pelist_coarse is not subset of pelist")
487  endif
488  enddo
489  do n = 1, npes_fine
490  pes_fine(n) = nest_domain%domain_fine%list(n-1)%pe
491  if( .NOT. any(nest_domain%pelist(:) == pes_fine(n)) ) then
492  call mpp_error(fatal, "mpp_define_nest_domains.inc: pelist_fine is not subset of pelist")
493  endif
494  enddo
495 
496  if (associated(nest_domain%pelist_fine)) deallocate(nest_domain%pelist_fine) !< Check if allocated
497  allocate(nest_domain%pelist_fine(npes_fine))
498  if (associated(nest_domain%pelist_coarse)) deallocate(nest_domain%pelist_coarse) !< Check if allocated
499  allocate(nest_domain%pelist_coarse(npes_coarse))
500  nest_domain%pelist_fine = pes_fine
501  nest_domain%pelist_coarse = pes_coarse
502  if( nest_domain%is_fine_pe .neqv. any(pes_fine(:) == mpp_pe()) ) then
503  call mpp_error(fatal, .neqv."mpp_define_nest_domains.inc: nest_domain%is_fine_pe ANY(pes_fine(:) == mpp_pe())")
504  endif
505  if( nest_domain%is_coarse_pe .neqv. any(pes_coarse(:) == mpp_pe()) ) then
506  call mpp_error(fatal, "mpp_define_nest_domains.inc: "// &
507  .neqv."nest_domain%is_coarse_pe ANY(pes_coarse(:) == mpp_pe())")
508  endif
509 
510  !--- figure out the the corresponding nested region.
511  !--- on coarse grid pe, it might overlap multiple fine regon.
512  !--- on fine grid pe, it always only overlap at most 1 coarse region.
513  my_num_nest= 0
514  my_nest_id(:) = 0
515  if( nest_domain%is_fine_pe ) then
516  !--- figure out the nest number on current pe
517  do n = 1, nest_domain%num_nest
518  if(nest_domain%domain_fine%tile_id(1) == nest_domain%tile_fine(n)) then
519  my_num_nest = my_num_nest + 1
520  my_nest_id(my_num_nest) = n
521  exit
522  end if
523  end do
524  if(my_num_nest .NE. 1) then
525  print*, "num_nest=", my_num_nest, nest_domain%domain_fine%tile_id(1), nest_domain%tile_fine(1)
526  call mpp_error(fatal, .ne."mpp_define_nest_domains.inc: my_num_nest 1 on fine pelist")
527  endif
528  else if( nest_domain%is_coarse_pe ) then
529  my_num_nest = nest_domain%num_nest
530  do n = 1, nest_domain%num_nest
531  my_nest_id(n) = n
532  enddo
533  endif
534 
535  nest_domain%my_num_nest = my_num_nest
536  if(my_num_nest>0) then
537  allocate(nest_domain%my_nest_id(my_num_nest))
538  nest_domain%my_nest_id(:) = my_nest_id(1:my_num_nest)
539  endif
540 
541  !--- We are assuming the fine grid is fully overlapped with coarse grid.
542  if( nest_domain%is_fine_pe ) then
543  if( nest_domain%iend_fine(my_nest_id(1))-nest_domain%istart_fine(my_nest_id(1))+1 &
544  .NE. nest_domain%domain_fine%x(1)%global%size .OR. &
545  nest_domain%jend_fine(my_nest_id(1))-nest_domain%jstart_fine(my_nest_id(1))+1 &
546  .NE. nest_domain%domain_fine%y(1)%global%size ) then
547  print*, "x size are", nest_domain%domain_fine%x(1)%global%size, &
548  nest_domain%istart_fine(my_nest_id(1)), nest_domain%iend_fine(my_nest_id(1))
549  print*, "y size are", nest_domain%domain_fine%y(1)%global%size, &
550  nest_domain%jstart_fine(my_nest_id(1)), nest_domain%jend_fine(my_nest_id(1))
551  call mpp_error(fatal, "mpp_define_nest_domains.inc: The fine global domain is not covered by coarse domain")
552  endif
553  endif
554 
555  ! only support concurrent run for fine and coarse domain, currently only check on coarse pe
556  if(nest_domain%is_coarse_pe) then
557 ! if( npes_fine + npes_coarse .NE. npes ) then
558 ! print*, "On pe =", mpp_pe(), npes_fine, npes_coarse, npes
559 ! call mpp_error(FATAL, "mpp_domains_define.inc: size(pelist_coarse)+size(pelist_fine) .NE. size(pelist)")
560 ! endif
561  endif
562 
563  !--- coarse grid and fine grid should be both symmetry or non-symmetry.
564  if(nest_domain%domain_coarse%symmetry .neqv. nest_domain%domain_fine%symmetry) then
565  print*,"symmetry is", nest_domain%domain_coarse%symmetry, nest_domain%domain_fine%symmetry, mpp_pe()
566  call mpp_error(fatal, .neqv..NOT."mpp_domains_define.inc: domain_coarse%symmetry domain_fine%symmetry")
567  endif
568 
569  nest_domain%x_refine = x_refine
570  nest_domain%y_refine = y_refine
571 
572  if (associated(nest_domain%C2F_T)) deallocate(nest_domain%C2F_T) !< Check if allocated
573  if (associated(nest_domain%C2F_C)) deallocate(nest_domain%C2F_C) !< Check if allocated
574  if (associated(nest_domain%C2F_E)) deallocate(nest_domain%C2F_E) !< Check if allocated
575  if (associated(nest_domain%C2F_N)) deallocate(nest_domain%C2F_N) !< Check if allocated
576  allocate( nest_domain%C2F_T, nest_domain%C2F_C, nest_domain%C2F_E, nest_domain%C2F_N )
577  nest_domain%C2F_T%next => null()
578  nest_domain%C2F_C%next => null()
579  nest_domain%C2F_N%next => null()
580  nest_domain%C2F_E%next => null()
581  if (associated(nest_domain%F2C_T)) deallocate(nest_domain%F2C_T) !< Check if allocated
582  if (associated(nest_domain%F2C_C)) deallocate(nest_domain%F2C_C) !< Check if allocated
583  if (associated(nest_domain%F2C_E)) deallocate(nest_domain%F2C_E) !< Check if allocated
584  if (associated(nest_domain%F2C_N)) deallocate(nest_domain%F2C_N) !< Check if allocated
585  allocate( nest_domain%F2C_T, nest_domain%F2C_C, nest_domain%F2C_E, nest_domain%F2C_N )
586 
587  call compute_overlap_fine_to_coarse(nest_domain, nest_domain%F2C_T, center, "F2C T-cell")
588  call compute_overlap_fine_to_coarse(nest_domain, nest_domain%F2C_E, east, "F2C E-cell")
589  call compute_overlap_fine_to_coarse(nest_domain, nest_domain%F2C_C, corner, "F2C C-cell")
590  call compute_overlap_fine_to_coarse(nest_domain, nest_domain%F2C_N, north, "F2C N-cell")
591 
592  call compute_overlap_coarse_to_fine(nest_domain, nest_domain%C2F_T, extra_halo, center, "C2F T-cell")
593  call compute_overlap_coarse_to_fine(nest_domain, nest_domain%C2F_E, extra_halo, east, "C2F E-cell")
594  call compute_overlap_coarse_to_fine(nest_domain, nest_domain%C2F_C, extra_halo, corner, "C2F C-cell")
595  call compute_overlap_coarse_to_fine(nest_domain, nest_domain%C2F_N, extra_halo, north, "C2F N-cell")
596 
597  deallocate(pes_fine, pes_coarse)
598 
599 
600 end subroutine define_nest_level_type
601 
602 
603 !###############################################################################
604 subroutine compute_overlap_coarse_to_fine(nest_domain, overlap, extra_halo, position, name)
605  type(nest_level_type), intent(inout) :: nest_domain
606  type(nestspec), intent(inout) :: overlap
607  integer, intent(in ) :: extra_halo
608  integer, intent(in ) :: position
609  character(len=*), intent(in ) :: name
610 
611  type(domain2d), pointer :: domain_fine =>null()
612  type(domain2d), pointer :: domain_coarse=>null()
613  type(overlap_type), allocatable :: overlapList(:)
614  logical :: is_first
615  integer :: tile_fine, tile_coarse
616  integer :: istart_fine, iend_fine, jstart_fine, jend_fine
617  integer :: istart_coarse, iend_coarse, jstart_coarse, jend_coarse
618  integer :: whalo, ehalo, shalo, nhalo
619  integer :: npes, npes_fine, npes_coarse, n, m
620  integer :: isg_fine, ieg_fine, jsg_fine, jeg_fine
621  integer :: isc_coarse, iec_coarse, jsc_coarse, jec_coarse
622  integer :: is_coarse, ie_coarse, js_coarse, je_coarse
623  integer :: is_coarse2, ie_coarse2, js_coarse2, je_coarse2
624  integer :: rotate
625  integer :: is_convert2(2), ie_convert2(2), js_convert2(2), je_convert2(2), rotate2(2)
626  integer :: isc_fine, iec_fine, jsc_fine, jec_fine
627  integer :: isd_fine, ied_fine, jsd_fine, jed_fine
628  integer :: x_refine, y_refine, ishift, jshift
629  integer :: nsend, nrecv, dir, l, nn
630  integer :: nconvert
631  integer, allocatable :: isl_coarse(:), iel_coarse(:), jsl_coarse(:), jel_coarse(:)
632  integer, allocatable :: isl_fine(:), iel_fine(:), jsl_fine(:), jel_fine(:)
633  integer, allocatable :: isgl_fine(:), iegl_fine(:), jsgl_fine(:), jegl_fine(:)
634  integer :: outunit
635 
636 
637  outunit = stdout()
638  domain_fine => nest_domain%domain_fine
639  domain_coarse => nest_domain%domain_coarse
640  call mpp_get_domain_shift (domain_coarse, ishift, jshift, position)
641  npes = mpp_npes()
642  npes_fine = size(nest_domain%pelist_fine(:))
643  npes_coarse = size(nest_domain%pelist_coarse(:))
644 
645  allocate(isl_coarse(npes_coarse), iel_coarse(npes_coarse))
646  allocate(jsl_coarse(npes_coarse), jel_coarse(npes_coarse))
647  allocate(isl_fine(npes_fine ), iel_fine(npes_fine ))
648  allocate(jsl_fine(npes_fine ), jel_fine(npes_fine ))
649  allocate(isgl_fine(npes_fine ), iegl_fine(npes_fine ))
650  allocate(jsgl_fine(npes_fine ), jegl_fine(npes_fine ))
651 
652  call mpp_get_global_domain (domain_fine, xbegin=isg_fine, xend=ieg_fine, &
653  ybegin=jsg_fine, yend=jeg_fine)
654  call mpp_get_compute_domain (domain_coarse, xbegin=isc_coarse, xend=iec_coarse, &
655  ybegin=jsc_coarse, yend=jec_coarse)
656  call mpp_get_compute_domain (domain_fine, xbegin=isc_fine, xend=iec_fine, &
657  ybegin=jsc_fine, yend=jec_fine)
658  call mpp_get_compute_domains(domain_coarse, xbegin=isl_coarse, xend=iel_coarse, &
659  ybegin=jsl_coarse, yend=jel_coarse)
660  call mpp_get_compute_domains(domain_fine, xbegin=isl_fine, xend=iel_fine, &
661  ybegin=jsl_fine, yend=jel_fine)
662  call mpp_get_global_domains(domain_fine, xbegin=isgl_fine, xend=iegl_fine, &
663  ybegin=jsgl_fine, yend=jegl_fine)
664 
665  if( nest_domain%is_coarse_pe ) then
666  allocate(overlaplist(npes_fine))
667  overlap%xbegin = isc_coarse - domain_coarse%whalo
668  overlap%xend = iec_coarse + domain_coarse%ehalo + ishift
669  overlap%ybegin = jsc_coarse - domain_coarse%shalo
670  overlap%yend = jec_coarse + domain_coarse%nhalo + jshift
671  else
672  allocate(overlaplist(npes_coarse))
673  overlap%xbegin = isc_fine - domain_fine%whalo
674  overlap%xend = iec_fine + domain_fine%ehalo + ishift
675  overlap%ybegin = jsc_fine - domain_fine%shalo
676  overlap%yend = jec_fine + domain_fine%nhalo + jshift
677  endif
678 
679  overlap%extra_halo = extra_halo
680  x_refine = nest_domain%x_refine
681  y_refine = nest_domain%y_refine
682  whalo = domain_fine%whalo + extra_halo
683  ehalo = domain_fine%ehalo + extra_halo
684  shalo = domain_fine%shalo + extra_halo
685  nhalo = domain_fine%nhalo + extra_halo
686 
687  isd_fine = isc_fine - whalo
688  ied_fine = iec_fine + ehalo
689  jsd_fine = jsc_fine - shalo
690  jed_fine = jec_fine + nhalo
691 
692  overlap%nsend = 0
693  overlap%nrecv = 0
694  call init_index_type(overlap%west)
695  call init_index_type(overlap%east)
696  call init_index_type(overlap%south)
697  call init_index_type(overlap%north)
698  nsend = 0
699  nrecv = 0
700 
701  do nn = 1, nest_domain%num_nest
702 
703  tile_fine = nest_domain%tile_fine(nn)
704  tile_coarse = nest_domain%tile_coarse(nn)
705  istart_fine = nest_domain%istart_fine(nn)
706  iend_fine = nest_domain%iend_fine(nn)
707  jstart_fine = nest_domain%jstart_fine(nn)
708  jend_fine = nest_domain%jend_fine(nn)
709  istart_coarse = nest_domain%istart_coarse(nn)
710  iend_coarse = nest_domain%iend_coarse(nn)
711  jstart_coarse = nest_domain%jstart_coarse(nn)
712  jend_coarse = nest_domain%jend_coarse(nn)
713 
714  !--- first compute the halo region and corresponding index in coarse grid.
715  if( nest_domain%is_fine_pe .and. domain_fine%tile_id(1) == tile_fine) then
716  if( ieg_fine == iec_fine ) then ! east halo
717  is_coarse = iend_coarse
718  ie_coarse = iend_coarse + ehalo
719  js_coarse = jstart_coarse + ( jsc_fine - jsg_fine )/y_refine
720  je_coarse = jstart_coarse + ( jec_fine - jsg_fine )/y_refine
721  js_coarse = js_coarse - shalo
722  je_coarse = je_coarse + nhalo
723 
724  overlap%east%is_me = iec_fine + 1
725  overlap%east%ie_me = ied_fine
726  overlap%east%js_me = jsd_fine
727  overlap%east%je_me = jed_fine
728  overlap%east%is_you = is_coarse
729  overlap%east%ie_you = ie_coarse
730  overlap%east%js_you = js_coarse
731  overlap%east%je_you = je_coarse
732  endif
733 
734  if( jsg_fine == jsc_fine ) then ! south
735  is_coarse = istart_coarse + ( isc_fine - isg_fine )/x_refine
736  ie_coarse = istart_coarse + ( iec_fine - isg_fine )/x_refine
737  is_coarse = is_coarse - whalo
738  ie_coarse = ie_coarse + ehalo
739  js_coarse = jstart_coarse - shalo
740  je_coarse = jstart_coarse
741  overlap%south%is_me = isd_fine
742  overlap%south%ie_me = ied_fine
743  overlap%south%js_me = jsd_fine
744  overlap%south%je_me = jsc_fine-1
745  overlap%south%is_you = is_coarse
746  overlap%south%ie_you = ie_coarse
747  overlap%south%js_you = js_coarse
748  overlap%south%je_you = je_coarse
749  endif
750 
751  if( isg_fine == isc_fine ) then ! west
752  is_coarse = istart_coarse - whalo
753  ie_coarse = istart_coarse
754  js_coarse = jstart_coarse + ( jsc_fine - jsg_fine )/y_refine
755  je_coarse = jstart_coarse + ( jec_fine - jsg_fine )/y_refine
756  js_coarse = js_coarse - shalo
757  je_coarse = je_coarse + nhalo
758  overlap%west%is_me = isd_fine
759  overlap%west%ie_me = isc_fine-1
760  overlap%west%js_me = jsd_fine
761  overlap%west%je_me = jed_fine
762  overlap%west%is_you = is_coarse
763  overlap%west%ie_you = ie_coarse
764  overlap%west%js_you = js_coarse
765  overlap%west%je_you = je_coarse
766  endif
767 
768  if( jeg_fine == jec_fine ) then ! north
769  is_coarse = istart_coarse + ( isc_fine - isg_fine )/x_refine
770  ie_coarse = istart_coarse + ( iec_fine - isg_fine )/x_refine
771  is_coarse = is_coarse - whalo
772  ie_coarse = ie_coarse + ehalo
773  js_coarse = jend_coarse
774  je_coarse = jend_coarse + nhalo
775  overlap%north%is_me = isd_fine
776  overlap%north%ie_me = ied_fine
777  overlap%north%js_me = jec_fine+1
778  overlap%north%je_me = jed_fine
779  overlap%north%is_you = is_coarse
780  overlap%north%ie_you = ie_coarse
781  overlap%north%js_you = js_coarse
782  overlap%north%je_you = je_coarse
783  endif
784 
785  !-------------------------------------------------------------------------
786  !
787  ! Receiving
788  !
789  !-------------------------------------------------------------------------
790  !--- loop through coarse pelist
791  do n = 1, npes_coarse
792  is_first = .true.
793  do m = 1, 4
794  select case (m)
795  case (1) !--- east halo receiving
796  dir = 1
797  is_coarse = overlap%east%is_you
798  ie_coarse = overlap%east%ie_you
799  js_coarse = overlap%east%js_you
800  je_coarse = overlap%east%je_you
801  case (2) !--- south halo receiving
802  dir = 3
803  is_coarse = overlap%south%is_you
804  ie_coarse = overlap%south%ie_you
805  js_coarse = overlap%south%js_you
806  je_coarse = overlap%south%je_you
807  case (3) !--- west halo receiving
808  dir = 5
809  is_coarse = overlap%west%is_you
810  ie_coarse = overlap%west%ie_you
811  js_coarse = overlap%west%js_you
812  je_coarse = overlap%west%je_you
813  case (4) !--- north halo receiving
814  dir = 7
815  is_coarse = overlap%north%is_you
816  ie_coarse = overlap%north%ie_you
817  js_coarse = overlap%north%js_you
818  je_coarse = overlap%north%je_you
819  end select
820  if( je_coarse .GE. js_coarse .AND. ie_coarse .GE. is_coarse ) then
821  ! convert coarse grid index to the nested grid coarse grid index.
822  nconvert = convert_index_to_nest(domain_coarse, 0, 0, tile_coarse, istart_coarse, iend_coarse, &
823  jstart_coarse, jend_coarse, domain_coarse%ntiles, domain_coarse%list(n-1)%tile_id(1), &
824  isl_coarse(n), iel_coarse(n), jsl_coarse(n), jel_coarse(n), &
825  is_convert2, ie_convert2, js_convert2, je_convert2, rotate2)
826  do l = 1, nconvert
827  is_coarse2 = max( is_coarse, is_convert2(l) )
828  ie_coarse2 = min( ie_coarse, ie_convert2(l) )
829  js_coarse2 = max( js_coarse, js_convert2(l) )
830  je_coarse2 = min( je_coarse, je_convert2(l) )
831  if( ie_coarse2 .GE. is_coarse2 .AND. je_coarse2 .GE. js_coarse2 ) then
832  select case (m)
833  case (1) !--- east halo
834  is_coarse2 = is_coarse2+ishift
835  ie_coarse2 = ie_coarse2+ishift
836  if(je_coarse2 == overlap%east%je_you) je_coarse2 = je_coarse2+jshift
837  case (2) !--- south halo
838  if(ie_coarse2 == overlap%south%ie_you) ie_coarse2 = ie_coarse2+ishift
839  case (3) !--- west halo
840  if(je_coarse2 == overlap%west%je_you) je_coarse2 = je_coarse2+jshift
841  case (4) !--- north halo
842  if(ie_coarse2 == overlap%north%ie_you) ie_coarse2 = ie_coarse2+ishift
843  js_coarse2 = js_coarse2+jshift
844  je_coarse2 = je_coarse2+jshift
845  end select
846 
847  if(is_first) then
848  nrecv = nrecv + 1
849  call allocate_nest_overlap(overlaplist(nrecv), maxoverlap)
850  is_first = .false.
851  endif
852  rotate = -rotate2(l)
853  call insert_nest_overlap(overlaplist(nrecv), nest_domain%pelist_coarse(n), &
854  is_coarse2, ie_coarse2, js_coarse2, je_coarse2 , dir, rotate2(l))
855  endif
856  enddo
857  endif
858  enddo
859  enddo
860 
861  endif
862  !-----------------------------------------------------------------------
863  !
864  ! Sending
865  !
866  !-----------------------------------------------------------------------
867 
868  if( nest_domain%is_coarse_pe ) then
869  do n = 1, npes_fine
870  if( domain_fine%list(n-1)%tile_id(1) .NE. tile_fine ) cycle
871  is_first = .true.
872  isg_fine = isgl_fine(n)
873  ieg_fine = iegl_fine(n)
874  jsg_fine = jsgl_fine(n)
875  jeg_fine = jegl_fine(n)
876 
877  !--- to_pe's east
878  if( ieg_fine == iel_fine(n) ) then
879  dir = 1
880  is_coarse = iend_coarse
881  ie_coarse = iend_coarse + ehalo
882  js_coarse = jstart_coarse + ( jsl_fine(n) - jsg_fine )/y_refine
883  je_coarse = jstart_coarse + ( jel_fine(n) - jsg_fine )/y_refine
884  js_coarse = js_coarse - shalo
885  je_coarse = je_coarse + nhalo
886  !--- convert the index to coarse grid index.
887  nconvert = convert_index_to_coarse(domain_coarse, 0, 0, tile_coarse, istart_coarse, iend_coarse, &
888  & jstart_coarse, jend_coarse, domain_coarse%ntiles, domain_coarse%tile_id(1), is_coarse, ie_coarse,&
889  & js_coarse, je_coarse, is_convert2, ie_convert2, js_convert2, je_convert2, rotate2)
890  do l = 1, nconvert
891  is_coarse = max(isc_coarse, is_convert2(l))
892  ie_coarse = min(iec_coarse, ie_convert2(l))
893  js_coarse = max(jsc_coarse, js_convert2(l))
894  je_coarse = min(jec_coarse, je_convert2(l))
895  if( ie_coarse .GE. is_coarse .AND. je_coarse .GE. js_coarse ) then
896  if(rotate2(l)==zero) then
897  is_coarse = is_coarse+ishift
898  ie_coarse = ie_coarse+ishift
899  if( je_coarse == je_convert2(l) ) je_coarse = je_coarse+jshift
900  else if(rotate2(l) == minus_ninety) then
901  js_coarse = js_coarse+ishift
902  je_coarse = je_coarse+ishift
903  if(is_coarse==is_convert2(l)) is_coarse = is_coarse-jshift
904  is_coarse = is_coarse+jshift
905  ie_coarse = ie_coarse+jshift
906  else if(rotate2(l) == ninety) then
907  if(ie_coarse==ie_convert2(l)) ie_coarse = ie_coarse+jshift
908  endif
909 
910  if(is_first) then
911  nsend = nsend + 1
912  call allocate_nest_overlap(overlaplist(nsend), maxoverlap)
913  is_first = .false.
914  endif
915  rotate = -rotate2(l)
916  call insert_nest_overlap(overlaplist(nsend), nest_domain%pelist_fine(n), &
917  is_coarse, ie_coarse, js_coarse, je_coarse , dir, rotate)
918  endif
919  enddo
920  endif
921 
922  !--- to_pe's south
923  if( jsg_fine == jsl_fine(n) ) then
924  dir = 3
925  is_coarse = istart_coarse + ( isl_fine(n) - isg_fine )/x_refine
926  ie_coarse = istart_coarse + ( iel_fine(n) - isg_fine )/x_refine
927  is_coarse = is_coarse - shalo
928  ie_coarse = ie_coarse + nhalo
929  js_coarse = jstart_coarse - shalo
930  je_coarse = jstart_coarse
931  !--- convert the index to coarse grid index.
932  nconvert=convert_index_to_coarse(domain_coarse, 0, 0, tile_coarse, istart_coarse, iend_coarse, &
933  & jstart_coarse, jend_coarse, domain_coarse%ntiles, domain_coarse%tile_id(1), is_coarse, &
934  & ie_coarse, js_coarse, je_coarse, is_convert2, ie_convert2, js_convert2, je_convert2, rotate2)
935  do l = 1, nconvert
936  is_coarse = max(isc_coarse, is_convert2(l))
937  ie_coarse = min(iec_coarse, ie_convert2(l))
938  js_coarse = max(jsc_coarse, js_convert2(l))
939  je_coarse = min(jec_coarse, je_convert2(l))
940 
941  if( ie_coarse .GE. is_coarse .AND. je_coarse .GE. js_coarse ) then
942  if(rotate2(l)==zero .AND. ie_coarse==ie_convert2(l)) then
943  ie_coarse = ie_coarse+ishift
944  else if( rotate2(l) .NE. zero .AND. je_coarse == je_convert2(l) ) then
945  je_coarse = je_coarse+ishift
946  endif
947  if(is_first) then
948  nsend = nsend + 1
949  call allocate_nest_overlap(overlaplist(nsend), maxoverlap)
950  is_first = .false.
951  endif
952  rotate = -rotate2(l)
953  call insert_nest_overlap(overlaplist(nsend), nest_domain%pelist_fine(n), &
954  is_coarse, ie_coarse, js_coarse, je_coarse , dir, rotate)
955  endif
956  enddo
957  endif
958 
959  !--- to_pe's west
960  if( isg_fine == isl_fine(n) ) then
961  dir = 5
962  is_coarse = istart_coarse - whalo
963  ie_coarse = istart_coarse
964  js_coarse = jstart_coarse + ( jsl_fine(n) - jsg_fine )/y_refine
965  je_coarse = jstart_coarse + ( jel_fine(n) - jsg_fine )/y_refine
966  js_coarse = js_coarse - shalo
967  je_coarse = je_coarse + nhalo
968  !--- convert the index to coarse grid index.
969  nconvert=convert_index_to_coarse(domain_coarse, 0, 0, tile_coarse, istart_coarse, iend_coarse, &
970  & jstart_coarse, jend_coarse, domain_coarse%ntiles, domain_coarse%tile_id(1), is_coarse, &
971  & ie_coarse, js_coarse, je_coarse, is_convert2, ie_convert2, js_convert2, je_convert2, rotate2)
972  do l = 1, nconvert
973  is_coarse = max(isc_coarse, is_convert2(l))
974  ie_coarse = min(iec_coarse, ie_convert2(l))
975  js_coarse = max(jsc_coarse, js_convert2(l))
976  je_coarse = min(jec_coarse, je_convert2(l))
977  if( ie_coarse .GE. is_coarse .AND. je_coarse .GE. js_coarse ) then
978  if(rotate2(l)==zero .and. je_coarse == je_convert2(l) ) then
979  je_coarse = je_coarse+jshift
980  else if(rotate2(l) .NE. zero .and. ie_coarse == ie_convert2(l) ) then
981  ie_coarse = ie_coarse+jshift
982  endif
983  if(is_first) then
984  nsend = nsend + 1
985  call allocate_nest_overlap(overlaplist(nsend), maxoverlap)
986  is_first = .false.
987  endif
988  rotate = -rotate2(l)
989  call insert_nest_overlap(overlaplist(nsend), nest_domain%pelist_fine(n), &
990  is_coarse, ie_coarse, js_coarse, je_coarse , dir, rotate)
991  endif
992  enddo
993  endif
994 
995  !--- to_pe's north
996  if( jeg_fine == jel_fine(n) ) then
997  dir = 7
998  is_coarse = istart_coarse + ( isl_fine(n) - isg_fine )/x_refine
999  ie_coarse = istart_coarse + ( iel_fine(n) - isg_fine )/x_refine
1000  is_coarse = is_coarse - shalo
1001  ie_coarse = ie_coarse + nhalo
1002  js_coarse = jend_coarse
1003  je_coarse = jend_coarse + nhalo
1004  !--- convert the index to coarse grid index.
1005  nconvert=convert_index_to_coarse(domain_coarse, 0, 0, tile_coarse, istart_coarse, iend_coarse, &
1006  & jstart_coarse, jend_coarse, domain_coarse%ntiles, domain_coarse%tile_id(1), is_coarse, &
1007  & ie_coarse, js_coarse, je_coarse, is_convert2, ie_convert2, js_convert2, je_convert2, rotate2)
1008  do l = 1, nconvert
1009  is_coarse = max(isc_coarse, is_convert2(l))
1010  ie_coarse = min(iec_coarse, ie_convert2(l))
1011  js_coarse = max(jsc_coarse, js_convert2(l))
1012  je_coarse = min(jec_coarse, je_convert2(l))
1013  if( ie_coarse .GE. is_coarse .AND. je_coarse .GE. js_coarse ) then
1014  if(rotate2(l)==zero) then
1015  if(ie_coarse==ie_convert2(l)) ie_coarse = ie_coarse+ishift
1016  js_coarse = js_coarse+jshift
1017  je_coarse = je_coarse+jshift
1018  else if(rotate2(l) == ninety) then
1019  if(js_coarse==js_convert2(l)) js_coarse = js_coarse-ishift
1020  js_coarse = js_coarse+ishift
1021  je_coarse = je_coarse+ishift
1022  is_coarse = is_coarse+jshift
1023  ie_coarse = ie_coarse+jshift
1024  else if(rotate2(l) == minus_ninety ) then
1025  if(je_coarse==je_convert2(l)) je_coarse = je_coarse+ishift
1026  endif
1027  if(is_first) then
1028  nsend = nsend + 1
1029  call allocate_nest_overlap(overlaplist(nsend), maxoverlap)
1030  is_first = .false.
1031  endif
1032  rotate = -rotate2(l)
1033  call insert_nest_overlap(overlaplist(nsend), nest_domain%pelist_fine(n), &
1034  is_coarse, ie_coarse, js_coarse, je_coarse , dir, rotate)
1035  endif
1036  enddo
1037  endif
1038  enddo
1039  endif
1040  enddo
1041 
1042  !--- copy the overlapping into nest_domain data.
1043  overlap%nrecv = nrecv
1044  if( nrecv > 0 ) then
1045  if (associated(overlap%recv)) deallocate(overlap%recv) !< Check if allocated
1046  allocate(overlap%recv(nrecv))
1047  do n = 1, nrecv
1048  call copy_nest_overlap( overlap%recv(n), overlaplist(n) )
1049 ! call print_nest_overlap(overlap%recv(n), "C2F RECV")
1050  call deallocate_nest_overlap( overlaplist(n) )
1051  enddo
1052  endif
1053 
1054  overlap%nsend = nsend
1055  if( nsend > 0 ) then
1056  if (associated(overlap%send)) deallocate(overlap%send) !< Check if allocated
1057  allocate(overlap%send(nsend))
1058  do n = 1, nsend
1059  call copy_nest_overlap( overlap%send(n), overlaplist(n) )
1060 ! call print_nest_overlap(overlap%send(n), "C2F SEND")
1061  call deallocate_nest_overlap( overlaplist(n) )
1062  enddo
1063  endif
1064  if(allocated(overlaplist))deallocate(overlaplist)
1065 
1066 
1067  deallocate(isl_coarse, iel_coarse, jsl_coarse, jel_coarse)
1068  deallocate(isl_fine, iel_fine, jsl_fine, jel_fine)
1069  deallocate(isgl_fine, iegl_fine, jsgl_fine, jegl_fine)
1070 
1071  !--- add shift value accoring grid position
1072  if( nest_domain%is_fine_pe ) then
1073  if( ieg_fine == iec_fine ) then ! east halo
1074  overlap%east%is_me = overlap%east%is_me + ishift
1075  overlap%east%ie_me = overlap%east%ie_me + ishift
1076  overlap%east%je_me = overlap%east%je_me + jshift
1077  overlap%east%is_you = overlap%east%is_you + ishift
1078  overlap%east%ie_you = overlap%east%ie_you + ishift
1079  overlap%east%je_you = overlap%east%je_you + jshift
1080  endif
1081 
1082  if( jsg_fine == jsc_fine ) then ! south
1083  overlap%south%ie_me = overlap%south%ie_me + ishift
1084  overlap%south%ie_you = overlap%south%ie_you + ishift
1085  endif
1086 
1087  if( isg_fine == isc_fine ) then ! west
1088  overlap%west%je_me = overlap%west%je_me + jshift
1089  overlap%west%je_you = overlap%west%je_you + jshift
1090  endif
1091 
1092  if( jeg_fine == jec_fine ) then ! north
1093  overlap%north%ie_me = overlap%north%ie_me + ishift
1094  overlap%north%js_me = overlap%north%js_me + jshift
1095  overlap%north%je_me = overlap%north%je_me + jshift
1096  overlap%north%ie_you = overlap%north%ie_you + ishift
1097  overlap%north%js_you = overlap%north%js_you + jshift
1098  overlap%north%je_you = overlap%north%je_you + jshift
1099  endif
1100  endif
1101 
1102  if(debug_message_passing) call debug_message_size(overlap, name)
1103 
1104 
1105 end subroutine compute_overlap_coarse_to_fine
1106 
1107 !###############################################################################
1108 !> This routine will compute the send and recv information between overlapped nesting
1109 !! region. The data is assumed on T-cell center.
1110 subroutine compute_overlap_fine_to_coarse(nest_domain, overlap, position, name)
1111  type(nest_level_type), intent(inout) :: nest_domain
1112  type(nestspec), intent(inout) :: overlap
1113  integer, intent(in ) :: position
1114  character(len=*), intent(in ) :: name
1115 
1116  !--- local variables
1117 
1118  type(domain2d), pointer :: domain_fine =>null()
1119  type(domain2d), pointer :: domain_coarse=>null()
1120  type(overlap_type), allocatable :: overlapList(:)
1121  integer :: tile_fine, tile_coarse
1122  integer :: istart_fine, iend_fine, jstart_fine, jend_fine
1123  integer :: istart_coarse, iend_coarse, jstart_coarse, jend_coarse
1124  integer :: npes_fine, npes_coarse, n
1125  integer :: isg_fine, ieg_fine, jsg_fine, jeg_fine
1126  integer :: isc_coarse, iec_coarse, jsc_coarse, jec_coarse
1127  integer :: is_coarse, ie_coarse, js_coarse, je_coarse
1128  integer :: isc_fine, iec_fine, jsc_fine, jec_fine
1129  integer :: is_you, ie_you, js_you, je_you
1130  integer :: x_refine, y_refine
1131  integer :: nsend, nrecv, dir
1132  integer, allocatable :: isl_coarse(:), iel_coarse(:), jsl_coarse(:), jel_coarse(:)
1133  integer, allocatable :: isl_fine(:), iel_fine(:), jsl_fine(:), jel_fine(:)
1134  integer :: is_convert2(2), ie_convert2(2), js_convert2(2), je_convert2(2), rotate2(2)
1135  integer :: is2, ie2, js2, je2, nconvert
1136  integer :: xbegin_c, xend_c, ybegin_c, yend_c
1137  integer :: ishift, jshift, l, is3, ie3, js3, je3, nn
1138 
1139  domain_fine => nest_domain%domain_fine
1140  domain_coarse => nest_domain%domain_coarse
1141  npes_fine = size(nest_domain%pelist_fine(:))
1142  npes_coarse = size(nest_domain%pelist_coarse(:))
1143 
1144  allocate(isl_coarse(npes_coarse), iel_coarse(npes_coarse) )
1145  allocate(jsl_coarse(npes_coarse), jel_coarse(npes_coarse) )
1146  allocate(isl_fine(npes_fine), iel_fine(npes_fine) )
1147  allocate(jsl_fine(npes_fine), jel_fine(npes_fine) )
1148  call mpp_get_domain_shift (domain_coarse, ishift, jshift, position)
1149 
1150  call mpp_get_compute_domain (domain_coarse, xbegin=isc_coarse, xend=iec_coarse, ybegin=jsc_coarse, yend=jec_coarse)
1151  call mpp_get_compute_domain (domain_fine, xbegin=isc_fine, xend=iec_fine, ybegin=jsc_fine, yend=jec_fine)
1152  call mpp_get_compute_domains(domain_coarse, xbegin=isl_coarse, xend=iel_coarse, ybegin=jsl_coarse, yend=jel_coarse)
1153  call mpp_get_compute_domains(domain_fine, xbegin=isl_fine, xend=iel_fine, ybegin=jsl_fine, yend=jel_fine)
1154  call mpp_get_global_domain (domain_fine, xbegin=isg_fine, xend=ieg_fine, ybegin=jsg_fine, yend=jeg_fine)
1155  overlap%center%is_you = 0; overlap%center%ie_you = -1
1156  overlap%center%js_you = 0; overlap%center%je_you = -1
1157 
1158  overlap%nsend = 0
1159  overlap%nrecv = 0
1160  call init_index_type(overlap%center)
1161 
1162  if( nest_domain%is_fine_pe ) then
1163  overlap%xbegin = 0; overlap%xend = -1
1164  overlap%ybegin = 0; overlap%yend = -1
1165  else
1166  overlap%xbegin = isc_coarse - domain_coarse%whalo
1167  overlap%xend = iec_coarse + domain_coarse%ehalo + ishift
1168  overlap%ybegin = jsc_coarse - domain_coarse%shalo
1169  overlap%yend = jec_coarse + domain_coarse%nhalo + jshift
1170  overlap%xsize_c = overlap%xend - overlap%xbegin + 1
1171  overlap%ysize_c = overlap%yend - overlap%ybegin + 1
1172  overlap%xbegin_f = 0
1173  overlap%xend_f = -1
1174  overlap%ybegin_f = 0
1175  overlap%yend_f = -1
1176  overlap%xbegin_c = 0
1177  overlap%xend_c = -1
1178  overlap%ybegin_c = 0
1179  overlap%yend_c = -1
1180  endif
1181 
1182  if(nest_domain%is_fine_pe) then
1183  nsend = 0
1184  allocate(overlaplist(npes_coarse))
1185  do nn = 1, nest_domain%num_nest
1186  tile_fine = nest_domain%tile_fine(nn)
1187  tile_coarse = nest_domain%tile_coarse(nn)
1188  istart_fine = nest_domain%istart_fine(nn)
1189  iend_fine = nest_domain%iend_fine(nn)
1190  jstart_fine = nest_domain%jstart_fine(nn)
1191  jend_fine = nest_domain%jend_fine(nn)
1192  istart_coarse = nest_domain%istart_coarse(nn)
1193  iend_coarse = nest_domain%iend_coarse(nn)
1194  jstart_coarse = nest_domain%jstart_coarse(nn)
1195  jend_coarse = nest_domain%jend_coarse(nn)
1196  x_refine = nest_domain%x_refine
1197  y_refine = nest_domain%y_refine
1198 
1199  !--- set up the data range for fine and coarse grid.
1200  !--- on coarse grid pelist, xbegin_f, ybegin_f, xend_f, yend_f is dummy value
1201  !--- on fine grid pelist, xbegin_c, xend_c, ybegin_c, yend_c are the coarse grid index that
1202  !--- the fine grid overlapped with.
1203  !--- One coarse grid box might overlap with multiple fine grid processor. We use
1204  !--- the west/south/southwest processor to store the coarse grid data.
1205  if(tile_fine .NE. domain_fine%tile_id(1)) cycle
1206  is_coarse = istart_coarse + (isc_fine-istart_fine)/x_refine
1207  ie_coarse = istart_coarse + (iec_fine-istart_fine)/x_refine
1208  if(mod(isc_fine-istart_fine, x_refine) .NE. 0 ) is_coarse = is_coarse + 1
1209  js_coarse = jstart_coarse + (jsc_fine-jstart_fine)/y_refine
1210  je_coarse = jstart_coarse + (jec_fine-jstart_fine)/y_refine
1211  if(mod(jsc_fine-jstart_fine, y_refine) .NE. 0 ) js_coarse = js_coarse + 1
1212  overlap%xbegin_c = is_coarse
1213  overlap%xend_c = ie_coarse
1214  overlap%ybegin_c = js_coarse
1215  overlap%yend_c = je_coarse
1216  overlap%xbegin_f = istart_fine + (overlap%xbegin_c-istart_coarse)*x_refine
1217  overlap%xend_f = istart_fine + (overlap%xend_c-istart_coarse+1)*x_refine - 1
1218  overlap%ybegin_f = jstart_fine + (overlap%ybegin_c-jstart_coarse)*y_refine
1219  overlap%yend_f = jstart_fine + (overlap%yend_c-jstart_coarse+1)*y_refine - 1
1220  xbegin_c = overlap%xbegin_c
1221  xend_c = overlap%xend_c
1222  ybegin_c = overlap%ybegin_c
1223  yend_c = overlap%yend_c
1224  ! if(iec_fine == ieg_fine) then
1225  overlap%xend_c = overlap%xend_c + ishift
1226  overlap%xend_f = overlap%xend_f + ishift
1227  ! endif
1228  ! if(jec_fine == jeg_fine) then
1229  overlap%yend_c = overlap%yend_c + jshift
1230  overlap%yend_f = overlap%yend_f + jshift
1231  ! endif
1232 
1233  overlap%xsize_c = overlap%xend_c - overlap%xbegin_c + 1
1234  overlap%ysize_c = overlap%yend_c - overlap%ybegin_c + 1
1235 
1236  !-----------------------------------------------------------------------------------------
1237  !
1238  ! Sending From fine to coarse.
1239  ! compute the send information from fine grid to coarse grid. This will only need to send
1240  ! the internal of fine grid to coarse grid.
1241  !-----------------------------------------------------------------------------------------
1242  do n = 1, npes_coarse
1243  nconvert = convert_index_to_nest(domain_coarse, ishift, jshift, tile_coarse, istart_coarse, iend_coarse, &
1244  jstart_coarse, jend_coarse, domain_coarse%ntiles, domain_coarse%list(n-1)%tile_id(1), &
1245  isl_coarse(n), iel_coarse(n), jsl_coarse(n), jel_coarse(n), &
1246  is_convert2, ie_convert2, js_convert2, je_convert2, rotate2)
1247  is2 = xbegin_c; ie2 = xend_c
1248  js2 = ybegin_c; je2 = yend_c
1249  is3 = is2; js3 = js2
1250  do l = 1, nconvert
1251  if(rotate2(l) == ninety .OR. rotate2(l) == minus_ninety) then
1252  ie3 = ie2 + jshift
1253  je3 = je2 + ishift
1254  else
1255  ie3 = ie2 + ishift
1256  je3 = je2 + jshift
1257  endif
1258  is_coarse = max( is3, is_convert2(l) )
1259  ie_coarse = min( ie3, ie_convert2(l) )
1260  js_coarse = max( js3, js_convert2(l) )
1261  je_coarse = min( je3, je_convert2(l) )
1262  if(ie_coarse .GE. is_coarse .AND. je_coarse .GE. js_coarse ) then
1263  dir = 0
1264  nsend = nsend + 1
1265  call allocate_nest_overlap(overlaplist(nsend), maxoverlap)
1266  call insert_nest_overlap(overlaplist(nsend), nest_domain%pelist_coarse(n), &
1267  is_coarse, ie_coarse, js_coarse, je_coarse, dir, rotate2(l))
1268  endif
1269  enddo
1270  enddo
1271  enddo
1272  overlap%nsend = nsend
1273  if(nsend > 0) then
1274  if (associated(overlap%send)) deallocate(overlap%send) !< Check if allocated
1275  allocate(overlap%send(nsend))
1276  do n = 1, nsend
1277  call copy_nest_overlap(overlap%send(n), overlaplist(n) )
1278 ! call print_nest_overlap(overlap%send(n), "SEND")
1279  call deallocate_nest_overlap(overlaplist(n))
1280  enddo
1281  endif
1282  if(allocated(overlaplist))deallocate(overlaplist)
1283  endif
1284  !--------------------------------------------------------------------------------
1285  ! compute the recv information from fine grid to coarse grid. This will only need to send
1286  ! the internal of fine grid to coarse grid.
1287  !--------------------------------------------------------------------------------
1288 
1289  if( nest_domain%is_coarse_pe ) then
1290  nrecv = 0
1291  allocate(overlaplist(npes_fine))
1292  do nn = 1, nest_domain%num_nest
1293  tile_fine = nest_domain%tile_fine(nn)
1294  tile_coarse = nest_domain%tile_coarse(nn)
1295  istart_fine = nest_domain%istart_fine(nn)
1296  iend_fine = nest_domain%iend_fine(nn)
1297  jstart_fine = nest_domain%jstart_fine(nn)
1298  jend_fine = nest_domain%jend_fine(nn)
1299  istart_coarse = nest_domain%istart_coarse(nn)
1300  iend_coarse = nest_domain%iend_coarse(nn)
1301  jstart_coarse = nest_domain%jstart_coarse(nn)
1302  jend_coarse = nest_domain%jend_coarse(nn)
1303  x_refine = nest_domain%x_refine
1304  y_refine = nest_domain%y_refine
1305 
1306  dir = 0
1307  do n = 1, npes_fine
1308  if(tile_fine .NE. domain_fine%list(n-1)%tile_id(1)) cycle
1309  is_you = istart_coarse + (isl_fine(n)-istart_fine)/x_refine
1310  ie_you = istart_coarse + (iel_fine(n)-istart_fine)/x_refine
1311  if(mod(isl_fine(n)-istart_fine, x_refine) .NE. 0 ) is_you = is_you + 1
1312  js_you = jstart_coarse + (jsl_fine(n)-jstart_fine)/y_refine
1313  je_you = jstart_coarse + (jel_fine(n)-jstart_fine)/y_refine
1314  if(mod(jsl_fine(n)-jstart_fine, y_refine) .NE. 0 ) js_you = js_you + 1
1315  nconvert=convert_index_to_coarse(domain_coarse, ishift, jshift, tile_coarse, istart_coarse, iend_coarse, &
1316  & jstart_coarse, jend_coarse, domain_coarse%ntiles, domain_coarse%tile_id(1), is_you, ie_you, &
1317  & js_you, je_you, is_convert2, ie_convert2, js_convert2, je_convert2, rotate2)
1318  do l = 1, nconvert
1319  is2 = max(is_convert2(l), isc_coarse)
1320  ie2 = min(ie_convert2(l), iec_coarse+ishift)
1321  js2 = max(js_convert2(l), jsc_coarse)
1322  je2 = min(je_convert2(l), jec_coarse+jshift)
1323 
1324  if( ie2 .GE. is2 .AND. je2 .GE. js2 ) then
1325  nrecv = nrecv + 1
1326  call allocate_nest_overlap(overlaplist(nrecv), maxoverlap)
1327  call insert_nest_overlap(overlaplist(nrecv), nest_domain%pelist_fine(n), &
1328  is2, ie2, js2, je2, dir, rotate2(l))
1329  endif
1330  enddo
1331  enddo
1332  enddo
1333  overlap%nrecv = nrecv
1334  if(nrecv > 0) then
1335  allocate(overlap%recv(nrecv))
1336  do n = 1, nrecv
1337  call copy_nest_overlap(overlap%recv(n), overlaplist(n) )
1338 ! call print_nest_overlap(overlap%recv(n), "RECV")
1339  call deallocate_nest_overlap( overlaplist(n) )
1340  enddo
1341  endif
1342  if(allocated(overlaplist))deallocate(overlaplist)
1343 
1344  endif
1345 
1346  if(debug_message_passing) call debug_message_size(overlap, name)
1347 
1348  deallocate(isl_coarse, iel_coarse, jsl_coarse, jel_coarse)
1349  deallocate(isl_fine, iel_fine, jsl_fine, jel_fine)
1350 
1351 end subroutine compute_overlap_fine_to_coarse
1352 
1353 function find_index(array, index_data, start_pos)
1354  integer, intent(in) :: array(:)
1355  integer, intent(in) :: index_data
1356  integer, intent(in) :: start_pos
1357  integer :: find_index
1358  integer :: i
1359 
1360  find_index = 0
1361  do i = start_pos, size(array)
1362  if(array(i) == index_data) then
1363  find_index = i
1364  exit
1365  endif
1366  enddo
1367  if(find_index == 0) then
1368  print*, "start_pos = ", start_pos, index_data, array
1369  call mpp_error(fatal, "mpp_define_nest_domains.inc: can not find data in array")
1370  endif
1371 
1372 end function find_index
1373 
1374 subroutine debug_message_size(overlap, name)
1375  type(nestspec), intent(in) :: overlap
1376  character(len=*), intent(in) :: name
1377  integer, allocatable :: msg1(:), msg2(:), msg3(:), pelist(:)
1378  integer :: m, n, l, npes, msgsize
1379  integer :: is, ie, js, je, outunit
1380 
1381  outunit = stdout()
1382  npes = mpp_npes()
1383  allocate(msg1(npes), msg2(npes), msg3(npes) )
1384  allocate(pelist(npes))
1385  call mpp_get_current_pelist(pelist)
1386  msg1 = 0
1387  msg2 = 0
1388  msg3 = 0
1389  l = 0
1390  do m = 1, overlap%nrecv
1391  msgsize = 0
1392  do n = 1, overlap%recv(m)%count
1393  is = overlap%recv(m)%is(n); ie = overlap%recv(m)%ie(n)
1394  js = overlap%recv(m)%js(n); je = overlap%recv(m)%je(n)
1395  msgsize = msgsize + (ie-is+1)*(je-js+1)
1396  end do
1397  l = find_index(pelist, overlap%recv(m)%pe, l+1)
1398  msg2(l) = msgsize
1399  enddo
1400  l = 0
1401  do m = 1, overlap%nsend
1402  msgsize = 0
1403  do n = 1, overlap%send(m)%count
1404  is = overlap%send(m)%is(n); ie = overlap%send(m)%ie(n)
1405  js = overlap%send(m)%js(n); je = overlap%send(m)%je(n)
1406  msgsize = msgsize + (ie-is+1)*(je-js+1)
1407  end do
1408  l = find_index(pelist, overlap%send(m)%pe, l+1)
1409  msg3(l) = msgsize
1410  enddo
1411 
1412  call mpp_alltoall(msg3, 1, msg1, 1)
1413 
1414  do m = 1, npes
1415  if(msg1(m) .NE. msg2(m)) then
1416  print*, "debug_message_size: My pe = ", mpp_pe(), ",name =", trim(name),", from pe=", &
1417  pelist(m), ":send size = ", msg1(m), ", recv size = ", msg2(m)
1418  call mpp_error(fatal, "debug_message_size: mismatch on send and recv size")
1419  endif
1420  enddo
1421  write(outunit,*)"NOTE from compute_overlap_fine_to_coarse: "// &
1422  "message sizes are matched between send and recv for "//trim(name)
1423  deallocate(msg1, msg2, msg3, pelist)
1424 
1425 end subroutine debug_message_size
1426 
1427 !###############################################################################
1428 
1429 subroutine init_index_type (indexData )
1430  type(index_type), intent(inout) :: indexData
1431 
1432  indexdata%is_me = 0
1433  indexdata%ie_me = -1
1434  indexdata%js_me = 0
1435  indexdata%je_me = -1
1436  indexdata%is_you = 0
1437  indexdata%ie_you = -1
1438  indexdata%js_you = 0
1439  indexdata%je_you = -1
1440 
1441 end subroutine init_index_type
1442 
1443 subroutine allocate_nest_overlap(overlap, count)
1444  type(overlap_type), intent(inout) :: overlap
1445  integer, intent(in ) :: count
1446 
1447  overlap%count = 0
1448  overlap%pe = null_pe
1449  if( ASSOCIATED(overlap%is) ) call mpp_error(fatal, &
1450  "mpp_define_nest_domains.inc: overlap is already been allocated")
1451 
1452  allocate(overlap%is (count) )
1453  allocate(overlap%ie (count) )
1454  allocate(overlap%js (count) )
1455  allocate(overlap%je (count) )
1456  allocate(overlap%dir (count) )
1457  allocate(overlap%rotation (count) )
1458  allocate(overlap%msgsize (count) )
1459 
1460 end subroutine allocate_nest_overlap
1461 
1462 !##############################################################################
1463 subroutine deallocate_nest_overlap(overlap)
1464  type(overlap_type), intent(inout) :: overlap
1465 
1466  overlap%count = 0
1467  overlap%pe = null_pe
1468  deallocate(overlap%is)
1469  deallocate(overlap%ie)
1470  deallocate(overlap%js)
1471  deallocate(overlap%je)
1472  deallocate(overlap%dir)
1473  deallocate(overlap%rotation)
1474  deallocate(overlap%msgsize)
1475 
1476 end subroutine deallocate_nest_overlap
1477 
1478 !##############################################################################
1479 subroutine insert_nest_overlap(overlap, pe, is, ie, js, je, dir, rotation)
1480  type(overlap_type), intent(inout) :: overlap
1481  integer, intent(in ) :: pe
1482  integer, intent(in ) :: is, ie, js, je
1483  integer, intent(in ) :: dir, rotation
1484  integer :: count
1485 
1486  if( overlap%count == 0 ) then
1487  overlap%pe = pe
1488  else
1489  if(overlap%pe .NE. pe) call mpp_error(fatal, &
1490  "mpp_define_nest_domains.inc: mismatch on pe")
1491  endif
1492  overlap%count = overlap%count+1
1493  count = overlap%count
1494  if(count > size(overlap%is(:))) call mpp_error(fatal, &
1495  "mpp_define_nest_domains.inc: overlap%count > size(overlap%is), contact developer")
1496  overlap%is (count) = is
1497  overlap%ie (count) = ie
1498  overlap%js (count) = js
1499  overlap%je (count) = je
1500  overlap%dir (count) = dir
1501  overlap%rotation (count) = rotation
1502  overlap%msgsize (count) = (ie-is+1)*(je-js+1)
1503 
1504 end subroutine insert_nest_overlap
1505 
1506 subroutine print_nest_overlap(overlap, msg)
1507  type(overlap_type), intent(in) :: overlap
1508  character(len=*), intent(in) :: msg
1509 
1510  integer :: i
1511  write(1000+mpp_pe(),*) trim(msg), ",pe=",overlap%pe, overlap%count
1512  do i = 1, overlap%count
1513  write(1000+mpp_pe(),*) trim(msg), ",index=", overlap%is(i), overlap%ie(i),overlap%js(i),overlap%je(i)
1514  write(1000+mpp_pe(),*) trim(msg), ",rotation=", overlap%dir(i), overlap%rotation(i), overlap%msgsize(i)
1515  enddo
1516  flush(1000+mpp_pe())
1517 
1518 end subroutine print_nest_overlap
1519 
1520 !#########################################################
1521 subroutine copy_nest_overlap(overlap_out, overlap_in)
1522  type(overlap_type), intent(inout) :: overlap_out
1523  type(overlap_type), intent(in) :: overlap_in
1524 
1525  if(overlap_in%count == 0) call mpp_error(fatal, &
1526  "mpp_define_nest_domains.inc: overlap_in%count is 0")
1527 
1528  if(associated(overlap_out%is)) call mpp_error(fatal, &
1529  "mpp_define_nest_domains.inc: overlap_out is already been allocated")
1530 
1531  call allocate_nest_overlap(overlap_out, overlap_in%count)
1532  overlap_out%count = overlap_in%count
1533  overlap_out%pe = overlap_in%pe
1534 
1535  overlap_out%is(:) = overlap_in%is(1:overlap_in%count)
1536  overlap_out%ie(:) = overlap_in%ie(1:overlap_in%count)
1537  overlap_out%js(:) = overlap_in%js(1:overlap_in%count)
1538  overlap_out%je(:) = overlap_in%je(1:overlap_in%count)
1539  overlap_out%is(:) = overlap_in%is(1:overlap_in%count)
1540  overlap_out%dir(:) = overlap_in%dir(1:overlap_in%count)
1541  overlap_out%rotation(:) = overlap_in%rotation(1:overlap_in%count)
1542  overlap_out%msgsize(:) = overlap_in%msgsize(1:overlap_in%count)
1543 
1544 
1545 end subroutine copy_nest_overlap
1546 
1547 
1548 !#######################################################################
1549  ! this routine found the domain has the same halo size with the input
1550  ! whalo, ehalo,
1551 function search_c2f_nest_overlap(nest_domain, nest_level, extra_halo, position)
1552  type(nest_domain_type), intent(inout) :: nest_domain
1553  integer, intent(in) :: extra_halo
1554  integer, intent(in) :: position, nest_level
1555  type(nestspec), pointer :: search_C2F_nest_overlap
1556  type(nestspec), pointer :: update_ref
1557  character(len=128) :: name
1558 
1559  if(nest_level < 1 .OR. nest_level > nest_domain%num_level) call mpp_error(fatal, &
1560  "mpp_define_nest_domains.inc(search_C2F_nest_overlap): nest_level should be between 1 and nest_domain%num_level")
1561 
1562  select case(position)
1563  case (center)
1564  name = trim(nest_domain%name)//" T-cell"
1565  update_ref => nest_domain%nest(nest_level)%C2F_T
1566  case (corner)
1567  update_ref => nest_domain%nest(nest_level)%C2F_C
1568  case (north)
1569  update_ref => nest_domain%nest(nest_level)%C2F_N
1570  case (east)
1571  update_ref => nest_domain%nest(nest_level)%C2F_E
1572  case default
1573  call mpp_error(fatal, &
1574  & "mpp_define_nest_domains.inc(search_C2F_nest_overlap): position should be CENTER|CORNER|EAST|NORTH")
1575  end select
1576 
1577  search_c2f_nest_overlap => update_ref
1578 
1579  do
1580  if(extra_halo == search_c2f_nest_overlap%extra_halo) then
1581  exit ! found domain
1582  endif
1583  !--- if not found, switch to next
1584  if(.NOT. ASSOCIATED(search_c2f_nest_overlap%next)) then
1585  allocate(search_c2f_nest_overlap%next)
1586  search_c2f_nest_overlap => search_c2f_nest_overlap%next
1587  call compute_overlap_coarse_to_fine(nest_domain%nest(nest_level), search_c2f_nest_overlap, &
1588  & extra_halo, position, name)
1589  exit
1590  else
1591  search_c2f_nest_overlap => search_c2f_nest_overlap%next
1592  end if
1593 
1594  end do
1595 
1596  update_ref => null()
1597 
1598  end function search_c2f_nest_overlap
1599 
1600 !#######################################################################
1601  ! this routine found the domain has the same halo size with the input
1602  ! whalo, ehalo,
1603  function search_f2c_nest_overlap(nest_domain, nest_level, position)
1604  type(nest_domain_type), intent(inout) :: nest_domain
1605  integer, intent(in) :: position, nest_level
1606  type(nestspec), pointer :: search_F2C_nest_overlap
1607 
1608  if(nest_level < 1 .OR. nest_level > nest_domain%num_level) call mpp_error(fatal, &
1609  "mpp_define_nest_domains.inc(search_F2C_nest_overlap): nest_level should be between 1 and nest_domain%num_level")
1610 
1611  select case(position)
1612  case (center)
1613  search_f2c_nest_overlap => nest_domain%nest(nest_level)%F2C_T
1614  case (corner)
1615  search_f2c_nest_overlap => nest_domain%nest(nest_level)%F2C_C
1616  case (north)
1617  search_f2c_nest_overlap => nest_domain%nest(nest_level)%F2C_N
1618  case (east)
1619  search_f2c_nest_overlap => nest_domain%nest(nest_level)%F2C_E
1620  case default
1621  call mpp_error(fatal, &
1622  & "mpp_define_nest_domains.inc(search_F2C_nest_overlap): position should be CENTER|CORNER|EAST|NORTH")
1623  end select
1624 
1625  end function search_f2c_nest_overlap
1626 
1627  !################################################################
1628  !> @brief Get the index of the data passed from coarse grid to fine grid.
1629  !!
1630  !> Get the index of the data passed from coarse grid to fine grid.
1631  !!
1632  !! <br>Example usage:
1633  !! @code{.F90}
1634  !! call mpp_get_C2F_index(nest_domain, is_fine, ie_fine, js_fine, je_fine,
1635  !! is_coarse, ie_coarse, js_coarse, je_coarse, dir,
1636  !! nest_level, position)
1637  !! @endcode
1638  subroutine mpp_get_c2f_index(nest_domain, is_fine, ie_fine, js_fine, je_fine, &
1639  is_coarse, ie_coarse, js_coarse, je_coarse, dir, nest_level, position)
1640 
1641  type(nest_domain_type), intent(in ) :: nest_domain !< holds the information to pass data
1642  !! between fine and coarse grids
1643  integer, intent(out) :: is_fine, ie_fine, js_fine, je_fine !< index in the fine
1644  !! grid of the nested region
1645  integer, intent(out) :: is_coarse, ie_coarse, js_coarse, je_coarse !< index in the coarse
1646  !! grid of the nested region
1647  integer, intent(in ) :: dir, nest_level !< direction of the halo update.
1648  !! Its value should be WEST, EAST, SOUTH or NORTH.;
1649  !! level of the nest (> 1 implies a telescoping nest)
1650  integer, optional, intent(in ) :: position !< Cell position. It value should be CENTER,
1651  !! EAST, CORNER, or NORTH.
1652 
1653  integer :: update_position
1654  type(nestspec), pointer :: update => null()
1655 
1656  if(nest_level < 1 .OR. nest_level > nest_domain%num_level) call mpp_error(fatal, &
1657  "mpp_define_nest_domains.inc(mpp_get_C2F_index): nest_level should be between 1 and nest_domain%num_level")
1658 
1659  update_position = center
1660  if(present(position)) update_position = position
1661 
1662  select case(update_position)
1663  case (center)
1664  update => nest_domain%nest(nest_level)%C2F_T
1665  case (east)
1666  update => nest_domain%nest(nest_level)%C2F_E
1667  case (corner)
1668  update => nest_domain%nest(nest_level)%C2F_C
1669  case (north)
1670  update => nest_domain%nest(nest_level)%C2F_N
1671  case default
1672  call mpp_error(fatal, "mpp_define_nest_domains.inc(mpp_get_C2F_index): invalid option argument position")
1673  end select
1674 
1675  select case(dir)
1676  case(west)
1677  is_fine = update%west%is_me
1678  ie_fine = update%west%ie_me
1679  js_fine = update%west%js_me
1680  je_fine = update%west%je_me
1681  is_coarse = update%west%is_you
1682  ie_coarse = update%west%ie_you
1683  js_coarse = update%west%js_you
1684  je_coarse = update%west%je_you
1685  case(east)
1686  is_fine = update%east%is_me
1687  ie_fine = update%east%ie_me
1688  js_fine = update%east%js_me
1689  je_fine = update%east%je_me
1690  is_coarse = update%east%is_you
1691  ie_coarse = update%east%ie_you
1692  js_coarse = update%east%js_you
1693  je_coarse = update%east%je_you
1694  case(south)
1695  is_fine = update%south%is_me
1696  ie_fine = update%south%ie_me
1697  js_fine = update%south%js_me
1698  je_fine = update%south%je_me
1699  is_coarse = update%south%is_you
1700  ie_coarse = update%south%ie_you
1701  js_coarse = update%south%js_you
1702  je_coarse = update%south%je_you
1703  case(north)
1704  is_fine = update%north%is_me
1705  ie_fine = update%north%ie_me
1706  js_fine = update%north%js_me
1707  je_fine = update%north%je_me
1708  is_coarse = update%north%is_you
1709  ie_coarse = update%north%ie_you
1710  js_coarse = update%north%js_you
1711  je_coarse = update%north%je_you
1712  case default
1713  call mpp_error(fatal, "mpp_define_nest_domains.inc: invalid value for argument dir")
1714  end select
1715 
1716 
1717  end subroutine mpp_get_c2f_index
1718 
1719  subroutine mpp_get_f2c_index_fine(nest_domain, is_coarse, ie_coarse, js_coarse, je_coarse, &
1720  is_fine, ie_fine, js_fine, je_fine, nest_level, position)
1721 
1722  type(nest_domain_type), intent(in ) :: nest_domain !< Holds the information to pass data
1723  !! between fine and coarse grid.
1724  integer, intent(out) :: is_fine, ie_fine, js_fine, je_fine !< index in the fine
1725  !! grid of the nested region
1726  integer, intent(out) :: is_coarse, ie_coarse, js_coarse, je_coarse !< index in
1727  !! the coarse grid of the nested region
1728  integer, intent(in) :: nest_level !< level of the nest (> 1 implies a telescoping nest)
1729  integer, optional, intent(in ) :: position !< Cell position. It value should be CENTER,
1730  !! EAST, CORNER, or NORTH.
1731 
1732  integer :: update_position
1733  type(nestspec), pointer :: update => null()
1734 
1735  if(nest_level < 1 .OR. nest_level > nest_domain%num_level) call mpp_error(fatal, &
1736  "mpp_define_nest_domains.inc(mpp_get_F2C_index): nest_level should be between 1 and nest_domain%num_level")
1737 
1738  if(.not. nest_domain%nest(nest_level)%on_level) call mpp_error(fatal, &
1739  "mpp_define_nest_domains.inc(mpp_get_F2C_index_fine): nest_domain%nest(nest_level)%on_level is false")
1740 
1741  update_position = center
1742  if(present(position)) update_position = position
1743 
1744  select case(update_position)
1745  case (center)
1746  update => nest_domain%nest(nest_level)%F2C_T
1747  case (east)
1748  update => nest_domain%nest(nest_level)%F2C_E
1749  case (corner)
1750  update => nest_domain%nest(nest_level)%F2C_C
1751  case (north)
1752  update => nest_domain%nest(nest_level)%F2C_N
1753  case default
1754  call mpp_error(fatal, "mpp_define_nest_domains.inc(mpp_get_F2C_index): invalid option argument position")
1755  end select
1756  is_fine = update%xbegin_f
1757  ie_fine = update%xend_f
1758  js_fine = update%ybegin_f
1759  je_fine = update%yend_f
1760  is_coarse = update%xbegin_c
1761  ie_coarse = update%xend_c
1762  js_coarse = update%ybegin_c
1763  je_coarse = update%yend_c
1764 
1765  end subroutine mpp_get_f2c_index_fine
1766 
1767  !################################################################
1768  subroutine mpp_get_f2c_index_coarse(nest_domain, is_coarse, ie_coarse, js_coarse, je_coarse, nest_level, position)
1769 
1770  type(nest_domain_type), intent(in ) :: nest_domain !< Holds the information to pass data
1771  !! between fine and coarse grid.
1772  integer, intent(out) :: is_coarse, ie_coarse, js_coarse, je_coarse !< index in the fine
1773  !! grid of the nested region
1774  integer, intent(in ) :: nest_level !< level of the nest (> 1 implies a telescoping nest)
1775  integer, optional, intent(in ) :: position !< Cell position. It value should be CENTER,
1776  !! EAST, CORNER, or NORTH.
1777 
1778  integer :: update_position
1779  type(nestspec), pointer :: update => null()
1780 
1781  if(nest_level < 1 .OR. nest_level > nest_domain%num_level) call mpp_error(fatal, &
1782  & "mpp_define_nest_domains.inc(mpp_get_F2C_index_coarse):"// &
1783  & "nest_level should be between 1 and nest_domain%num_level")
1784 
1785 
1786  update_position = center
1787  if(present(position)) update_position = position
1788 
1789  select case(update_position)
1790  case (center)
1791  update => nest_domain%nest(nest_level)%F2C_T
1792  case (east)
1793  update => nest_domain%nest(nest_level)%F2C_E
1794  case (corner)
1795  update => nest_domain%nest(nest_level)%F2C_C
1796  case (north)
1797  update => nest_domain%nest(nest_level)%F2C_N
1798  case default
1799  call mpp_error(fatal, &
1800  & "mpp_define_nest_domains.inc(mpp_get_F2C_index_coarse): invalid option argument position")
1801  end select
1802  is_coarse = update%xbegin_c
1803  ie_coarse = update%xend_c
1804  js_coarse = update%ybegin_c
1805  je_coarse = update%yend_c
1806 
1807  end subroutine mpp_get_f2c_index_coarse
1808 
1809  subroutine get_coarse_index(rotate, is, ie, js, je, iadd, jadd, is_c, ie_c, js_c, je_c)
1810  integer, intent(in) :: rotate, is, ie, js, je, iadd, jadd
1811  integer, intent(out) :: is_c, ie_c, js_c, je_c
1812 
1813  if(rotate == 0) then
1814  is_c = is; ie_c = ie
1815  js_c = js; je_c = je
1816  else
1817  is_c = js; ie_c = je
1818  js_c = is; je_c = ie
1819  endif
1820  is_c = is_c + iadd; ie_c = ie_c + iadd
1821  js_c = js_c + jadd; je_c = je_c + jadd
1822 
1823  end subroutine get_coarse_index
1824 
1825  !--- this routine will get number of nest.
1826  subroutine get_nnest(domain, num_nest, tile_coarse, istart_coarse, iend_coarse, jstart_coarse, jend_coarse, &
1827  x_refine, y_refine, nnest, t_coarse, ncross_coarse, rotate_coarse, &
1828  is_coarse, ie_coarse, js_coarse, je_coarse, is_fine, ie_fine, js_fine, je_fine)
1829  type(domain2d), intent(in) :: domain
1830  integer, intent(in) :: num_nest, istart_coarse(:), iend_coarse(:), jstart_coarse(:), jend_coarse(:)
1831  integer, intent(in) :: tile_coarse(:)
1832  integer, intent(in) :: x_refine, y_refine
1833  integer, intent(out) :: nnest, is_coarse(:), ie_coarse(:), js_coarse(:), je_coarse(:)
1834  integer, intent(out) :: is_fine(:), ie_fine(:), js_fine(:), je_fine(:)
1835  integer, intent(out) :: t_coarse(:), ncross_coarse(:), rotate_coarse(:)
1836  integer :: is, ie, js, je, tile, isg, ieg, jsg, jeg
1837  integer :: ncross, rotate, i1, i2
1838  integer :: is_c, ie_c, js_c, je_c
1839  integer :: n, iadd, jadd
1840 
1841 
1842  call mpp_get_global_domain(domain, isg, ieg, jsg, jeg)
1843  nnest = 0
1844  do n = 1, num_nest
1845  is = istart_coarse(n); ie = iend_coarse(n)
1846  js = jstart_coarse(n); je = jend_coarse(n)
1847  tile = tile_coarse(n)
1848  iadd = 0 ; jadd = 0
1849  ncross = 0
1850  rotate = 0
1851  do while ( ie .GE. is .AND. je .GE. js)
1852  nnest = nnest+1
1853  t_coarse(nnest) = tile
1854  ncross_coarse(nnest) = ncross
1855  rotate_coarse(nnest) = rotate
1856  !--- rotate should be 0, 90 or -90.
1857  if(rotate .NE. 0 .AND. rotate .NE. 90 .AND. rotate .NE. -90) then
1858  call mpp_error(fatal, "get_nnest: roate should be 0, 90 or -90")
1859  endif
1860  if( ieg .GE. ie .AND. jeg .GE. je) then
1861  is_coarse(nnest) = is; ie_coarse(nnest) = ie
1862  js_coarse(nnest) = js; je_coarse(nnest) = je
1863  call get_coarse_index(rotate, is_coarse(nnest), ie_coarse(nnest), js_coarse(nnest), je_coarse(nnest), &
1864  iadd, jadd, is_c, ie_c, js_c, je_c)
1865  is = ie + 1; js = je + 1
1866  else if( ieg .GE. ie ) then ! jeg < je, will cross the north edge
1867  is_coarse(nnest) = is; ie_coarse(nnest) = ie
1868  js_coarse(nnest) = js; je_coarse(nnest) = jeg
1869  call get_coarse_index(rotate, is_coarse(nnest), ie_coarse(nnest), js_coarse(nnest), je_coarse(nnest), &
1870  iadd, jadd, is_c, ie_c, js_c, je_c)
1871  if(rotate ==0) then
1872  jadd = jadd + jeg
1873  else
1874  iadd = iadd + ieg
1875  endif
1876  js = 1; je = je-jeg
1877  ncross = ncross+1
1878  if(mod(tile,2) ==0) then ! tile 2 4 6
1879  tile = tile + 1
1880  if(tile>6) tile=tile-6
1881  else ! rotate 90 degree
1882  tile = tile + 2
1883  if(tile>6) tile=tile-6
1884  i1 = is; i2 = ie
1885  is = js; ie = je
1886  js = i1; je = i2
1887  rotate = rotate + 90
1888  endif
1889 
1890 
1891  else if( jeg .GE. je ) then ! ieg < ie, will cross the east edge
1892  is_coarse(nnest) = is; ie_coarse(nnest) = ieg
1893  js_coarse(nnest) = js; je_coarse(nnest) = je
1894  call get_coarse_index(rotate, is_coarse(nnest), ie_coarse(nnest), js_coarse(nnest), je_coarse(nnest), &
1895  iadd, jadd, is_c, ie_c, js_c, je_c)
1896  if(rotate ==0) then
1897  iadd = iadd + ieg
1898  else
1899  jadd = jadd + jeg
1900  endif
1901  is = 1; ie = ie-ieg
1902  ncross = ncross+1
1903  if(mod(tile,2) ==0) then ! rotate -90
1904  tile = tile + 2
1905  if(tile>6) tile=tile-6
1906  i1 = is; i2 = ie
1907  is = js; ie = je
1908  js = i1; je = i2
1909  rotate = rotate - 90
1910  else
1911  tile = tile + 1
1912  if(tile>6) tile=tile-6
1913  endif
1914  else
1915  call mpp_error(fatal, "get_nnest: do not support cross the corner")
1916  endif
1917 
1918  !--- is_c:ie_c,js_c:je_c must be inside istart_coarse(n):iend_coarse(n), jstart_coarse(n):jend_coarse(n)
1919  if(is_c < istart_coarse(n)) call mpp_error(fatal, "get_nnest: is_c < istart_coarse")
1920  if(ie_c > iend_coarse(n)) call mpp_error(fatal, "get_nnest: ie_c > iend_coarse")
1921  if(js_c < jstart_coarse(n)) call mpp_error(fatal, "get_nnest: js_c < jstart_coarse")
1922  if(je_c > jend_coarse(n)) call mpp_error(fatal, "get_nnest: je_c > jend_coarse")
1923  is_fine(nnest) = (is_c - istart_coarse(n)) * x_refine + 1
1924  ie_fine(nnest) = (ie_c - istart_coarse(n)+1) * x_refine
1925  js_fine(nnest) = (js_c - jstart_coarse(n)) * y_refine + 1
1926  je_fine(nnest) = (je_c - jstart_coarse(n)+1) * y_refine
1927 
1928  !--- it should not cross the edge more than 3 times.
1929  if(ncross > 3) call mpp_error(fatal, "get_nnest: nncross > 3")
1930  enddo
1931  enddo
1932 
1933 
1934  end subroutine get_nnest
1935 
1936 
1937  !> This routine will convert the global coarse grid index to nest grid index.
1938  function convert_index_to_nest(domain, ishift, jshift, tile_coarse, istart_coarse, iend_coarse, jstart_coarse, &
1939  & jend_coarse, ntiles_coarse, tile_in, is_in, ie_in, js_in, je_in, is_out, ie_out,&
1940  & js_out, je_out, rotate_out)
1941  type(domain2d), intent(in) :: domain
1942  integer, intent(in) :: ishift, jshift
1943  integer, intent(in) :: istart_coarse, iend_coarse, jstart_coarse, jend_coarse
1944  integer, intent(in) :: tile_coarse
1945  integer, intent(in) :: ntiles_coarse, tile_in, is_in, ie_in, js_in, je_in
1946  integer, intent(out) :: is_out(:), ie_out(:), js_out(:), je_out(:), rotate_out(:)
1947  integer :: convert_index_to_nest
1948  integer :: is, ie, js, je, tile, isg, ieg, jsg, jeg
1949  integer :: ncross, rotate, nout, diff, ntiles
1950 
1951  ntiles = ntiles_coarse
1952  call mpp_get_global_domain(domain, isg, ieg, jsg, jeg)
1953  is = istart_coarse; ie = iend_coarse
1954  js = jstart_coarse; je = jend_coarse
1955  tile = tile_coarse
1956 
1957  if(size(is_out(:)) < 2) call mpp_error(fatal,"convert_index_to_nest: size(is_out(:)) < 2")
1958  if(size(ie_out(:)) < 2) call mpp_error(fatal,"convert_index_to_nest: size(ie_out(:)) < 2")
1959  if(size(js_out(:)) < 2) call mpp_error(fatal,"convert_index_to_nest: size(js_out(:)) < 2")
1960  if(size(je_out(:)) < 2) call mpp_error(fatal,"convert_index_to_nest: size(je_out(:)) < 2")
1961  if(size(rotate_out(:)) < 2) call mpp_error(fatal,"convert_index_to_nest: size(rotate_out(:)) < 2")
1962  if( ie > ieg .AND. je > jeg) then
1963  call mpp_error(fatal, "convert_index_to_nest: do not support cross the corner, contact developer")
1964  endif
1965  if( is > ieg .or. js > jeg) call mpp_error(fatal,.or."convert_index_to_nest: is > ieg js > jeg")
1966 
1967 
1968  nout = 0
1969 
1970  if(tile == tile_in) then
1971  nout = nout+1
1972  rotate_out(nout) = zero
1973  is_out(nout) = is_in; ie_out(nout) = ie_in + ishift
1974  js_out(nout) = js_in; je_out(nout) = je_in + jshift
1975  endif
1976 
1977  diff = tile_in - tile
1978  if(diff < 0) diff = diff + ntiles
1979  ncross = -1
1980  if( ie > ieg ) then
1981  select case(diff)
1982  case (0)
1983  rotate = zero
1984  ncross = 4
1985  case (1)
1986  if(mod(tile,2) ==1) then ! tile 1 3 5
1987  rotate = zero
1988  ncross = 1
1989  endif
1990  case (2)
1991  if(mod(tile,2) ==0) then ! tile 2 4 6
1992  rotate = minus_ninety
1993  ncross = 1
1994  endif
1995  case (3)
1996  rotate = minus_ninety
1997  ncross = 2
1998  case (4)
1999  if(mod(tile,2) ==1) then ! tile 1 3 5
2000  rotate = minus_ninety
2001  ncross = 3
2002  endif
2003  case (5)
2004  if(mod(tile,2) ==0) then ! tile 2 4 6
2005  rotate = zero
2006  ncross = 3
2007  endif
2008  case default
2009  call mpp_error(fatal,"convert_index_to_nest: invalid value of diff")
2010  end select
2011 
2012  if(ncross > 0) then
2013  nout =nout+1
2014  rotate_out(nout) = rotate
2015  if(rotate_out(nout) == zero) then
2016  js_out(nout) = js_in
2017  je_out(nout) = je_in + jshift
2018  is_out(nout) = is_in+ncross*ieg
2019  ie_out(nout) = ie_in+ncross*ieg + ishift
2020  else if(rotate_out(nout) == minus_ninety) then
2021  js_out(nout) = ieg-ie_in + 1
2022  je_out(nout) = ieg-is_in + 1 + ishift
2023  is_out(nout) = js_in+ncross*jeg
2024  ie_out(nout) = je_in+ncross*jeg + jshift
2025  endif
2026  endif
2027  else if( je > jeg ) then
2028  select case(diff)
2029  case (0)
2030  rotate = zero
2031  ncross = 4
2032  case (1)
2033  if(mod(tile,2) ==0) then ! tile 2 4 6
2034  rotate = zero
2035  ncross = 1
2036  endif
2037  case (2)
2038  if(mod(tile,2) ==1) then ! tile 1 3 5
2039  rotate = ninety
2040  ncross = 1
2041  endif
2042  case (3)
2043  rotate = ninety
2044  ncross = 2
2045  case (4)
2046  if(mod(tile,2) ==0) then ! tile 2 4 6
2047  rotate = ninety
2048  ncross = 3
2049  endif
2050  case (5)
2051  if(mod(tile,2) ==1) then ! tile 1 3 5
2052  rotate = zero
2053  ncross = 3
2054  endif
2055  end select
2056 
2057  if(ncross > 0) then
2058  nout =nout+1
2059  rotate_out(nout) = rotate
2060 
2061  if(rotate_out(nout) == zero) then
2062  js_out(nout) = js_in
2063  je_out(nout) = je_in + jshift
2064  is_out(nout) = is_in+ncross*ieg
2065  ie_out(nout) = ie_in+ncross*ieg + ishift
2066  else if(rotate_out(nout) == ninety) then
2067  is_out(nout) = ieg-je_in + 1
2068  ie_out(nout) = ieg-js_in+1 + jshift
2069  js_out(nout) = is_in+ncross*jeg
2070  je_out(nout) = ie_in+ncross*jeg + ishift
2071  endif
2072  endif
2073  endif
2074 
2075  convert_index_to_nest = nout
2076 
2077  end function convert_index_to_nest
2078 
2079  function convert_index_to_coarse(domain, ishift, jshift, tile_coarse, istart_coarse, iend_coarse, jstart_coarse, &
2080  & jend_coarse, ntiles_coarse, tile_in, is_in, ie_in, js_in, je_in, is_out, ie_out,&
2081  & js_out, je_out, rotate_out)
2082  type(domain2d), intent(in) :: domain
2083  integer, intent(in) :: ishift, jshift
2084  integer, intent(in) :: istart_coarse, iend_coarse, jstart_coarse, jend_coarse
2085  integer, intent(in) :: tile_coarse
2086  integer, intent(in) :: ntiles_coarse, tile_in, is_in, ie_in, js_in, je_in
2087  integer, intent(out) :: is_out(:), ie_out(:), js_out(:), je_out(:), rotate_out(:)
2088  integer :: convert_index_to_coarse
2089  integer :: is, ie, js, je, isg, ieg, jsg, jeg
2090  integer :: ncross, rotate, ntiles, nout, diff, tile
2091 
2092  ntiles = ntiles_coarse
2093  call mpp_get_global_domain(domain, isg, ieg, jsg, jeg)
2094  is = istart_coarse; ie = iend_coarse
2095  js = jstart_coarse; je = jend_coarse
2096  tile = tile_coarse
2097 
2098  if(size(is_out(:)) < 2) call mpp_error(fatal,"convert_index_to_coarse: size(is_out(:)) < 2")
2099  if(size(ie_out(:)) < 2) call mpp_error(fatal,"convert_index_to_coarse: size(ie_out(:)) < 2")
2100  if(size(js_out(:)) < 2) call mpp_error(fatal,"convert_index_to_coarse: size(js_out(:)) < 2")
2101  if(size(je_out(:)) < 2) call mpp_error(fatal,"convert_index_to_coarse: size(je_out(:)) < 2")
2102  if(size(rotate_out(:)) < 2) call mpp_error(fatal,"convert_index_to_coarse: size(rotate_out(:)) < 2")
2103  if( ie > ieg .AND. je > jeg) then
2104  call mpp_error(fatal, "convert_index_to_coarse: do not support cross the corner, contact developer")
2105  endif
2106  if( is > ieg .or. js > jeg) call mpp_error(fatal,.or."convert_index_to_coarse: is > ieg js > jeg")
2107 
2108  nout = 0
2109 
2110  if(tile_coarse == tile_in) then
2111  nout = nout+1
2112  rotate_out(nout) = zero
2113  is_out(nout) = is_in; ie_out(nout) = ie_in + ishift
2114  js_out(nout) = js_in; je_out(nout) = je_in + jshift
2115  endif
2116 
2117  diff = tile_in - tile
2118  if(diff < 0) diff = diff + ntiles
2119  ncross = -1
2120  if( ie > ieg ) then
2121  select case(diff)
2122  case (0)
2123  rotate = zero
2124  ncross = 4
2125  case (1)
2126  if(mod(tile,2) ==1) then ! tile 1 3 5
2127  rotate = zero
2128  ncross = 1
2129  endif
2130  case (2)
2131  if(mod(tile,2) ==0) then ! tile 2 4 6
2132  rotate = minus_ninety
2133  ncross = 1
2134  endif
2135  case (3)
2136  rotate = minus_ninety
2137  ncross = 2
2138  case (4)
2139  if(mod(tile,2) ==1) then ! tile 1 3 5
2140  rotate = minus_ninety
2141  ncross = 3
2142  endif
2143  case (5)
2144  if(mod(tile,2) ==0) then ! tile 2 4 6
2145  rotate = zero
2146  ncross = 3
2147  endif
2148  case default
2149  call mpp_error(fatal,"convert_index_to_coarse: invalid value of diff")
2150  end select
2151 
2152  if(ncross > 0) then
2153  nout =nout+1
2154  rotate_out(nout) = rotate
2155  if(rotate_out(nout) == zero) then
2156  js_out(nout) = js_in
2157  je_out(nout) = je_in + jshift
2158  is_out(nout) = is_in-ncross*ieg
2159  ie_out(nout) = ie_in-ncross*ieg + ishift
2160  else if(rotate_out(nout) == minus_ninety) then
2161  is_out(nout) = ieg-je_in + 1
2162  ie_out(nout) = ieg-js_in + 1 + ishift
2163  js_out(nout) = is_in-ncross*jeg
2164  je_out(nout) = ie_in-ncross*jeg + jshift
2165  endif
2166  endif
2167  else if( je > jeg ) then
2168  select case(diff)
2169  case (0)
2170  rotate = zero
2171  ncross = 4
2172  case (1)
2173  if(mod(tile,2) ==0) then ! tile 2 4 6
2174  rotate = zero
2175  ncross = 1
2176  endif
2177  case (2)
2178  if(mod(tile,2) ==1) then ! tile 1 3 5
2179  rotate = ninety
2180  ncross = 1
2181  endif
2182  case (3)
2183  rotate = ninety
2184  ncross = 2
2185  case (4)
2186  if(mod(tile,2) ==0) then ! tile 2 4 6
2187  rotate = ninety
2188  ncross = 3
2189  endif
2190  case (5)
2191  if(mod(tile,2) ==1) then ! tile 1 3 5
2192  rotate = zero
2193  ncross = 3
2194  endif
2195  end select
2196 
2197  if(ncross > 0) then
2198  nout =nout+1
2199  rotate_out(nout) = rotate
2200 
2201  if(rotate_out(nout) == zero) then
2202  js_out(nout) = js_in
2203  je_out(nout) = je_in + jshift
2204  is_out(nout) = is_in-ncross*ieg
2205  ie_out(nout) = ie_in-ncross*ieg + ishift
2206  else if(rotate_out(nout) == ninety) then
2207  is_out(nout) = js_in - ncross*jeg
2208  ie_out(nout) = je_in - ncross*jeg + ishift
2209  js_out(nout) = jeg - ie_in + 1
2210  je_out(nout) = jeg - is_in + 1 + jshift
2211  endif
2212  endif
2213  endif
2214 
2215  convert_index_to_coarse = nout
2216 
2217 
2218  end function convert_index_to_coarse
2219 
2220 
2221  subroutine convert_index_back(domain, ishift, jshift, rotate, is_in, ie_in, js_in, je_in, is_out, ie_out, &
2222  & js_out, je_out)
2223  type(domain2d), intent(in) :: domain
2224  integer, intent(in) :: ishift, jshift
2225  integer, intent(in) :: is_in, ie_in, js_in, je_in, rotate
2226  integer, intent(out) :: is_out, ie_out, js_out, je_out
2227  integer :: isg, ieg, jsg, jeg
2228  integer :: ncross
2229 
2230  call mpp_get_global_domain(domain, isg, ieg, jsg, jeg)
2231  ncross = 0
2232  if( je_in > jeg+jshift .and. ie_in > ieg+ishift ) then
2233  call mpp_error(fatal,.and."convert_index_back: je_in > jeg ie_in > ieg")
2234  else if (je_in > jeg+jshift) then
2235  ncross = je_in/jeg
2236  select case(rotate)
2237  case(0)
2238  is_out = is_in
2239  ie_out = ie_in
2240  js_out = js_in - ncross*jeg
2241  je_out = je_in - ncross*jeg
2242  case(90)
2243  is_out = js_in - ncross*jeg
2244  ie_out = je_in - ncross*jeg
2245  js_out = jeg - ie_in + 1
2246  je_out = jeg - is_in + 1
2247  case default
2248  call mpp_error(fatal, "convert_index_back: rotate should be 0 or 90 when je_in>jeg")
2249  end select
2250  else if (ie_in > ieg+ishift) then
2251  ncross = ie_in/ieg
2252  select case(rotate)
2253  case(0)
2254  is_out = is_in - ncross*ieg
2255  ie_out = ie_in - ncross*ieg
2256  js_out = js_in
2257  je_out = je_in
2258  case(-90)
2259  js_out = is_in - ncross*ieg
2260  je_out = ie_in - ncross*ieg
2261  is_out = ieg - je_in + 1
2262  ie_out = ieg - js_in + 1
2263  case default
2264  call mpp_error(fatal, "convert_index_back: rotate should be 0 or -90 when ie_in>ieg")
2265  end select
2266  else
2267  is_out = is_in
2268  ie_out = ie_in
2269  js_out = js_in
2270  je_out = je_in
2271  endif
2272 
2273  end subroutine convert_index_back
2274 
2275 
2276 
2277  function get_nest_vector_recv(nest_domain, update_x, update_y, ind_x, ind_y, start_pos, pelist)
2278  type(nest_level_type), intent(in) :: nest_domain
2279  type(nestspec), intent(in) :: update_x, update_y
2280  integer, intent(out) :: ind_x(:), ind_y(:)
2281  integer, intent(out) :: start_pos(:)
2282  integer, intent(out) :: pelist(:)
2283  integer :: get_nest_vector_recv
2284  integer :: nlist, nrecv_x, nrecv_y, ntot, n
2285  integer :: ix, iy, rank_x, rank_y, cur_pos
2286  integer :: nrecv
2287 
2288  nlist = size(nest_domain%pelist)
2289  nrecv_x = update_x%nrecv
2290  nrecv_y = update_y%nrecv
2291 
2292  ntot = nrecv_x + nrecv_y
2293 
2294  n = 1
2295  ix = 1
2296  iy = 1
2297  ind_x = -1
2298  ind_y = -1
2299  nrecv = 0
2300  cur_pos = 0
2301  do while (n<=ntot)
2302  if ( ix <= nrecv_x ) then
2303  rank_x = update_x%recv(ix)%pe
2304  else
2305  rank_x = -1
2306  endif
2307  if ( iy <= nrecv_y ) then
2308  rank_y = update_y%recv(iy)%pe
2309  else
2310  rank_y = -1
2311  endif
2312  nrecv = nrecv + 1
2313  start_pos(nrecv) = cur_pos
2314  if ( (rank_x == rank_y) .and. ( (rank_x >= 0) .and. (rank_y >= 0) ) ) then
2315  n = n+2
2316  ind_x(nrecv) = ix
2317  ind_y(nrecv) = iy
2318  cur_pos = cur_pos + update_x%recv(ix)%totsize + update_y%recv(iy)%totsize
2319  pelist(nrecv) = update_x%recv(ix)%pe
2320  ix = ix + 1
2321  iy = iy + 1
2322  else if ( rank_x < rank_y ) then
2323  n = n+1
2324  if ( rank_x < 0 ) then
2325  ind_x(nrecv) = -1
2326  ind_y(nrecv) = iy
2327  cur_pos = cur_pos + update_y%recv(iy)%totsize
2328  pelist(nrecv) = update_y%recv(iy)%pe
2329  iy = iy + 1
2330  else
2331  ind_x(nrecv) = ix
2332  ind_y(nrecv) = -1
2333  cur_pos = cur_pos + update_x%recv(ix)%totsize
2334  pelist(nrecv) = update_x%recv(ix)%pe
2335  ix = ix + 1
2336  endif
2337  else if ( rank_y < rank_x ) then
2338  n = n+1
2339  if ( rank_y < 0 ) then
2340  ind_x(nrecv) = ix
2341  ind_y(nrecv) = -1
2342  cur_pos = cur_pos + update_x%recv(ix)%totsize
2343  pelist(nrecv) = update_x%recv(ix)%pe
2344  ix = ix + 1
2345  else
2346  ind_x(nrecv) = -1
2347  ind_y(nrecv) = iy
2348  cur_pos = cur_pos + update_y%recv(iy)%totsize
2349  pelist(nrecv) = update_y%recv(iy)%pe
2350  iy = iy + 1
2351  endif
2352  endif
2353  end do
2354 
2355  get_nest_vector_recv = nrecv
2356 
2357 
2358  end function get_nest_vector_recv
2359 
2360 
2361  function get_nest_vector_send(nest_domain, update_x, update_y, ind_x, ind_y, start_pos, pelist)
2362  type(nest_level_type), intent(in) :: nest_domain
2363  type(nestspec), intent(in) :: update_x, update_y
2364  integer, intent(out) :: ind_x(:), ind_y(:)
2365  integer, intent(out) :: start_pos(:)
2366  integer, intent(out) :: pelist(:)
2367  integer :: get_nest_vector_send
2368  integer :: nlist, nsend_x, nsend_y, ntot, n
2369  integer :: ix, iy, rank_x, rank_y, cur_pos
2370  integer :: nsend
2371 
2372  nlist = size(nest_domain%pelist_fine(:)) + size(nest_domain%pelist_coarse(:))
2373  nsend_x = update_x%nsend
2374  nsend_y = update_y%nsend
2375 
2376  ntot = nsend_x + nsend_y
2377 
2378  n = 1
2379  ix = 1
2380  iy = 1
2381  ind_x = -1
2382  ind_y = -1
2383  nsend = 0
2384  cur_pos = 0
2385  do while (n<=ntot)
2386  if ( ix <= nsend_x ) then
2387  rank_x = update_x%send(ix)%pe
2388  else
2389  rank_x = -1
2390  endif
2391  if ( iy <= nsend_y ) then
2392  rank_y = update_y%send(iy)%pe
2393  else
2394  rank_y = -1
2395  endif
2396  nsend = nsend + 1
2397  start_pos(nsend) = cur_pos
2398  if ( (rank_x == rank_y) .and. ( (rank_x >= 0) .and. (rank_y >= 0) ) ) then
2399  n = n+2
2400  ind_x(nsend) = ix
2401  ind_y(nsend) = iy
2402  cur_pos = cur_pos + update_x%send(ix)%totsize + update_y%send(iy)%totsize
2403  pelist(nsend) = update_x%send(ix)%pe
2404  ix = ix + 1
2405  iy = iy + 1
2406  else if ( rank_x < rank_y ) then
2407  n = n+1
2408  if ( rank_x < 0 ) then
2409  ind_x(nsend) = -1
2410  ind_y(nsend) = iy
2411  cur_pos = cur_pos + update_y%send(iy)%totsize
2412  pelist(nsend) = update_y%send(iy)%pe
2413  iy = iy + 1
2414  else
2415  ind_x(nsend) = ix
2416  ind_y(nsend) = -1
2417  cur_pos = cur_pos + update_x%send(ix)%totsize
2418  pelist(nsend) = update_x%send(ix)%pe
2419  ix = ix + 1
2420  endif
2421  else if ( rank_y < rank_x ) then
2422  n = n+1
2423  if ( rank_y < 0 ) then
2424  ind_x(nsend) = ix
2425  ind_y(nsend) = -1
2426  cur_pos = cur_pos + update_x%send(ix)%totsize
2427  pelist(nsend) = update_x%send(ix)%pe
2428  ix = ix + 1
2429  else
2430  ind_x(nsend) = -1
2431  ind_y(nsend) = iy
2432  cur_pos = cur_pos + update_y%send(iy)%totsize
2433  pelist(nsend) = update_y%send(iy)%pe
2434  iy = iy + 1
2435  endif
2436  endif
2437  end do
2438 
2439  get_nest_vector_send = nsend
2440 
2441 
2442  end function get_nest_vector_send
2443 
2444  subroutine check_data_size_1d(module, str1, size1, str2, size2)
2445  character(len=*), intent(in) :: module, str1, str2
2446  integer, intent(in) :: size1, size2
2447 
2448 
2449  if(size2 > 0 .AND. size1 .NE. size2 ) then
2450  print '(a, 3I5)', trim(module), mpp_pe(), size1, size2
2451  call mpp_error(fatal, trim(module)//": mismatch between size of "//trim(str1)//" and "//trim(str2))
2452  endif
2453 
2454  end subroutine check_data_size_1d
2455 
2456 
2457  subroutine check_data_size_2d(module, str1, isize1, jsize1, str2, isize2, jsize2)
2458  character(len=*), intent(in) :: module, str1, str2
2459  integer, intent(in) :: isize1, jsize1, isize2, jsize2
2460 
2461 
2462  if(isize2 > 0 .AND. jsize2 > 0 .AND. (isize1 .NE. isize2 .OR. jsize1 .NE. jsize2) ) then
2463  print '(a, 5I5)', trim(module), mpp_pe(), isize1, jsize1, isize2, jsize2
2464  call mpp_error(fatal, trim(module)//": mismatch between size of "//trim(str1)//" and "//trim(str2))
2465  endif
2466 
2467  end subroutine check_data_size_2d
2468 
2469  function mpp_get_nest_coarse_domain(nest_domain, nest_level)
2470  type(nest_domain_type), intent(in) :: nest_domain
2471  integer, intent(in) :: nest_level
2472  type(domain2d), pointer :: mpp_get_nest_coarse_domain
2473 
2474  if(nest_level < 1 .OR. nest_level > nest_domain%num_level) call mpp_error(fatal, &
2475  & "mpp_define_nest_domains.inc(mpp_get_nest_coarse_domain):"// &
2476  & "nest_level should be between 1 and nest_domain%num_level")
2477 
2478  if(.not. nest_domain%nest(nest_level)%on_level) call mpp_error(fatal, &
2479  "mpp_define_nest_domains.inc(mpp_get_nest_coarse_domain): nest_domain%nest(nest_level)%on_level is false")
2480  mpp_get_nest_coarse_domain => nest_domain%nest(nest_level)%domain_coarse
2481 
2482  end function mpp_get_nest_coarse_domain
2483 
2484  function mpp_get_nest_fine_domain(nest_domain, nest_level)
2485  type(nest_domain_type), intent(in) :: nest_domain
2486  integer, intent(in) :: nest_level
2487  type(domain2d), pointer :: mpp_get_nest_fine_domain
2488 
2489  if(nest_level < 1 .OR. nest_level > nest_domain%num_level) call mpp_error(fatal, &
2490  & "mpp_define_nest_domains.inc(mpp_get_nest_fine_domain):"// &
2491  & "nest_level should be between 1 and nest_domain%num_level")
2492 
2493  if(.not. nest_domain%nest(nest_level)%on_level) call mpp_error(fatal, &
2494  "mpp_define_nest_domains.inc(mpp_get_nest_fine_domain): nest_domain%nest(nest_level)%on_level is false")
2495  mpp_get_nest_fine_domain => nest_domain%nest(nest_level)%domain_fine
2496 
2497  end function mpp_get_nest_fine_domain
2498 
2499  function mpp_get_nest_npes(nest_domain, nest_level)
2500  type(nest_domain_type), intent(in) :: nest_domain
2501  integer, intent(in) :: nest_level
2502  integer :: mpp_get_nest_npes
2503 
2504  if(nest_level < 1 .OR. nest_level > nest_domain%num_level) call mpp_error(fatal, &
2505  "mpp_define_nest_domains.inc(mpp_get_nest_npes): nest_level should be between 1 and nest_domain%num_level")
2506 
2507  mpp_get_nest_npes = size(nest_domain%nest(nest_level)%pelist(:))
2508 
2509  end function mpp_get_nest_npes
2510 
2511  subroutine mpp_get_nest_pelist(nest_domain, nest_level, pelist)
2512  type(nest_domain_type), intent(in) :: nest_domain
2513  integer, intent(in) :: nest_level
2514  integer, intent(out) :: pelist(:)
2515  if(nest_level < 1 .OR. nest_level > nest_domain%num_level) call mpp_error(fatal, &
2516  "mpp_define_nest_domains.inc(mpp_get_nest_pelist): nest_level should be between 1 and nest_domain%num_level")
2517 
2518  if(size(pelist) .NE. size(nest_domain%nest(nest_level)%pelist)) call mpp_error(fatal, &
2519  .NE."mpp_define_nest_domains.inc(mpp_get_nest_pelist): size(pelist) size(nest_domain%nest(nest_level)%pelist)")
2520 
2521  pelist = nest_domain%nest(nest_level)%pelist
2522 
2523  end subroutine mpp_get_nest_pelist
2524 
2525  function mpp_get_nest_fine_npes(nest_domain, nest_level)
2526  type(nest_domain_type), intent(in) :: nest_domain
2527  integer, intent(in) :: nest_level
2528  integer :: mpp_get_nest_fine_npes
2529 
2530  if(nest_level < 1 .OR. nest_level > nest_domain%num_level) call mpp_error(fatal, &
2531  "mpp_define_nest_domains.inc(mpp_get_nest_fine_npes): nest_level should be between 1 and nest_domain%num_level")
2532 
2533  mpp_get_nest_fine_npes = size(nest_domain%nest(nest_level)%pelist_fine(:))
2534 
2535  end function mpp_get_nest_fine_npes
2536 
2537  subroutine mpp_get_nest_fine_pelist(nest_domain, nest_level, pelist)
2538  type(nest_domain_type), intent(in) :: nest_domain
2539  integer, intent(in) :: nest_level
2540  integer, intent(out) :: pelist(:)
2541  if(nest_level < 1 .OR. nest_level > nest_domain%num_level) call mpp_error(fatal, &
2542  & "mpp_define_nest_domains.inc(mpp_get_nest_fine_pelist):"// &
2543  & "nest_level should be between 1 and nest_domain%num_level")
2544 
2545  if(size(pelist) .NE. size(nest_domain%nest(nest_level)%pelist_fine)) call mpp_error(fatal, &
2546  "mpp_define_nest_domains.inc(mpp_get_nest_fine_pelist): size(pelist) "// &
2547  .NE." size(nest_domain%nest(nest_level)%pelist)")
2548 
2549  pelist = nest_domain%nest(nest_level)%pelist_fine
2550 
2551  end subroutine mpp_get_nest_fine_pelist
2552 
2553 
2554 
2555  function mpp_is_nest_fine(nest_domain, nest_level)
2556  type(nest_domain_type), intent(in) :: nest_domain
2557  integer, intent(in) :: nest_level
2558  logical :: mpp_is_nest_fine
2559 
2560  if(nest_level < 1 .OR. nest_level > nest_domain%num_level) call mpp_error(fatal, &
2561  "mpp_define_nest_domains.inc(mpp_is_nest_fine): nest_level should be between 1 and nest_domain%num_level")
2562 
2563  if(.not. nest_domain%nest(nest_level)%on_level) call mpp_error(fatal, &
2564  "mpp_define_nest_domains.inc(mpp_is_nest_fine): nest_domain%nest(nest_level)%on_level is false")
2565 
2566  mpp_is_nest_fine = nest_domain%nest(nest_level)%is_fine_pe
2567 
2568  end function mpp_is_nest_fine
2569 
2570  function mpp_is_nest_coarse(nest_domain, nest_level)
2571  type(nest_domain_type), intent(in) :: nest_domain
2572  integer, intent(in) :: nest_level
2573  logical :: mpp_is_nest_coarse
2574 
2575  if(nest_level < 1 .OR. nest_level > nest_domain%num_level) call mpp_error(fatal, &
2576  "mpp_define_nest_domains.inc(mpp_is_nest_coarse): nest_level should be between 1 and nest_domain%num_level")
2577 
2578  if(.not. nest_domain%nest(nest_level)%on_level) call mpp_error(fatal, &
2579  "mpp_define_nest_domains.inc(mpp_is_nest_coarse): nest_domain%nest(nest_level)%on_level is false")
2580 
2581  mpp_is_nest_coarse = nest_domain%nest(nest_level)%is_coarse_pe
2582 
2583  end function mpp_is_nest_coarse
2584 !> @}
subroutine mpp_define_nest_domains(nest_domain, domain, num_nest, nest_level, tile_fine, tile_coarse, istart_coarse, icount_coarse, jstart_coarse, jcount_coarse, npes_nest_tile, x_refine, y_refine, extra_halo, name)
Set up a domain to pass data between aligned coarse and fine grid of nested model.
subroutine mpp_shift_nest_domains(nest_domain, domain, delta_i_coarse, delta_j_coarse, extra_halo)
Based on mpp_define_nest_domains, but just resets positioning of nest Modifies the parent/coarse star...
subroutine compute_overlap_fine_to_coarse(nest_domain, overlap, position, name)
This routine will compute the send and recv information between overlapped nesting region....
integer function convert_index_to_nest(domain, ishift, jshift, tile_coarse, istart_coarse, iend_coarse, jstart_coarse, jend_coarse, ntiles_coarse, tile_in, is_in, ie_in, js_in, je_in, is_out, ie_out, js_out, je_out, rotate_out)
This routine will convert the global coarse grid index to nest grid index.
subroutine mpp_get_f2c_index_coarse(nest_domain, is_coarse, ie_coarse, js_coarse, je_coarse, nest_level, position)
subroutine mpp_get_f2c_index_fine(nest_domain, is_coarse, ie_coarse, js_coarse, je_coarse, is_fine, ie_fine, js_fine, je_fine, nest_level, position)
subroutine define_nest_level_type(nest_domain, x_refine, y_refine, extra_halo)
subroutine mpp_get_domain_shift(domain, ishift, jshift, position)
Returns the shift value in x and y-direction according to domain position..
subroutine mpp_get_c2f_index(nest_domain, is_fine, ie_fine, js_fine, je_fine, is_coarse, ie_coarse, js_coarse, je_coarse, dir, nest_level, position)
Get the index of the data passed from coarse grid to fine grid.
integer function stdout()
This function returns the current standard fortran unit numbers for output.
Definition: mpp_util.inc:43
subroutine mpp_set_current_pelist(pelist, no_sync)
Set context pelist.
Definition: mpp_util.inc:490
subroutine mpp_declare_pelist(pelist, name, commID)
Declare a pelist.
Definition: mpp_util.inc:461
integer function mpp_npes()
Returns processor count for current pelist.
Definition: mpp_util.inc:421
integer function mpp_pe()
Returns processor ID.
Definition: mpp_util.inc:407