Salome HOME
OverlapDEC: many improvements/fixes:
[modules/med.git] / src / ParaMEDMEMTest / ParaMEDMEMTest_OverlapDEC.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 "ParaMEDMEMTest.hxx"
21 #include <cppunit/TestAssert.h>
22
23 #include "CommInterface.hxx"
24 #include "ProcessorGroup.hxx"
25 #include "MPIProcessorGroup.hxx"
26 #include "Topology.hxx"
27 #include "OverlapDEC.hxx"
28 #include "ParaMESH.hxx"
29 #include "ParaFIELD.hxx"
30 #include "ComponentTopology.hxx"
31
32 #include "MEDCouplingUMesh.hxx"
33
34 #include <set>
35
36 using namespace std;
37
38 #include "MEDCouplingAutoRefCountObjectPtr.hxx"
39 #include "MEDLoader.hxx"
40 #include "MEDLoaderBase.hxx"
41 #include "MEDCouplingFieldDouble.hxx"
42 #include "MEDCouplingMemArray.hxx"
43 #include "MEDCouplingRemapper.hxx"
44
45 using namespace ParaMEDMEM;
46
47 typedef  MEDCouplingAutoRefCountObjectPtr<MEDCouplingUMesh> MUMesh;
48 typedef  MEDCouplingAutoRefCountObjectPtr<MEDCouplingFieldDouble> MFDouble;
49
50 //void ParaMEDMEMTest::testOverlapDEC_LMEC_seq()
51 //{
52 //  //  T_SC_Trio_src.med  -- "SupportOf_"
53 //  //  T_SC_Trio_dst.med  -- "SupportOf_T_SC_Trio"
54 //  //  h_TH_Trio_src.med  -- "SupportOf_"
55 //  //  h_TH_Trio_dst.med  -- "SupportOf_h_TH_Trio"
56 //  string rep("/export/home/adrien/support/antoine_LMEC/");
57 //  string src_mesh_nam(rep + string("T_SC_Trio_src.med"));
58 //  string tgt_mesh_nam(rep + string("T_SC_Trio_dst.med"));
59 ////  string src_mesh_nam(rep + string("h_TH_Trio_src.med"));
60 ////  string tgt_mesh_nam(rep + string("h_TH_Trio_dst.med"));
61 //  MUMesh src_mesh=MEDLoader::ReadUMeshFromFile(src_mesh_nam,"SupportOf_",0);
62 //  MUMesh tgt_mesh=MEDLoader::ReadUMeshFromFile(tgt_mesh_nam,"SupportOf_T_SC_Trio",0);
63 ////  MUMesh tgt_mesh=MEDLoader::ReadUMeshFromFile(tgt_mesh_nam,"SupportOf_h_TH_Trio",0);
64 //
65 //  MFDouble srcField = MEDCouplingFieldDouble::New(ON_CELLS, ONE_TIME);
66 //  srcField->setMesh(src_mesh);
67 //  DataArrayDouble * dad = DataArrayDouble::New(); dad->alloc(src_mesh->getNumberOfCells(),1);
68 //  dad->fillWithValue(1.0);
69 //  srcField->setArray(dad);
70 //  srcField->setNature(ConservativeVolumic);
71 //
72 //  MEDCouplingRemapper remap;
73 //  remap.setOrientation(2); // always consider surface intersections as absolute areas.
74 //  remap.prepare(src_mesh, tgt_mesh, "P0P0");
75 //  MFDouble tgtField = remap.transferField(srcField, 1.0e+300);
76 //  tgtField->setName("result");
77 //  string out_nam(rep + string("adrien.med"));
78 //  MEDLoader::WriteField(out_nam,tgtField, true);
79 //  cout << "wrote: " << out_nam << "\n";
80 //  double integ1 = 0.0, integ2 = 0.0;
81 //  srcField->integral(true, &integ1);
82 //  tgtField->integral(true, &integ2);
83 ////  tgtField->reprQuickOverview(cout);
84 //  CPPUNIT_ASSERT_DOUBLES_EQUAL(integ1,integ2,1e-8);
85 //
86 //  dad->decrRef();
87 //}
88 //
89 //void ParaMEDMEMTest::testOverlapDEC_LMEC_para()
90 //{
91 //  using namespace ParaMEDMEM;
92 //
93 //  int size;
94 //  int rank;
95 //  MPI_Comm_size(MPI_COMM_WORLD,&size);
96 //  MPI_Comm_rank(MPI_COMM_WORLD,&rank);
97 //
98 //  if (size != 1) return ;
99 //
100 //  int nproc = 1;
101 //  std::set<int> procs;
102 //
103 //  for (int i=0; i<nproc; i++)
104 //    procs.insert(i);
105 //
106 //  CommInterface interface;
107 //  OverlapDEC dec(procs);
108 //
109 //  ParaMESH* parameshS=0;
110 //  ParaMESH* parameshT=0;
111 //  ParaFIELD* parafieldS=0;
112 //  ParaFIELD* parafieldT=0;
113 //  MFDouble srcField;
114 //
115 //  // **** FILE LOADING
116 //  //  T_SC_Trio_src.med  -- "SupportOf_"
117 //  //  T_SC_Trio_dst.med  -- "SupportOf_T_SC_Trio"
118 //  //  h_TH_Trio_src.med  -- "SupportOf_"
119 //  //  h_TH_Trio_dst.med  -- "SupportOf_h_TH_Trio"
120 //  string rep("/export/home/adrien/support/antoine_LMEC/");
121 //  string src_mesh_nam(rep + string("T_SC_Trio_src.med"));
122 //  string tgt_mesh_nam(rep + string("T_SC_Trio_dst.med"));
123 //
124 //
125 //  MPI_Barrier(MPI_COMM_WORLD);
126 //  if(rank==0)
127 //    {
128 //    //  string src_mesh_nam(rep + string("h_TH_Trio_src.med"));
129 //    //  string tgt_mesh_nam(rep + string("h_TH_Trio_dst.med"));
130 //      MUMesh src_mesh=MEDLoader::ReadUMeshFromFile(src_mesh_nam,"SupportOf_",0);
131 //      MUMesh tgt_mesh=MEDLoader::ReadUMeshFromFile(tgt_mesh_nam,"SupportOf_T_SC_Trio",0);
132 //    //  MUMesh tgt_mesh=MEDLoader::ReadUMeshFromFile(tgt_mesh_nam,"SupportOf_h_TH_Trio",0);
133 //
134 //      // **** SOURCE
135 //      srcField = MEDCouplingFieldDouble::New(ON_CELLS, ONE_TIME);
136 //      srcField->setMesh(src_mesh);
137 //      DataArrayDouble * dad = DataArrayDouble::New(); dad->alloc(src_mesh->getNumberOfCells(),1);
138 //      dad->fillWithValue(1.0);
139 //      srcField->setArray(dad);
140 //      srcField->setNature(ConservativeVolumic);
141 //
142 //      ComponentTopology comptopo;
143 //      parameshS = new ParaMESH(src_mesh,*dec.getGroup(),"source mesh");
144 //      parafieldS = new ParaFIELD(ON_CELLS,ONE_TIME,parameshS,comptopo);
145 //      parafieldS->getField()->setNature(ConservativeVolumic);//IntegralGlobConstraint
146 //      parafieldS->getField()->setArray(dad);
147 //
148 //      // **** TARGET
149 //      parameshT=new ParaMESH(tgt_mesh,*dec.getGroup(),"target mesh");
150 //      parafieldT=new ParaFIELD(ON_CELLS,ONE_TIME,parameshT,comptopo);
151 //      parafieldT->getField()->setNature(ConservativeVolumic);//IntegralGlobConstraint
152 //      parafieldT->getField()->getArray()->fillWithValue(1.0e300);
153 ////      valsT[0]=7.;
154 //    }
155 //  dec.setOrientation(2);
156 //  dec.attachSourceLocalField(parafieldS);
157 //  dec.attachTargetLocalField(parafieldT);
158 //  dec.synchronize();
159 //  dec.sendRecvData(true);
160 //  //
161 //  if(rank==0)
162 //    {
163 //      double integ1 = 0.0, integ2 = 0.0;
164 //      MEDCouplingFieldDouble * tgtField;
165 //
166 //      srcField->integral(true, &integ1);
167 //      tgtField = parafieldT->getField();
168 ////      tgtField->reprQuickOverview(cout);
169 //      tgtField->integral(true, &integ2);
170 //      tgtField->setName("result");
171 //      string out_nam(rep + string("adrien_para.med"));
172 //      MEDLoader::WriteField(out_nam,tgtField, true);
173 //      cout << "wrote: " << out_nam << "\n";
174 //      CPPUNIT_ASSERT_DOUBLES_EQUAL(integ1,integ2,1e-8);
175 //    }
176 //  delete parafieldS;
177 //  delete parafieldT;
178 //  delete parameshS;
179 //  delete parameshT;
180 //
181 //  MPI_Barrier(MPI_COMM_WORLD);
182 //}
183
184 void prepareData1(int rank, ProcessorGroup * grp, NatureOfField nature,
185                                   MEDCouplingUMesh *& meshS, MEDCouplingUMesh *& meshT,
186                                   ParaMESH*& parameshS, ParaMESH*& parameshT,
187                                   ParaFIELD*& parafieldS, ParaFIELD*& parafieldT)
188 {
189   if(rank==0)
190     {
191       const double coordsS[10]={0.,0.,0.5,0.,1.,0.,0.,0.5,0.5,0.5};
192       const double coordsT[6]={0.,0.,1.,0.,1.,1.};
193       meshS=MEDCouplingUMesh::New();
194       meshS->setMeshDimension(2);
195       DataArrayDouble *myCoords=DataArrayDouble::New();
196       myCoords->alloc(5,2);
197       std::copy(coordsS,coordsS+10,myCoords->getPointer());
198       meshS->setCoords(myCoords);
199       myCoords->decrRef();
200       int connS[7]={0,3,4,1, 1,4,2};
201       meshS->allocateCells(2);
202       meshS->insertNextCell(INTERP_KERNEL::NORM_QUAD4,4,connS);
203       meshS->insertNextCell(INTERP_KERNEL::NORM_TRI3,3,connS+4);
204       meshS->finishInsertingCells();
205       ComponentTopology comptopo;
206       parameshS=new ParaMESH(meshS, *grp,"source mesh");
207       parafieldS=new ParaFIELD(ON_CELLS,NO_TIME,parameshS,comptopo);
208       parafieldS->getField()->setNature(nature);
209       double *valsS=parafieldS->getField()->getArray()->getPointer();
210       valsS[0]=7.; valsS[1]=8.;
211       //
212       meshT=MEDCouplingUMesh::New();
213       meshT->setMeshDimension(2);
214       myCoords=DataArrayDouble::New();
215       myCoords->alloc(3,2);
216       std::copy(coordsT,coordsT+6,myCoords->getPointer());
217       meshT->setCoords(myCoords);
218       myCoords->decrRef();
219       int connT[3]={0,2,1};
220       meshT->allocateCells(1);
221       meshT->insertNextCell(INTERP_KERNEL::NORM_TRI3,3,connT);
222       meshT->finishInsertingCells();
223       parameshT=new ParaMESH(meshT,*grp,"target mesh");
224       parafieldT=new ParaFIELD(ON_CELLS,NO_TIME,parameshT,comptopo);
225       parafieldT->getField()->setNature(nature);
226       double *valsT=parafieldT->getField()->getArray()->getPointer();
227       valsT[0]=7.;
228     }
229   //
230   if(rank==1)
231     {
232       const double coordsS[10]={1.,0.,0.5,0.5,1.,0.5,0.5,1.,1.,1.};
233       const double coordsT[6]={0.,0.,0.5,0.5,0.,1.};
234       meshS=MEDCouplingUMesh::New();
235       meshS->setMeshDimension(2);
236       DataArrayDouble *myCoords=DataArrayDouble::New();
237       myCoords->alloc(5,2);
238       std::copy(coordsS,coordsS+10,myCoords->getPointer());
239       meshS->setCoords(myCoords);
240       myCoords->decrRef();
241       int connS[7]={0,1,2, 1,3,4,2};
242       meshS->allocateCells(2);
243       meshS->insertNextCell(INTERP_KERNEL::NORM_TRI3,3,connS);
244       meshS->insertNextCell(INTERP_KERNEL::NORM_QUAD4,4,connS+3);
245       meshS->finishInsertingCells();
246       ComponentTopology comptopo;
247       parameshS=new ParaMESH(meshS,*grp,"source mesh");
248       parafieldS=new ParaFIELD(ON_CELLS,NO_TIME,parameshS,comptopo);
249       parafieldS->getField()->setNature(nature);
250       double *valsS=parafieldS->getField()->getArray()->getPointer();
251       valsS[0]=9.;
252       valsS[1]=11.;
253       //
254       meshT=MEDCouplingUMesh::New();
255       meshT->setMeshDimension(2);
256       myCoords=DataArrayDouble::New();
257       myCoords->alloc(3,2);
258       std::copy(coordsT,coordsT+6,myCoords->getPointer());
259       meshT->setCoords(myCoords);
260       myCoords->decrRef();
261       int connT[3]={0,2,1};
262       meshT->allocateCells(1);
263       meshT->insertNextCell(INTERP_KERNEL::NORM_TRI3,3,connT);
264       meshT->finishInsertingCells();
265       parameshT=new ParaMESH(meshT,*grp,"target mesh");
266       parafieldT=new ParaFIELD(ON_CELLS,NO_TIME,parameshT,comptopo);
267       parafieldT->getField()->setNature(nature);
268       double *valsT=parafieldT->getField()->getArray()->getPointer();
269       valsT[0]=8.;
270     }
271   //
272   if(rank==2)
273     {
274       const double coordsS[8]={0.,0.5, 0.5,0.5, 0.,1., 0.5,1.};
275       const double coordsT[6]={0.5,0.5,0.,1.,1.,1.};
276       meshS=MEDCouplingUMesh::New();
277       meshS->setMeshDimension(2);
278       DataArrayDouble *myCoords=DataArrayDouble::New();
279       myCoords->alloc(4,2);
280       std::copy(coordsS,coordsS+8,myCoords->getPointer());
281       meshS->setCoords(myCoords);
282       myCoords->decrRef();
283       int connS[4]={0,2,3,1};
284       meshS->allocateCells(1);
285       meshS->insertNextCell(INTERP_KERNEL::NORM_QUAD4,4,connS);
286       meshS->finishInsertingCells();
287       ComponentTopology comptopo;
288       parameshS=new ParaMESH(meshS,*grp,"source mesh");
289       parafieldS=new ParaFIELD(ON_CELLS,NO_TIME,parameshS,comptopo);
290       parafieldS->getField()->setNature(nature);
291       double *valsS=parafieldS->getField()->getArray()->getPointer();
292       valsS[0]=10.;
293       //
294       meshT=MEDCouplingUMesh::New();
295       meshT->setMeshDimension(2);
296       myCoords=DataArrayDouble::New();
297       myCoords->alloc(3,2);
298       std::copy(coordsT,coordsT+6,myCoords->getPointer());
299       meshT->setCoords(myCoords);
300       myCoords->decrRef();
301       int connT[3]={0,1,2};
302       meshT->allocateCells(1);
303       meshT->insertNextCell(INTERP_KERNEL::NORM_TRI3,3,connT);
304       meshT->finishInsertingCells();
305       parameshT=new ParaMESH(meshT,*grp,"target mesh");
306       parafieldT=new ParaFIELD(ON_CELLS,NO_TIME,parameshT,comptopo);
307       parafieldT->getField()->setNature(nature);
308       double *valsT=parafieldT->getField()->getArray()->getPointer();
309       valsT[0]=9.;
310     }
311 }
312
313 /*! Test case from the official doc of the OverlapDEC.
314  *  WARNING: bounding boxes are tweaked here to make the case more interesting (i.e. to avoid an all to all exchange
315  *  between all procs).
316  */
317 void ParaMEDMEMTest::testOverlapDEC1()
318 {
319   int size, rank;
320   MPI_Comm_size(MPI_COMM_WORLD,&size);
321   MPI_Comm_rank(MPI_COMM_WORLD,&rank);
322   //  char hostname[256];
323   //  printf("(%d) PID %d on localhost ready for attach\n", rank, getpid());
324   //  fflush(stdout);
325
326 //    if (rank == 1)
327 //      {
328 //        int i=1, j=0;
329 //        while (i!=0)
330 //          j=2;
331 //      }
332
333   if (size != 3) return ;
334   int nproc = 3;
335   std::set<int> procs;
336   for (int i=0; i<nproc; i++)
337     procs.insert(i);
338
339   CommInterface interface;
340   OverlapDEC dec(procs);
341   ProcessorGroup * grp = dec.getGroup();
342   MEDCouplingUMesh* meshS=0, *meshT=0;
343   ParaMESH* parameshS=0, *parameshT=0;
344   ParaFIELD* parafieldS=0, *parafieldT=0;
345
346   MPI_Barrier(MPI_COMM_WORLD);
347   prepareData1(rank, grp, ConservativeVolumic, meshS, meshT, parameshS, parameshT, parafieldS, parafieldT);
348
349   /*!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
350    *  HACK ON BOUNDING BOX TO MAKE THIS CASE SIMPLE AND USABLE IN DEBUG
351    * Bounding boxes are slightly smaller than should be, thus localizing the work to be done
352    * and avoiding every proc talking to everyone else.
353    * Obviously this is NOT a good idea to do this in production code :-)
354    * !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
355    */
356   dec.setBoundingBoxAdjustmentAbs(-1.0e-12);
357
358   dec.attachSourceLocalField(parafieldS);
359   dec.attachTargetLocalField(parafieldT);
360   dec.synchronize();
361   dec.sendRecvData(true);
362   //
363   MPI_Barrier(MPI_COMM_WORLD);
364   if(rank==0)
365     {
366       CPPUNIT_ASSERT_DOUBLES_EQUAL(8.75,parafieldT->getField()->getArray()->getIJ(0,0),1e-12);
367     }
368   if(rank==1)
369     {
370       CPPUNIT_ASSERT_DOUBLES_EQUAL(8.5,parafieldT->getField()->getArray()->getIJ(0,0),1e-12);
371     }
372   if(rank==2)
373     {
374       CPPUNIT_ASSERT_DOUBLES_EQUAL(10.5,parafieldT->getField()->getArray()->getIJ(0,0),1e-12);
375     }
376
377   delete parafieldS;
378   delete parafieldT;
379   delete parameshS;
380   delete parameshT;
381   meshS->decrRef();
382   meshT->decrRef();
383
384   MPI_Barrier(MPI_COMM_WORLD);
385 }
386
387 /*!
388  * Same as testOverlapDEC1() but with regular bounding boxes. If you're looking for a nice debug case,
389  * testOverlapDEC1() is more appropriate.
390  */
391 void ParaMEDMEMTest::testOverlapDEC2()
392 {
393   int size, rank;
394   MPI_Comm_size(MPI_COMM_WORLD,&size);
395   MPI_Comm_rank(MPI_COMM_WORLD,&rank);
396
397   if (size != 3) return ;
398   int nproc = 3;
399   std::set<int> procs;
400   for (int i=0; i<nproc; i++)
401     procs.insert(i);
402
403 //      if (rank == 0)
404 //        {
405 //          int i=1, j=0;
406 //          while (i!=0)
407 //            j=2;
408 //        }
409
410
411   CommInterface interface;
412   OverlapDEC dec(procs);
413   ProcessorGroup * grp = dec.getGroup();
414   MEDCouplingUMesh* meshS=0, *meshT=0;
415   ParaMESH* parameshS=0, *parameshT=0;
416   ParaFIELD* parafieldS=0, *parafieldT=0;
417
418   MPI_Barrier(MPI_COMM_WORLD);
419   prepareData1(rank, grp, ConservativeVolumic,meshS, meshT, parameshS, parameshT, parafieldS, parafieldT);
420
421   /* Normal bounding boxes here: */
422   dec.setBoundingBoxAdjustmentAbs(+1.0e-12);
423
424   dec.attachSourceLocalField(parafieldS);
425   dec.attachTargetLocalField(parafieldT);
426   dec.synchronize();
427   dec.sendRecvData(true);
428   //
429   if(rank==0)
430     {
431       CPPUNIT_ASSERT_DOUBLES_EQUAL(8.75,parafieldT->getField()->getArray()->getIJ(0,0),1e-12);
432     }
433   if(rank==1)
434     {
435       CPPUNIT_ASSERT_DOUBLES_EQUAL(8.5,parafieldT->getField()->getArray()->getIJ(0,0),1e-12);
436     }
437   if(rank==2)
438     {
439       CPPUNIT_ASSERT_DOUBLES_EQUAL(10.5,parafieldT->getField()->getArray()->getIJ(0,0),1e-12);
440     }
441   delete parafieldS;
442   delete parafieldT;
443   delete parameshS;
444   delete parameshT;
445   meshS->decrRef();
446   meshT->decrRef();
447
448   MPI_Barrier(MPI_COMM_WORLD);
449 }
450
451 void prepareData2_buildOneSquare(MEDCouplingUMesh* & meshS_0, MEDCouplingUMesh* & meshT_0)
452 {
453   const double coords[10] = {0.0,0.0,  0.0,1.0,  1.0,1.0,  1.0,0.0, 0.5,0.5};
454   meshS_0 = MEDCouplingUMesh::New("source", 2);
455   DataArrayDouble *myCoords=DataArrayDouble::New();
456   myCoords->alloc(5,2);
457   std::copy(coords,coords+10,myCoords->getPointer());
458   meshS_0->setCoords(myCoords);  myCoords->decrRef();
459   int connS[4]={0,1,2,3};
460   meshS_0->allocateCells(2);
461   meshS_0->insertNextCell(INTERP_KERNEL::NORM_QUAD4,4,connS);
462   //
463   meshT_0 = MEDCouplingUMesh::New("target", 2);
464   myCoords=DataArrayDouble::New();
465   myCoords->alloc(5,2);
466   std::copy(coords,coords+10,myCoords->getPointer());
467   meshT_0->setCoords(myCoords);
468   myCoords->decrRef();
469   int connT[12]={0,1,4,  1,2,4,  2,3,4,  3,0,4};
470   meshT_0->allocateCells(4);
471   meshT_0->insertNextCell(INTERP_KERNEL::NORM_TRI3,3,connT);
472   meshT_0->insertNextCell(INTERP_KERNEL::NORM_TRI3,3,connT+3);
473   meshT_0->insertNextCell(INTERP_KERNEL::NORM_TRI3,3,connT+6);
474   meshT_0->insertNextCell(INTERP_KERNEL::NORM_TRI3,3,connT+9);
475 }
476
477 /**
478  * Prepare five (detached) QUAD4 disposed like this:
479  *   (0)  (1)  (2)
480  *   (3)  (4)
481  *
482  * On the target side the global mesh is identical except that each QUAD4 is split in 4 TRI3 (along the diagonals).
483  * This is a case for two procs:
484  *    - proc #0 has source squares 0,1,2 and target squares 0,3 (well, sets of TRI3s actually)
485  *    - proc #1 has source squares 3,4 and target squares 1,2,4
486  */
487 void prepareData2(int rank, ProcessorGroup * grp, NatureOfField nature,
488                   MEDCouplingUMesh *& meshS, MEDCouplingUMesh *& meshT,
489                   ParaMESH*& parameshS, ParaMESH*& parameshT,
490                   ParaFIELD*& parafieldS, ParaFIELD*& parafieldT)
491 {
492   MEDCouplingUMesh *meshS_0 = 0, *meshT_0 = 0;
493   prepareData2_buildOneSquare(meshS_0, meshT_0);
494
495   if(rank==0)
496     {
497       const double tr1[] = {1.5, 0.0};
498       MEDCouplingUMesh *meshS_1 = static_cast<MEDCouplingUMesh*>(meshS_0->deepCpy());
499       meshS_1->translate(tr1);
500       const double tr2[] = {3.0, 0.0};
501       MEDCouplingUMesh *meshS_2 = static_cast<MEDCouplingUMesh*>(meshS_0->deepCpy());
502       meshS_2->translate(tr2);
503
504       std::vector<const MEDCouplingUMesh*> vec;
505       vec.push_back(meshS_0);vec.push_back(meshS_1);vec.push_back(meshS_2);
506       meshS = MEDCouplingUMesh::MergeUMeshes(vec);
507       meshS_1->decrRef(); meshS_2->decrRef();
508
509       ComponentTopology comptopo;
510       parameshS=new ParaMESH(meshS, *grp,"source mesh");
511       parafieldS=new ParaFIELD(ON_CELLS,ONE_TIME,parameshS,comptopo);
512       parafieldS->getField()->setNature(nature);
513       double *valsS=parafieldS->getField()->getArray()->getPointer();
514       valsS[0]=1.; valsS[1]=2.; valsS[2]=3.;
515
516       //
517       const double tr3[] = {0.0, -1.5};
518       MEDCouplingUMesh *meshT_3 = static_cast<MEDCouplingUMesh*>(meshT_0->deepCpy());
519       meshT_3->translate(tr3);
520       vec.clear();
521       vec.push_back(meshT_0);vec.push_back(meshT_3);
522       meshT = MEDCouplingUMesh::MergeUMeshes(vec);
523       meshT_3->decrRef();
524
525       parameshT=new ParaMESH(meshT,*grp,"target mesh");
526       parafieldT=new ParaFIELD(ON_CELLS,ONE_TIME,parameshT,comptopo);
527       parafieldT->getField()->setNature(nature);
528     }
529   //
530   if(rank==1)
531     {
532       const double tr3[] = {0.0, -1.5};
533       MEDCouplingUMesh *meshS_3 = static_cast<MEDCouplingUMesh*>(meshS_0->deepCpy());
534       meshS_3->translate(tr3);
535       const double tr4[] = {1.5, -1.5};
536       MEDCouplingUMesh *meshS_4 = static_cast<MEDCouplingUMesh*>(meshS_0->deepCpy());
537       meshS_4->translate(tr4);
538
539       std::vector<const MEDCouplingUMesh*> vec;
540       vec.push_back(meshS_3);vec.push_back(meshS_4);
541       meshS = MEDCouplingUMesh::MergeUMeshes(vec);
542       meshS_3->decrRef(); meshS_4->decrRef();
543
544       ComponentTopology comptopo;
545       parameshS=new ParaMESH(meshS, *grp,"source mesh");
546       parafieldS=new ParaFIELD(ON_CELLS,ONE_TIME,parameshS,comptopo);
547       parafieldS->getField()->setNature(nature);
548       double *valsS=parafieldS->getField()->getArray()->getPointer();
549       valsS[0]=4.; valsS[1]=5.;
550
551       //
552       const double tr5[] = {1.5, 0.0};
553       MEDCouplingUMesh *meshT_1 = static_cast<MEDCouplingUMesh*>(meshT_0->deepCpy());
554       meshT_1->translate(tr5);
555       const double tr6[] = {3.0, 0.0};
556       MEDCouplingUMesh *meshT_2 = static_cast<MEDCouplingUMesh*>(meshT_0->deepCpy());
557       meshT_2->translate(tr6);
558       const double tr7[] = {1.5, -1.5};
559       MEDCouplingUMesh *meshT_4 = static_cast<MEDCouplingUMesh*>(meshT_0->deepCpy());
560       meshT_4->translate(tr7);
561
562       vec.clear();
563       vec.push_back(meshT_1);vec.push_back(meshT_2);vec.push_back(meshT_4);
564       meshT = MEDCouplingUMesh::MergeUMeshes(vec);
565       meshT_1->decrRef(); meshT_2->decrRef(); meshT_4->decrRef();
566
567       parameshT=new ParaMESH(meshT,*grp,"target mesh");
568       parafieldT=new ParaFIELD(ON_CELLS,ONE_TIME,parameshT,comptopo);
569       parafieldT->getField()->setNature(nature);
570     }
571   meshS_0->decrRef();
572   meshT_0->decrRef();
573 }
574
575 /*! Test focused on the mapping of cell IDs.
576  * (i.e. when only part of the source/target mesh is transmitted)
577  */
578 void ParaMEDMEMTest::testOverlapDEC3()
579 {
580   int size, rank;
581   MPI_Comm_size(MPI_COMM_WORLD,&size);
582   MPI_Comm_rank(MPI_COMM_WORLD,&rank);
583
584   int nproc = 2;
585   if (size != nproc) return ;
586   std::set<int> procs;
587   for (int i=0; i<nproc; i++)
588     procs.insert(i);
589
590   CommInterface interface;
591   OverlapDEC dec(procs);
592   ProcessorGroup * grp = dec.getGroup();
593   MEDCouplingUMesh* meshS=0, *meshT=0;
594   ParaMESH* parameshS=0, *parameshT=0;
595   ParaFIELD* parafieldS=0, *parafieldT=0;
596
597   prepareData2(rank, grp, ConservativeVolumic, meshS, meshT, parameshS, parameshT, parafieldS, parafieldT);
598 //  if (rank == 1)
599 //    {
600 //      int i=1, j=0;
601 //      while (i!=0)
602 //        j=2;
603 //    }
604
605   dec.attachSourceLocalField(parafieldS);
606   dec.attachTargetLocalField(parafieldT);
607   dec.synchronize();
608   dec.sendRecvData(true);
609   //
610   MEDCouplingFieldDouble * resField = parafieldT->getField();
611   if(rank==0)
612     {
613       CPPUNIT_ASSERT_EQUAL(8, resField->getNumberOfTuples());
614       for(int i=0;i<4;i++)
615         CPPUNIT_ASSERT_DOUBLES_EQUAL(1.0,resField->getArray()->getIJ(i,0),1e-12);
616       for(int i=4;i<8;i++)
617         CPPUNIT_ASSERT_DOUBLES_EQUAL(4.0,resField->getArray()->getIJ(i,0),1e-12);
618     }
619   if(rank==1)
620     {
621       CPPUNIT_ASSERT_EQUAL(12, resField->getNumberOfTuples());
622       for(int i=0;i<4;i++)
623         CPPUNIT_ASSERT_DOUBLES_EQUAL(2.0,resField->getArray()->getIJ(i,0),1e-12);
624       for(int i=4;i<8;i++)
625         CPPUNIT_ASSERT_DOUBLES_EQUAL(3.0,resField->getArray()->getIJ(i,0),1e-12);
626       for(int i=8;i<12;i++)
627         CPPUNIT_ASSERT_DOUBLES_EQUAL(5.0,resField->getArray()->getIJ(i,0),1e-12);
628     }
629   delete parafieldS;
630   delete parafieldT;
631   delete parameshS;
632   delete parameshT;
633   meshS->decrRef();
634   meshT->decrRef();
635
636   MPI_Barrier(MPI_COMM_WORLD);
637 }
638
639 /*!
640  * Test default values and multi-component fields
641  */
642 void ParaMEDMEMTest::testOverlapDEC4()
643 {
644 }