Salome HOME
Synchronize adm files
[modules/med.git] / src / ParaMEDMEMTest / ParaMEDMEMTest_Gauthier1.cxx
1 // Copyright (C) 2007-2014  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 "ParaMEDMEMTest.hxx"
21 #include <cppunit/TestAssert.h>
22
23 #include "CommInterface.hxx"
24 #include "ProcessorGroup.hxx"
25 #include "MPIProcessorGroup.hxx"
26 #include "DEC.hxx"
27 #include "InterpKernelDEC.hxx"
28 #include "MEDCouplingUMesh.hxx"
29 #include "MEDCouplingFieldDouble.hxx"
30 #include "ParaMESH.hxx"
31 #include "ParaFIELD.hxx"
32 #include "ComponentTopology.hxx"
33 #include "BlockTopology.hxx"
34
35 #include <set>
36 #include <time.h>
37 #include <iostream>
38 #include <assert.h>
39 #include <string>
40 #include <math.h>
41
42 using namespace std;
43 using namespace ParaMEDMEM;
44 using namespace ICoCo;
45
46 void afficheGauthier1(const ParaFIELD& field, const double *vals, int lgth)
47 {
48   const DataArrayDouble *valsOfField(field.getField()->getArray());
49   CPPUNIT_ASSERT_EQUAL(lgth,valsOfField->getNumberOfTuples());
50   for (int ele=0;ele<valsOfField->getNumberOfTuples();ele++)
51     CPPUNIT_ASSERT_DOUBLES_EQUAL(vals[ele],valsOfField->getIJ(ele,0),1e-12);
52 }
53
54 MEDCouplingUMesh *init_quadGauthier1(int is_master)
55 {
56   MEDCouplingAutoRefCountObjectPtr<MEDCouplingUMesh> m(MEDCouplingUMesh::New("champ_quad",2));
57   MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> coo(DataArrayDouble::New());
58   if(is_master)
59     {
60       const double dataCoo[24]={0,0,0,1,0,0,0,0,1,1,0,1,0,1,0,1,1,0,0,1,1,1,1,1};
61       coo->alloc(8,3);
62       std::copy(dataCoo,dataCoo+24,coo->getPointer());
63       const int conn[8]={0,1,3,2,4,5,7,6};
64       m->allocateCells(2);
65       m->insertNextCell(INTERP_KERNEL::NORM_QUAD4,4,conn);
66       m->insertNextCell(INTERP_KERNEL::NORM_QUAD4,4,conn+4);
67     }
68   else
69     {
70       coo->alloc(0,3);
71       m->allocateCells(0);
72     }
73   m->setCoords(coo);
74   return m.retn();
75 }
76
77 MEDCouplingUMesh *init_triangleGauthier1(int is_master)
78 {
79   MEDCouplingAutoRefCountObjectPtr<MEDCouplingUMesh> m(MEDCouplingUMesh::New("champ_triangle",2));
80   MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> coo(DataArrayDouble::New());
81   if(is_master)
82     {
83       const double dataCoo[24]={0,0,0,1,0,0,0,0,1,1,0,1,0,1,0,1,1,0,0,1,1,1,1,1};
84       coo->alloc(8,3);
85       std::copy(dataCoo,dataCoo+24,coo->getPointer());
86       const int conn[12]={0,1,2,1,2,3,4,5,7,4,6,7};
87       m->allocateCells(2);
88       for(int i=0;i<4;i++)
89         m->insertNextCell(INTERP_KERNEL::NORM_TRI3,3,conn+3*i);
90     }
91   else
92     {
93       coo->alloc(0,3);
94       m->allocateCells(0);
95     }
96   m->setCoords(coo);
97   return m.retn();
98 }
99
100
101 void ParaMEDMEMTest::testGauthier1()
102 {
103   int num_cas=0;
104   int rank, size;
105   MPI_Comm_rank(MPI_COMM_WORLD,&rank);
106   MPI_Comm_size(MPI_COMM_WORLD,&size);
107   
108   int is_master=0;
109
110   CommInterface comm;
111   set<int> emetteur_ids;
112   set<int> recepteur_ids;
113   emetteur_ids.insert(0);
114   if(size!=4)
115     return;
116   recepteur_ids.insert(1);
117   if (size >2) 
118     recepteur_ids.insert(2);
119   if (size >2) 
120     emetteur_ids.insert(3);
121   if ((rank==0)||(rank==1)) 
122     is_master=1;
123   
124   MPIProcessorGroup recepteur_group(comm,recepteur_ids);
125   MPIProcessorGroup emetteur_group(comm,emetteur_ids);
126
127   string cas;
128   if (recepteur_group.containsMyRank())
129     {
130       cas="recepteur";
131       //freopen("recpeteur.out","w",stdout);
132       //freopen("recepteur.err","w",stderr);
133     }
134   else
135     {
136       cas="emetteur";
137       // freopen("emetteur.out","w",stdout);
138       //freopen("emetteur.err","w",stderr);
139     }
140   double expected[8][4]={
141     {1.,1.,1.,1.},
142     {40., 40., 1., 1.},
143     {1.,1.,1e200,1e200},
144     {40.,1.,1e200,1e200},
145     {1.,1.,1.,1.},
146     {40.,1.,1.,1.},
147     {1.,1.,1e200,1e200},
148     {20.5,1.,1e200,1e200}
149   };
150   int expectedLgth[8]={4,4,2,2,4,4,2,2};
151   
152   for (int send=0;send<2;send++)
153     for (int rec=0;rec<2;rec++)
154       {
155         InterpKernelDEC dec_emetteur(emetteur_group, recepteur_group);
156         ParaMEDMEM::ParaFIELD *champ_emetteur(0),*champ_recepteur(0);
157         ParaMEDMEM::ParaMESH *paramesh(0);
158         MEDCouplingAutoRefCountObjectPtr<MEDCouplingUMesh> mesh;
159         dec_emetteur.setOrientation(2);
160         if (send==0)
161           {
162             mesh=init_quadGauthier1(is_master);
163           }
164         else
165           {
166             mesh=init_triangleGauthier1(is_master);
167           }
168         paramesh=new ParaMEDMEM::ParaMESH(mesh,recepteur_group.containsMyRank()?recepteur_group:emetteur_group,"emetteur mesh");
169         ParaMEDMEM::ComponentTopology comptopo;
170         champ_emetteur=new ParaMEDMEM::ParaFIELD(ON_CELLS,ONE_TIME,paramesh,comptopo);
171         champ_emetteur->getField()->setNature(ConservativeVolumic);
172         champ_emetteur->setOwnSupport(true);
173         if (rec==0)
174           {
175             mesh=init_triangleGauthier1(is_master);
176           }
177         else
178           {
179             mesh=init_quadGauthier1(is_master);
180           }
181         paramesh=new ParaMEDMEM::ParaMESH(mesh,recepteur_group.containsMyRank()?recepteur_group:emetteur_group,"recepteur mesh");
182         champ_recepteur=new ParaMEDMEM::ParaFIELD(ON_CELLS,ONE_TIME,paramesh,comptopo);
183         champ_recepteur->getField()->setNature(ConservativeVolumic);
184         champ_recepteur->setOwnSupport(true);
185         if (cas=="emetteur") 
186           {
187             champ_emetteur->getField()->getArray()->fillWithValue(1.);
188           }
189   
190   
191         MPI_Barrier(MPI_COMM_WORLD);
192
193         //clock_t clock0= clock ();
194         int compti=0;
195
196         bool init=true; // first time step ??
197         bool stop=false;
198         //boucle sur les pas de quads
199         while (!stop) {
200   
201           compti++;
202           //clock_t clocki= clock ();
203           //cout << compti << " CLOCK " << (clocki-clock0)*1.e-6 << endl; 
204           for (int non_unif=0;non_unif<2;non_unif++)
205             {
206               if (cas=="emetteur") 
207                 {
208                   if (non_unif)
209                     if(rank!=3)
210                       champ_emetteur->getField()->getArray()->setIJ(0,0,40);
211                 }
212               //bool ok=false; // Is the time interval successfully solved ?
213     
214               // Loop on the time interval tries
215               if(1) {
216       
217
218                 if (cas=="emetteur")
219                   dec_emetteur.attachLocalField(champ_emetteur);
220                 else
221                   dec_emetteur.attachLocalField(champ_recepteur);
222
223
224                 if(init) dec_emetteur.synchronize();
225                 init=false;
226
227                 if (cas=="emetteur") {
228                   //    affiche(champ_emetteur);
229                   dec_emetteur.sendData();
230                 }
231                 else if (cas=="recepteur")
232                   {
233                     dec_emetteur.recvData();
234                     if (is_master)
235                       afficheGauthier1(*champ_recepteur,expected[num_cas],expectedLgth[num_cas]);
236                   }
237                 else
238                   throw 0;
239                 MPI_Barrier(MPI_COMM_WORLD);
240               }
241               stop=true;
242               num_cas++;
243             }
244         }
245         delete champ_emetteur;
246         delete champ_recepteur;
247       }
248 }
249
250 void ParaMEDMEMTest::testGauthier2()
251 {
252   double valuesExpected1[2]={0.,0.};
253   double valuesExpected2[2]={0.95,0.970625};
254   
255   double valuesExpected30[]={0., 0., 0.05, 0., 0., 0.15, 0., 0., 0.25, 0., 0., 0.35, 0., 0., 0.45, 0., 0., 0.55, 0., 0., 0.65, 0., 0., 0.75, 0., 0., 0.85, 0., 0., 0.95};
256   double valuesExpected31[]={0.,  0.,  0.029375,  0.,  0.,  0.029375,  0.,  0.,  0.1,  0.,  0.,  0.1,  0.,  0.,  0.2,  0.,  0.,  0.2,  0.,  0.,  0.3,  0.,  0.,  0.3,  0.,  0.,  0.4,  0.,  0.,  0.4,  0.,  0.,  0.5,  0.,  0.,  0.5,  0.,  0.,  0.6,  0.,  0.,  0.6,  0.,  0.,  0.7,  0.,  0.,  0.7,  0.,  0.,  0.8,  0.,  0.,  0.8,  0.,  0.,  0.9,  0.,  0.,  0.9,  0.,  0.,  0.970625,  0.,  0.,  0.970625 };
257
258   double *valuesExpected3[2]={valuesExpected30,valuesExpected31};
259
260   int rank, size;
261   MPI_Comm_rank(MPI_COMM_WORLD,&rank);
262   MPI_Comm_size(MPI_COMM_WORLD,&size);
263   if (size <2)
264     return ;
265   CommInterface comm;
266   set<int> Genepi_ids;
267   set<int> entree_chaude_ids;
268   Genepi_ids.insert(0);
269   for (int i=1;i<size;i++)
270     entree_chaude_ids.insert(i);
271   for (int type=0;type<2;type++)
272     {
273       MPIProcessorGroup entree_chaude_group(comm,entree_chaude_ids);
274       MPIProcessorGroup Genepi_group(comm,Genepi_ids);
275
276       ParaMEDMEM::ParaFIELD *vitesse(0);
277       InterpKernelDEC dec_vit_in_chaude(entree_chaude_group, Genepi_group);
278
279       if ( entree_chaude_group.containsMyRank())
280         {
281           MEDCouplingAutoRefCountObjectPtr<MEDCouplingUMesh> mesh(MEDCouplingUMesh::New("mesh",2));
282           MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> arr(DataArrayDouble::New()); arr->alloc(63,3);
283           const double cooData[189]={0.,0.,0.,0.5,0.,0.,0.5,0.05,0.,0.,0.1,0.,0.5,0.1,0.,0.5,0.15,0.,0.,0.2,0.,0.5,0.2,0.,0.5,0.25,0.,0.,0.3,0.,0.5,0.3,0.,0.5,0.35,0.,0.,0.4,0.,0.5,0.4,0.,0.5,0.45,0.,0.,0.5,0.,0.5,0.5,0.,0.5,0.55,0.,0.,0.6,0.,0.5,0.6,0.,0.5,0.65,0.,0.,0.7,0.,0.5,0.7,0.,0.5,0.75,0.,0.,0.8,0.,0.5,0.8,0.,0.5,0.85,0.,0.,0.9,0.,0.5,0.9,0.,0.5,0.95,0.,1.,0.,0.,1.,0.1,0.,1.,0.2,0.,1.,0.3,0.,1.,0.4,0.,1.,0.5,0.,1.,0.6,0.,1.,0.7,0.,1.,0.8,0.,1.,0.9,0.,1.,0.05,0.,1.,0.15,0.,1.,0.25,0.,1.,0.35,0.,1.,0.45,0.,1.,0.55,0.,1.,0.65,0.,1.,0.75,0.,1.,0.85,0.,1.,0.95,0.,1.,1.,0.,0.,1.,0.,0.5,1.,0.,0.,0.05,0.,0.,0.15,0.,0.,0.25,0.,0.,0.35,0.,0.,0.45,0.,0.,0.55,0.,0.,0.65,0.,0.,0.75,0.,0.,0.85,0.,0.,0.95,0.};
284           std::copy(cooData,cooData+189,arr->getPointer());
285           mesh->setCoords(arr);
286           mesh->allocateCells(80);
287           const int conn[240]={0,1,2,3,4,5,6,7,8,9,10,11,12,13,14,15,16,17,18,19,20,21,22,23,24,25,26,27,28,29,30,2,1,31,5,4,32,8,7,33,11,10,34,14,13,35,17,16,36,20,19,37,23,22,38,26,25,39,29,28,30,40,2,31,41,5,32,42,8,33,43,11,34,44,14,35,45,17,36,46,20,37,47,23,38,48,26,39,49,29,31,2,40,32,5,41,33,8,42,34,11,43,35,14,44,36,17,45,37,20,46,38,23,47,39,26,48,50,29,49,3,2,4,6,5,7,9,8,10,12,11,13,15,14,16,18,17,19,21,20,22,24,23,25,27,26,28,51,29,52,31,4,2,32,7,5,33,10,8,34,13,11,35,16,14,36,19,17,37,22,20,38,25,23,39,28,26,50,52,29,0,2,53,3,5,54,6,8,55,9,11,56,12,14,57,15,17,58,18,20,59,21,23,60,24,26,61,27,29,62,3,53,2,6,54,5,9,55,8,12,56,11,15,57,14,18,58,17,21,59,20,24,60,23,27,61,26,51,62,29};
288           for(int i=0;i<80;i++)
289             mesh->insertNextCell(INTERP_KERNEL::NORM_TRI3,3,conn+3*i);
290           MEDCouplingAutoRefCountObjectPtr<MEDCouplingFieldDouble> f(MEDCouplingFieldDouble::New(ON_NODES,ONE_TIME));
291           const double valsOfField[189]={0.,0.,0.,0.,0.,0.,0.,0.,0.05,0.,0.,0.1,0.,0.,0.1,0.,0.,0.15,0.,0.,0.2,0.,0.,0.2,0.,0.,0.25,0.,0.,0.3,0.,0.,0.3,0.,0.,0.35,0.,0.,0.4,0.,0.,0.4,0.,0.,0.45,0.,0.,0.5,0.,0.,0.5,0.,0.,0.55,0.,0.,0.6,0.,0.,0.6,0.,0.,0.65,0.,0.,0.7,0.,0.,0.7,0.,0.,0.75,0.,0.,0.8,0.,0.,0.8,0.,0.,0.85,0.,0.,0.9,0.,0.,0.9,0.,0.,0.95,0.,0.,0.,0.,0.,0.1,0.,0.,0.2,0.,0.,0.3,0.,0.,0.4,0.,0.,0.5,0.,0.,0.6,0.,0.,0.7,0.,0.,0.8,0.,0.,0.9,0.,0.,0.05,0.,0.,0.15,0.,0.,0.25,0.,0.,0.35,0.,0.,0.45,0.,0.,0.55,0.,0.,0.65,0.,0.,0.75,0.,0.,0.85,0.,0.,0.95,0.,0.,1.,0.,0.,1.,0.,0.,1.,0.,0.,0.05,0.,0.,0.15,0.,0.,0.25,0.,0.,0.35,0.,0.,0.45,0.,0.,0.55,0.,0.,0.65,0.,0.,0.75,0.,0.,0.85,0.,0.,0.95};
292           f->setMesh(mesh); f->setName("VITESSE_P1_OUT");
293           arr=DataArrayDouble::New(); arr->alloc(63,3);
294           std::copy(valsOfField,valsOfField+189,arr->getPointer());
295           f->setArray(arr); f->setNature(ConservativeVolumic);
296           ParaMEDMEM::ParaMESH *paramesh(new ParaMEDMEM::ParaMESH(mesh,entree_chaude_group,"emetteur mesh"));
297           vitesse=new ParaMEDMEM::ParaFIELD(f,paramesh,entree_chaude_group);
298           vitesse->setOwnSupport(true);
299           dec_vit_in_chaude.setMethod("P1");
300         }
301       else
302         {
303           MEDCouplingAutoRefCountObjectPtr<MEDCouplingUMesh> mesh(MEDCouplingUMesh::New("mesh",2));
304           MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> arr(DataArrayDouble::New()); arr->alloc(22,3);
305           const double cooData[66]={0,0,0,1,0,0,0,0.1,0,1,0.1,0,0,0.2,0,1,0.2,0,0,0.3,0,1,0.3,0,0,0.4,0,1,0.4,0,0,0.5,0,1,0.5,0,0,0.6,0,1,0.6,0,0,0.7,0,1,0.7,0,0,0.8,0,1,0.8,0,0,0.9,0,1,0.9,0,0,1,0,1,1,0};
306           std::copy(cooData,cooData+66,arr->getPointer());
307           mesh->setCoords(arr);
308           mesh->allocateCells(10);
309           const int conn[40]={0,1,3,2,2,3,5,4,4,5,7,6,6,7,9,8,8,9,11,10,10,11,13,12,12,13,15,14,14,15,17,16,16,17,19,18,18,19,21,20};
310           for(int i=0;i<10;i++)
311             mesh->insertNextCell(INTERP_KERNEL::NORM_QUAD4,4,conn+4*i);
312           MEDCouplingAutoRefCountObjectPtr<MEDCouplingFieldDouble> f(MEDCouplingFieldDouble::New(type==0?ON_CELLS:ON_NODES,ONE_TIME));
313           f->setMesh(mesh); f->setName("vitesse_in_chaude");
314           arr=DataArrayDouble::New(); arr->alloc(f->getNumberOfTuplesExpected()*3); arr->fillWithZero(); arr->rearrange(3);
315           f->setArray(arr); f->setNature(ConservativeVolumic);
316           ParaMEDMEM::ParaMESH *paramesh(new ParaMEDMEM::ParaMESH(mesh,Genepi_group,"recepteur mesh"));
317           vitesse=new ParaMEDMEM::ParaFIELD(f,paramesh,Genepi_group);
318           vitesse->setOwnSupport(true);
319           dec_vit_in_chaude.setMethod(f->getDiscretization()->getRepr());
320         }
321
322       dec_vit_in_chaude.attachLocalField(vitesse);
323       
324       dec_vit_in_chaude.synchronize();
325   
326   
327       // Envois - receptions
328       if (entree_chaude_group.containsMyRank())
329         {
330           dec_vit_in_chaude.sendData();
331         }
332       else
333         {
334           dec_vit_in_chaude.recvData(); 
335         }
336       if ( !entree_chaude_group.containsMyRank() )
337         {
338           double pmin=1e38, pmax=-1e38;
339           const double *p(vitesse->getField()->getArray()->begin());
340           for(std::size_t i=0;i<vitesse->getField()->getArray()->getNbOfElems();i++,p++)
341             {
342               if (*p<pmin) pmin=*p;
343               if (*p>pmax) pmax=*p;
344             }
345           CPPUNIT_ASSERT_DOUBLES_EQUAL(valuesExpected1[type],pmin,1e-12);
346           CPPUNIT_ASSERT_DOUBLES_EQUAL(valuesExpected2[type],pmax,1e-12);
347       
348           int nbCompo(vitesse->getField()->getNumberOfComponents());
349           p=vitesse->getField()->getArray()->begin();
350           for(int i=0;i<vitesse->getField()->getNumberOfTuples();i++)
351             for(int c=0;c<nbCompo;c++,p++)
352               CPPUNIT_ASSERT_DOUBLES_EQUAL(valuesExpected3[type][i*nbCompo+c],*p,1e-12);
353         }
354       delete vitesse;
355     }
356 }
357
358 /*!
359  * Non regression test testing copy constructor of InterpKernelDEC. 
360  */
361 void ParaMEDMEMTest::testGauthier3()
362 {
363   int num_cas=0;
364   int rank, size;
365   MPI_Comm_rank(MPI_COMM_WORLD,&rank);
366   MPI_Comm_size(MPI_COMM_WORLD,&size);
367   
368   int is_master=0;
369
370   CommInterface comm;
371   set<int> emetteur_ids;
372   set<int> recepteur_ids;
373   emetteur_ids.insert(0);
374   if(size!=4)
375     return;
376   recepteur_ids.insert(1);
377   if (size >2) 
378     recepteur_ids.insert(2);
379   if (size >2) 
380     emetteur_ids.insert(3);
381   if ((rank==0)||(rank==1)) 
382     is_master=1;
383   
384   MPIProcessorGroup recepteur_group(comm,recepteur_ids);
385   MPIProcessorGroup emetteur_group(comm,emetteur_ids);
386
387   string cas;
388   if (recepteur_group.containsMyRank())
389     {
390       cas="recepteur";
391       //freopen("recpeteur.out","w",stdout);
392       //freopen("recepteur.err","w",stderr);
393     }
394   else
395     {
396       cas="emetteur";
397       // freopen("emetteur.out","w",stdout);
398       //freopen("emetteur.err","w",stderr);
399     }
400   double expected[8][4]={
401     {1.,1.,1.,1.},
402     {40., 40., 1., 1.},
403     {1.,1.,1e200,1e200},
404     {40.,1.,1e200,1e200},
405     {1.,1.,1.,1.},
406     {40.,1.,1.,1.},
407     {1.,1.,1e200,1e200},
408     {20.5,1.,1e200,1e200}
409   };
410   int expectedLgth[8]={4,4,2,2,4,4,2,2};
411   
412   for (int send=0;send<2;send++)
413     for (int rec=0;rec<2;rec++)
414       {
415         std::vector<InterpKernelDEC> decu(1);
416         decu[0]=InterpKernelDEC(emetteur_group,recepteur_group);
417         InterpKernelDEC& dec_emetteur=decu[0];
418         ParaMEDMEM::ParaFIELD *champ_emetteur(0),*champ_recepteur(0);
419         ParaMEDMEM::ParaMESH *paramesh(0);
420         MEDCouplingAutoRefCountObjectPtr<MEDCouplingUMesh> mesh;
421         dec_emetteur.setOrientation(2);
422         if (send==0)
423           {
424             mesh=init_quadGauthier1(is_master);
425           }
426         else
427           {
428             mesh=init_triangleGauthier1(is_master);
429           }
430         paramesh=new ParaMEDMEM::ParaMESH(mesh,recepteur_group.containsMyRank()?recepteur_group:emetteur_group,"emetteur mesh");
431         ParaMEDMEM::ComponentTopology comptopo;
432         champ_emetteur=new ParaMEDMEM::ParaFIELD(ON_CELLS,ONE_TIME,paramesh,comptopo);
433         champ_emetteur->getField()->setNature(ConservativeVolumic);
434         champ_emetteur->setOwnSupport(true);
435         if (rec==0)
436           {
437             mesh=init_triangleGauthier1(is_master);
438           }
439         else
440           {
441             mesh=init_quadGauthier1(is_master);
442           }
443         paramesh=new ParaMEDMEM::ParaMESH(mesh,recepteur_group.containsMyRank()?recepteur_group:emetteur_group,"recepteur mesh");
444         champ_recepteur=new ParaMEDMEM::ParaFIELD(ON_CELLS,ONE_TIME,paramesh,comptopo);
445         champ_recepteur->getField()->setNature(ConservativeVolumic);
446         champ_recepteur->setOwnSupport(true);
447         if (cas=="emetteur") 
448           {
449             champ_emetteur->getField()->getArray()->fillWithValue(1.);
450           }
451   
452   
453         MPI_Barrier(MPI_COMM_WORLD);
454
455         //clock_t clock0= clock ();
456         int compti=0;
457
458         bool init=true; // first time step ??
459         bool stop=false;
460         //boucle sur les pas de quads
461         while (!stop) {
462   
463           compti++;
464           //clock_t clocki= clock ();
465           //cout << compti << " CLOCK " << (clocki-clock0)*1.e-6 << endl; 
466           for (int non_unif=0;non_unif<2;non_unif++)
467             {
468               if (cas=="emetteur") 
469                 {
470                   if (non_unif)
471                     if(rank!=3)
472                       champ_emetteur->getField()->getArray()->setIJ(0,0,40);
473                 }
474               //bool ok=false; // Is the time interval successfully solved ?
475     
476               // Loop on the time interval tries
477               if(1) {
478       
479
480                 if (cas=="emetteur")
481                   dec_emetteur.attachLocalField(champ_emetteur);
482                 else
483                   dec_emetteur.attachLocalField(champ_recepteur);
484
485
486                 if(init) dec_emetteur.synchronize();
487                 init=false;
488
489                 if (cas=="emetteur") {
490                   //    affiche(champ_emetteur);
491                   dec_emetteur.sendData();
492                 }
493                 else if (cas=="recepteur")
494                   {
495                     dec_emetteur.recvData();
496                     if (is_master)
497                       afficheGauthier1(*champ_recepteur,expected[num_cas],expectedLgth[num_cas]);
498                   }
499                 else
500                   throw 0;
501                 MPI_Barrier(MPI_COMM_WORLD);
502               }
503               stop=true;
504               num_cas++;
505             }
506         }
507         delete champ_emetteur;
508         delete champ_recepteur;
509       }
510 }
511
512 /*!
513  * This test is the parallel version of MEDCouplingBasicsTest.test3D1DOnP1P0_1 test.
514  */
515 void ParaMEDMEMTest::testGauthier4()
516 {
517   //
518   const double sourceCoords[19*3]={0.5,0.5,0.1,0.5,0.5,1.2,0.5,0.5,1.6,0.5,0.5,1.8,0.5,0.5,2.43,0.5,0.5,2.55,0.5,0.5,4.1,0.5,0.5,4.4,0.5,0.5,4.9,0.5,0.5,5.1,0.5,0.5,7.6,0.5,0.5,7.7,0.5,0.5,8.2,0.5,0.5,8.4,0.5,0.5,8.6,0.5,0.5,8.8,0.5,0.5,9.2,0.5,0.5,9.6,0.5,0.5,11.5};
519   const int sourceConn[18*2]={0,1,1,2,2,3,3,4,4,5,5,6,6,7,7,8,8,9,9,10,10,11,11,12,12,13,13,14,14,15,15,16,16,17,17,18};
520   const double sourceVals[19]={0.49,2.8899999999999997,7.29,13.69,22.09,32.49,44.89,59.29,75.69,94.09, 114.49,136.89,161.29,187.69,216.09,246.49,278.89,313.29,349.69};
521   const double targetCoords0[20*3]={0.,0.,0.,1.,0.,0.,0.,1.,0.,1.,1.,0.,0.,0.,1.,1.,0.,1.,0.,1.,1.,1.,1.,1.,0.,0.,2.,1.,0.,2.,0.,1.,2.,1.,1.,2.,0.,0.,3.,1.,0.,3.,0.,1.,3.,1.,1.,3.,0.,0.,4.,1.,0.,4.,0.,1.,4.,1.,1.,4.};
522   const int targetConn0[8*4]={1,0,2,3,5,4,6,7,5,4,6,7,9,8,10,11,9,8,10,11,13,12,14,15,13,12,14,15,17,16,18,19};
523   const double targetCoords1[28*3]={0.,0.,4.,1.,0.,4.,0.,1.,4.,1.,1.,4.,0.,0.,5.,1.,0.,5.,0.,1.,5.,1.,1.,5.,0.,0.,6.,1.,0.,6.,0.,1.,6.,1.,1.,6.,0.,0.,7.,1.,0.,7.,0.,1.,7.,1.,1.,7.,0.,0.,8.,1.,0.,8.,0.,1.,8.,1.,1.,8.,0.,0.,9.,1.,0.,9.,0.,1.,9.,1.,1.,9.,0.,0.,10.,1.,0.,10.,0.,1.,10.,1.,1.,10.};
524   const int targetConn1[8*6]={1,0,2,3,5,4,6,7,5,4,6,7,9,8,10,11,9,8,10,11,13,12,14,15,13,12,14,15,17,16,18,19,17,16,18,19,21,20,22,23,21,20,22,23,25,24,26,27};
525   //
526   int size;
527   int rank;
528   MPI_Comm_size(MPI_COMM_WORLD,&size);
529   MPI_Comm_rank(MPI_COMM_WORLD,&rank);
530   //
531   if(size!=3)
532     return ;
533   int nproc_source = 1;
534   set<int> self_procs;
535   set<int> procs_source;
536   set<int> procs_target;
537   
538   for (int i=0; i<nproc_source; i++)
539     procs_source.insert(i);
540   for (int i=nproc_source; i<size; i++)
541     procs_target.insert(i);
542   self_procs.insert(rank);
543   //
544   ParaMEDMEM::MEDCouplingUMesh *mesh=0;
545   ParaMEDMEM::ParaMESH *paramesh=0;
546   ParaMEDMEM::ParaFIELD* parafield=0;
547   //
548   ParaMEDMEM::CommInterface interface;
549   //
550   ProcessorGroup* self_group = new ParaMEDMEM::MPIProcessorGroup(interface,self_procs);
551   ProcessorGroup* target_group = new ParaMEDMEM::MPIProcessorGroup(interface,procs_target);
552   ProcessorGroup* source_group = new ParaMEDMEM::MPIProcessorGroup(interface,procs_source);
553   //
554   MPI_Barrier(MPI_COMM_WORLD);
555   if(source_group->containsMyRank())
556     {
557       std::ostringstream stream; stream << "sourcemesh2D proc " << rank;
558       mesh=MEDCouplingUMesh::New(stream.str().c_str(),1);
559       mesh->allocateCells();
560       for(int i=0;i<18;i++)
561         mesh->insertNextCell(INTERP_KERNEL::NORM_SEG2,2,sourceConn+2*i);
562       mesh->finishInsertingCells();
563       DataArrayDouble *myCoords=DataArrayDouble::New();
564       myCoords->alloc(19,3);
565       std::copy(sourceCoords,sourceCoords+19*3,myCoords->getPointer());
566       mesh->setCoords(myCoords);
567       myCoords->decrRef();
568       paramesh=new ParaMESH(mesh,*source_group,"source mesh");
569       ParaMEDMEM::ComponentTopology comptopo;
570       parafield = new ParaFIELD(ON_NODES,NO_TIME,paramesh,comptopo);
571       double *value=parafield->getField()->getArray()->getPointer();
572       std::copy(sourceVals,sourceVals+19,value);
573     }
574   else
575     {
576       if(rank==1)
577         {
578           std::ostringstream stream; stream << "targetmesh2D proc " << rank-nproc_source;
579           mesh=MEDCouplingUMesh::New(stream.str().c_str(),3);
580           mesh->allocateCells();
581           for(int i=0;i<4;i++)
582             mesh->insertNextCell(INTERP_KERNEL::NORM_HEXA8,8,targetConn0+8*i);
583           mesh->finishInsertingCells();
584           DataArrayDouble *myCoords=DataArrayDouble::New();
585           myCoords->alloc(20,3);
586           std::copy(targetCoords0,targetCoords0+20*3,myCoords->getPointer());
587           mesh->setCoords(myCoords);
588           myCoords->decrRef();
589           paramesh=new ParaMESH (mesh,*target_group,"target mesh");
590           ParaMEDMEM::ComponentTopology comptopo;
591           parafield = new ParaFIELD(ON_CELLS,NO_TIME,paramesh, comptopo);
592         }
593       else if(rank==2)
594         {
595           std::ostringstream stream; stream << "targetmesh2D proc " << rank-nproc_source;
596           mesh=MEDCouplingUMesh::New(stream.str().c_str(),3);
597           mesh->allocateCells();
598           for(int i=0;i<6;i++)
599             mesh->insertNextCell(INTERP_KERNEL::NORM_HEXA8,8,targetConn1+8*i);
600           mesh->finishInsertingCells();
601           DataArrayDouble *myCoords=DataArrayDouble::New();
602           myCoords->alloc(28,3);
603           std::copy(targetCoords1,targetCoords1+28*3,myCoords->getPointer());
604           mesh->setCoords(myCoords);
605           myCoords->decrRef();
606           paramesh=new ParaMESH (mesh,*target_group,"target mesh");
607           ParaMEDMEM::ComponentTopology comptopo;
608           parafield = new ParaFIELD(ON_CELLS,NO_TIME,paramesh, comptopo);
609         }
610     }
611   //test 1 - primaire -> secondaire
612   ParaMEDMEM::InterpKernelDEC dec(*source_group,*target_group);
613   dec.setIntersectionType(INTERP_KERNEL::PointLocator);
614   parafield->getField()->setNature(ConservativeVolumic);//very important
615   if (source_group->containsMyRank())
616     { 
617       dec.setMethod("P1");
618       dec.attachLocalField(parafield);
619       dec.synchronize();
620       dec.setForcedRenormalization(false);
621       dec.sendData();
622     }
623   else
624     {
625       dec.setMethod("P0");
626       dec.attachLocalField(parafield);
627       dec.synchronize();
628       dec.setForcedRenormalization(false);
629       dec.recvData();
630       const double *res(parafield->getField()->getArray()->getConstPointer());
631       if(rank==1)
632         {
633           const double expected0[4]={0.49,7.956666666666667,27.29,0.};
634           for(int i=0;i<4;i++)
635             CPPUNIT_ASSERT_DOUBLES_EQUAL(expected0[i],res[i],1e-13);
636         }
637       else
638         {
639           const double expected1[6]={59.95666666666667,94.09,0.,125.69,202.89,296.09};
640           for(int i=0;i<6;i++)
641             CPPUNIT_ASSERT_DOUBLES_EQUAL(expected1[i],res[i],1e-13);
642         }
643     }
644   MPI_Barrier(MPI_COMM_WORLD);
645   if (source_group->containsMyRank())
646     {
647       dec.recvData();
648       const double expected2[19]={0.49,7.956666666666667,7.956666666666667,7.956666666666667,27.29,27.29,59.95666666666667,59.95666666666667,59.95666666666667,94.09,125.69,125.69,202.89,202.89,202.89,202.89,296.09,296.09,0.};
649       const double *res(parafield->getField()->getArray()->getConstPointer());
650       for(int i=0;i<19;i++)
651         CPPUNIT_ASSERT_DOUBLES_EQUAL(expected2[i],res[i],1e-13);
652     }
653   else
654     {
655       dec.sendData();
656     }
657   delete parafield;
658   mesh->decrRef();
659   delete paramesh;
660   delete self_group;
661   delete target_group;
662   delete source_group;
663   //
664   MPI_Barrier(MPI_COMM_WORLD);
665 }