FV3 Bundle
test_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 #ifdef test_mpp_pset
20 #include <fms_platform.h>
21 program test
22  use mpp_mod, only: mpp_init, mpp_exit, mpp_pe, mpp_npes, stderr, stdout, &
23  mpp_clock_id, mpp_clock_begin, mpp_clock_end
27  implicit none
28 !test program demonstrates how to create PSETs
29 ! how to distribute allocatable arrays
30 ! how to distribute automatic arrays
31  integer, parameter :: n=96 !divisible by lots of numbers
32  real, allocatable, dimension(:,:,:) :: a, b, cc
33  real :: c(n,n,n)
34 #ifdef use_CRI_pointers
35  pointer( ptr_c, c )
36 #endif
37  integer(POINTER_KIND) :: ptr !useless declaration, but it will compile
38  integer :: i, j, k, ks, ke
39 !MPP
40  integer :: pe, npes
41 !MPP_PSET
42  type(mpp_pset_type) :: pset
43  logical :: root
44 !clocks
45  integer :: id_full, id_alloc, id_auto
46  integer :: out_unit, errunit
47 
48  call mpp_init()
49  pe = mpp_pe()
50  npes = mpp_npes()
51  out_unit = stdout()
52  errunit = stderr()
53  write( out_unit,'(a,i6)' )'Starting MPP_PSET unit test, npes=', npes
54  call mpp_pset_create( npes, pset )
55  root = mpp_pset_root(pset)
56  id_full = mpp_clock_id( 'Full array' )
57  id_alloc = mpp_clock_id( 'Allocatable array, PSETs' )
58  id_auto = mpp_clock_id( 'Automatic array, PSETs' )
59 !allocate a and b
60  allocate( a(n,n,n) )
61  allocate( b(n,n,n) )
62 !allocate shared array c
63  if( root )then
64  allocate( cc(n,n,n) )
65 #ifdef use_CRI_pointers
66  ptr = loc(cc)
67 #endif
68  end if
69  call mpp_pset_broadcast_ptr( pset, ptr )
70 #ifdef use_CRI_pointers
71  ptr_c = ptr
72 #endif
73 !initialize a and b
74  call random_number(a)
75  call mpp_clock_begin(id_full)
76  do k = 1,n
77  do j = 1,n
78  do i = 1,n
79  b(i,j,k) = 2*a(i,j,k)
80  end do
81  end do
82  end do
83  call mpp_clock_end(id_full)
84 !divide up among PSETs
85  call mpp_pset_segment_array( pset, 1, n, ks, ke )
86  write( errunit,'(a,4i6)' )'pe, n, ks, ke=', pe, n, ks, ke
87  call mpp_clock_begin(id_alloc)
88  do k = ks,ke
89  do j = 1,n
90  do i = 1,n
91  c(i,j,k) = 2*a(i,j,k)
92  end do
93  end do
94  end do
95  call mpp_pset_sync(pset)
96  call mpp_clock_end(id_alloc)
97  write( errunit,'(a,i6,2es23.15)' )'b, c should be equal: pe b c=', &
98  pe, sum(b), sum(c)
99  call mpp_pset_print_chksum( pset, 'test_alloc', c(:,:,ks:ke) )
100  call test_auto(n)
101  call mpp_pset_delete(pset)
102  call mpp_exit()
103 
104 contains
105 
106  subroutine test_auto(m)
107 !same test as above, on auto array d
108 !this is how you create shared auto arrays
109  integer, intent(in) :: m
110  real :: d(m,m,m)
111  integer :: js, je
112 #ifdef use_CRI_pointers
113  pointer( pd, d )
114  call mpp_pset_stack_push( pset, pd, size(d) )
115 #endif
116  call mpp_pset_segment_array( pset, 1, m, js, je )
117  call mpp_clock_begin(id_auto)
118  do k = 1,m
119  do j = js,je
120  do i = 1,m
121  d(i,j,k) = 2*a(i,j,k)
122  end do
123  end do
124  end do
125  call mpp_pset_sync(pset)
126  call mpp_clock_end(id_auto)
127  write( errunit,'(a,i6,2es23.15)' )'b, d should be equal: pe b d=', &
128  pe, sum(b), sum(d)
129  call mpp_pset_print_chksum( pset, 'test_auto ', d(:,js:je,:) )
130  end subroutine test_auto
131 
132 end program test
133 #else
135 end module
136 #endif
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
Definition: mpp.F90:39
subroutine, public mpp_pset_create(npset, pset, stacksize, pelist, commID)
Definition: mpp_pset.F90:125
subroutine, public mpp_pset_delete(pset)
Definition: mpp_pset.F90:236
logical function, public mpp_pset_root(pset)
Definition: mpp_pset.F90:601
subroutine, public mpp_pset_sync(pset)
Definition: mpp_pset.F90:332