Salome HOME
Update copyrights
[tools/medcoupling.git] / src / MEDPartitioner / MEDPARTITIONER_UtilsPara.cxx
1 // Copyright (C) 2007-2019  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 "MEDPARTITIONER_Utils.hxx"
21
22 #include "MEDLoader.hxx"
23 #include "MEDLoaderBase.hxx"
24 #include "MEDFileUtilities.hxx"
25 #include "CellModel.hxx"
26 #include "MEDCouplingUMesh.hxx"
27 #include "MEDCouplingFieldDouble.hxx"
28 #include "InterpKernelException.hxx"
29 #include "MCAuto.hxx"
30 #include "InterpKernelAutoPtr.hxx"
31
32 #include <fstream>
33 #include <iostream>
34 #include <iomanip>
35 #include <sstream>
36 #include <string>
37
38 #ifdef HAVE_MPI
39 #include <mpi.h>
40 #endif
41
42 using namespace MEDPARTITIONER;
43
44 /*!
45  * not optimized but suffisant
46  * return empty vector if i am not target
47  */
48 std::vector<std::string> MEDPARTITIONER::SendAndReceiveVectorOfString(const std::vector<std::string>& vec, const int source, const int target)
49 {
50   int rank=MyGlobals::_Rank;
51
52   MPI_Status status;
53   int tag = 111001;
54   if (rank == source)
55     {
56       std::string str=SerializeFromVectorOfString(vec);
57       int size=str.length();
58       MPI_Send( &size, 1, MPI_INT, target, tag, MPI_COMM_WORLD );
59       MPI_Send( (void*)str.data(), str.length(), MPI_CHAR, target, tag+100, MPI_COMM_WORLD );
60     }
61   
62   int recSize=0;
63   if (rank == target)
64     {
65       MPI_Recv(&recSize, 1, MPI_INT, source, tag, MPI_COMM_WORLD, &status);
66       std::string recData(recSize,'x');
67       MPI_Recv((void*)recData.data(), recSize, MPI_CHAR, source, tag+100, MPI_COMM_WORLD, &status);
68       return DeserializeToVectorOfString(recData); //not empty one for target proc
69     }
70   std::vector<std::string> res;
71   return res; //empty one for other proc
72 }
73
74 /*!
75  * strings NO need all same size!!!!
76  */
77 std::vector<std::string> MEDPARTITIONER::AllgathervVectorOfString(const std::vector<std::string>& vec)
78 {
79   if (MyGlobals::_World_Size==1) //nothing to do
80     return vec;
81
82   int world_size=MyGlobals::_World_Size;
83   std::string str=SerializeFromVectorOfString(vec);
84   
85   std::vector<int> indexes(world_size);
86   int size=str.length();
87   MPI_Allgather(&size, 1, MPI_INT, 
88                 &indexes[0], 1, MPI_INT, MPI_COMM_WORLD);
89   
90   //calcul of displacement
91   std::vector< int > disp(1,0);
92   for (int i=0; i<world_size; i++) disp.push_back( disp.back() + indexes[i] );
93   
94   std::string recData(disp.back(),'x');
95   MPI_Allgatherv((void*)str.data(), str.length(), MPI_CHAR,
96                  (void*)recData.data(), &indexes[0], &disp[0], MPI_CHAR,
97                  MPI_COMM_WORLD);
98   
99   //really extraordinary verbose for debug
100   std::vector<std::string> deserial=DeserializeToVectorOfString(recData);
101   if (MyGlobals::_Verbose>1000) 
102     {
103       std::cout << "proc "<<MyGlobals::_Rank<<" : receive '" << recData << "'" << std::endl;
104       std::cout << "deserialize is : a vector of size " << deserial.size() << std::endl;
105       std::cout << ReprVectorOfString(deserial) << std::endl;
106     }
107   return deserial;
108 }
109
110 /*!
111   Sends content of \a vec to processor \a target. To be used with \a RecvDoubleVec method.
112   \param vec vector to be sent
113   \param target processor id of the target
114 */
115 void MEDPARTITIONER::SendDoubleVec(const std::vector<double>& vec, const int target)
116 {
117   int tag = 111002;
118   int size=vec.size();
119   if (MyGlobals::_Verbose>1000) 
120     std::cout << "proc " << MyGlobals::_Rank << " : --> SendDoubleVec " << size << std::endl;
121 #ifdef HAVE_MPI
122   MPI_Send(&size, 1, MPI_INT, target, tag, MPI_COMM_WORLD);
123   MPI_Send(const_cast<double*>(&vec[0]), size, MPI_DOUBLE, target, tag+100, MPI_COMM_WORLD);
124 #endif
125 }
126
127 /*! Receives messages from proc \a source to fill vector<int> vec.
128   To be used with \a SendDoubleVec method.
129
130   \param vec vector that is filled
131   \param source processor id of the incoming messages
132 */
133 std::vector<double>* MEDPARTITIONER::RecvDoubleVec(const int source)
134 {
135   int tag = 111002;
136   int size;
137 #ifdef HAVE_MPI
138   MPI_Status status;  
139   MPI_Recv(&size, 1, MPI_INT, source, tag, MPI_COMM_WORLD, &status);
140   if (MyGlobals::_Verbose>1000) 
141     std::cout << "proc " << MyGlobals::_Rank << " : <-- RecvDoubleVec " << size << std::endl;
142   std::vector<double>* vec=new std::vector<double>;
143   vec->resize(size);
144   MPI_Recv(&vec[0], size, MPI_DOUBLE, source, tag+100, MPI_COMM_WORLD, &status);
145 #endif
146   return vec;
147 }
148
149 void MEDPARTITIONER::RecvDoubleVec(std::vector<double>& vec, const int source)
150 {
151   int tag = 111002;
152   int size;
153 #ifdef HAVE_MPI
154   MPI_Status status;  
155   MPI_Recv(&size, 1, MPI_INT, source, tag, MPI_COMM_WORLD, &status);
156   if (MyGlobals::_Verbose>1000)
157     std::cout<< "proc " << MyGlobals::_Rank << " : <-- RecvDoubleVec " << size << std::endl;;
158   vec.resize(size);
159   MPI_Recv(&vec[0], size, MPI_DOUBLE, source, tag+100, MPI_COMM_WORLD, &status);
160 #endif
161 }
162 /*!
163   Sends content of \a vec to processor \a target. To be used with \a RecvIntVec method.
164   \param vec vector to be sent
165   \param target processor id of the target
166 */
167 void MEDPARTITIONER::SendIntVec(const std::vector<int>& vec, const int target)
168 {
169   int tag = 111003;
170   int size=vec.size();
171   if (MyGlobals::_Verbose>1000)
172     std::cout << "proc " << MyGlobals::_Rank << " : --> SendIntVec " << size << std::endl;
173 #ifdef HAVE_MPI
174   MPI_Send(&size, 1, MPI_INT, target, tag, MPI_COMM_WORLD);
175   MPI_Send(const_cast<int*>(&vec[0]), size,MPI_INT, target, tag+100, MPI_COMM_WORLD);
176 #endif
177 }
178
179 /*! Receives messages from proc \a source to fill vector<int> vec.
180   To be used with \a SendIntVec method.
181   \param vec vector that is filled
182   \param source processor id of the incoming messages
183 */
184 std::vector<int> *MEDPARTITIONER::RecvIntVec(const int source)
185 {
186   int tag = 111003;
187   int size;
188 #ifdef HAVE_MPI
189   MPI_Status status;  
190   MPI_Recv(&size, 1, MPI_INT, source, tag, MPI_COMM_WORLD, &status);
191   if (MyGlobals::_Verbose>1000)
192     std::cout << "proc " << MyGlobals::_Rank << " : <-- RecvIntVec " << size << std::endl;
193   std::vector<int> *vec=new std::vector<int>;
194   vec->resize(size);
195   MPI_Recv(&vec[0], size, MPI_INT, source, tag+100, MPI_COMM_WORLD, &status);
196 #endif
197   return vec;
198 }
199
200 void MEDPARTITIONER::RecvIntVec(std::vector<int>& vec, const int source)
201 {
202   int tag = 111003;
203   int size;
204 #ifdef HAVE_MPI
205   MPI_Status status;  
206   MPI_Recv(&size, 1, MPI_INT, source, tag, MPI_COMM_WORLD, &status);
207   if (MyGlobals::_Verbose>1000)
208     std::cout << "proc " << MyGlobals::_Rank << " : <-- RecvIntVec " << size << std::endl;
209   vec.resize(size);
210   MPI_Recv(&vec[0], size, MPI_INT, source, tag+100, MPI_COMM_WORLD,&status);
211 #endif
212 }
213
214 /*!
215   Sends content of \a dataArrayInt to processor \a target. 
216   To be used with \a RecvDataArrayInt method.
217   \param da dataArray to be sent
218   \param target processor id of the target
219 */
220 void MEDPARTITIONER::SendDataArrayInt(const MEDCoupling::DataArrayInt *da, const int target)
221 {
222   if (da==0)
223     throw INTERP_KERNEL::Exception("Problem send DataArrayInt* NULL");
224   int tag = 111004;
225   int size[3];
226   size[0]=da->getNbOfElems();
227   size[1]=da->getNumberOfTuples();
228   size[2]=da->getNumberOfComponents();
229   if (MyGlobals::_Verbose>1000) 
230     std::cout << "proc " << MyGlobals::_Rank << " : --> SendDataArrayInt " << size[0] << std::endl;
231 #ifdef HAVE_MPI
232   MPI_Send(&size, 3, MPI_INT, target, tag, MPI_COMM_WORLD);
233   const int *p=da->getConstPointer();
234   MPI_Send(const_cast<int*>(&p[0]), size[0] ,MPI_INT, target, tag+100, MPI_COMM_WORLD);
235 #endif
236 }
237
238 /*! Receives messages from proc \a source to fill dataArrayInt da.
239   To be used with \a SendIntVec method.
240   \param da dataArrayInt that is filled
241   \param source processor id of the incoming messages
242 */
243 MEDCoupling::DataArrayInt *MEDPARTITIONER::RecvDataArrayInt(const int source)
244 {
245   int tag = 111004;
246   int size[3];
247 #ifdef HAVE_MPI
248   MPI_Status status;
249   MPI_Recv(size, 3, MPI_INT, source, tag, MPI_COMM_WORLD, &status);
250   if (MyGlobals::_Verbose>1000)
251     std::cout << "proc " << MyGlobals::_Rank << " : <-- RecvDataArrayInt " << size[0] << std::endl;
252   if (size[0]!=(size[1]*size[2]))
253     throw INTERP_KERNEL::Exception("Problem in RecvDataArrayInt incoherent sizes");
254   MEDCoupling::DataArrayInt* da=MEDCoupling::DataArrayInt::New();
255   da->alloc(size[1],size[2]);
256   int *p=da->getPointer();
257   MPI_Recv(const_cast<int*>(&p[0]), size[0], MPI_INT, source, tag+100, MPI_COMM_WORLD, &status);
258 #endif
259   return da;
260 }
261
262 /*!
263   Sends content of \a dataArrayInt to processor \a target. 
264   To be used with \a RecvDataArrayDouble method.
265   \param da dataArray to be sent
266   \param target processor id of the target
267 */
268 void MEDPARTITIONER::SendDataArrayDouble(const MEDCoupling::DataArrayDouble *da, const int target)
269 {
270   if (da==0)
271     throw INTERP_KERNEL::Exception("Problem send DataArrayDouble* NULL");
272   int tag = 111005;
273   int size[3];
274   size[0]=da->getNbOfElems();
275   size[1]=da->getNumberOfTuples();
276   size[2]=da->getNumberOfComponents();
277   if (MyGlobals::_Verbose>1000) 
278     std::cout << "proc " << MyGlobals::_Rank << " : --> SendDataArrayDouble " << size[0] << std::endl;
279 #ifdef HAVE_MPI
280   MPI_Send(&size, 3, MPI_INT, target, tag, MPI_COMM_WORLD);
281   const double *p=da->getConstPointer();
282   MPI_Send(const_cast<double*>(&p[0]), size[0] ,MPI_DOUBLE, target, tag+100, MPI_COMM_WORLD);
283 #endif
284 }
285
286 /*! Receives messages from proc \a source to fill dataArrayDouble da.
287   To be used with \a SendDoubleVec method.
288   \param da dataArrayDouble that is filled
289   \param source processor id of the incoming messages
290 */
291 MEDCoupling::DataArrayDouble* MEDPARTITIONER::RecvDataArrayDouble(const int source)
292 {
293   int tag = 111005;
294   int size[3];
295 #ifdef HAVE_MPI
296   MPI_Status status;
297   MPI_Recv(size, 3, MPI_INT, source, tag, MPI_COMM_WORLD, &status);
298   if (MyGlobals::_Verbose>1000)
299     std::cout << "proc " << MyGlobals::_Rank << " : <-- RecvDataArrayDouble " << size[0] << std::endl;
300   if (size[0]!=(size[1]*size[2]))
301     throw INTERP_KERNEL::Exception("Problem in RecvDataArrayDouble incoherent sizes");
302   MEDCoupling::DataArrayDouble* da=MEDCoupling::DataArrayDouble::New();
303   da->alloc(size[1],size[2]);
304   double *p=da->getPointer();
305   MPI_Recv(const_cast<double*>(&p[0]), size[0], MPI_DOUBLE, source, tag+100, MPI_COMM_WORLD, &status);
306 #endif
307   return da;
308 }
309
310 void MEDPARTITIONER::TestVectorOfStringMpi()
311 {
312   int rank=MyGlobals::_Rank;
313   int world_size=MyGlobals::_World_Size;
314   std::vector<std::string> myVector;
315   std::ostringstream oss;
316   oss << "hello from " << std::setw(5) << rank << " " << std::string(rank+1,'n') <<
317     " next is an empty one";
318   myVector.push_back(oss.str());
319   myVector.push_back("");
320   myVector.push_back("next is an singleton");
321   myVector.push_back("1");
322   
323   if (rank==0)
324     {
325       std::string s0=SerializeFromVectorOfString(myVector);
326       std::vector<std::string> res=DeserializeToVectorOfString(s0);
327       if (res.size()!=myVector.size()) 
328         throw INTERP_KERNEL::Exception("Problem in (de)serialise VectorOfString incoherent sizes");
329       for (std::size_t i=0; i<myVector.size(); i++)
330         if (res[i]!=myVector[i])
331           throw INTERP_KERNEL::Exception("Problem in (de)serialise VectorOfString incoherent elements");
332     }
333
334   for (int i=0; i<world_size; i++)
335     {
336       for (int j=0; j<world_size; j++)
337         {
338           std::vector<std::string> res=SendAndReceiveVectorOfString(myVector, i, j);
339           if ((rank==j) && MyGlobals::_Verbose>20)
340             std::cout << "proc " << rank << " : receive \n" << ReprVectorOfString(res) << std::endl;
341           if (rank==j)
342             {
343               if (res.size()!=myVector.size()) 
344                 throw INTERP_KERNEL::Exception("Problem in SendAndReceiveVectorOfString incoherent sizes");
345               for (std::size_t ii=1; ii<myVector.size(); ii++) //first is different
346                 if (res[i]!=myVector[ii])
347                   throw INTERP_KERNEL::Exception("Problem in SendAndReceiveVectorOfString incoherent elements");
348             }
349           else 
350             {
351               if (res.size()!=0) 
352                 throw INTERP_KERNEL::Exception("Problem in SendAndReceiveVectorOfString size have to be 0");
353             }
354         }
355     }
356   std::vector<std::string> res=AllgathervVectorOfString(myVector);
357   //sometimes for test
358   res=AllgathervVectorOfString(myVector);
359   res=AllgathervVectorOfString(myVector);
360   if (rank==0 && MyGlobals::_Verbose>20)
361     std::cout << "proc " << rank << " : receive \n" << ReprVectorOfString(res) << std::endl;
362   if (res.size()!=myVector.size()*world_size) 
363     throw INTERP_KERNEL::Exception("Problem in AllgathervVectorOfString incoherent sizes");
364   int jj=-1;
365   for (int j=0; j<world_size; j++)
366     {
367       for (int i=0; i<(int)myVector.size(); i++)
368         {
369           jj=jj+1;
370           if (i==0)
371             continue; //first is different
372           if (res[jj]!=myVector[i])
373             throw INTERP_KERNEL::Exception("Problem in AllgathervVectorOfString incoherent elements");
374         }
375     }
376   if (MyGlobals::_Verbose)
377     std::cout << "proc " << rank << " : OK TestVectorOfStringMpi END" << std::endl;
378 }
379
380 void MEDPARTITIONER::TestMapOfStringIntMpi()
381 {
382   int rank=MyGlobals::_Rank;
383   std::map<std::string,int> myMap;
384   myMap["one"]=1;
385   myMap["two"]=22;  //a bug
386   myMap["three"]=3;
387   myMap["two"]=2; //last speaking override
388   
389   if (rank==0)
390     {
391       std::vector<std::string> v2=VectorizeFromMapOfStringInt(myMap);
392       std::map<std::string,int> m3=DevectorizeToMapOfStringInt(v2);
393       if (ReprMapOfStringInt(m3)!=ReprMapOfStringInt(myMap))
394         throw INTERP_KERNEL::Exception("Problem in (de)vectorize MapOfStringInt");
395     }
396     
397   std::vector<std::string> v2=AllgathervVectorOfString(VectorizeFromMapOfStringInt(myMap));
398   if (rank==0 && MyGlobals::_Verbose>20)
399     {
400       std::cout << "v2 is : a vector of size " << v2.size() << std::endl;
401       std::cout << ReprVectorOfString(v2) << std::endl;
402       std::map<std::string,int> m2=DevectorizeToMapOfStringInt(v2);
403       std::cout << "m2 is : a map of size " << m2.size() << std::endl;
404       std::cout << ReprMapOfStringInt(m2) << std::endl;
405     }
406   if (MyGlobals::_Verbose)
407     std::cout << "proc " << rank << " : OK TestMapOfStringIntMpi END" << std::endl;
408 }
409
410 void MEDPARTITIONER::TestMapOfStringVectorOfStringMpi()
411 {
412   int rank=MyGlobals::_Rank;
413   std::vector<std::string> myVector;
414   std::ostringstream oss;
415   oss << "hello from " << std::setw(5) << MyGlobals::_Rank << " " << std::string(rank+1,'n') << " next is an empty one";
416   myVector.push_back(oss.str());
417   myVector.push_back("");
418   myVector.push_back("next is an singleton");
419   myVector.push_back("1");
420   
421   if (rank==0)
422     {
423       std::map< std::string,std::vector<std::string> > m2;
424       m2["first key"]=myVector;
425       m2["second key"]=myVector;
426       std::vector<std::string> v2=VectorizeFromMapOfStringVectorOfString(m2);
427       std::map< std::string,std::vector<std::string> > m3=DevectorizeToMapOfStringVectorOfString(v2);
428       if (rank==0 && MyGlobals::_Verbose>20)
429         std::cout << "m2 is : a MapOfStringVectorOfString of size " << m2.size() << std::endl;
430       std::cout << ReprMapOfStringVectorOfString(m2) << std::endl;
431       std::cout << "v2 is : a vector of size " << v2.size() << std::endl;
432       std::cout << ReprVectorOfString(v2) << std::endl;
433       std::cout << "m3 is : a map of size "<<m3.size() << std::endl;
434       std::cout << ReprMapOfStringVectorOfString(m3) << std::endl;
435       if (ReprMapOfStringVectorOfString(m3)!=ReprMapOfStringVectorOfString(m2))
436         throw INTERP_KERNEL::Exception("Problem in (de)vectorize MapOfStringVectorOfString");
437     }
438     
439   std::map< std::string,std::vector<std::string> > m4;
440   m4["1rst key"]=myVector;
441   m4["2snd key"]=myVector;
442   std::vector<std::string> v4=AllgathervVectorOfString(VectorizeFromMapOfStringVectorOfString(m4));
443   if (rank==0 && MyGlobals::_Verbose>20)
444     {
445       std::map< std::string,std::vector<std::string> > m5=DevectorizeToMapOfStringVectorOfString(v4);
446       std::map< std::string,std::vector<std::string> > m6=DeleteDuplicatesInMapOfStringVectorOfString(m5);
447       std::cout<< "m5 is : a map of size "<<m5.size() << std::endl;
448       std::cout<< ReprMapOfStringVectorOfString(m5) << std::endl;
449       std::cout<< "m6 is : a map from m5 with deleteDuplicates of size " << m6.size() << std::endl;
450       std::cout<< ReprMapOfStringVectorOfString(m6) << std::endl;
451     }
452   if (MyGlobals::_Verbose)
453     std::cout<<"proc " << rank << " : OK TestMapOfStringVectorOfStringMpi END" << std::endl;
454 }
455
456 void MEDPARTITIONER::TestDataArrayMpi()
457 {
458   int rank=MyGlobals::_Rank;
459   //int
460   {
461     MEDCoupling::DataArrayInt* send=MEDCoupling::DataArrayInt::New();
462     MEDCoupling::DataArrayInt* recv=0;
463     int nbOfTuples=5;
464     int numberOfComponents=3;
465     send->alloc(nbOfTuples,numberOfComponents);
466     std::vector<int> vals;
467     for (int j=0; j<nbOfTuples; j++)
468       for (int i=0; i<numberOfComponents; i++) vals.push_back((j+1)*10+i+1);
469     std::copy(vals.begin(),vals.end(),send->getPointer());
470     if (rank==0)
471       SendDataArrayInt(send, 1);
472     if (rank==1)
473       recv=RecvDataArrayInt(0);
474     if (rank==1 && MyGlobals::_Verbose>20)
475       {
476         std::cout << send->repr() << std::endl;
477         std::cout << recv->repr() << std::endl;
478       }
479     if (rank==1)
480       {
481         if (send->repr()!=recv->repr())
482           throw INTERP_KERNEL::Exception("Problem in send&recv DataArrayInt");
483       }
484     send->decrRef();
485     if (rank==1)
486       recv->decrRef();
487   }
488   //double
489   {
490     MEDCoupling::DataArrayDouble* send=MEDCoupling::DataArrayDouble::New();
491     MEDCoupling::DataArrayDouble* recv=0;
492     int nbOfTuples=5;
493     int numberOfComponents=3;
494     send->alloc(nbOfTuples,numberOfComponents);
495     std::vector<double> vals;
496     for (int j=0; j<nbOfTuples; j++)
497       for (int i=0; i<numberOfComponents; i++) vals.push_back(double(j+1)+double(i+1)/10);
498     std::copy(vals.begin(),vals.end(),send->getPointer());
499     if (rank==0) SendDataArrayDouble(send, 1);
500     if (rank==1) recv=RecvDataArrayDouble(0);
501     if (rank==1 && MyGlobals::_Verbose>20)
502       {
503         std::cout << send->repr() << std::endl;
504         std::cout << recv->repr() << std::endl;
505       }
506     if (rank==1)
507       {
508         if (send->repr()!=recv->repr())
509           throw INTERP_KERNEL::Exception("Problem in send&recv DataArrayDouble");
510       }
511     send->decrRef();
512     if (rank==1) recv->decrRef();
513   }
514   
515   if (MyGlobals::_Verbose)
516     std::cout << "proc " << rank << " : OK TestDataArrayMpi END" << std::endl;
517 }
518
519 void MEDPARTITIONER::TestPersistantMpi0To1(int taille, int nb)
520 {
521   double temps_debut=MPI_Wtime();
522   int rank=MyGlobals::_Rank;
523   std::vector<int> x, y;
524   int tag=111111;
525   MPI_Request requete0, requete1;
526   MPI_Status statut;
527   int ok=0;
528   std::string res;
529   if (rank==0)
530     {
531       x.resize(taille);
532       MPI_Ssend_init(&x[0], taille, MPI_INT, 1, tag, MPI_COMM_WORLD , &requete0);
533       for(int k=0; k<nb; k++)
534         {
535           for (int i=0; i<taille; ++i) x[i]=k;
536           //Envoi d’un gros message --> cela peut prendre du temps
537           MPI_Start(&requete0);
538           //Traitement sequentiel independant de "x"
539           //...
540           MPI_Wait(&requete0, &statut);
541           //Traitement sequentiel impliquant une modification de "x" en memoire
542           //x=...
543         }
544       MPI_Request_free(&requete0);
545     }
546   else if (rank == 1)
547     {
548       y.resize(taille);
549       MPI_Recv_init(&y[0], taille,  MPI_INT, 0, tag, MPI_COMM_WORLD , &requete1);
550       for(int k=0; k<nb; k++)
551         {
552           //Pre-traitement sequentiel
553           //...
554           for (int i=0; i<taille; ++i) y[i]=(-1);
555           //Reception du gros message --> cela peut prendre du temps
556           MPI_Start(&requete1);
557           //Traitement sequentiel independant de "y"
558           //...
559           MPI_Wait(&requete1, &statut);
560           //Traitement sequentiel dependant de "y"
561           //...=f(y)
562           int nbb=0;
563           for (int i=0; i<taille; ++i)
564             if (y[i]==k)
565               nbb++;
566           if (nbb==taille)
567             ok++;
568           if (MyGlobals::_Verbose>9)
569             {
570               res="0K";
571               if (nbb!=taille)
572                 res="KO";
573               std::cout << res << k << " ";
574             }
575         }
576       res="0K";
577       if (ok!=nb)
578         res="BAD";
579       if (MyGlobals::_Verbose>1) 
580         std::cout << "result " << res << " time(sec) " << MPI_Wtime()-temps_debut << std::endl;
581       MPI_Request_free(&requete1);
582     }
583   //end_time=(MPI_WTIME()-start_time);
584 }
585
586 void MEDPARTITIONER::TestPersistantMpiRing(int taille, int nb)
587 {
588   double temps_debut=MPI_Wtime();
589   int befo, next, rank, wsize, tagbefo, tagnext;
590   rank=MyGlobals::_Rank;
591   wsize=MyGlobals::_World_Size;
592   befo=rank-1; if (befo<0) befo=wsize-1;
593   next=rank+1; if (next>=wsize) next=0;
594   std::vector<int> x, y;
595   tagbefo=111111+befo;
596   tagnext=111111+rank;
597   MPI_Request requete0, requete1;
598   MPI_Status statut1, statut2;
599   int ok=0;
600   std::string res;
601   //cout<<"ini|"<<rank<<'|'<<befo<<'|'<<next<<' ';
602   {
603     x.resize(taille);
604     y.resize(taille);
605     MPI_Ssend_init(&x[0], taille, MPI_INT, next, tagnext, MPI_COMM_WORLD , &requete0);
606     MPI_Recv_init(&y[0], taille,  MPI_INT, befo, tagbefo, MPI_COMM_WORLD , &requete1);
607     for(int k=0; k<nb; k++)
608       {
609         for (int i=0; i<taille; ++i) x[i]=k+rank;
610         //Envoi d’un gros message --> cela peut prendre du temps
611         MPI_Start(&requete0);
612         //Reception du gros message --> cela peut prendre du temps
613         for (int i=0; i<taille; ++i) y[i]=(-1);
614         MPI_Start(&requete1);
615         //Traitement sequentiel independant de "x"
616         //...
617         //Traitement sequentiel independant de "y"
618         //...
619         MPI_Wait(&requete1, &statut1);
620         //Traitement sequentiel dependant de "y"
621         //...=f(y)
622         int nbb=0;
623         for (int i=0; i<taille; ++i)
624           if (y[i]==k+befo)
625             nbb++;
626         if (nbb==taille)
627           ok++;
628         if (MyGlobals::_Verbose>9)
629           {
630             res="0K"+IntToStr(rank);
631             if (nbb!=taille)
632               res="KO"+IntToStr(rank);
633             std::cout << res << k << " ";
634           }
635         MPI_Wait(&requete0, &statut2);
636         //Traitement sequentiel impliquant une modification de "x" en memoire
637         //x=...
638       }
639     res="0K"; if (ok!=nb) res="MAUVAIS";
640     temps_debut=MPI_Wtime()-temps_debut;
641     MPI_Request_free(&requete1);
642     MPI_Request_free(&requete0);
643   }
644   //end_time=(MPI_WTIME()-start_time);
645   if (MyGlobals::_Verbose>1) 
646     std::cout << "result on proc " << rank << " " << res << " time(sec) " << temps_debut << std::endl;
647 }
648
649 void MEDPARTITIONER::TestPersistantMpiRingOnCommSplit(int size, int nb)
650 {
651   double temps_debut=MPI_Wtime();
652   int rank=MyGlobals::_Rank;
653   MPI_Comm newcomm;
654   int color=1;
655   int rankMax=4;
656   if (rank>=rankMax)
657     color=MPI_UNDEFINED;
658   //MPI_Comm_dup (MPI_COMM_WORLD, &newcomm) ;
659   MPI_Comm_split(MPI_COMM_WORLD, color, rank, &newcomm);
660   
661   int befo, next, wsize, tagbefo, tagnext;
662   wsize=rankMax;
663   if (wsize>MyGlobals::_World_Size)
664     wsize=MyGlobals::_World_Size;
665   befo=rank-1;
666   if (befo<0)
667     befo=wsize-1;
668   next=rank+1;
669   if (next>=wsize)
670     next=0;
671   std::vector<int> x, y;
672   tagbefo=111111+befo;
673   tagnext=111111+rank;
674   MPI_Request requete0, requete1;
675   MPI_Status statut1, statut2;
676   int ok=0;
677   std::string res;
678   
679   if (color==1)
680     {
681       x.resize(size);
682       y.resize(size);
683       MPI_Ssend_init(&x[0], size, MPI_INT, next, tagnext, newcomm , &requete0);
684       MPI_Recv_init(&y[0], size,  MPI_INT, befo, tagbefo, newcomm , &requete1);
685       for(int k=0; k<nb; k++)
686         {
687           for (int i=0; i<size; ++i)
688             x[i]=k+rank;
689           //Send of big message --> time consuming
690           MPI_Start(&requete0);
691           //Reception of big message --> time consuming
692           for (int i=0; i<size; ++i)
693             y[i]=-1;
694           MPI_Start(&requete1);
695           //Traitement sequentiel independant de "x"
696           //...
697           //Traitement sequentiel independant de "y"
698           //...
699           //cout<<"dsr|"<<rank<<' ';
700           MPI_Wait(&requete1, &statut1);
701           //Traitement sequentiel dependant de "y"
702           //...=f(y)
703           int nbb=0;
704           for (int i=0; i<size; ++i)
705             if (y[i]==k+befo)
706               nbb++;
707           if (nbb==size)
708             ok++;
709           if (MyGlobals::_Verbose>9)
710             {
711               res="0K"+IntToStr(rank);
712               if (nbb!=size)
713                 res="KO"+IntToStr(rank);
714               std::cout << res << k << " ";
715             }
716           MPI_Wait(&requete0, &statut2);
717           //Traitement sequentiel impliquant une modification de "x" en memoire
718           //x=...
719         }
720       res="0K";
721       if (ok!=nb)
722         res="MAUVAIS";
723       temps_debut=MPI_Wtime()-temps_debut;
724       MPI_Request_free(&requete1);
725       MPI_Request_free(&requete0);
726     }
727   //MPI_Barrier(MPI_COMM_WORLD);
728   if (color==1)
729     MPI_Comm_free(&newcomm);
730   if (MyGlobals::_Verbose>1) 
731     std::cout << "resultat proc " << rank <<" " << res << " time(sec) " << temps_debut << std::endl;
732 }