FV3 Bundle
test_mpp.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 #ifdef test_mpp
20 #ifdef SYSTEM_CLOCK
21 #undef SYSTEM_CLOCK
22 #endif
23 
24 program test !test various aspects of mpp_mod
25 #include <fms_platform.h>
26 
27 #ifdef sgi_mipspro
28  use shmem_interface
29 #endif
30 
31  use mpp_mod, only : mpp_init, mpp_exit, mpp_pe, mpp_npes, mpp_root_pe, stdout
32  use mpp_mod, only : mpp_clock_id, mpp_clock_begin, mpp_clock_end, mpp_sync, mpp_malloc
33  use mpp_mod, only : mpp_declare_pelist, mpp_set_current_pelist, mpp_set_stack_size
35  use mpp_mod, only : mpp_gather, mpp_error, fatal, mpp_sync_self
36  use mpp_io_mod, only: mpp_io_init, mpp_flush
37 #ifdef use_MPI_GSM
38  use mpp_mod, only : mpp_gsm_malloc, mpp_gsm_free
39 #endif
40 
41  implicit none
42 
43  integer, parameter :: n=1048576
44  real, allocatable, dimension(:) :: a, b, c
45 #ifdef use_MPI_GSM
46  real :: d(n)
47  pointer(locd, d)
48 #else
49  real, allocatable, dimension(:) :: d
50  integer(LONG_KIND) :: locd
51 #endif
52  integer :: tick, tick0, ticks_per_sec, id
53  integer :: pe, npes, root, i, j, k, l, m, n2, istat
54  integer :: out_unit
55  real :: dt
56 
57  call mpp_init()
58  call mpp_io_init()
59  call mpp_set_stack_size(3145746)
60  pe = mpp_pe()
61  npes = mpp_npes()
62  root = mpp_root_pe()
63 
64  out_unit = stdout()
65  call test_gather(npes,pe,root,out_unit)
66  call test_gatherv(npes,pe,root,out_unit)
67  call test_gather2dv(npes,pe,root,out_unit)
68 
69  if(.false.) then
70 
71  ! first test broadcast
72  call test_broadcast()
73 
74  call system_clock( count_rate=ticks_per_sec )
75 
76  allocate( a(n), b(n) )
77  id = mpp_clock_id( 'Random number' )
78  call mpp_clock_begin(id)
79  call random_number(a)
80  call mpp_clock_end (id)
81  !---------------------------------------------------------------------!
82  ! time transmit, compare against shmem_put and get !
83  !---------------------------------------------------------------------!
84  if( pe.EQ.root )then
85  print *, 'Time mpp_transmit for various lengths...'
86 #ifdef SGICRAY
87  print *, 'For comparison, times for shmem_get and shmem_put are also provided.'
88 #endif
89  print *
90  end if
91  id = mpp_clock_id( 'mpp_transmit' )
92  call mpp_clock_begin(id)
93  !timing is done for cyclical pass (more useful than ping-pong etc)
94  l = n
95  do while( l.GT.0 )
96  !--- mpp_transmit -------------------------------------------------
97  call mpp_sync()
98  call system_clock(tick0)
99  do i = 1,npes
100  call mpp_transmit( put_data=a(1), plen=l, to_pe=modulo(pe+npes-i,npes), &
101  get_data=b(1), glen=l, from_pe=modulo(pe+i,npes) )
102  call mpp_sync_self()
103  ! call mpp_sync_self( (/modulo(pe+npes-i,npes)/) )
104  end do
105  call mpp_sync()
106  call system_clock(tick)
107  dt = real(tick-tick0)/(npes*ticks_per_sec)
108  dt = max( dt, epsilon(dt) )
109  if( pe.EQ.root )write( out_unit,'(/a,i8,f13.6,f8.2)' )'MPP_TRANSMIT length, time, bw(Mb/s)=', l, dt, l*8e-6/dt
110 !#ifdef SGICRAY
111 ! !--- shmem_put ----------------------------------------------------
112 ! call mpp_sync()
113 ! call SYSTEM_CLOCK(tick0)
114 ! do i = 1,npes
115 ! call shmem_real_put( b, a, l, modulo(pe+1,npes) )
116 ! end do
117 ! call mpp_sync()
118 ! call SYSTEM_CLOCK(tick)
119 ! dt = real(tick-tick0)/(npes*ticks_per_sec)
120 ! dt = max( dt, epsilon(dt) )
121 ! if( pe.EQ.root )write( stdout(),'( a,i8,f13.6,f8.2)' )'SHMEM_PUT length, time, bw(Mb/s)=', l, dt, l*8e-6/dt
122 ! !--- shmem_get ----------------------------------------------------
123 ! call mpp_sync()
124 ! call SYSTEM_CLOCK(tick0)
125 ! do i = 1,npes
126 ! call shmem_real_get( b, a, l, modulo(pe+1,npes) )
127 ! end do
128 ! call SYSTEM_CLOCK(tick)
129 ! dt = real(tick-tick0)/(npes*ticks_per_sec)
130 ! dt = max( dt, epsilon(dt) )
131 ! if( pe.EQ.root )write( stdout(),'( a,i8,f13.6,f8.2)' )'SHMEM_GET length, time, bw(Mb/s)=', l, dt, l*8e-6/dt
132 !#endif
133  l = l/2
134  end do
135  !---------------------------------------------------------------------!
136  ! test mpp_sum !
137  !---------------------------------------------------------------------!
138  if( pe.EQ.root )then
139  print '(/a)', 'Time mpp_sum...'
140  end if
141  a = real(pe+1)
142  call mpp_sync()
143  call system_clock(tick0)
144  call mpp_sum(a(1:1000),1000)
145  call system_clock(tick)
146  dt = real(tick-tick0)/ticks_per_sec
147  dt = max( dt, epsilon(dt) )
148  if( pe.EQ.root )write( out_unit,'(a,2i6,f9.1,i8,f13.6,f8.2/)' ) &
149  'mpp_sum: pe, npes, sum(pe+1), length, time, bw(Mb/s)=', pe, npes, a(1), n, dt, n*8e-6/dt
150  call mpp_clock_end(id)
151  !---------------------------------------------------------------------!
152  ! test mpp_max !
153  !---------------------------------------------------------------------!
154  if( pe.EQ.root )then
155  print *
156  print *, 'Test mpp_max...'
157  end if
158  a = real(pe+1)
159  print *, 'pe, pe+1 =', pe, a(1)
160  call mpp_max( a(1) )
161  print *, 'pe, max(pe+1)=', pe, a(1)
162  !pelist check
163  call mpp_sync()
164  call flush(out_unit,istat)
165  if( npes.GE.2 )then
166  if( pe.EQ.root )print *, 'Test of pelists: bcast, sum and max using PEs 0...npes-2 (excluding last PE)'
167  call mpp_declare_pelist( (/(i,i=0,npes-2)/) )
168  a = real(pe+1)
169  if( pe.NE.npes-1 )call mpp_broadcast( a, n, npes-2, (/(i,i=0,npes-2)/) )
170  print *, 'bcast(npes-1) from 0 to npes-2=', pe, a(1)
171  a = real(pe+1)
172  if( pe.NE.npes-1 )then
173  call mpp_set_current_pelist( (/(i,i=0,npes-2)/) )
174  id = mpp_clock_id( 'Partial mpp_sum' )
175  call mpp_clock_begin(id)
176  call mpp_sum( a(1:1000), 1000, (/(i,i=0,npes-2)/) )
177  call mpp_clock_end (id)
178  end if
179  if( pe.EQ.root )print *, 'sum(pe+1) from 0 to npes-2=', a(1)
180  a = real(pe+1)
181  if( pe.NE.npes-1 )call mpp_max( a(1), (/(i,i=0,npes-2)/) )
182  if( pe.EQ.root )print *, 'max(pe+1) from 0 to npes-2=', a(1)
183  end if
184  call mpp_set_current_pelist()
185 
186 #ifdef use_CRI_pointers
187  !---------------------------------------------------------------------!
188  ! test mpp_chksum !
189  !---------------------------------------------------------------------!
190  if( modulo(n,npes).EQ.0 )then !only set up for even division
191  n2 = 1024
192  a = 0.d0
193  if( pe.EQ.root )call random_number(a(1:n2))
194 ! if( pe.EQ.root )call random_number(a)
195  call mpp_sync()
196  call mpp_transmit( put_data=a(1), plen=n2, to_pe=all_pes, &
197  get_data=a(1), glen=n2, from_pe=root )
198  call mpp_sync_self ()
199 ! call mpp_transmit( put_data=a(1), plen=n, to_pe=ALL_PES, &
200 ! get_data=a(1), glen=n, from_pe=root )
201  m= n2/npes
202 ! m= n/npes
203  allocate( c(m) )
204  c = a(pe*m+1:pe*m+m)
205 
206  if( pe.EQ.root )then
207  print *
208  print *, 'Test mpp_chksum...'
209  print *, 'This test shows that a whole array and a distributed array give identical checksums.'
210  end if
211  print *, 'chksum(a(1:1024))=', mpp_chksum(a(1:n2),(/pe/))
212  print *, 'chksum(c(1:1024))=', mpp_chksum(c)
213 ! print *, 'chksum(a)=', mpp_chksum(a,(/pe/))
214 ! print *, 'chksum(c)=', mpp_chksum(c)
215  end if
216 !test of pointer sharing
217 #ifdef use_MPI_GSM
218  call mpp_gsm_malloc( locd, sizeof(d) )
219 #else
220  if( pe.EQ.root )then
221  allocate( d(n) )
222  locd = loc(d)
223  end if
224  call mpp_broadcast(locd,root)
225 #endif
226  if( pe.EQ.root )then
227  call random_number(d)
228  end if
229  call mpp_sync()
230 ! call test_shared_pointers(locd,n)
231 
232 #ifdef use_MPI_GSM
233  call mpp_gsm_free( locd )
234 #else
235  if( pe.EQ.root )then
236  deallocate( d )
237  end if
238 #endif
239 #endif
240  endif ! if(.false.)
241 
242  call mpp_exit()
243 
244 contains
245 
246  !***********************************************
247  !currently only test the mpp_broadcast_char
248  subroutine test_broadcast()
249  integer, parameter :: ARRAYSIZE = 3
250  integer, parameter :: STRINGSIZE = 256
251  character(len=STRINGSIZE), dimension(ARRAYSIZE) :: textA, textB
252  integer :: n
253 
254  texta(1) = "This is line 1 "
255  texta(2) = "Here comes the line 2 "
256  texta(3) = "Finally is line 3 "
257  do n = 1, arraysize
258  textb(n) = texta(n)
259  enddo
260 
261  if(mpp_pe() .NE. mpp_root_pe()) then
262  do n =1, arraysize
263  texta(n) = ""
264  enddo
265  endif
266 
267  !--- comparing textA and textB. textA and textB are supposed to be different on pe other than root_pe
268  if(mpp_pe() == mpp_root_pe()) then
269  do n = 1, arraysize
270  if(texta(n) .NE. textb(n)) call mpp_error(fatal, "test_broadcast: on root_pe, textA should equal textB")
271  enddo
272  else
273  do n = 1, arraysize
274  if(texta(n) == textb(n)) call mpp_error(fatal, "test_broadcast: on root_pe, textA should not equal textB")
275  enddo
276  endif
277  call mpp_broadcast(texta, stringsize, mpp_root_pe())
278  !--- after broadcast, textA and textB should be the same
279  do n = 1, arraysize
280  if(texta(n) .NE. textb(n)) call mpp_error(fatal, "test_broadcast: after broadcast, textA should equal textB")
281  enddo
282 
283  write(out_unit,*) "==> NOTE from test_broadcast: The test is succesful"
284 
285  end subroutine test_broadcast
286 
287  subroutine test_gather(npes,pe,root,out_unit)
288  integer, intent(in) :: npes,pe,root,out_unit
289 
290  integer :: pelist(npes)
291  integer :: i
292  real :: rdata(npes)
293  real :: val
294 
295  if(npes < 3)then
296  write(out_unit,*) "Minimum of 3 ranks required. Not testing gather; too few ranks."
297  return
298  endif
299  write(out_unit,*)
300 
301  val = pe
302  rdata = -1.0
303  do i=1,npes
304  pelist(i) = i-1
305  enddo
306 
307  call mpp_gather((/val/),rdata)
308  if(pe == root)then
309  do i=1,npes
310  if(int(rdata(i)) /= pelist(i))then
311  write(6,*) "Gathered data ",int(rdata(i)), " NE reference ",pelist(i), "at i=",i
312  call mpp_error(fatal, "Test gather uniform vector with global pelist failed")
313  endif
314  enddo
315  endif
316 
317  call mpp_sync()
318  write(out_unit,*) "Test gather uniform vector with global pelist successful"
319 
320  rdata = -1.0
321  if(any(pe == pelist(2:npes)))call mpp_gather((/val/),rdata(2:npes),pelist(2:npes))
322  if(pe == pelist(2))then
323  do i=2,npes
324  if(int(rdata(i)) /= pelist(i))then
325  write(6,*) "Gathered data ",int(rdata(i)), " NE reference ",pelist(i), "at i=",i
326  call mpp_error(fatal, "Test gather uniform vector with reduced pelist failed")
327  endif
328  enddo
329  endif
330  call mpp_sync()
331  write(out_unit,*) "Test gather uniform vector with reduced pelist successful"
332 
333  end subroutine test_gather
334 
335 
336  subroutine test_gatherv(npes,pe,root,out_unit)
337  implicit none
338  integer, intent(in) :: npes,pe,root,out_unit
339 
340  integer :: pelist(npes),rsize(npes)
341  integer :: i,j,k,dsize,ssize
342  real,allocatable :: sdata(:), rdata(:), ref(:)
343 
344  if(npes < 3)then
345  write(out_unit,*) "Minimum of 3 ranks required. Not testing gatherV; too few ranks."
346  return
347  elseif(npes > 9999)then
348  write(out_unit,*) "Maximum of 9999 ranks supported. Not testing gatherV; too many ranks."
349  return
350  endif
351  write(out_unit,*)
352 
353  ssize = pe+1
354  allocate(sdata(ssize))
355  do i=1,ssize
356  sdata(i) = pe + 0.0001*i
357  enddo
358  do i=1,npes
359  pelist(i) = i-1
360  rsize(i) = i
361  enddo
362 
363  dsize = sum(rsize)
364  allocate(rdata(dsize),ref(dsize))
365  rdata = -1.0
366  k=1
367  do j=1,npes
368  do i=1,rsize(j)
369  ref(k) = pelist(j) + 0.0001*i
370  k = k+1
371  enddo;enddo
372 
373  call mpp_gather(sdata,ssize,rdata,rsize)
374 
375  if(pe == root)then
376  k = 1
377  do j=1,npes
378  do i=1,rsize(j)
379  if(rdata(k) /= ref(k))then
380  write(6,*) "Gathered data ",rdata(k), " NE reference ",ref(k), "at k=",k
381  call mpp_error(fatal, "Test gatherV global pelist failed")
382  endif
383  k = k+1
384  enddo;enddo
385  endif
386 
387  call mpp_sync()
388  write(out_unit,*) "Test gatherV with global pelist successful"
389 
390  rdata = -1.0
391  ref(1) = -1.0
392 
393  if(any(pe == pelist(2:npes)))call mpp_gather(sdata,ssize,rdata(2:),rsize(2:),pelist(2:npes))
394 
395  if(pe == pelist(2))then
396  k = 1
397  do j=1,npes
398  do i=1,rsize(j)
399  if(rdata(k) /= ref(k))then
400  write(6,*) "Gathered data ",rdata(k), " NE reference ",ref(k), "at k=",k
401  call mpp_error(fatal, "Test gatherV with reduced pelist failed")
402  endif
403  k = k+1
404  enddo;enddo
405  endif
406  call mpp_sync()
407 
408  write(out_unit,*) "Test gatherV with reduced pelist successful"
409  deallocate(sdata,rdata,ref)
410  end subroutine test_gatherv
411 
412 subroutine test_gather2dv(npes,pe,root,out_unit)
413  implicit none
414  integer, intent(in) :: npes,pe,root,out_unit
415 
416  integer :: pelist(npes),rsize(npes)
417  integer :: pelist2(npes),rsize2(npes)
418  integer :: i,j,k,l,nz,ssize,nelems
419  real,allocatable,dimension(:,:) :: data, cdata, sbuff,rbuff
420  real,allocatable :: ref(:,:)
421  integer, parameter :: KSIZE=10
422 
423  real :: sbuff1D(size(sbuff))
424  real :: rbuff1D(size(rbuff))
425  pointer(sptr,sbuff1d); pointer(rptr,rbuff1d)
426 
427 
428  if(npes < 3)then
429  write(out_unit,*) "Minimum of 3 ranks required. Not testing gather2DV; too few ranks."
430  return
431  elseif(npes > 9999)then
432  write(out_unit,*) "Maximum of 9999 ranks supported. Not testing gather2DV; too many ranks."
433  return
434  endif
435  write(out_unit,*)
436 
437  ssize = pe+1
438  allocate(data(ssize,ksize))
439  do k=1,ksize; do i=1,ssize
440  data(i,k) = 10000.0*k + pe + 0.0001*i
441  enddo; enddo
442  do i=1,npes
443  pelist(i) = i-1
444  rsize(i) = i
445  enddo
446 
447  nz = ksize
448  nelems = sum(rsize(:))
449 
450  allocate(rbuff(nz,nelems)); rbuff = -1.0
451  allocate(ref(nelems,nz),cdata(nelems,nz))
452  ref = 0.0; cdata = 0.0
453  if(pe == root)then
454  do k=1,ksize
455  l=1
456  do j=1,npes
457  do i=1,rsize(j)
458  ref(l,k) = 10000.0*k + pelist(j) + 0.0001*i
459  l = l+1
460  enddo; enddo;enddo
461  endif
462  allocate(sbuff(nz,ssize))
463  ! this matrix inversion makes for easy gather to the IO root
464  ! and a clear, concise unpack
465  do j=1,ssize
466  do i=1,nz
467  sbuff(i,j) = data(j,i)
468  enddo; enddo
469 
470  ! Note that the gatherV implied here is asymmetric; only root needs to know the vector of recv size
471  sptr = loc(sbuff); rptr = loc(rbuff)
472  call mpp_gather(sbuff1d,size(sbuff),rbuff1d,nz*rsize(:))
473 
474  if(pe == root)then
475  do j=1,nz
476  do i=1,nelems
477  cdata(i,j) = rbuff(j,i)
478  enddo; enddo
479  do j=1,nz
480  do i=1,nelems
481  if(cdata(i,j) /= ref(i,j))then
482  write(6,*) "Gathered data ",cdata(i,j), " NE reference ",ref(i,j), "at i,j=",i,j
483  call mpp_error(fatal, "Test gather2DV global pelist failed")
484  endif
485  enddo;enddo
486  endif
487 
488  call mpp_sync()
489  write(out_unit,*) "Test gather2DV with global pelist successful"
490 
491  do i=1,npes
492  pelist2(i) = pelist(npes-i+1)
493  rsize2(i) = rsize(npes-i+1)
494  enddo
495 
496  rbuff = -1.0
497  ref = 0.0; cdata = 0.0
498  if(pe == pelist2(1))then
499  do k=1,ksize
500  l=1
501  do j=1,npes
502  do i=1,rsize2(j)
503  ref(l,k) = 10000.0*k + pelist2(j) + 0.0001*i
504  l = l+1
505  enddo; enddo;enddo
506  endif
507 
508  call mpp_gather(sbuff1d,size(sbuff),rbuff1d,nz*rsize2(:),pelist2)
509 
510  if(pe == pelist2(1))then
511  do j=1,nz
512  do i=1,nelems
513  cdata(i,j) = rbuff(j,i)
514  enddo; enddo
515  do j=1,nz
516  do i=1,nelems
517  if(cdata(i,j) /= ref(i,j))then
518  write(6,*) "Gathered data ",cdata(i,j), " NE reference ",ref(i,j), "at i,j=",i,j
519  call mpp_error(fatal, "Test gather2DV with reversed pelist failed")
520  endif
521  enddo;enddo
522  endif
523  call mpp_sync()
524  write(out_unit,*) "Test gather2DV with reversed pelist successful"
525  deallocate(data,sbuff,rbuff,cdata,ref)
526  end subroutine test_gather2dv
527 
528  subroutine test_shared_pointers(locd,n)
529  integer(LONG_KIND), intent(in) :: locd
530  integer :: n
531  real :: dd(n)
532  pointer( p, dd )
533 
534  p = locd
535  print *, 'TEST_SHARED_POINTERS: pe, locd=', pe, locd
536 ! print *, 'TEST_SHARED_POINTERS: pe, chksum(d)=', pe, mpp_chksum(dd,(/pe/))
537  print *, 'TEST_SHARED_POINTERS: pe, sum(d)=', pe, sum(dd)
538  return
539  end subroutine test_shared_pointers
540 end program test
541 
542 #else
544 end module
545 
546 #endif /* test_mpp */
Definition: mpp.F90:39
only root needs to know the vector of recv size nz do nelems cdata(i, j)
************************************************************************GNU Lesser General Public License **This file is part of the GFDL Flexible Modeling System(FMS). ! *! *FMS is free software without even the implied warranty of MERCHANTABILITY or *FITNESS FOR A PARTICULAR PURPOSE See the GNU General Public License *for more details **You should have received a copy of the GNU Lesser General Public *License along with FMS If see< http:! ***********************************************************************!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! !! MPP_TRANSMIT !! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! subroutine MPP_TRANSMIT_(put_data, put_len, to_pe, get_data, get_len, from_pe, block, tag, recv_request, send_request)!a message-passing routine intended to be reminiscent equally of both MPI and SHMEM!put_data and get_data are contiguous MPP_TYPE_ arrays!at each call, your put_data array is put to to_pe 's get_data! your get_data array is got from from_pe 's put_data!i.e we assume that typically(e.g updating halo regions) each PE performs a put _and_ a get!special PE designations:! NULL_PE:to disable a put or a get(e.g at boundaries)! ANY_PE:if remote PE for the put or get is to be unspecific! ALL_PES:broadcast and collect operations(collect not yet implemented)!ideally we would not pass length, but this f77-style call performs better(arrays passed by address, not descriptor)!further, this permits< length > contiguous words from an array of any rank to be passed(avoiding f90 rank conformance check)!caller is responsible for completion checks(mpp_sync_self) before and after integer, intent(in) ::put_len, to_pe, get_len, from_pe MPP_TYPE_, intent(in) ::put_data(*) MPP_TYPE_, intent(out) ::get_data(*) logical, intent(in), optional ::block integer, intent(in), optional ::tag integer, intent(out), optional ::recv_request, send_request logical ::block_comm integer ::i MPP_TYPE_, allocatable, save ::local_data(:) !local copy used by non-parallel code(no SHMEM or MPI) integer ::comm_tag integer ::rsize if(.NOT.module_is_initialized) call mpp_error(FATAL, 'MPP_TRANSMIT:You must first call mpp_init.') if(to_pe.EQ.NULL_PE .AND. from_pe.EQ.NULL_PE) return block_comm=.true. if(PRESENT(block)) block_comm=block if(debug) then call SYSTEM_CLOCK(tick) write(stdout_unit,'(a, i18, a, i6, a, 2i6, 2i8)')&'T=', tick, ' PE=', pe, ' MPP_TRANSMIT begin:to_pe, from_pe, put_len, get_len=', to_pe, from_pe, put_len, get_len end if comm_tag=DEFAULT_TAG if(present(tag)) comm_tag=tag!do put first and then get if(to_pe.GE.0 .AND. to_pe.LT.npes) then!use non-blocking sends if(debug .and.(current_clock.NE.0)) call SYSTEM_CLOCK(start_tick)!z1l:truly non-blocking send.! if(request(to_pe).NE.MPI_REQUEST_NULL) then !only one message from pe-> to_pe in queue *PE waiting for to_pe ! call error else get_len get_data(i)
#define max(a, b)
Definition: mosaic_util.h:33
************************************************************************GNU Lesser General Public License **This file is part of the GFDL Flexible Modeling System(FMS). ! *! *FMS is free software without even the implied warranty of MERCHANTABILITY or *FITNESS FOR A PARTICULAR PURPOSE See the GNU General Public License *for more details **You should have received a copy of the GNU Lesser General Public *License along with FMS If see< http:! ***********************************************************************subroutine MPP_WRITE_COMPRESSED_1D_(unit, field, domain, data, nelems_io, tstamp, default_data) integer, intent(in) ::unit type(fieldtype), intent(inout) ::field type(domain2D), intent(inout) ::domain MPP_TYPE_, intent(inout) ::data(:) integer, intent(in) ::nelems_io(:) ! number of compressed elements real, intent(in), optional ::tstamp MPP_TYPE_, intent(in), optional ::default_data MPP_TYPE_ ::data2D(size(data, 1), 1) pointer(ptr, data2D) ptr=LOC(data) call mpp_write_compressed(unit, field, domain, data2D, nelems_io, tstamp, default_data) return end subroutine MPP_WRITE_COMPRESSED_1D_ subroutine MPP_WRITE_COMPRESSED_3D_(unit, field, domain, data, nelems_io, tstamp, default_data) integer, intent(in) ::unit type(fieldtype), intent(inout) ::field type(domain2D), intent(inout) ::domain MPP_TYPE_, intent(inout) ::data(:,:,:) integer, intent(in) ::nelems_io(:) ! number of compressed elements real, intent(in), optional ::tstamp MPP_TYPE_, intent(in), optional ::default_data MPP_TYPE_ ::data2D(size(data, 1), size(data, 2) *size(data, 3)) pointer(ptr, data2D) ptr=LOC(data) call mpp_write_compressed(unit, field, domain, data2D, nelems_io, tstamp, default_data) return end subroutine MPP_WRITE_COMPRESSED_3D_ subroutine MPP_WRITE_COMPRESSED_2D_(unit, field, domain, data, nelems_io, tstamp, default_data) integer, intent(in) ::unit type(fieldtype), intent(inout) ::field type(domain2D), intent(inout) ::domain MPP_TYPE_, intent(inout) ::data(:,:) integer, intent(in) ::nelems_io(:) ! number of compressed elements from each ! member of the io_domain. It MUST have the ! same order as the io_domain pelist. real, intent(in), optional ::tstamp MPP_TYPE_, intent(in), optional ::default_data!cdata is used to store the io-domain compressed data MPP_TYPE_, allocatable, dimension(:,:) ::cdata MPP_TYPE_, allocatable, dimension(:,:) ::sbuff, rbuff MPP_TYPE_ ::fill MPP_TYPE_ ::sbuff1D(size(data)) MPP_TYPE_ ::rbuff1D(size(data, 2) *sum(nelems_io(:))) pointer(sptr, sbuff1D);pointer(rptr, rbuff1D) integer, allocatable ::pelist(:) integer, allocatable ::nz_gather(:) integer ::i, j, nz, nelems, mynelems, idx, npes type(domain2d), pointer ::io_domain=> pelist concise unpack do mynelems do nz sbuff(i, j)