Salome HOME
Initiating medtool
[modules/med.git] / src / ParaMEDMEMTest / test_MPI_Access_Time.cxx
1 // Copyright (C) 2007-2015  CEA/DEN, EDF R&D
2 //
3 // This library is free software; you can redistribute it and/or
4 // modify it under the terms of the GNU Lesser General Public
5 // License as published by the Free Software Foundation; either
6 // version 2.1 of the License, or (at your option) any later version.
7 //
8 // This library is distributed in the hope that it will be useful,
9 // but WITHOUT ANY WARRANTY; without even the implied warranty of
10 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
11 // Lesser General Public License for more details.
12 //
13 // You should have received a copy of the GNU Lesser General Public
14 // License along with this library; if not, write to the Free Software
15 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA
16 //
17 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
18 //
19
20 #include <string>
21 #include <vector>
22 #include <map>
23 #include <iostream>
24 #include <mpi.h>
25
26 #include "MPIAccessTest.hxx"
27 #include <cppunit/TestAssert.h>
28
29 //#include "CommInterface.hxx"
30 //#include "ProcessorGroup.hxx"
31 //#include "MPIProcessorGroup.hxx"
32 #include "MPIAccess.hxx"
33
34 // use this define to enable lines, execution of which leads to Segmentation Fault
35 #define ENABLE_FAULTS
36
37 // use this define to enable CPPUNIT asserts and fails, showing bugs
38 #define ENABLE_FORCED_FAILURES
39
40 using namespace std;
41 using namespace ParaMEDMEM;
42
43 void MPIAccessTest::test_MPI_Access_Time() {
44
45   debugStream << "test_MPI_Access_Time" << endl ;
46
47   //  MPI_Init(&argc, &argv) ; 
48
49   int size ;
50   int myrank ;
51   MPI_Comm_size(MPI_COMM_WORLD,&size) ;
52   MPI_Comm_rank(MPI_COMM_WORLD,&myrank) ;
53
54   if ( size < 2 ) {
55     ostringstream strstream ;
56     strstream << "test_MPI_Access_Time must be runned with 2 procs" << endl ;
57     cerr << strstream.str() << endl ;
58     //CPPUNIT_FAIL( strstream.str() ) ;
59     return;
60   }
61
62   debugStream << "test_MPI_Access_Time" << myrank << endl ;
63
64   ParaMEDMEM::CommInterface interface ;
65
66   ParaMEDMEM::MPIProcessorGroup* group = new ParaMEDMEM::MPIProcessorGroup(interface) ;
67
68   ParaMEDMEM::MPIAccess mpi_access( group ) ;
69
70 #define maxreq 10
71
72   if ( myrank >= 2 ) {
73     debugStream << "test_MPI_Access_Time_0 rank" << myrank << " --> mpi_access->Barrier" << endl ;
74     mpi_access.barrier() ;
75     debugStream << "test_MPI_Access_Time_0 rank" << myrank << " <-- mpi_access->Barrier" << endl ;
76     delete group ;
77     debugStream << "test_MPI_Access_Time" << myrank << " OK" << endl ;
78     return ;
79   }
80
81   int target = 1 - myrank ;
82   int SendTimeRequestId[maxreq] ;
83   int RecvTimeRequestId[maxreq] ;
84   int SendRequestId[maxreq] ;
85   int RecvRequestId[maxreq] ;
86   int sts ;
87   int sendbuf[maxreq] ;
88   int recvbuf[maxreq] ;
89   int i = 0 ;
90   ParaMEDMEM::TimeMessage aSendTimeMsg[maxreq] ;
91   ParaMEDMEM::TimeMessage aRecvTimeMsg[maxreq] ;
92   double t ;
93   double dt = 1. ;
94   double maxt = 10. ;
95   for ( t = 0 ; t < maxt ; t = t+dt ) {
96     if ( myrank == 0 ) {
97       aSendTimeMsg[i].time = t ;
98       aSendTimeMsg[i].deltatime = dt ;
99       //aSendTimeMsg[i].maxtime = maxt ;
100       //sts = mpi_access.ISend( &aSendTimeMsg , mpi_access.timeExtent() ,
101       sts = mpi_access.ISend( &aSendTimeMsg[i] , 1 ,
102                               mpi_access.timeType() , target ,
103                               SendTimeRequestId[i]) ;
104       debugStream << "test" << myrank << " ISend RequestId " << SendTimeRequestId[i]
105            << " tag " << mpi_access.sendMPITag(target) << endl ;
106       sendbuf[i] = i ;
107       sts = mpi_access.ISend(&sendbuf[i],1,MPI_INT,target, SendRequestId[i]) ;
108       debugStream << "test" << myrank << " ISend RequestId " << SendRequestId[i]
109            << " tag " << mpi_access.sendMPITag(target) << endl ;
110     }
111     else {
112       //sts = mpi_access.IRecv( &aRecvTimeMsg , mpi_access.timeExtent() ,
113       sts = mpi_access.IRecv( &aRecvTimeMsg[i] , 1 ,
114                               mpi_access.timeType() , target ,
115                               RecvTimeRequestId[i]) ;
116       debugStream << "test" << myrank << " IRecv RequestId " << RecvTimeRequestId[i]
117            << " tag " << mpi_access.recvMPITag(target) << endl ;
118       sts = mpi_access.IRecv(&recvbuf[i],1,MPI_INT,target, RecvRequestId[i]) ;
119       debugStream << "test" << myrank << " IRecv RequestId " << RecvRequestId[i]
120            << " tag " << mpi_access.recvMPITag(target) << endl ;
121     }
122     int j ;
123     for (j = 0 ; j <= i ; j++) {
124       int flag ;
125       if ( myrank == 0 ) {
126         mpi_access.test( SendTimeRequestId[j], flag ) ;
127       }
128       else {
129         mpi_access.test( RecvTimeRequestId[j], flag ) ;
130       }
131       if ( flag ) {
132         int target,source, tag, error, outcount ;
133         if ( myrank == 0 ) {
134           mpi_access.status( SendTimeRequestId[j], target, tag, error, outcount,
135                              true ) ;
136           debugStream << "test" << myrank << " Test(Send TimeRequestId " << SendTimeRequestId[j]
137                << ") : target " << target << " tag " << tag << " error " << error
138                << " flag " << flag << aSendTimeMsg[j] << endl ;
139         }
140         else {
141           mpi_access.status( RecvTimeRequestId[j], source, tag, error, outcount,
142                              true ) ;
143           debugStream << "test" << myrank << " Test(Recv TimeRequestId "
144                << RecvTimeRequestId[j] << ") : source " << source << " tag " << tag
145                << " error " << error << " outcount " << outcount
146                << " flag " << flag << aRecvTimeMsg[j] << endl ;
147           if ( (outcount != 1) | (aRecvTimeMsg[j].time != j) ) {
148             ostringstream strstream ;
149             strstream << "==========================================================="
150                       << endl << "test" << myrank << " outcount " << outcount << " KO"
151                       << " RecvTimeRequestId " << RecvTimeRequestId[j] << endl
152                       << "==========================================================="
153                       << endl ;
154             debugStream << strstream.str() << endl ;
155             CPPUNIT_FAIL( strstream.str() ) ;
156           }
157           else {
158             debugStream << "==========================================================="
159                  << endl << "test" << myrank << " outcount " << outcount
160                  << " RecvTimeRequestId " << RecvTimeRequestId[j] << " OK" << endl
161                  << "==========================================================="
162                  << endl ;
163           }
164         }
165       }
166       if ( myrank == 0 ) {
167         mpi_access.test( SendRequestId[j], flag ) ;
168       }
169       else {
170         mpi_access.test( RecvRequestId[j], flag ) ;
171       }
172       if ( flag ) {
173         int target,source, tag, error, outcount ;
174         if ( myrank == 0 ) {
175           mpi_access.status( SendRequestId[j], target, tag, error, outcount,
176                              true ) ;
177           debugStream << "test" << myrank << " Test(Send RequestId " << SendRequestId[j]
178                << ") : target " << target << " tag " << tag << " error " << error
179                << " flag " << flag << endl ;
180         }
181         else {
182           mpi_access.status( RecvRequestId[j], source, tag, error, outcount,
183                              true ) ;
184           debugStream << "test" << myrank << " Test(Recv RequestId "
185                << RecvRequestId[j] << ") : source " << source << " tag " << tag
186                << " error " << error << " outcount " << outcount
187                << " flag " << flag << endl ;
188           if ( (outcount != 1) | (recvbuf[j] != j) ) {
189             ostringstream strstream ;
190             strstream << "==========================================================="
191                       << endl << "test" << myrank << " outcount "
192                       << outcount << " recvbuf " << recvbuf[j] << " KO" << endl
193                       << "==========================================================="
194                       << endl ;
195             debugStream << strstream.str() << endl ;
196             CPPUNIT_FAIL( strstream.str() ) ;
197           }
198           else {
199             debugStream << "==========================================================="
200                  << endl << "test" << myrank << " outcount " << outcount
201                  << " RequestId " << RecvRequestId[j] << " OK" << endl
202                  << "==========================================================="
203                  << endl ;
204           }
205         }
206       }
207     }
208     char msgerr[MPI_MAX_ERROR_STRING] ;
209     int lenerr ;
210     mpi_access.errorString(sts, msgerr, &lenerr) ;
211     debugStream << "test" << myrank << " lenerr " << lenerr << " "
212          << msgerr << endl ;
213
214     if ( sts != MPI_SUCCESS ) {
215       ostringstream strstream ;
216       strstream << "==========================================================="
217                 << "test" << myrank << " KO"
218                 << "==========================================================="
219                 << endl ;
220       debugStream << strstream.str() << endl ;
221       CPPUNIT_FAIL( strstream.str() ) ;
222     }
223     i = i + 1 ;
224   }
225
226   if(MPI_ACCESS_VERBOSE) mpi_access.check() ;
227   if ( myrank == 0 ) {
228     mpi_access.waitAll(maxreq, SendTimeRequestId) ;
229     mpi_access.deleteRequests(maxreq, SendTimeRequestId) ;
230     mpi_access.waitAll(maxreq, SendRequestId) ;
231     mpi_access.deleteRequests(maxreq, SendRequestId) ;
232   }
233   else {
234     mpi_access.waitAll(maxreq, RecvTimeRequestId) ;
235     mpi_access.deleteRequests(maxreq, RecvTimeRequestId) ;
236     mpi_access.waitAll(maxreq, RecvRequestId) ;
237     mpi_access.deleteRequests(maxreq, RecvRequestId) ;
238   }
239   if(MPI_ACCESS_VERBOSE) mpi_access.check() ;
240
241   if ( myrank == 0 ) {
242     int sendrequests[2*maxreq] ;
243     int sendreqsize = mpi_access.sendRequestIds( target , 2*maxreq , sendrequests ) ;
244     if ( sendreqsize != 0 ) {
245       ostringstream strstream ;
246       strstream << "=========================================================" << endl
247                 << "test" << myrank << " sendreqsize " << sendreqsize << " KO" << endl
248                 << "=========================================================" << endl ;
249       debugStream << strstream.str() << endl ;
250       CPPUNIT_FAIL( strstream.str() ) ;
251     }
252     else {
253       debugStream << "=========================================================" << endl
254            << "test" << myrank << " sendreqsize " << sendreqsize << " OK" << endl
255            << "=========================================================" << endl ;
256     }
257   }
258   else {
259     int recvrequests[2*maxreq] ;
260     int recvreqsize = mpi_access.sendRequestIds( target , 2*maxreq , recvrequests ) ;
261     if ( recvreqsize != 0 ) {
262       ostringstream strstream ;
263       strstream << "=========================================================" << endl
264                 << "test" << myrank << " recvreqsize " << recvreqsize << " KO" << endl
265                 << "=========================================================" << endl ;
266       debugStream << strstream.str() << endl ;
267       CPPUNIT_FAIL( strstream.str() ) ;
268     }
269     else {
270       debugStream << "=========================================================" << endl
271            << "test" << myrank << " recvreqsize " << recvreqsize << " OK" << endl
272            << "=========================================================" << endl ;
273     }
274   }
275
276   debugStream << "test_MPI_Access_Time_0 rank" << myrank << " --> mpi_access->Barrier" << endl ;
277   mpi_access.barrier() ;
278   debugStream << "test_MPI_Access_Time_0 rank" << myrank << " <-- mpi_access->Barrier" << endl ;
279
280   delete group ;
281
282   //  MPI_Finalize();
283
284   debugStream << "test_MPI_Access_Time" << myrank << " OK" << endl ;
285
286   return ;
287 }
288
289
290
291