]> SALOME platform Git repositories - tools/medcoupling.git/blob - src/ParaMEDMEM/Test/test_MPI_Access_IProbe.cxx
Salome HOME
Merge from BR_V5_DEV 16Feb09
[tools/medcoupling.git] / src / ParaMEDMEM / Test / test_MPI_Access_IProbe.cxx
1 //  Copyright (C) 2007-2008  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.
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 #include <time.h>
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_IProbe() {
44
45   cout << "test_MPI_Access_IProbe" << 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_IProbe must be runned with 2 procs" << endl ;
57     cout << strstream.str() << endl ;
58     CPPUNIT_FAIL( strstream.str() ) ;
59   }
60
61   cout << "test_MPI_Access_IProbe" << myrank << endl ;
62
63   ParaMEDMEM::CommInterface interface ;
64
65   ParaMEDMEM::MPIProcessorGroup* group = new ParaMEDMEM::MPIProcessorGroup(interface) ;
66
67   ParaMEDMEM::MPIAccess mpi_access( group ) ;
68
69   if ( myrank >= 2 ) {
70     mpi_access.barrier() ;
71     delete group ;
72     return ;
73   }
74
75   int target = 1 - myrank ;
76   int sendbuf[10] ;
77   int RequestId[10] ;
78   int sts ;
79   int i ;
80   for ( i = 0 ; i < 10 ; i++ ) {
81      if ( myrank == 0 ) {
82        sendbuf[i] = i ;
83        sts = mpi_access.ISend(&sendbuf[i],1,MPI_INT,target, RequestId[i]) ;
84        cout << "test" << myrank << " Send RequestId " << RequestId[i]
85             << endl ;
86      }
87      else {
88        int flag = false ;
89        while ( !flag ) {
90             int source, tag, outcount ;
91             MPI_Datatype datatype ;
92             sts = mpi_access.IProbe(target, source, tag, datatype, outcount, flag ) ;
93             if ( flag ) {
94               cout << "test" << myrank << " " << i << " IProbe target " << target
95                    << " source " << source << " tag " << tag
96                    << " outcount " << outcount << " flag " << flag << endl ;
97             }
98             else {
99               cout << "test" << myrank << " IProbe flag " << flag << endl ;
100               sleep( 1 ) ;
101             }
102             if ( flag ) {
103               int recvbuf ;
104               sts = mpi_access.recv(&recvbuf,outcount,datatype,source, RequestId[i],
105                                     &outcount) ;
106               if ( (outcount != 1) | (recvbuf != i) ) {
107                 ostringstream strstream ;
108                 strstream << "==========================================================="
109                           << endl << "test" << myrank << " outcount " << outcount
110                           << " recvbuf " << recvbuf << " KO" << endl
111                           << "==========================================================="
112                           << endl ;
113                 cout << strstream.str() << endl ;
114                 CPPUNIT_FAIL( strstream.str() ) ;
115               }
116               cout << "==========================================================="
117                    << endl << "test" << myrank << " outcount " << outcount
118                    << " recvbuf " << recvbuf << " OK" << endl
119                    << "==========================================================="
120                    << endl ;
121             }
122        }
123      }
124      char msgerr[MPI_MAX_ERROR_STRING] ;
125      int lenerr ;
126      mpi_access.errorString(sts, msgerr, &lenerr) ;
127      cout << "test" << myrank << " lenerr " << lenerr << " "
128           << msgerr << endl ;
129
130      if ( sts != MPI_SUCCESS ) {
131        ostringstream strstream ;
132        strstream << "==========================================================="
133                  << "test" << myrank << " KO"
134                  << "==========================================================="
135                  << endl ;
136        cout << strstream.str() << endl ;
137        CPPUNIT_FAIL( strstream.str() ) ;
138      }
139      mpi_access.check() ;
140   }
141   int flag ;
142   mpi_access.testAll(10,RequestId,flag) ;
143   mpi_access.waitAll(10,RequestId) ;
144   mpi_access.deleteRequests(10,RequestId) ;
145   mpi_access.testAll(10,RequestId,flag) ;
146   if ( !flag ) {
147     ostringstream strstream ;
148     strstream << "test" << myrank << " flag " << flag << " KO" << endl ;
149     cout << strstream.str() << endl ;
150     CPPUNIT_FAIL( strstream.str() ) ;
151   }
152   mpi_access.check() ;
153
154   mpi_access.barrier() ;
155
156   delete group ;
157
158 //  MPI_Finalize();
159
160   cout << "test" << myrank << " OK" << endl ;
161
162   return ;
163 }
164
165
166
167