FV3 Bundle
mpp_pset.F90
Go to the documentation of this file.
1 !***********************************************************************
2 !* GNU Lesser General Public License
3 !*
4 !* This file is part of the GFDL Flexible Modeling System (FMS).
5 !*
6 !* FMS is free software: you can redistribute it and/or modify it under
7 !* the terms of the GNU Lesser General Public License as published by
8 !* the Free Software Foundation, either version 3 of the License, or (at
9 !* your option) any later version.
10 !*
11 !* FMS is distributed in the hope that it will be useful, but WITHOUT
12 !* ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
13 !* FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
14 !* for more details.
15 !*
16 !* You should have received a copy of the GNU Lesser General Public
17 !* License along with FMS. If not, see <http://www.gnu.org/licenses/>.
18 !***********************************************************************
19 ! module within MPP for handling PSETs:
20 ! PSET: Persistent Shared-memory Execution Thread
21 !
22 ! AUTHOR: V. Balaji (v.balaji@noaa.gov)
23 ! DATE: 2006-01-15
24 #ifdef test_mpp_pset
25 !PSET_DEBUG is always turned on in the test program
26 #define PSET_DEBUG
27 #endif
28 
30 #include <fms_platform.h>
31  use mpp_mod, only: mpp_pe, mpp_npes, mpp_root_pe, mpp_send, mpp_recv, &
32  mpp_sync, mpp_error, fatal, warning, stdout, stderr, mpp_chksum, &
33  mpp_declare_pelist, mpp_get_current_pelist, mpp_set_current_pelist, &
34  mpp_init, comm_tag_1, comm_tag_2, comm_tag_3, mpp_sync_self
35  implicit none
36  private
37 
38 !private variables
39  integer :: pe
40  integer :: commid !MPI communicator, copy here from pset
41  logical :: verbose=.false.
42  logical :: module_is_initialized=.false.
43  character(len=256) :: text
44 #ifdef use_SGI_GSM
45 #include <mpp/shmem.fh>
46  integer :: psync(shmem_barrier_sync_size)
47  pointer( p_psync, psync ) !used by SHPALLOC
48 #endif
49 !generic interfaces
51  module procedure mpp_pset_broadcast_ptr_scalar
52  module procedure mpp_pset_broadcast_ptr_array
53  end interface
54  interface mpp_send_ptr
55  module procedure mpp_send_ptr_scalar
56  module procedure mpp_send_ptr_array
57  end interface
58  interface mpp_recv_ptr
59  module procedure mpp_recv_ptr_scalar
60  module procedure mpp_recv_ptr_array
61  end interface
63  module procedure mpp_pset_print_chksum_1d
64  module procedure mpp_pset_print_chksum_2d
65  module procedure mpp_pset_print_chksum_3d
66  module procedure mpp_pset_print_chksum_4d
67  end interface
68 !public type
69  type :: mpp_pset_type
70  private
71 #ifdef IBM_FIX
72  sequence
73 #endif
74  integer :: npset !number of PSETs
75  integer :: next_in_pset, prev_in_pset !next and prev PE in PSET (cyclic)
76  integer :: root_in_pset !PE designated to be the root within PSET
77  logical :: root !true if you are the root PSET
78  integer :: pos !position of current PE within pset
79 !stack is allocated by root
80 !it is then mapped to mpp_pset_stack by mpp_pset_broadcast_ptr
81  real, _allocatable :: stack(:) _null
82  integer, _allocatable :: pelist(:) _null !base PElist
83  integer, _allocatable :: root_pelist(:) _null !a PElist of all the roots
84  integer, _allocatable :: pset(:) _null !PSET IDs
85  integer(POINTER_KIND) :: p_stack
86  integer :: lstack, maxstack, hiwm !current stack length, max, hiWM
87  integer :: commid
88  character(len=32) :: name
89  logical :: initialized=.false.
90  end type mpp_pset_type
91 !public types
92  public :: mpp_pset_type
93 !public variables
94 !public member functions
100 
101 
102 contains
103  subroutine mpp_pset_init
104 #ifdef use_SGI_GSM
105  integer :: err
106 #ifdef sgi_mipspro
107  character(len=8) :: value !won't be read
108  integer :: lenname, lenval!won't be read
109 #endif
110  if( module_is_initialized )return
111 !this part needs to be called _all_ PEs
112  call shmem_barrier_all()
113  call shpalloc( p_psync, shmem_barrier_sync_size, err, -1 )
114  call shmem_barrier_all()
115 #ifdef sgi_mipspro
116  call pxfgetenv( 'SMA_GLOBAL_ALLOC', 0, value, lenval, err )
117  if( err.NE.0 )call mpp_error( fatal, &
118  'The environment variable SMA_GLOBAL_ALLOC must be set on Irix.' )
119 #endif
120 #endif
121  module_is_initialized = .true.
122  end subroutine mpp_pset_init
123 
124  subroutine mpp_pset_create(npset,pset,stacksize,pelist, commID)
125 !create PSETs
126 ! called by all PEs in parent pelist
127 ! mpset must be exact divisor of npes
128  integer, intent(in) :: npset !number of PSETs per set
129  type(mpp_pset_type), intent(inout) :: pset
130  integer, intent(in), optional :: stacksize
131  integer, intent(in), optional :: pelist(:)
132  integer, intent(in), optional :: commid
133 
134  integer :: npes, my_commid
135  integer :: i, j, k, out_unit, errunit
136  integer, allocatable :: my_pelist(:), root_pelist(:)
137 
138  call mpp_init()
139  call mpp_pset_init()
140 
141 #ifdef PSET_DEBUG
142  verbose=.true.
143 #endif
144  out_unit = stdout()
145  errunit = stderr()
146  pe = mpp_pe()
147  if(present(pelist)) then
148  npes = size(pelist(:))
149  else
150  npes = mpp_npes()
151  endif
152  if( mod(npes,npset).NE.0 )then
153  write( text,'(a,2i6)' ) &
154  'MPP_PSET_CREATE: PSET size (npset) must divide npes exactly:'// &
155  ' npset, npes=', npset, npes
156  call mpp_error( fatal, text )
157  end if
158 
159  !configure out root_pelist
160  allocate(my_pelist(0:npes-1) )
161  allocate(root_pelist(0:npes/npset-1) )
162  if(present(pelist)) then
163  if(.not. present(commid)) call mpp_error(fatal, &
164  'MPP_PSET_CREATE: when pelist is present, commID should also be present')
165  my_pelist = pelist
166  my_commid = commid
167  else
168  call mpp_get_current_pelist(my_pelist, commid = my_commid)
169  endif
170  do i = 0,npes/npset-1
171  root_pelist(i) = my_pelist(npset*i)
172  enddo
173  write( out_unit,'(a,i6)' )'MPP_PSET_CREATE creating PSETs... npset=', npset
174  if(any(my_pelist == pe) ) then
175  if( pset%initialized )call mpp_error( fatal, &
176  'MPP_PSET_CREATE: PSET already initialized!' )
177  pset%npset = npset
178  allocate( pset%pelist(0:npes-1) )
179  allocate( pset%root_pelist(0:npes/npset-1) )
180  pset%commID = my_commid
181  pset%pelist = my_pelist
182 !create the root PElist
183  pset%root_pelist = root_pelist
184  allocate( pset%pset(0:npset-1) )
185  do i = 0,npes/npset-1
186  k = npset*i
187 !designate the root PE, next PE, prev PE
188  do j = 0,npset-1
189  if( pe.EQ.pset%pelist(k+j) )then
190  pset%pset(:) = pset%pelist(k:k+npset-1)
191  pset%pos = j
192  pset%root_in_pset = pset%root_pelist(i)
193  if( j.EQ.0 )then
194  pset%prev_in_pset = pset%pelist(k+npset-1)
195  else
196  pset%prev_in_pset = pset%pelist(k+j-1)
197  end if
198  if( j.EQ.npset-1 )then
199  pset%next_in_pset = pset%pelist(k)
200  else
201  pset%next_in_pset = pset%pelist(k+j+1)
202  end if
203  end if
204  end do
205  end do
206 
207  pset%root = pe.EQ.pset%root_in_pset
208 
209 !stack
210  pset%hiWM = 0 !initialize hi-water-mark
211  pset%maxstack = 1000000 !default
212  if( PRESENT(stacksize) )pset%maxstack = stacksize
213  write( out_unit,'(a,i8)' ) &
214  'MPP_PSET_CREATE: setting stacksize=', pset%maxstack
215  if( pset%root )then
216  allocate( pset%stack(pset%maxstack) )
217 #ifdef use_CRI_pointers
218  pset%p_stack = loc(pset%stack)
219 #endif
220  end if
221  pset%initialized = .true. !must be called before using pset
222  call mpp_pset_broadcast_ptr(pset,pset%p_stack)
223  endif
224 
225  call mpp_declare_pelist(root_pelist)
226 
227  if( verbose )then
228  write( errunit,'(a,4i6)' )'MPP_PSET_CREATE: pe, root, next, prev=', &
229  pe, pset%root_in_pset, pset%next_in_pset, pset%prev_in_pset
230  write( errunit,* )'PE ', pe, ' pset=', pset%pset(:)
231  write( out_unit,* )'root pelist=', pset%root_pelist(:)
232  end if
233  end subroutine mpp_pset_create
234 
235  subroutine mpp_pset_delete(pset)
236  type(mpp_pset_type), intent(inout) :: pset
237  integer :: out_unit
238 
239  out_unit = stdout()
240  if( .NOT.pset%initialized )call mpp_error( fatal, &
241  'MPP_PSET_DELETE: called with uninitialized PSET.' )
242 !deallocate arrays...
243  deallocate( pset%pelist )
244  deallocate( pset%root_pelist )
245  deallocate( pset%pset )
246  if( pset%root )deallocate( pset%stack )
247  write( out_unit, '(a,i10)' ) &
248  'Deleting PSETs... stack high-water-mark=', pset%hiWM
249 !... and set status flag
250  pset%initialized = .false.
251  end subroutine mpp_pset_delete
252 
253  subroutine mpp_send_ptr_scalar( ptr, pe )
254  integer(POINTER_KIND), intent(in) :: ptr
255  integer, intent(in) :: pe
256 
257 !currently only wraps mpp_send
258 !on some architectures, mangling might occur
259  call mpp_send( ptr, pe, tag=comm_tag_1 )
260  end subroutine mpp_send_ptr_scalar
261 
262  subroutine mpp_send_ptr_array( ptr, pe )
263  integer(POINTER_KIND), intent(in) :: ptr(:)
264  integer, intent(in) :: pe
265 
266 !currently only wraps mpp_send
267 !on some architectures, mangling might occur
268  call mpp_send( ptr, size(ptr), pe, tag=comm_tag_2 )
269  end subroutine mpp_send_ptr_array
270 
271  subroutine mpp_recv_ptr_scalar( ptr, pe )
272  integer(POINTER_KIND), intent(inout) :: ptr
273  integer, intent(in) :: pe
274 
275  call mpp_recv( ptr, pe, tag=comm_tag_1 )
276  call mpp_translate_remote_ptr( ptr, pe )
277  return
278  end subroutine mpp_recv_ptr_scalar
279 
280  subroutine mpp_recv_ptr_array( ptr, pe )
281  integer(POINTER_KIND), intent(inout) :: ptr(:)
282  integer, intent(in) :: pe
283  integer :: i
284 
285  call mpp_recv( ptr, size(ptr), pe, tag=comm_tag_2 )
286  do i = 1, size(ptr)
287  call mpp_translate_remote_ptr( ptr(i), pe )
288  end do
289  return
290  end subroutine mpp_recv_ptr_array
291 
292  subroutine mpp_translate_remote_ptr( ptr, pe )
293 !modifies the received pointer to correct numerical address
294  integer(POINTER_KIND), intent(inout) :: ptr
295  integer, intent(in) :: pe
296 #ifdef use_SGI_GSM
297 !from the MPI_SGI_GLOBALPTR manpage
298 ! POINTER(global_ptr, global_addr)
299 ! INTEGER rem_rank, comm, ierror
300 ! INTEGER(KIND=MPI_ADDRESS_KIND) rem_addr, size, global_addr
301 !
302 ! CALL MPI_SGI_GLOBALPTR(rem_addr, size, rem_rank, comm, global_ptr, ierror)
303  real :: dummy
304  pointer( p, dummy )
305  integer :: ierror
306 !length goes in the second argument to MPI_SGI_GLOBALPTR
307 ! according to Kim Mcmahon, this is only used to ensure the requested array
308 ! length is within the valid memory-mapped region. We do not have access to
309 ! the actual array length, so we are only going to set it to 1. This might
310 ! unexpectedly fail on some large model.
311  integer(POINTER_KIND) :: length=1
312 #ifdef sgi_mipspro
313  return !no translation needed on sgi_mipspro if SMA_GLOBAL_ALLOC is set
314 #endif
315 #ifdef use_libMPI
316 !the MPI communicator was stored in pset%commID
317 !since this routine doesn't take a pset argument, we let the caller store
318 !it in the module global variable commID (see broadcast_ptr and check_ptr)
319  p = ptr
320  call mpi_sgi_globalptr( dummy, length, pe, commid, ptr, ierror )
321  if( ierror.EQ.-1 )call mpp_error( fatal, &
322  'MPP_TRANSLATE_REMOTE_PTR: unknown MPI_SGI_GLOBALPTR error.' )
323 #else
324  call mpp_error( fatal, &
325  'MPP_TRANSLATE_REMOTE_PTR now only works under -Duse_libMPI' )
326 #endif
327 #endif
328  return
329  end subroutine mpp_translate_remote_ptr
330 
331  subroutine mpp_pset_sync(pset)
332 !this is a replacement for mpp_sync, doing syncs across
333 !shared arrays without calling mpp_sync
334  type(mpp_pset_type), intent(in) :: pset
335 
336  if( .NOT.pset%initialized )call mpp_error( fatal, &
337  'MPP_PSET_SYNC: called with uninitialized PSET.' )
338 #ifdef use_SGI_GSM
339 !assumes npset contiguous PEs starting with root_in_pset
340  call shmem_barrier( pset%root_in_pset, 0, pset%npset, psync )
341 #else
342 !currently does mpp_sync!!! slow!!!
343 !try and make a lightweight pset sync
344  call mpp_sync
345 #endif
346  end subroutine mpp_pset_sync
347 
348  subroutine mpp_pset_broadcast(pset,a)
349 !broadcast value on the root to its sub-threads
350  type(mpp_pset_type), intent(in) :: pset
351  real, intent(inout) :: a
352  integer :: i
353 
354  if( .NOT.pset%initialized )call mpp_error( fatal, &
355  'MPP_PSET_BROADCAST: called with uninitialized PSET.' )
356  if( pset%root )then
357  do i = 1,pset%npset-1
358  call mpp_send( a, pset%pset(i), tag=comm_tag_3 )
359  end do
360  else
361  call mpp_recv( a, pset%root_in_pset, tag=comm_tag_3 )
362  end if
363  call mpp_pset_sync(pset)
364  end subroutine mpp_pset_broadcast
365 
366  subroutine mpp_pset_broadcast_ptr_scalar(pset,ptr)
367 !create a shared array by broadcasting pointer
368 !root allocates memory and passes pointer in
369 !on return all other PSETs will have the pointer to a shared object
370  type(mpp_pset_type), intent(in) :: pset
371  integer(POINTER_KIND), intent(inout) :: ptr
372  integer :: i
373 
374  if( .NOT.pset%initialized )call mpp_error( fatal, &
375  'MPP_PSET_BROADCAST_PTR: called with uninitialized PSET.' )
376  commid = pset%commID !pass to mpp_translate_remote_ptr
377  if( pset%root )then
378  do i = 1,pset%npset-1
379  call mpp_send_ptr( ptr, pset%pset(i) )
380  end do
381  else
382  call mpp_recv_ptr( ptr, pset%root_in_pset )
383  end if
384  call mpp_sync_self()
385  end subroutine mpp_pset_broadcast_ptr_scalar
386 
387  subroutine mpp_pset_broadcast_ptr_array(pset,ptr)
388 !create a shared array by broadcasting pointer
389 !root allocates memory and passes pointer in
390 !on return all other PSETs will have the pointer to a shared object
391  type(mpp_pset_type), intent(in) :: pset
392  integer(POINTER_KIND), intent(inout) :: ptr(:)
393  integer :: i
394 
395  if( .NOT.pset%initialized )call mpp_error( fatal, &
396  'MPP_PSET_BROADCAST_PTR: called with uninitialized PSET.' )
397  commid = pset%commID !pass to mpp_translate_remote_ptr
398  if( pset%root )then
399  do i = 1,pset%npset-1
400  call mpp_send_ptr( ptr, pset%pset(i) )
401  end do
402  else
403  call mpp_recv_ptr( ptr, pset%root_in_pset )
404  end if
405  call mpp_sync_self()
406 
407  end subroutine mpp_pset_broadcast_ptr_array
408 
409  subroutine mpp_pset_check_ptr(pset,ptr)
410 !checks if the supplied pointer is indeed shared
411  type(mpp_pset_type), intent(in) :: pset
412 #ifdef use_CRI_pointers
413  real :: dummy
414  pointer( ptr, dummy )
415 #else
416  integer(POINTER_KIND), intent(in) :: ptr
417 #endif
418 #ifdef PSET_DEBUG
419  integer(POINTER_KIND) :: p
420  integer :: i
421  if( .NOT.pset%initialized )call mpp_error( fatal, &
422  'MPP_PSET_CHECK_PTR: called with uninitialized PSET.' )
423  commid = pset%commID !pass to mpp_translate_remote_ptr
424 !check if this is a shared pointer
425  p = ptr
426  if( pset%root )then
427  do i = 1,pset%npset-1
428  call mpp_send_ptr( p, pset%pset(i) )
429  end do
430  else
431  call mpp_recv_ptr( p, pset%root_in_pset )
432  end if
433  call mpp_pset_sync(pset)
434  if( p.NE.ptr )call mpp_error( fatal, &
435  'MPP_PSET_CHECK_PTR: pointers do not match!' )
436 #else
437 !do nothing if the debug CPP flag isn't on
438 #endif
439  end subroutine mpp_pset_check_ptr
440 
441  subroutine mpp_pset_segment_array( pset, ls, le, lsp, lep )
442 !given input indices ls, le, returns indices lsp, lep
443 !so that segments span the range ls:le with no overlaps.
444 !attempts load balance: also some PSETs might get lsp>lep
445 !so that do-loops will be null
446  type(mpp_pset_type), intent(in) :: pset
447  integer, intent(in) :: ls, le
448  integer, intent(out) :: lsp, lep
449  integer :: i
450 
451  if( .NOT.pset%initialized )call mpp_error( fatal, &
452  'MPP_PSET_SEGMENT_ARRAY: called with uninitialized PSET.' )
453 #ifdef PSET_DEBUG
454  if( le-ls+1.LT.pset%npset )then
455  write( text,'(3(a,i6))' ) &
456  'MPP_PSET_ARRAY_SEGMENT: parallel range (', ls, ',', le, &
457  ') is smaller than the number of threads:', pset%npset
458  call mpp_error( warning, text )
459  end if
460 #endif
461  lep = ls-1 !initialize so that lsp is correct on first pass
462  do i = 0,pset%pos
463  lsp = lep + 1
464  lep = lsp + ceiling( REAL(le-lsp+1)/(pset%npset-i) ) - 1
465  end do
466  end subroutine mpp_pset_segment_array
467 
468  subroutine mpp_pset_stack_push( pset, ptr, len )
469 !mpp_malloc specialized for shared arrays
470 !len is the length of the required array
471 !lstack is the stack already in play
472 !user should zero lstack (call mpp_pset_stack_reset) when the stack is to be cleared
473  type(mpp_pset_type), intent(inout) :: pset
474  integer, intent(in) :: len
475 #ifdef use_CRI_pointers
476  real :: dummy
477  pointer( ptr, dummy )
478  real :: stack(pset%maxstack)
479  pointer( p, stack )
480 
481  if( .NOT.pset%initialized )call mpp_error( fatal, &
482  'MPP_PSET_STACK_PUSH: called with uninitialized PSET.' )
483  if( pset%lstack+len.GT.pset%maxstack )then
484  write( text, '(a,3i12)' ) &
485  'MPP_PSET_STACK_PUSH: mpp_pset_stack overflow: '// &
486  .GT.'len+lstackmaxstack. len, lstack, maxstack=', &
487  len, pset%lstack, pset%maxstack
488  call mpp_error( fatal, text )
489  end if
490  p = pset%p_stack !point stack to shared stack pointer
491  ptr = loc( stack(pset%lstack+1) )
492  call mpp_pset_check_ptr(pset,ptr) !make sure ptr is the same across PSETs
493  pset%lstack = pset%lstack + len
494  pset%hiWM = max( pset%hiWM, pset%lstack )
495 #else
496  integer(POINTER_KIND), intent(out) :: ptr
497  call mpp_error( fatal, &
498  'MPP_PSET_STACK_PUSH only works with Cray pointers.' )
499 #endif
500  end subroutine mpp_pset_stack_push
501 
502  subroutine mpp_pset_stack_reset(pset)
503  type(mpp_pset_type), intent(inout) :: pset
504 !reset stack... will reuse any temporary arrays! USE WITH CARE
505 !next few lines are to zero stack contents...
506 !but it's better noone tries to use uninitialized stack variables!
507 ! integer :: l1, l2
508 ! real :: mpp_pset_stack(maxstack)
509 ! pointer( p_mpp_pset_stack, mpp_pset_stack )
510 ! p_mpp_pset_stack = ptr_mpp_pset_stack
511 ! call mpp_pset_array_segment( 1, lstack, l1, l2 )
512 ! mpp_pset_stack(l1:l2) = 0.
513  if( .NOT.pset%initialized )call mpp_error( fatal, &
514  'MPP_PSET_STACK_RESET: called with uninitialized PSET.' )
515  pset%lstack = 0
516  end subroutine mpp_pset_stack_reset
517 
518  subroutine mpp_pset_print_chksum_1d(pset, caller, array)
519 !print a checksum of an array
520 !pass the whole domain seen by root PSET
521 !add lines to check on shared array?
522  type(mpp_pset_type), intent(in) :: pset
523  character(len=*), intent(in) :: caller
524  real, intent(in) :: array(:)
525  integer :: errunit
526 
527 #ifdef PSET_DEBUG
528  logical :: do_print
529  integer(LONG_KIND) :: chksum
530 
531  if( .NOT.pset%initialized )call mpp_error( fatal, &
532  'MPP_PSET_PRINT_CHKSUM: called with uninitialized PSET.' )
533  errunit = stderr()
534 
535  if( pset%root )then
536  do_print = pe.EQ.mpp_root_pe() !set to T to print from all PEs
537  call mpp_set_current_pelist(pset%root_pelist)
538  chksum = mpp_chksum( array )
539  if( do_print ) &
540  write( errunit, '(a,z18)' )trim(caller)//' chksum=', chksum
541  end if
542  call mpp_set_current_pelist(pset%pelist)
543 #endif
544  return
545  end subroutine mpp_pset_print_chksum_1d
546 
547  subroutine mpp_pset_print_chksum_2d(pset, caller, array)
548  type(mpp_pset_type), intent(in) :: pset
549  character(len=*), intent(in) :: caller
550  real, intent(in) :: array(:,:)
551  real :: array1D( size(array) )
552 #ifdef use_CRI_pointers
553  pointer( p, array1d )
554  p = loc(array)
555 #else
556  array1d = transfer( array, array1d )
557 #endif
558  call mpp_pset_print_chksum(pset, caller, array1d)
559  end subroutine mpp_pset_print_chksum_2d
560 
561  subroutine mpp_pset_print_chksum_3d(pset, caller, array)
562  type(mpp_pset_type), intent(in) :: pset
563  character(len=*), intent(in) :: caller
564  real, intent(in) :: array(:,:,:) !overload for other ranks
565  real :: array1D( size(array) )
566 #ifdef use_CRI_pointers
567  pointer( p, array1d )
568  p = loc(array)
569 #else
570  array1d = transfer( array, array1d )
571 #endif
572  call mpp_pset_print_chksum(pset, caller, array1d)
573  end subroutine mpp_pset_print_chksum_3d
574 
575  subroutine mpp_pset_print_chksum_4d(pset, caller, array)
576  type(mpp_pset_type), intent(in) :: pset
577  character(len=*), intent(in) :: caller
578  real, intent(in) :: array(:,:,:,:)
579  real :: array1D( size(array) )
580 #ifdef use_CRI_pointers
581  pointer( p, array1d )
582  p = loc(array)
583 #else
584  array1d = transfer( array, array1d )
585 #endif
586  call mpp_pset_print_chksum(pset, caller, array1d)
587  end subroutine mpp_pset_print_chksum_4d
588 
589  subroutine mpp_pset_print_stack_chksum( pset, caller )
590  type(mpp_pset_type), intent(in) :: pset
591  character(len=*), intent(in) :: caller
592 
593  if( .NOT.pset%initialized )call mpp_error( fatal, &
594  'MPP_PSET_PRINT_STACK_CHKSUM: called with uninitialized PSET.' )
595  call mpp_pset_print_chksum( pset, trim(caller)//' stack', &
596  pset%stack(1:pset%lstack) )
597  end subroutine mpp_pset_print_stack_chksum
598 
599 !accessor functions
600  function mpp_pset_root(pset)
601  logical :: mpp_pset_root
602  type(mpp_pset_type), intent(in) :: pset
603 
604  if( .NOT.pset%initialized )call mpp_error( fatal, &
605  'MPP_PSET_ROOT: called with uninitialized PSET.' )
606  mpp_pset_root = pset%root
607  end function mpp_pset_root
608 
609  function mpp_pset_numroots(pset)
610 !necessary to export root_pelist: caller needs to pre-allocate
611  integer :: mpp_pset_numroots
612  type(mpp_pset_type), intent(in) :: pset
613 
614  if( .NOT.pset%initialized )call mpp_error( fatal, &
615  'MPP_PSET_NUMROOTS: called with uninitialized PSET.' )
616  mpp_pset_numroots = size(pset%root_pelist)
617  end function mpp_pset_numroots
618 
619  subroutine mpp_pset_get_root_pelist(pset,pelist,commID)
620  type(mpp_pset_type), intent(in) :: pset
621  integer, intent(out) :: pelist(:)
622  integer, intent(out), optional :: commid
623 
624  if( .NOT.pset%initialized )call mpp_error( fatal, &
625  'MPP_PSET_GET_ROOT_PELIST: called with uninitialized PSET.' )
626  if( size(pelist).NE.size(pset%root_pelist) )then
627  write( text,'(a,2i6)' ) &
628  'pelist argument has wrong size: requested, actual=', &
629  size(pelist), size(pset%root_pelist)
630  call mpp_error( fatal, 'MPP_PSET_GET_ROOT_PELIST: '//text )
631  end if
632  pelist(:) = pset%root_pelist(:)
633  if( PRESENT(commid) )then
634 #ifdef use_libMPI
635  commid = pset%commID
636 #else
637  call mpp_error( warning, &
638  'MPP_PSET_GET_ROOT_PELIST: commID is only defined under -Duse_libMPI.' )
639 #endif
640  end if
641  end subroutine mpp_pset_get_root_pelist
642 
643 end module mpp_pset_mod
644 
subroutine, public mpp_pset_stack_reset(pset)
Definition: mpp_pset.F90:503
subroutine mpp_translate_remote_ptr(ptr, pe)
Definition: mpp_pset.F90:293
subroutine mpp_send_ptr_array(ptr, pe)
Definition: mpp_pset.F90:263
subroutine, public mpp_pset_print_stack_chksum(pset, caller)
Definition: mpp_pset.F90:590
logical module_is_initialized
Definition: mpp_pset.F90:42
subroutine, public mpp_pset_segment_array(pset, ls, le, lsp, lep)
Definition: mpp_pset.F90:442
subroutine, public mpp_pset_stack_push(pset, ptr, len)
Definition: mpp_pset.F90:469
subroutine, public mpp_pset_check_ptr(pset, ptr)
Definition: mpp_pset.F90:410
character(len=256) text
Definition: mpp_pset.F90:43
subroutine mpp_send_ptr_scalar(ptr, pe)
Definition: mpp_pset.F90:254
subroutine mpp_pset_broadcast_ptr_scalar(pset, ptr)
Definition: mpp_pset.F90:367
logical verbose
Definition: mpp_pset.F90:41
Definition: mpp.F90:39
subroutine, public mpp_pset_create(npset, pset, stacksize, pelist, commID)
Definition: mpp_pset.F90:125
subroutine, public mpp_pset_init
Definition: mpp_pset.F90:104
subroutine mpp_pset_print_chksum_1d(pset, caller, array)
Definition: mpp_pset.F90:519
subroutine, public mpp_pset_delete(pset)
Definition: mpp_pset.F90:236
subroutine mpp_recv_ptr_array(ptr, pe)
Definition: mpp_pset.F90:281
logical function, public mpp_pset_root(pset)
Definition: mpp_pset.F90:601
subroutine mpp_pset_print_chksum_2d(pset, caller, array)
Definition: mpp_pset.F90:548
integer function, public mpp_pset_numroots(pset)
Definition: mpp_pset.F90:610
subroutine, public mpp_pset_get_root_pelist(pset, pelist, commID)
Definition: mpp_pset.F90:620
#define max(a, b)
Definition: mosaic_util.h:33
subroutine, public mpp_pset_sync(pset)
Definition: mpp_pset.F90:332
integer commid
Definition: mpp_pset.F90:40
subroutine, public mpp_pset_broadcast(pset, a)
Definition: mpp_pset.F90:349
subroutine mpp_recv_ptr_scalar(ptr, pe)
Definition: mpp_pset.F90:272
subroutine mpp_pset_print_chksum_3d(pset, caller, array)
Definition: mpp_pset.F90:562
subroutine mpp_pset_print_chksum_4d(pset, caller, array)
Definition: mpp_pset.F90:576
subroutine mpp_pset_broadcast_ptr_array(pset, ptr)
Definition: mpp_pset.F90:388