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_wgtf.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       integer     src_wgts(2), dest_wgts(2)
21
22       integer     world_rank, world_size;
23       integer     idx, nbr_sep
24
25       comm_topo = MPI_UNDEFINED
26       call MPI_Topo_test(dgraph_comm, comm_topo, ierr);
27       if (comm_topo .ne. MPI_DIST_GRAPH) then
28           validate_dgraph = .false.
29           write(6,*) "dgraph_comm is NOT of type MPI_DIST_GRAPH."
30           return
31       endif
32
33       call MPI_Dist_graph_neighbors_count(dgraph_comm,
34      &                                    src_sz, dest_sz, wgt_flag,
35      &                                    ierr)
36       if (ierr .ne. MPI_SUCCESS) then
37           validate_dgraph = .false.
38           write(6,*) "MPI_Dist_graph_neighbors_count() fails!"
39           return
40       endif
41       if (.not. wgt_flag) then
42           validate_dgraph = .false.
43           write(6,*) "dgraph_comm is created with MPI_UNWEIGHTED!"
44           return
45       endif
46
47       if (src_sz .ne. 2 .or. dest_sz .ne. 2) then
48           validate_dgraph = .false.
49           write(6,*) "source or destination edge array is not size 2."
50           write(6,"('src_sz = ',I3,', dest_sz = ',I3)") src_sz, dest_sz
51           return
52       endif
53
54       call MPI_Dist_graph_neighbors(dgraph_comm,
55      &                              src_sz, srcs, src_wgts,
56      &                              dest_sz, dests, dest_wgts,
57      &                              ierr)
58       if (ierr .ne. MPI_SUCCESS) then
59           validate_dgraph = .false.
60           write(6,*) "MPI_Dist_graph_neighbors() fails!"
61           return
62       endif
63
64 C     Check if the neighbors returned from MPI are really
65 C     the nearest neighbors that within a ring.
66       call MPI_Comm_rank(MPI_COMM_WORLD, world_rank, ierr)
67       call MPI_Comm_size(MPI_COMM_WORLD, world_size, ierr)
68
69       do idx = 1, src_sz
70           nbr_sep = iabs(srcs(idx) - world_rank)
71           if (nbr_sep .ne. 1 .and. nbr_sep .ne. (world_size-1)) then
72               validate_dgraph = .false.
73               write(6,"('srcs[',I3,']=',I3,
74      &                  ' is NOT a neighbor of my rank',I3)")
75      &              idx, srcs(idx), world_rank
76               return
77           endif
78       enddo
79       if (src_wgts(1) .ne. src_wgts(2)) then
80           validate_dgraph = .false.
81           write(6,"('src_wgts = [',I3,',',I3,']')")
82      &          src_wgts(1), src_wgts(2)
83           return
84       endif
85       do idx = 1, dest_sz
86           nbr_sep = iabs(dests(idx) - world_rank)
87           if (nbr_sep .ne. 1 .and. nbr_sep .ne. (world_size-1)) then
88               validate_dgraph = .false.
89               write(6,"('dests[',I3,']=',I3,
90      &                  ' is NOT a neighbor of my rank',I3)")
91      &              idx, dests(idx), world_rank
92               return
93           endif
94       enddo
95       if (dest_wgts(1) .ne. dest_wgts(2)) then
96           validate_dgraph = .false.
97           write(6,"('dest_wgts = [',I3,',',I3,']')")
98      &          dest_wgts(1), dest_wgts(2)
99           return
100       endif
101
102       validate_dgraph = .true.
103       return
104       end
105
106       integer function ring_rank(world_size, in_rank)
107       implicit none
108       integer world_size, in_rank
109       if (in_rank .ge. 0 .and. in_rank .lt. world_size) then
110           ring_rank = in_rank
111           return
112       endif
113       if (in_rank .lt. 0 ) then
114           ring_rank = in_rank + world_size
115           return
116       endif
117       if (in_rank .ge. world_size) then
118           ring_rank = in_rank - world_size
119           return
120       endif
121       ring_rank = -99999
122       return
123       end
124
125
126
127       program dgraph_unwgt
128       implicit none
129       include 'mpif.h'
130
131       integer    ring_rank
132       external   ring_rank
133       logical    validate_dgraph
134       external   validate_dgraph
135       integer    errs, ierr
136
137       integer    dgraph_comm
138       integer    world_size, world_rank
139       integer    src_sz, dest_sz
140       integer    degs(1)
141       integer    srcs(2), dests(2)
142       integer    src_wgts(2), dest_wgts(2)
143
144       errs = 0
145       call MTEST_Init(ierr)
146       call MPI_Comm_rank(MPI_COMM_WORLD, world_rank, ierr)
147       call MPI_Comm_size(MPI_COMM_WORLD, world_size, ierr)
148
149       srcs(1)      = world_rank
150       degs(1)      = 2;
151       dests(1)     = ring_rank(world_size, world_rank-1)
152       dests(2)     = ring_rank(world_size, world_rank+1)
153       dest_wgts(1) = 1
154       dest_wgts(2) = 1
155       call MPI_Dist_graph_create(MPI_COMM_WORLD, 1, srcs, degs, dests,
156      &                           dest_wgts, MPI_INFO_NULL,
157      &                          .true., dgraph_comm, ierr)
158       if (ierr .ne. MPI_SUCCESS) then
159           write(6,*) "MPI_Dist_graph_create() fails!"
160           call MPI_Abort(MPI_COMM_WORLD, 1, ierr)
161           stop
162       endif
163       if (.not. validate_dgraph(dgraph_comm)) then
164           write(6,*) "MPI_Dist_graph_create() does not create "
165      &               //"a bidirectional ring graph!"
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       src_wgts(1)  = 1
175       src_wgts(2)  = 1
176       dest_sz      = 2
177       dests(1)     = ring_rank(world_size, world_rank-1)
178       dests(2)     = ring_rank(world_size, world_rank+1)
179       dest_wgts(1) = 1
180       dest_wgts(2) = 1
181       call MPI_Dist_graph_create_adjacent(MPI_COMM_WORLD,
182      &                                    src_sz, srcs, src_wgts,
183      &                                    dest_sz, dests, dest_wgts,
184      &                                    MPI_INFO_NULL, .true.,
185      &                                    dgraph_comm, ierr)
186       if (ierr .ne. MPI_SUCCESS) then
187           write(6,*) "MPI_Dist_graph_create_adjacent() fails!"
188           call MPI_Abort(MPI_COMM_WORLD, 1, ierr)
189           stop
190       endif
191       if (.not. validate_dgraph(dgraph_comm)) then
192           write(6,*) "MPI_Dist_graph_create_adjacent() does not create "
193      &               //"a bidirectional ring graph!"
194           call MPI_Abort(MPI_COMM_WORLD, 1, ierr)
195           stop
196       endif
197       call MPI_Comm_free(dgraph_comm, ierr)
198
199       call MTEST_Finalize(errs)
200       call MPI_Finalize(ierr)
201       end