FMS  2024.03
Flexible Modeling System
mpp_domains_comm.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 domain communications via @ref domaincommunicator2d
25 
26 !> @addtogroup mpp_domains_mod
27 !> @{
28  function mpp_redistribute_init_comm(domain_in,l_addrs_in, domain_out,l_addrs_out, &
29  isize_in,jsize_in,ksize_in,isize_out,jsize_out,ksize_out) RESULT(d_comm)
30  type(DomainCommunicator2D), pointer :: d_comm
31  type(domain2D),target, intent(in) :: domain_in
32  integer(i8_kind), intent(in) :: l_addrs_in(:)
33  type(domain2D),target, intent(in) :: domain_out
34  integer(i8_kind), intent(in) :: l_addrs_out(:)
35  integer, intent(in) :: isize_in
36  integer, intent(in) :: jsize_in
37  integer, intent(in) :: ksize_in
38  integer, intent(in) :: isize_out
39  integer, intent(in) :: jsize_out
40  integer, intent(in) :: ksize_out
41 
42  integer(i8_kind) :: domain_id
43  integer :: m, list
44  integer :: is, ie, js, je, ke, ioff, joff, list_size
45  integer :: isc, iec, jsc, jec, mytile
46  integer :: lsize,rsize,msgsize,to_pe,from_pe
47  integer, allocatable,dimension(:) :: isL, jsL
48  integer(i8_kind),allocatable,dimension(:,:) :: slist_addr
49  character(len=8) :: text
50 
51 
52  ! This test determines whether input fields are from allocated memory (LOC gets global
53  ! address) or "static" memory (need shmem_ptr). This probably needs to be generalized
54  ! to determine appropriate mechanism for each incoming address.
55 
56  ! "Concurrent" run mode may leave field_in or field_out unassociated if pe does not
57  ! contain in/out data. Use of STATIC option for ocean complicates this as ocean component
58  ! always defined. Field_out is always a boundary structure and so is always allocated or
59  ! not depending on whether it's used. If field out is defined (>0), then it is used otherwise
60  ! field in must be defined.
61 
62 !fix ke
63  ke = 0
64  if( domain_in%pe /= null_pe )ke = ksize_in
65  if( domain_out%pe /= null_pe )then
66  if( ke /= 0 .AND. ke /= ksize_out ) &
67  call mpp_error( fatal, 'MPP_REDISTRIBUTE_INIT_COMM: mismatch between field_in and field_out.' )
68  ke = ksize_out
69  end if
70  if( ke == 0 )call mpp_error( fatal, &
71  & 'MPP_REDISTRIBUTE_INIT_COMM: either domain_in or domain_out must be native.' )
72 !check sizes
73  if( domain_in%pe /= null_pe )then
74  if( isize_in /= domain_in%x(1)%domain_data%size .OR. jsize_in /= domain_in%y(1)%domain_data%size ) &
75  call mpp_error( fatal, 'MPP_REDISTRIBUTE_INIT_COMM: field_in must be on data domain of domain_in.' )
76  end if
77  if( domain_out%pe /= null_pe )then
78  if( isize_out /= domain_out%x(1)%domain_data%size .OR. jsize_out /= domain_out%y(1)%domain_data%size ) &
79  call mpp_error( fatal, 'MPP_REDISTRIBUTE_INIT_COMM: field_out must be on data domain of domain_out.' )
80  end if
81 
82 
83  ! Create unique domain identifier
84  list_size = size(l_addrs_in(:))
85  if(l_addrs_out(1) > 0)then
86  domain_id = set_domain_id(domain_out%id,ke+list_size)
87  else
88  domain_id = set_domain_id(domain_in%id,ke+list_size)
89  endif
90 
91  d_comm =>get_comm(domain_id,l_addrs_in(1),l_addrs_out(1))
92 
93  if(d_comm%initialized)return ! Found existing field/domain communicator
94 
95  d_comm%l_addr = l_addrs_in(1)
96  d_comm%domain_in =>domain_in
97  d_comm%Slist_size = size(domain_out%list(:))
98  d_comm%isize_in = isize_in
99  d_comm%jsize_in = jsize_in
100  d_comm%ke = ke
101 
102 !send
103  lsize = d_comm%Slist_size-1
104  allocate(d_comm%sendis(1,0:lsize), d_comm%sendie(1,0:lsize), &
105  d_comm%sendjs(1,0:lsize), d_comm%sendje(1,0:lsize), &
106  d_comm%S_msize(0:lsize),isl(0:lsize),jsl(0:lsize))
107  allocate(slist_addr(list_size,0:lsize))
108  allocate(d_comm%cto_pe(0:lsize), d_comm%S_do_buf(0:lsize))
109 
110  isl=0;jsl=0
111  slist_addr = -9999
112  d_comm%cto_pe=-1
113  d_comm%sendis=0; d_comm%sendie=0
114  d_comm%sendjs=0; d_comm%sendje=0;
115  d_comm%S_msize=0
116  d_comm%S_do_buf=.false.
117 
118  ioff = domain_in%x(1)%domain_data%begin
119  joff = domain_in%y(1)%domain_data%begin
120  mytile = domain_in%tile_id(1)
121 
122  call mpp_get_compute_domain( domain_in, isc, iec, jsc, jec )
123  do list = 0,lsize
124  m = mod( domain_out%pos+list+lsize+1, lsize+1 )
125  if( mytile .NE. domain_out%list(m)%tile_id(1) ) cycle
126  d_comm%cto_pe(list) = domain_out%list(m)%pe
127  to_pe = d_comm%cto_pe(list)
128  is = domain_out%list(m)%x(1)%compute%begin
129  ie = domain_out%list(m)%x(1)%compute%end
130  js = domain_out%list(m)%y(1)%compute%begin
131  je = domain_out%list(m)%y(1)%compute%end
132  is = max(is,isc); ie = min(ie,iec)
133  js = max(js,jsc); je = min(je,jec)
134  if( ie >= is .AND. je >= js )then
135  d_comm%S_do_buf(list) = .true.
136  d_comm%sendis(1,list)=is; d_comm%sendie(1,list)=ie
137  d_comm%sendjs(1,list)=js; d_comm%sendje(1,list)=je
138  d_comm%S_msize(list) = (ie-is+1)*(je-js+1)*ke
139  isl(list) = is-ioff+1; jsl(list) = js-joff+1
140  end if
141  end do
142 
143  call mpp_sync_self()
144 !recv
145  d_comm%domain_out =>domain_out
146  d_comm%Rlist_size = size(domain_in%list(:))
147  d_comm%isize_out = isize_out
148  d_comm%jsize_out = jsize_out
149 
150  rsize = d_comm%Rlist_size-1
151  allocate(d_comm%recvis(1,0:rsize), d_comm%recvie(1,0:rsize), &
152  d_comm%recvjs(1,0:rsize), d_comm%recvje(1,0:rsize), &
153  d_comm%R_msize(0:rsize))
154  allocate(d_comm%cfrom_pe(0:rsize), d_comm%R_do_buf(0:rsize))
155  allocate(d_comm%isizeR(0:rsize), d_comm%jsizeR(0:rsize))
156  allocate(d_comm%sendisR(1,0:rsize), d_comm%sendjsR(1,0:rsize))
157  allocate(d_comm%rem_addrl(list_size,0:rsize))
158  d_comm%rem_addrl=-9999
159  d_comm%cfrom_pe=-1
160  d_comm%recvis=0; d_comm%recvie=0
161  d_comm%recvjs=0; d_comm%recvje=0;
162  d_comm%R_msize=0
163  d_comm%R_do_buf=.false.
164  d_comm%isizeR=0; d_comm%jsizeR=0
165  d_comm%sendisR=0; d_comm%sendjsR=0
166 
167  mytile = domain_out%tile_id(1)
168  call mpp_get_compute_domain( domain_out, isc, iec, jsc, jec )
169  do list = 0,rsize
170  m = mod( domain_in%pos+rsize+1-list, rsize+1 )
171  if( mytile .NE. domain_in%list(m)%tile_id(1) ) cycle
172  d_comm%cfrom_pe(list) = domain_in%list(m)%pe
173  from_pe = d_comm%cfrom_pe(list)
174  is = domain_in%list(m)%x(1)%compute%begin
175  ie = domain_in%list(m)%x(1)%compute%end
176  js = domain_in%list(m)%y(1)%compute%begin
177  je = domain_in%list(m)%y(1)%compute%end
178  is = max(is,isc); ie = min(ie,iec)
179  js = max(js,jsc); je = min(je,jec)
180  if( ie >= is .AND. je >= js )then
181  d_comm%R_do_buf(list) = .true.
182  d_comm%recvis(1,list)=is; d_comm%recvie(1,list)=ie
183  d_comm%recvjs(1,list)=js; d_comm%recvje(1,list)=je
184  d_comm%R_msize(list) = (ie-is+1)*(je-js+1)*ke
185  end if
186  end do
187 
188  d_comm%isize_max = isize_in; call mpp_max(d_comm%isize_max)
189  d_comm%jsize_max = jsize_in; call mpp_max(d_comm%jsize_max)
190 
191  ! Handles case where S_msize and/or R_msize are 0 size array
192  msgsize = ( maxval( (/0,sum(d_comm%S_msize(:))/) ) + maxval( (/0,sum(d_comm%R_msize(:))/) ) ) * list_size
193  if(msgsize>0)then
194  mpp_domains_stack_hwm = max( mpp_domains_stack_hwm, msgsize )
195  if( mpp_domains_stack_hwm.GT.mpp_domains_stack_size )then
196  write( text,'(i8)' )mpp_domains_stack_hwm
197  call mpp_error( fatal, &
198  & 'MPP_REDISTRIBUTE_INIT_COMM: mpp_domains_stack overflow, call mpp_domains_set_stack_size(' &
199  & //trim(text)//') from all PEs.' )
200  end if
201  end if
202 
203  DEALLOCATE(slist_addr,isl,jsl)
204 
205  d_comm%initialized = .true.
206 
207  end function mpp_redistribute_init_comm
208 
209  !> initializes a @ref DomainCommunicator2D type for use in @ref mpp_global_field
210  function mpp_global_field_init_comm(domain,l_addr,isize_g,jsize_g,isize_l, &
211  jsize_l, ksize,l_addr2,flags, position) RESULT(d_comm)
212  type(domaincommunicator2d), pointer :: d_comm
213  type(domain2d),target, intent(in) :: domain
214  integer(i8_kind), intent(in) :: l_addr
215  integer, intent(in) :: isize_g
216  integer, intent(in) :: jsize_g
217  integer, intent(in) :: isize_l
218  integer, intent(in) :: jsize_l
219  integer, intent(in) :: ksize
220  integer(i8_kind),optional,intent(in) :: l_addr2
221  integer, optional, intent(in) :: flags
222  integer, optional, intent(in) :: position
223 
224  integer(i8_kind) :: domain_id
225  integer :: n, lpos, rpos, list, nlist, tile_id
226  integer :: update_flags
227  logical :: xonly, yonly
228  integer :: is, ie, js, je, ioff, joff, ishift, jshift
229  integer :: lsize,msgsize,from_pe
230  integer, allocatable,dimension(:) :: isl, jsl
231  integer(i8_kind),allocatable,dimension(:,:) :: slist_addr
232  integer(i8_kind),save ,dimension(2) :: rem_addr
233  character(len=8) :: text
234 
235  if( .NOT.module_is_initialized )call mpp_error( fatal, 'MPP_GLOBAL_FIELD: must first call mpp_domains_init.' )
236  update_flags=xupdate+yupdate; xonly = .false.; yonly = .false.
237  if( PRESENT(flags) )then
238  update_flags = flags
239  xonly = btest(flags,east)
240  yonly = btest(flags,south)
241  if( .NOT.xonly .AND. .NOT.yonly )call mpp_error( warning, &
242  'MPP_GLOBAL_FIELD: you must have flags=XUPDATE, YUPDATE or XUPDATE+YUPDATE' )
243  if(xonly .AND. yonly) then
244  xonly = .false.; yonly = .false.
245  endif
246  end if
247 
248  call mpp_get_domain_shift(domain, ishift, jshift, position=position)
249 
250  if( isize_g /= (domain%x(1)%global%size+ishift) .OR. jsize_g /= (domain%y(1)%global%size+jshift) ) &
251  call mpp_error( fatal, 'MPP_GLOBAL_FIELD_INIT_COMM: incoming arrays do not match domain.' )
252 
253  if( isize_l == (domain%x(1)%compute%size+ishift) .AND. jsize_l == (domain%y(1)%compute%size+jshift) )then
254 !local is on compute domain
255  ioff = -domain%x(1)%compute%begin + 1
256  joff = -domain%y(1)%compute%begin + 1
257  elseif( isize_l == (domain%x(1)%memory%size+ishift) .AND. jsize_l == (domain%y(1)%memory%size+jshift) )then
258 !local is on data domain
259  ioff = -domain%x(1)%domain_data%begin + 1
260  joff = -domain%y(1)%domain_data%begin + 1
261  else
262  call mpp_error(fatal, &
263  & 'MPP_GLOBAL_FIELD_INIT_COMM: incoming field array must match either compute domain or data domain.')
264  endif
265 
266  ! Create unique domain identifier
267  domain_id=set_domain_id(domain%id,ksize,update_flags, position=position)
268  d_comm =>get_comm(domain_id,l_addr,l_addr2)
269 
270  if(d_comm%initialized)return ! Found existing field/domain communicator
271  d_comm%domain => domain
272  d_comm%isize_in = isize_l; d_comm%isize_out = isize_g
273  d_comm%jsize_in = jsize_l; d_comm%jsize_out = jsize_g
274  d_comm%ke = ksize
275  d_comm%gf_ioff=ioff; d_comm%gf_joff=joff
276 
277 !fill off-domains (note loops begin at an offset of 1)
278  if( xonly )then
279  lsize = size(domain%x(1)%list(:))
280 !send
281  allocate(d_comm%cto_pe(0:lsize-1))
282  d_comm%cto_pe=-1
283  do list = 0,lsize-1
284  lpos = mod(domain%x(1)%pos+lsize-list,lsize)
285  d_comm%cto_pe(list) = domain%x(1)%list(lpos)%pe
286  end do
287 !recv
288  allocate(d_comm%cfrom_pe(0:lsize-1))
289  allocate(d_comm%recvis(1,0:lsize-1), d_comm%recvie(1,0:lsize-1), &
290  d_comm%recvjs(1,0:lsize-1), d_comm%recvje(1,0:lsize-1), &
291  d_comm%R_msize(0:lsize-1))
292  d_comm%cfrom_pe=-1
293  d_comm%recvis=0; d_comm%recvie=0
294  d_comm%recvjs=0; d_comm%recvje=0;
295  d_comm%R_msize=0
296  do list = 0,lsize-1
297  rpos = mod(domain%x(1)%pos+list,lsize)
298  from_pe = domain%x(1)%list(rpos)%pe
299  d_comm%cfrom_pe(list) = from_pe
300  is = domain%list(from_pe)%x(1)%compute%begin; ie = domain%list(from_pe)%x(1)%compute%end+ishift
301  js = domain%y(1)%compute%begin; je = domain%y(1)%compute%end+jshift
302  d_comm%recvis(1,list)=is; d_comm%recvie(1,list)=ie
303  d_comm%recvjs(1,list)=js; d_comm%recvje(1,list)=je
304  d_comm%R_msize(list) = (ie-is+1) * (je-js+1) * ksize
305  end do
306 
307  elseif( yonly )then
308  lsize = size(domain%y(1)%list(:))
309 !send
310  allocate(d_comm%cto_pe(0:lsize-1))
311  d_comm%cto_pe=-1
312  do list = 0,lsize
313  lpos = mod(domain%y(1)%pos+lsize-list,lsize)
314  d_comm%cto_pe(list) = domain%y(1)%list(lpos)%pe
315  end do
316 !recv
317  allocate(d_comm%cfrom_pe(0:lsize-1))
318  allocate(d_comm%recvis(1,0:lsize-1), d_comm%recvie(1,0:lsize-1), &
319  d_comm%recvjs(1,0:lsize-1), d_comm%recvje(1,0:lsize-1), &
320  d_comm%R_msize(0:lsize-1))
321  d_comm%cfrom_pe=-1
322  d_comm%recvis=0; d_comm%recvie=0
323  d_comm%recvjs=0; d_comm%recvje=0;
324  d_comm%R_msize=0
325  do list = 0,lsize-1
326  rpos = mod(domain%y(1)%pos+list,lsize)
327  from_pe = domain%y(1)%list(rpos)%pe
328  d_comm%cfrom_pe(list) = from_pe
329  is = domain%x(1)%compute%begin; ie = domain%x(1)%compute%end+ishift
330  js = domain%list(from_pe)%y(1)%compute%begin; je = domain%list(from_pe)%y(1)%compute%end+jshift
331  d_comm%recvis(1,list)=is; d_comm%recvie(1,list)=ie
332  d_comm%recvjs(1,list)=js; d_comm%recvje(1,list)=je
333  d_comm%R_msize(list) = (ie-is+1) * (je-js+1) * ksize
334  end do
335 
336  else
337  nlist = size(domain%list(:))
338  tile_id = domain%tile_id(1)
339 
340  lsize = 0
341  do list = 0,nlist-1
342  if( domain%list(list)%tile_id(1) .NE. tile_id ) cycle
343  lsize = lsize+1
344  end do
345 
346  !send
347  allocate(d_comm%cto_pe(0:lsize-1))
348  d_comm%cto_pe=-1
349  n = 0
350  do list = 0,nlist-1
351  lpos = mod(domain%pos+nlist-list,nlist)
352  if( domain%list(lpos)%tile_id(1) .NE. tile_id ) cycle
353  d_comm%cto_pe(n) = domain%list(lpos)%pe
354  n = n + 1
355  end do
356  !recv
357  allocate(d_comm%cfrom_pe(0:lsize-1))
358  allocate(d_comm%recvis(1,0:lsize-1), d_comm%recvie(1,0:lsize-1), &
359  d_comm%recvjs(1,0:lsize-1), d_comm%recvje(1,0:lsize-1), &
360  d_comm%R_msize(0:lsize-1))
361  d_comm%cfrom_pe=-1
362  d_comm%recvis=0; d_comm%recvie=0
363  d_comm%recvjs=0; d_comm%recvje=0;
364  d_comm%R_msize=0
365  n = 0
366  do list = 0,nlist-1
367  rpos = mod(domain%pos+list,nlist)
368  if( domain%list(rpos)%tile_id(1) .NE. tile_id ) cycle
369  d_comm%cfrom_pe(n) = domain%list(rpos)%pe
370  is = domain%list(rpos)%x(1)%compute%begin; ie = domain%list(rpos)%x(1)%compute%end+ishift
371  js = domain%list(rpos)%y(1)%compute%begin; je = domain%list(rpos)%y(1)%compute%end+jshift
372  d_comm%recvis(1,n)=is; d_comm%recvie(1,n)=ie
373  d_comm%recvjs(1,n)=js; d_comm%recvje(1,n)=je
374  d_comm%R_msize(n) = (je-js+1) * (ie-is+1) * ksize
375  n = n+1
376  end do
377 
378  endif
379 
380  d_comm%Slist_size = lsize
381  d_comm%Rlist_size = lsize
382 
383 !send
384  allocate(d_comm%sendis(1,0:lsize-1), d_comm%sendie(1,0:lsize-1), &
385  d_comm%sendjs(1,0:lsize-1), d_comm%sendje(1,0:lsize-1), &
386  d_comm%S_msize(0:lsize-1),isl(0:lsize-1),jsl(0:lsize-1))
387  allocate(slist_addr(2,0:lsize-1))
388  isl=0; jsl=0
389  slist_addr = -9999
390  d_comm%sendis=0; d_comm%sendie=0
391  d_comm%sendjs=0; d_comm%sendje=0;
392  d_comm%S_msize=0
393  do list = 0,lsize-1
394  is=domain%x(1)%compute%begin; ie=domain%x(1)%compute%end+ishift
395  js=domain%y(1)%compute%begin; je=domain%y(1)%compute%end+jshift
396  d_comm%sendis(1,list)=is; d_comm%sendie(1,list)=ie
397  d_comm%sendjs(1,list)=js; d_comm%sendje(1,list)=je
398  d_comm%S_msize(list) = (je-js+1) * (ie-is+1) * ksize
399  isl(list) = ioff+domain%x(1)%compute%begin; jsl(list) = joff+domain%y(1)%compute%begin
400  end do
401 
402 !recv
403  allocate(d_comm%isizeR(0:lsize-1), d_comm%jsizeR(0:lsize-1))
404  allocate(d_comm%sendisR(1,0:lsize-1), d_comm%sendjsR(1,0:lsize-1))
405  if(.not.PRESENT(l_addr2))then
406  allocate(d_comm%rem_addr(0:lsize-1))
407  d_comm%rem_addr=-9999
408  else
409  allocate(d_comm%rem_addrx(0:lsize-1),d_comm%rem_addry(0:lsize-1))
410  d_comm%rem_addrx=-9999; d_comm%rem_addry=-9999
411  endif
412  d_comm%isizeR=0; d_comm%jsizeR=0
413  d_comm%sendisR=0; d_comm%sendjsR=0
414  rem_addr = -9999
415 
416  ! Handles case where S_msize and/or R_msize are 0 size array
417  msgsize = maxval( (/0,sum(d_comm%S_msize(:))/) ) + maxval( (/0,sum(d_comm%R_msize(:))/) )
418  if(msgsize>0)then
419  mpp_domains_stack_hwm = max( mpp_domains_stack_hwm, msgsize )
420  if( mpp_domains_stack_hwm.GT.mpp_domains_stack_size )then
421  write( text,'(i8)' )mpp_domains_stack_hwm
422  call mpp_error( fatal, &
423  & 'MPP_GLOBAL_FIELD_INIT_COMM: mpp_domains_stack overflow, call mpp_domains_set_stack_size(' &
424  & //trim(text)//') from all PEs.' )
425  end if
426  end if
427 
428  DEALLOCATE(slist_addr,isl,jsl)
429 
430  d_comm%initialized = .true.
431 
432  end function mpp_global_field_init_comm
433 
434 
435  subroutine mpp_redistribute_free_comm(domain_in,l_addr,domain_out,l_addr2,ksize,lsize)
436  ! Since initialization of the d_comm type is expensive, freeing should be a rare
437  ! event. Thus no attempt is made to salvage freed d_comm's.
438  type(domain2d), intent(in) :: domain_in
439  integer(i8_kind), intent(in) :: l_addr
440  type(domain2d), intent(in) :: domain_out
441  integer(i8_kind), intent(in) :: l_addr2
442  integer, intent(in) :: ksize,lsize
443 
444  integer(i8_kind) :: domain_id
445 
446  if(l_addr2 > 0)then
447  domain_id = set_domain_id(domain_out%id,ksize+lsize)
448  else
449  domain_id = set_domain_id(domain_in%id,ksize+lsize)
450  endif
451  call free_comm(domain_id,l_addr,l_addr2)
452  end subroutine mpp_redistribute_free_comm
453 
454 
455  subroutine mpp_global_field_free_comm(domain,l_addr,ksize,l_addr2,flags)
456  ! Since initialization of the d_comm type is expensive, freeing should be a rare
457  ! event. Thus no attempt is made to salvage freed d_comm's.
458  type(domain2d), intent(in) :: domain
459  integer(i8_kind), intent(in) :: l_addr
460  integer, intent(in) :: ksize
461  integer(i8_kind),optional,intent(in) :: l_addr2
462  integer, optional,intent(in) :: flags
463 
464  integer :: update_flags
465  integer(i8_kind) :: domain_id
466 
467  update_flags=0; if(PRESENT(flags))update_flags=flags
468  domain_id=set_domain_id(domain%id,ksize,update_flags)
469  call free_comm(domain_id,l_addr,l_addr2)
470  end subroutine mpp_global_field_free_comm
471 
472 
473  subroutine free_comm(domain_id,l_addr,l_addr2)
474  ! Since initialization of the d_comm type is expensive, freeing should be a rare
475  ! event. Thus no attempt is made to salvage freed d_comm's.
476  integer(i8_kind), intent(in) :: domain_id
477  integer(i8_kind), intent(in) :: l_addr
478  integer(i8_kind),optional,intent(in) :: l_addr2
479 
480  integer(i8_kind) :: dc_key,a_key
481  integer :: dc_idx,a_idx,i_idx,insert,insert_a,insert_i
482  integer :: a2_idx,insert_a2
483 
484 
485  i_idx = find_key(domain_id,ids_sorted(1:n_ids),insert_i)
486  a_idx = find_key(l_addr,addrs_sorted(1:a_sort_len),insert_a)
487  a_key = int(addrs_idx(a_idx),kind(i8_kind))
488  if(PRESENT(l_addr2))then
489  a2_idx = find_key(l_addr2,addrs2_sorted(1:a2_sort_len),insert_a2)
490  a_key = a_key + addr2_base*int(addrs2_idx(a2_idx),kind(i8_kind))
491  endif
492  dc_key = domain_id_base*int(ids_idx(i_idx),kind(i8_kind)) + a_key
493  dc_idx = find_key(dc_key,dckey_sorted(1:dc_sort_len),insert)
494 
495  if(dc_idx < 0)then
496  call mpp_error(fatal,'FREE_COMM: attempt to remove nonexistent domains communicator key')
497  endif
498  call deallocate_comm(d_comm(dc_idx))
499  call pop_key(dckey_sorted,d_comm_idx,dc_sort_len,dc_idx)
500  call pop_key(addrs_sorted,addrs_idx,a_sort_len,a_idx)
501  if(PRESENT(l_addr2))call pop_key(addrs2_sorted,addrs2_idx,a2_sort_len,a2_idx)
502  end subroutine free_comm
503 
504 
505  function get_comm(domain_id,l_addr,l_addr2)
506  integer(i8_kind),intent(in) :: domain_id
507  integer(i8_kind),intent(in) :: l_addr
508  integer(i8_kind),intent(in),optional :: l_addr2
509  type(domaincommunicator2d), pointer :: get_comm
510 
511  integer(i8_kind) :: dc_key,a_key
512  integer :: i,dc_idx,a_idx,i_idx,insert,insert_a,insert_i
513  integer :: a2_idx,insert_a2
514 
515  if(.not.ALLOCATED(d_comm))ALLOCATE(d_comm(max_fields))
516  i_idx = find_key(domain_id,ids_sorted(1:n_ids),insert_i)
517  a_idx = find_key(l_addr,addrs_sorted(1:a_sort_len),insert_a)
518  a_key = int(addrs_idx(a_idx),kind(i8_kind))
519  if(PRESENT(l_addr2))then
520  a2_idx = find_key(l_addr2,addrs2_sorted(1:a2_sort_len),insert_a2)
521  a_key = a_key + addr2_base*int(addrs2_idx(a2_idx),kind(i8_kind))
522  endif
523  dc_key = domain_id_base*int(ids_idx(i_idx),kind(i8_kind)) + a_key
524  dc_idx = find_key(dc_key,dckey_sorted(1:dc_sort_len),insert)
525  if(dc_idx > 0)then
526  get_comm =>d_comm(d_comm_idx(dc_idx))
527  else
528  if(i_idx<0)then
529  if(n_ids == max_dom_ids)then
530  call mpp_error(fatal,'GET_COMM: Maximum number of domains exceeded')
531  endif
532  n_ids = n_ids+1
533  i_idx = push_key(ids_sorted,ids_idx,i_sort_len,insert_i,domain_id,n_ids)
534  endif
535  if(a_idx<0)then
536  if(n_addrs == max_addrs)then
537  call mpp_error(fatal,'GET_COMM: Maximum number of memory addresses exceeded')
538  endif
539  n_addrs = n_addrs + 1
540  a_idx = push_key(addrs_sorted,addrs_idx,a_sort_len,insert_a,l_addr,n_addrs)
541  endif
542  if(PRESENT(l_addr2))then
543  if(a2_idx<0)then
544  if(n_addrs2 == max_addrs2)then
545  call mpp_error(fatal,'GET_COMM: Maximum number of 2nd memory addresses exceeded')
546  endif
547  n_addrs2 = n_addrs2 + 1
548  a2_idx = push_key(addrs2_sorted,addrs2_idx,a2_sort_len,insert_a2,l_addr2,n_addrs2)
549  endif
550  endif
551  if(n_comm == max_fields)then
552  call mpp_error(fatal,'GET_COMM: Maximum number of fields exceeded')
553  endif
554  a_key = int(addrs_idx(a_idx),kind(8))
555  if(PRESENT(l_addr2))a_key = a_key + addr2_base*int(addrs2_idx(a2_idx),kind(8))
556  dc_key = domain_id_base*int(ids_idx(i_idx),kind(i8_kind)) + a_key
557  dc_idx = find_key(dc_key,dckey_sorted(1:dc_sort_len),insert)
558  if(dc_idx /= -1)call mpp_error(fatal,'GET_COMM: attempt to insert existing key')
559  n_comm = n_comm + 1
560  i = push_key(dckey_sorted,d_comm_idx,dc_sort_len,insert,dc_key,n_comm)
561  d_comm_idx(insert) = n_comm
562  if(PRESENT(l_addr2))then
563  d_comm(n_comm)%l_addrx = l_addr
564  d_comm(n_comm)%l_addry = l_addr2
565  else
566  d_comm(n_comm)%l_addr = l_addr
567  endif
568  get_comm =>d_comm(n_comm)
569  endif
570  end function get_comm
571 
572 
573  function push_key(sorted,idx,n_idx,insert,key,ival)
574  integer(i8_kind),intent(inout),dimension(:) :: sorted
575  integer, intent(inout),dimension(-1:) :: idx ! Start -1 to simplify first call logic in get_comm
576  integer, intent(inout) :: n_idx
577  integer, intent(in) :: insert
578  integer(i8_kind),intent(in) :: key
579  integer, intent(in) :: ival
580 
581  integer :: push_key,i
582 
583  do i=n_idx,insert,-1
584  sorted(i+1) = sorted(i)
585  idx(i+1) = idx(i)
586  end do
587  sorted(insert) = key
588  n_idx = n_idx + 1
589  idx(insert) = ival
590  push_key = insert
591  end function push_key
592 
593 
594  subroutine pop_key(sorted,idx,n_idx,key_idx)
595  integer(i8_kind),intent(inout),dimension(:) :: sorted
596  integer, intent(inout),dimension(-1:) :: idx ! Start -1 to simplify first call logic in get_comm
597  integer, intent(inout) :: n_idx
598  integer, intent(in) :: key_idx
599 
600  integer :: i
601 
602  do i=key_idx,n_idx-1
603  sorted(i) = sorted(i+1)
604  idx(i) = idx(i+1)
605  end do
606  sorted(n_idx) = -9999
607  idx(n_idx) = -9999
608  n_idx = n_idx - 1
609  end subroutine pop_key
610 
611 
612  function find_key(key,sorted,insert) RESULT(n)
613  ! The algorithm used here requires monotonic keys w/out repetition.
614  integer(i8_kind),intent(in) :: key ! new address to be found in list
615  integer(i8_kind),dimension(:),intent(in) :: sorted ! list of sorted local addrs
616  integer, intent(out) :: insert
617  integer :: n, n_max, n_min, n_key
618  logical :: not_found
619 
620  n_key = size(sorted(:))
621  insert = 1
622  n = -1 ! value not in list
623  if(n_key == 0)return ! first call
624 
625  if(key < sorted(1))then
626  insert = 1; return
627  elseif(key > sorted(n_key))then
628  insert = n_key+1; return
629  endif
630 
631  if(key == sorted(1))then
632  n = 1; return
633  elseif(key == sorted(n_key))then
634  n = n_key; return
635  endif
636 
637  not_found = .true.
638  n = n_key/2 + 1
639  n_min=1; n_max=n_key
640  do while(not_found)
641  if(key == sorted(n))then
642  not_found = .false.
643  elseif(key > sorted(n))then
644  if(key < sorted(n+1))then
645  insert = n+1; exit
646  endif
647  n_min = n
648  n = (n+1+n_max)/2
649  else
650  if(key > sorted(n-1))then
651  insert = n; exit
652  endif
653  n_max = n
654  n = (n+n_min)/2
655  endif
656  if(n==1 .or. n==n_key)exit
657  end do
658  if(not_found)n = -1 ! value not in list
659  end function find_key
660 
661 
662  subroutine deallocate_comm(d_comm)
663  type(domaincommunicator2d), intent(inout) :: d_comm
664 
665  d_comm%domain =>null()
666  d_comm%domain_in =>null()
667  d_comm%domain_out =>null()
668 
669  d_comm%initialized=.false.
670  d_comm%id=-9999
671  d_comm%l_addr =-9999
672  d_comm%l_addrx =-9999
673  d_comm%l_addry =-9999
674 
675  if( allocated(d_comm%sendis) ) DEALLOCATE(d_comm%sendis); !!d_comm%sendis =>NULL()
676  if( allocated(d_comm%sendie) ) DEALLOCATE(d_comm%sendie); !!d_comm%sendie =>NULL()
677  if( allocated(d_comm%sendjs) ) DEALLOCATE(d_comm%sendjs); !!d_comm%sendjs =>NULL()
678  if( allocated(d_comm%sendje) ) DEALLOCATE(d_comm%sendje); !!d_comm%sendje =>NULL()
679  if( allocated(d_comm%S_msize) ) DEALLOCATE(d_comm%S_msize); !!d_comm%S_msize =>NULL()
680  if( allocated(d_comm%S_do_buf) ) DEALLOCATE(d_comm%S_do_buf); !!d_comm%S_do_buf =>NULL()
681  if( allocated(d_comm%cto_pe) ) DEALLOCATE(d_comm%cto_pe); !!d_comm%cto_pe =>NULL()
682  if( allocated(d_comm%recvis) ) DEALLOCATE(d_comm%recvis); !!d_comm%recvis =>NULL()
683  if( allocated(d_comm%recvie) ) DEALLOCATE(d_comm%recvie); !!d_comm%recvie =>NULL()
684  if( allocated(d_comm%recvjs) ) DEALLOCATE(d_comm%recvjs); !!d_comm%recvjs =>NULL()
685  if( allocated(d_comm%recvje) ) DEALLOCATE(d_comm%recvje); !!d_comm%recvje =>NULL()
686  if( allocated(d_comm%R_msize) ) DEALLOCATE(d_comm%R_msize); !!d_comm%R_msize =>NULL()
687  if( allocated(d_comm%R_do_buf) ) DEALLOCATE(d_comm%R_do_buf); !!d_comm%R_do_buf =>NULL()
688  if( allocated(d_comm%cfrom_pe) ) DEALLOCATE(d_comm%cfrom_pe); !!d_comm%cfrom_pe =>NULL()
689  d_comm%Slist_size=0; d_comm%Rlist_size=0
690  d_comm%isize=0; d_comm%jsize=0; d_comm%ke=0
691  d_comm%isize_in=0; d_comm%jsize_in=0
692  d_comm%isize_out=0; d_comm%jsize_out=0
693  d_comm%isize_max=0; d_comm%jsize_max=0
694  d_comm%gf_ioff=0; d_comm%gf_joff=0
695  ! Remote data
696  if( allocated(d_comm%isizeR) ) DEALLOCATE(d_comm%isizeR); !!dd_comm%isizeR =>NULL()
697  if( allocated(d_comm%jsizeR) ) DEALLOCATE(d_comm%jsizeR); !!dd_comm%jsizeR =>NULL()
698  if( allocated(d_comm%sendisR) ) DEALLOCATE(d_comm%sendisR); !!dd_comm%sendisR =>NULL()
699  if( allocated(d_comm%sendjsR) ) DEALLOCATE(d_comm%sendjsR); !!dd_comm%sendjsR =>NULL()
700  if( allocated(d_comm%rem_addr) ) DEALLOCATE(d_comm%rem_addr); !!dd_comm%rem_addr =>NULL()
701  if( allocated(d_comm%rem_addrx) )DEALLOCATE(d_comm%rem_addrx); !!dd_comm%rem_addrx =>NULL()
702  if( allocated(d_comm%rem_addry) )DEALLOCATE(d_comm%rem_addry); !!dd_comm%rem_addry =>NULL()
703  if( allocated(d_comm%rem_addrl) )DEALLOCATE(d_comm%rem_addrl); !!dd_comm%rem_addrl =>NULL()
704  end subroutine deallocate_comm
705 
706 
707  function set_domain_id(d_id,ksize,flags,gtype, position, whalo, ehalo, shalo, nhalo)
708  integer(i8_kind), intent(in) :: d_id
709  integer , intent(in) :: ksize
710  integer , optional, intent(in) :: flags
711  integer , optional, intent(in) :: gtype
712  integer , optional, intent(in) :: position
713  integer , optional, intent(in) :: whalo, ehalo, shalo, nhalo
714 
715  integer(i8_kind) :: set_domain_id
716 
717  set_domain_id=d_id + ke_base*int(ksize,kind(d_id))
718  if(PRESENT(flags))set_domain_id=set_domain_id+int(flags,kind(d_id))
719  if(PRESENT(gtype))set_domain_id=set_domain_id+gt_base*int(gtype,kind(d_id)) ! Must be i8_kind arithmetic
720  !--- gtype is never been used to set id. we need to add position to calculate id to seperate
721  !--- BGRID and CGRID or scalar variable.
722  if(present(position)) set_domain_id=set_domain_id+gt_base*int(2**position, kind(d_id))
723  !z1l ???? the following calculation may need to be revised
724  if(present(whalo)) then
725  if(whalo>=0) then
726  set_domain_id=set_domain_id+gt_base*int(2**4*2**whalo, kind(d_id))
727  else
728  set_domain_id=set_domain_id-gt_base*int(2**4*2**(-whalo), kind(d_id))
729  endif
730  end if
731  if(present(ehalo)) then
732  if(ehalo>=0) then
733  set_domain_id=set_domain_id+gt_base*int(2**4*2**ehalo, kind(d_id))
734  else
735  set_domain_id=set_domain_id-gt_base*int(2**4*2**(-ehalo), kind(d_id))
736  endif
737  end if
738  if(present(shalo)) then
739  if(shalo>=0) then
740  set_domain_id=set_domain_id+gt_base*int(2**4*2**shalo, kind(d_id))
741  else
742  set_domain_id=set_domain_id-gt_base*int(2**4*2**(-shalo), kind(d_id))
743  endif
744  end if
745  if(present(nhalo)) then
746  if(nhalo>=0) then
747  set_domain_id=set_domain_id+gt_base*int(2**4*2**nhalo, kind(d_id))
748  else
749  set_domain_id=set_domain_id-gt_base*int(2**4*2**(-nhalo), kind(d_id))
750  endif
751  end if
752  end function set_domain_id
753 !> @}
type(domaincommunicator2d) function, pointer mpp_global_field_init_comm(domain, l_addr, isize_g, jsize_g, isize_l, jsize_l, ksize, l_addr2, flags, position)
initializes a DomainCommunicator2D type for use in mpp_global_field
subroutine mpp_get_domain_shift(domain, ishift, jshift, position)
Returns the shift value in x and y-direction according to domain position..
subroutine mpp_sync_self(pelist, check, request, msg_size, msg_type)
This is to check if current PE's outstanding puts are complete but we can't use shmem_fence because w...