source: CPL/oasis3-mct_5.0/lib/mct/mpi-serial/tests/ftest_internal.F90 @ 6328

Last change on this file since 6328 was 6328, checked in by aclsce, 17 months ago

First import of oasis3-mct_5.0 (from oasis git server, branch OASIS3-MCT_5.0)

File size: 8.9 KB
Line 
1      program test
2      use mpi
3      implicit none
4
5        call test_contiguous()
6        call test_vector()
7        call test_simple_hvector()
8        call test_simple_indexed()
9        call test_simple_bindexed()
10        call test_simple_hindexed()
11        call test_packed()
12        call test_multiple()
13        stop
14      end
15
16!!!!!!!!!!!!!!!!!!!
17! Contiguous type.  Simplest example.  Strings 5
18!  integers together and tests their equality after
19!  a send operation
20!!!!!!!!!!!!!!!!!!!
21
22      subroutine test_contiguous()
23      use mpi
24      integer ierr
25      integer datatype
26      integer a(5)
27      integer b(5)
28      integer i
29      data a/1,2,3,4,5/
30      data b/5 * 0/
31
32      print *, "Test Contiguous of 5 x MPI_INTEGER"
33      call mpi_type_contiguous(5, mpi_integer, datatype,ierr)
34
35      call mpi_type_commit(datatype, ierr)
36
37      call print_typemap(datatype,ierr)
38      call copy_data2(a,1,datatype, b,1,datatype, ierr)
39
40      do i=1,5
41        if (a(i) .ne. b(i)) then
42          print *,">>>FAILED: mpi_type_contiguous"
43          stop
44        end if
45      end do
46      print *, ">>>PASSED: mpi_type_contiguous"
47      end
48
49!!!!!!!!!!!!!!!!!!!!!!!!
50! Vector type.  collect a series of indices with
51! set stride from an array.
52!!!!!!!!!!!!!!!!!!!!!!!!
53
54      subroutine test_vector()
55      use mpi
56      integer ierr
57      integer datatype
58      integer a(10) != (1,2,3,4,5,6,7,8,9,0)
59      integer b(10)
60      integer check_index(6)
61      data a/1,2,3,4,5,6,7,8,9,10/
62      data b/10 * 0/
63      data check_index/1,2,4,5,7,8/
64      integer i
65
66      print *, "Test vector of MPI_INTEGER"
67
68      call mpi_type_vector(3, 2, 3, mpi_integer, datatype, ierr)
69      call mpi_type_commit(datatype, ierr)
70      call print_typemap(datatype,ierr)
71      call copy_data2(a,1,datatype,b,1,datatype,ierr)
72
73      do i=1,6
74        if (a(check_index(i)) .ne. b(check_index(i))) then
75          print *,">>>FAILED: mpi_type_vector"
76          stop
77        end if
78      end do
79      print *, ">>>PASSED: mpi_type_vector"
80      end
81
82!!!!!!!!!!!!!!!!!!!!!
83! Byte-addressed vector.
84! values calculated with mpi_type_extent(),
85! so basically we are doing the work here in the
86! test program instead of in the library
87!!!!!!!!!!!!!!!!!!!!!
88
89      subroutine test_simple_hvector()
90      use mpi
91        integer vector_type
92        integer (kind=mpi_address_kind) extent
93        integer i
94        integer a(10)
95        integer b(10)
96        integer index_test(6)
97        integer ierr
98
99        data a/1,2,3,4,5,6,7,8,9,10/, b/0,0,0,0,0,0,0,0,0,0/
100        data index_test/1,2,5,6,9,10/
101
102        print *, "Vector type of 3 groups of 2 MPI_INTEGER"
103        print *, "Stride of 4 (in bytes)"
104
105        call mpi_type_extent(mpi_integer, extent, ierr)
106        call mpi_type_hvector(3, 2, 4 * extent, mpi_integer, &
107             vector_type, ierr)
108        call mpi_type_commit(vector_type, ierr)
109        call print_typemap(vector_type,ierr)
110        call copy_data2(a,1,vector_type, b,1,vector_type,ierr)
111
112        do i=1,7
113          if (a(index_test(i)) .ne. (b(index_test(i)))) then
114            print *, ">>>FAILED: test_simple_hvector"
115            stop
116          end if
117        end do
118      print *, ">>>PASSED: test_simple_hvector"
119      end subroutine
120
121!!!!!!!!!!!!!!!!!!!!
122! indexed type.  test certain indices of an array
123!!!!!!!!!!!!!!!!!!!!
124
125      subroutine test_simple_indexed()
126      use mpi
127        integer i
128        complex a(15)
129        complex b(15)
130        integer index_test(6)
131        integer blens(3)
132        integer disps(3)
133        integer indexed_type
134        integer ierr
135
136        data a/1,2,3,4,5,6,7,8,9,10,11,12,13,14,15/
137        data b/0,0,0,0,0,0,0,0,0,0,0,0,0,0,0/
138        data index_test/1,6,7,11,12,13/
139        data blens/2,1,3/
140        data disps/5,0,10/
141        print *, "Indexed type"
142
143        call mpi_type_indexed(3, blens, disps, mpi_complex, &
144                         indexed_type, ierr)
145        call mpi_type_commit(indexed_type, ierr)
146        call print_typemap(indexed_type, ierr)
147        call copy_data2(a,1,indexed_type, b,1,indexed_type,ierr)
148
149        do i=1,6
150          if (a(index_test(i)) .ne. b(index_test(i))) then
151            print *, ">>>FAILED: test_simple_indexed"
152            stop
153          end if
154        end do
155        print *, ">>>PASSED: test_simple_indexed"
156      end subroutine
157
158!!!!!!!!!!!!!!!!
159! Block indexed.  All blocks have same length
160!!!!!!!!!!!!!!!!
161
162      subroutine test_simple_bindexed()
163      use mpi
164        integer i
165        integer disps(3)
166        integer a(10), b(10)
167        integer index_test(6)
168        integer indexed_type
169        integer ierr
170
171        data disps/0,4,7/
172        data a/1,2,3,4,5,6,7,8,9,10/
173        data b/0,0,0,0,0,0,0,0,0,0/
174        data index_test/1,2,5,6,8,9/
175        print *, "Block indexed type"
176
177        call mpi_type_indexed_block(3,2,disps,mpi_integer, &
178                        indexed_type, ierr)
179
180        call mpi_type_commit(indexed_type, ierr)
181        call print_typemap(indexed_type, ierr)
182        call copy_data2(a,1,indexed_type, b,1,indexed_type, ierr)
183
184        do i=1,6
185          if (a(index_test(i)) .ne. b(index_test(i))) then
186            print *, ">>>FAILED: test_simple_bindexed"
187            stop
188          end if
189        end do
190        print *, ">>>PASSED: test_simple_bindexed"
191      end subroutine
192
193!!!!!!!!!!!!!!!!
194! test_simple_indexed
195!  test equality of a byte-addressed
196!  type of integer array
197!  (disps calculated through mpi_type_extent()
198!!!!!!!!!!!!!!!
199      subroutine test_simple_hindexed()
200      use mpi
201        integer i
202        integer a(10), b(10)
203        integer index_test(6)
204        integer blens(3)
205        integer*8 disps(3)
206        integer indexed_type
207        integer*8 extent
208        integer ierr
209
210        data a/1,2,3,4,5,6,7,8,9,10/
211        data b/0,0,0,0,0,0,0,0,0,0/
212        data index_test/1,3,4,6,7,8/
213        data blens/2,1,3/
214
215        call mpi_type_extent(mpi_integer, extent, ierr)
216        disps(1) = 2*extent
217        disps(2) = 0
218        disps(3) = 5*extent
219
220
221        print *, "Byte addressed indexed type"
222        call mpi_type_hindexed(3,blens,disps, MPI_INTEGER, &
223             indexed_type,ierr)
224        call mpi_type_commit(indexed_type, ierr)
225        call print_typemap(indexed_type, ierr)
226        call copy_data2(a,1,indexed_type, b,1,indexed_type, ierr)
227
228        do i=1,6
229          if (a(index_test(i)) .ne. b(index_test(i))) then
230            print *, ">>>FAILED: test_simple_hindexed"
231            stop
232          end if
233        end do
234        print *, ">>>PASSED: test_simple_hindexed"
235      end subroutine
236
237!!!!!!!!!!!!!!!!!!!!!!!!!!!!
238! test_packed()
239! Creates a few variable pairs, assigns the first
240!  of each pair, then packs their values and unpacks
241!  them to the other set.
242!!!!!!!!!!!!!!!!!!!!!!!!!!!!
243      subroutine test_packed()
244      use mpi
245        integer size
246        integer x, y
247        real f, g
248        complex c, d
249        character*5 a, b
250        character buf(100), rbuf(100)
251        integer blens(3)
252        integer(kind=mpi_address_kind) disps(3)
253        integer pos
254
255
256        x = 10
257        f = 14.333
258        c = (100, 20)
259        a = "xyzab"
260
261        pos = 0
262        data blens/1,2,1/, disps/0,4,8/
263
264        print *, "Packed type "
265
266        call mpi_pack(x, 1, mpi_integer, buf, 100, pos, 0, ierr)
267        call mpi_pack(f, 1, mpi_real, buf, 100, pos, 0, ierr)
268        call mpi_pack(c, 1, mpi_complex, buf, 100, pos, 0, ierr)
269        call mpi_pack(a, 5, mpi_character, buf, 100, pos, 0, ierr)
270
271        call copy_data2(buf, pos, mpi_packed, rbuf, pos, &
272                        mpi_packed, ierr)
273
274        pos = 0;
275
276        call mpi_unpack(rbuf, 100, pos, y, 1, mpi_integer, 0, ierr)
277        call mpi_unpack(rbuf, 100, pos, g, 1, mpi_real, 0, ierr)
278        call mpi_unpack(rbuf, 100, pos, d, 1, mpi_complex, 0, ierr)
279        call mpi_unpack(rbuf, 100, pos, b, 5, mpi_character, &
280                        0, ierr)
281
282        if (x .ne. y .OR. f .ne. g &
283            .OR. c .ne. d .OR. a .ne. b) &
284            then
285          print *, ">>>FAILED: mpi_pack"
286          stop
287        end if
288
289        print *, ">>>PASSED: mpi_pack"
290
291      end subroutine
292
293!!!!!!!!!!!!!!!!!!!!!!!!!
294! Test an indexed send with a multiple receive
295!!!!!!!!!!!!!!!!!!!!!!!!!
296
297      subroutine test_multiple()
298      use mpi
299        integer i
300        complex a(15)
301        complex b(15)
302        integer index_test(6)
303        integer blens(3)
304        integer disps(3)
305        integer indexed_type
306        integer ierr
307
308        data a/1,2,3,4,5,6,7,8,9,10,11,12,13,14,15/
309        data b/0,0,0,0,0,0,0,0,0,0,0,0,0,0,0/
310        data index_test/1,6,7,11,12,13/
311        data blens/1,2,3/
312        data disps/0,5,10/
313        print *, "Indexed type"
314
315        call mpi_type_indexed(3, blens, disps, mpi_complex, &
316                         indexed_type, ierr)
317        call mpi_type_commit(indexed_type, ierr)
318        call copy_data2(a,1,indexed_type, b,6, mpi_complex, ierr)
319
320        do i=1,6
321          if (a(index_test(i)) .ne. b(i)) then
322            print *, ">>>FAILED: test_multiple"
323            stop
324          end if
325        end do
326        print *, ">>>PASSED: test_multiple"
327      end subroutine
328
Note: See TracBrowser for help on using the repository browser.