Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Kill trailing whitespaces in teshsuite/smpi/{isp,mpich3-test}.
[simgrid.git] / teshsuite / smpi / mpich3-test / f77 / topo / dgraph_unwgtf.f
1 C -*- Mode: Fortran; -*-
2 C
3 C  (C) 2011 by Argonne National Laboratory.
4 C      See COPYRIGHT in top-level directory.
5 C
6 C     This program is Fortran version of dgraph_unwgt.c
7 C     Specify a distributed graph of a bidirectional ring of the MPI_COMM_WORLD,
8 C     i.e. everyone only talks to left and right neighbors.
9
10       logical function validate_dgraph(dgraph_comm)
11       implicit none
12       include 'mpif.h'
13
14       integer     dgraph_comm
15       integer     comm_topo
16       integer     src_sz, dest_sz
17       integer     ierr;
18       logical     wgt_flag;
19       integer     srcs(2), dests(2)
20
21       integer     world_rank, world_size;
22       integer     idx, nbr_sep
23
24       comm_topo = MPI_UNDEFINED
25       call MPI_Topo_test(dgraph_comm, comm_topo, ierr);
26       if (comm_topo .ne. MPI_DIST_GRAPH) then
27           validate_dgraph = .false.
28           write(6,*) "dgraph_comm is NOT of type MPI_DIST_GRAPH."
29           return
30       endif
31
32       call MPI_Dist_graph_neighbors_count(dgraph_comm,
33      &                                    src_sz, dest_sz, wgt_flag,
34      &                                    ierr)
35       if (ierr .ne. MPI_SUCCESS) then
36           validate_dgraph = .false.
37           write(6,*) "MPI_Dist_graph_neighbors_count() fails!"
38           return
39       endif
40       if (wgt_flag) then
41           validate_dgraph = .false.
42           write(6,*) "dgraph_comm is NOT created with MPI_UNWEIGHTED."
43           return
44       endif
45
46       if (src_sz .ne. 2 .or. dest_sz .ne. 2) then
47           validate_dgraph = .false.
48           write(6,*) "source or destination edge array is not size 2."
49           write(6,"('src_sz = ',I3,', dest_sz = ',I3)") src_sz, dest_sz
50           return
51       endif
52
53       call MPI_Dist_graph_neighbors(dgraph_comm,
54      &                              src_sz, srcs, MPI_UNWEIGHTED,
55      &                              dest_sz, dests, MPI_UNWEIGHTED,
56      &                              ierr)
57       if (ierr .ne. MPI_SUCCESS) then
58           validate_dgraph = .false.
59           write(6,*) "MPI_Dist_graph_neighbors() fails!"
60           return
61       endif
62
63 C     Check if the neighbors returned from MPI are really
64 C     the nearest neighbors that within a ring.
65       call MPI_Comm_rank(MPI_COMM_WORLD, world_rank, ierr)
66       call MPI_Comm_size(MPI_COMM_WORLD, world_size, ierr)
67
68       do idx = 1, src_sz
69           nbr_sep = iabs(srcs(idx) - world_rank)
70           if (nbr_sep .ne. 1 .and. nbr_sep .ne. (world_size-1)) then
71               validate_dgraph = .false.
72               write(6,"('srcs[',I3,']=',I3,
73      &                  ' is NOT a neighbor of my rank',I3)")
74      &              idx, srcs(idx), world_rank
75               return
76           endif
77       enddo
78       do idx = 1, dest_sz
79           nbr_sep = iabs(dests(idx) - world_rank)
80           if (nbr_sep .ne. 1 .and. nbr_sep .ne. (world_size-1)) then
81               validate_dgraph = .false.
82               write(6,"('dests[',I3,']=',I3,
83      &                  ' is NOT a neighbor of my rank',I3)")
84      &              idx, dests(idx), world_rank
85               return
86           endif
87       enddo
88
89       validate_dgraph = .true.
90       return
91       end
92
93       integer function ring_rank(world_size, in_rank)
94       implicit none
95       integer world_size, in_rank
96       if (in_rank .ge. 0 .and. in_rank .lt. world_size) then
97           ring_rank = in_rank
98           return
99       endif
100       if (in_rank .lt. 0 ) then
101           ring_rank = in_rank + world_size
102           return
103       endif
104       if (in_rank .ge. world_size) then
105           ring_rank = in_rank - world_size
106           return
107       endif
108       ring_rank = -99999
109       return
110       end
111
112
113
114       program dgraph_unwgt
115       implicit none
116       include 'mpif.h'
117
118       integer    ring_rank
119       external   ring_rank
120       logical    validate_dgraph
121       external   validate_dgraph
122       integer    errs, ierr
123
124       integer    dgraph_comm
125       integer    world_size, world_rank
126       integer    src_sz, dest_sz
127       integer    degs(1)
128       integer    srcs(2), dests(2)
129
130       errs = 0
131       call MTEST_Init(ierr)
132       call MPI_Comm_rank(MPI_COMM_WORLD, world_rank, ierr)
133       call MPI_Comm_size(MPI_COMM_WORLD, world_size, ierr)
134
135       srcs(1) = world_rank
136       degs(1) = 2;
137       dests(1) = ring_rank(world_size, world_rank-1)
138       dests(2) = ring_rank(world_size, world_rank+1)
139       call MPI_Dist_graph_create(MPI_COMM_WORLD, 1, srcs, degs, dests,
140      &                           MPI_UNWEIGHTED, MPI_INFO_NULL,
141      &                          .true., dgraph_comm, ierr)
142       if (ierr .ne. MPI_SUCCESS) then
143           write(6,*) "MPI_Dist_graph_create() fails!"
144           call MPI_Abort(MPI_COMM_WORLD, 1, ierr)
145           stop
146       endif
147       if (.not. validate_dgraph(dgraph_comm)) then
148           write(6,*) "MPI_Dist_graph_create() does not create"
149      &               //"a bidirectional ring graph!"
150           call MPI_Abort(MPI_COMM_WORLD, 1, ierr)
151           stop
152       endif
153       call MPI_Comm_free(dgraph_comm, ierr)
154
155 C now create one with MPI_WEIGHTS_EMPTY
156 C NOTE that MPI_WEIGHTS_EMPTY was added in MPI-3 and does not
157 C appear before then.  Including this test means that this test cannot
158 C be compiled if the MPI version is less than 3 (see the testlist file)
159
160       degs(1) = 0;
161       call MPI_Dist_graph_create(MPI_COMM_WORLD, 1, srcs, degs, dests,
162      &                           MPI_WEIGHTS_EMPTY, MPI_INFO_NULL,
163      &                          .true., dgraph_comm, ierr)
164       if (ierr .ne. MPI_SUCCESS) then
165           write(6,*) "MPI_Dist_graph_create() fails!"
166           call MPI_Abort(MPI_COMM_WORLD, 1, ierr)
167           stop
168       endif
169       call MPI_Comm_free(dgraph_comm, ierr)
170
171       src_sz   = 2
172       srcs(1)  = ring_rank(world_size, world_rank-1)
173       srcs(2)  = ring_rank(world_size, world_rank+1)
174       dest_sz  = 2
175       dests(1) = ring_rank(world_size, world_rank-1)
176       dests(2) = ring_rank(world_size, world_rank+1)
177       call MPI_Dist_graph_create_adjacent(MPI_COMM_WORLD,
178      &                                    src_sz, srcs,
179      &                                    MPI_UNWEIGHTED,
180      &                                    dest_sz, dests,
181      &                                    MPI_UNWEIGHTED,
182      &                                    MPI_INFO_NULL, .true.,
183      &                                    dgraph_comm, ierr)
184       if (ierr .ne. MPI_SUCCESS) then
185           write(6,*) "MPI_Dist_graph_create_adjacent() fails!"
186           call MPI_Abort(MPI_COMM_WORLD, 1, ierr)
187           stop
188       endif
189       if (.not. validate_dgraph(dgraph_comm)) then
190           write(6,*) "MPI_Dist_graph_create_adjacent() does not create"
191      &               //"a bidirectional ring graph!"
192           call MPI_Abort(MPI_COMM_WORLD, 1, ierr)
193           stop
194       endif
195       call MPI_Comm_free(dgraph_comm, ierr)
196
197 C now create one with MPI_WEIGHTS_EMPTY
198       src_sz   = 0
199       dest_sz  = 0
200       call MPI_Dist_graph_create_adjacent(MPI_COMM_WORLD,
201      &                                    src_sz, srcs,
202      &                                    MPI_WEIGHTS_EMPTY,
203      &                                    dest_sz, dests,
204      &                                    MPI_WEIGHTS_EMPTY,
205      &                                    MPI_INFO_NULL, .true.,
206      &                                    dgraph_comm, ierr)
207       if (ierr .ne. MPI_SUCCESS) then
208           write(6,*) "MPI_Dist_graph_create_adjacent() fails!"
209           call MPI_Abort(MPI_COMM_WORLD, 1, ierr)
210           stop
211       endif
212       call MPI_Comm_free(dgraph_comm, ierr)
213
214       call MTEST_Finalize(errs)
215       call MPI_Finalize(ierr)
216       end