Salome HOME
bc904c487909e0b140b0aaad476b945785bbd1fb
[tools/medcoupling.git] / src / ParaMEDMEMTest / ParaMEDMEMTest_InterpKernelDEC.cxx
1 // Copyright (C) 2007-2020  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 "DEC.hxx"
28 #include "MxN_Mapping.hxx"
29 #include "InterpKernelDEC.hxx"
30 #include "ParaMESH.hxx"
31 #include "ParaFIELD.hxx"
32 #include "ComponentTopology.hxx"
33 #include "ICoCoMEDField.hxx"
34 #include "ParaMEDLoader.hxx"
35 #include "MEDLoader.hxx"
36 #include "TestInterpKernelUtils.hxx"
37
38  
39 #include <string>
40 #include <iterator>
41
42 // use this define to enable lines, execution of which leads to Segmentation Fault
43 #define ENABLE_FAULTS
44
45 // use this define to enable CPPUNIT asserts and fails, showing bugs
46 #define ENABLE_FORCED_FAILURES
47
48
49 using namespace std;
50 using namespace MEDCoupling;
51
52 void ParaMEDMEMTest::testInterpKernelDEC_2D()
53 {
54   testInterpKernelDEC_2D_("P0","P0");
55 }
56
57 void ParaMEDMEMTest::testInterpKernelDEC2_2D()
58 {
59   testInterpKernelDEC2_2D_("P0","P0");
60 }
61
62 void ParaMEDMEMTest::testInterpKernelDEC_3D()
63 {
64   testInterpKernelDEC_3D_("P0","P0");
65 }
66
67 void ParaMEDMEMTest::testInterpKernelDEC_2DP0P1()
68 {
69   //testInterpKernelDEC_2D_("P0","P1");
70 }
71
72 void ParaMEDMEMTest::testInterpKernelDEC_1D()
73 {
74   int size;
75   int rank;
76   MPI_Comm_size(MPI_COMM_WORLD,&size);
77   MPI_Comm_rank(MPI_COMM_WORLD,&rank);
78   //
79   if(size!=5)
80     return ;
81   int nproc_source = 3;
82   set<int> self_procs;
83   set<int> procs_source;
84   set<int> procs_target;
85   
86   for (int i=0; i<nproc_source; i++)
87     procs_source.insert(i);
88   for (int i=nproc_source; i<size; i++)
89     procs_target.insert(i);
90   self_procs.insert(rank);
91   //
92   MEDCoupling::MEDCouplingUMesh *mesh=0;
93   MEDCoupling::ParaMESH *paramesh=0;
94   MEDCoupling::ParaFIELD *parafieldP0=0;
95   //
96   MEDCoupling::CommInterface interface;
97   //
98   ProcessorGroup* self_group = new MEDCoupling::MPIProcessorGroup(interface,self_procs);
99   ProcessorGroup* target_group = new MEDCoupling::MPIProcessorGroup(interface,procs_target);
100   ProcessorGroup* source_group = new MEDCoupling::MPIProcessorGroup(interface,procs_source);
101   //
102   MPI_Barrier(MPI_COMM_WORLD);
103   if(source_group->containsMyRank())
104     {
105       if(rank==0)
106         {
107           double coords[4]={0.3,0.7, 0.9,1.0};
108           mcIdType conn[4]={0,1,2,3};
109           mesh=MEDCouplingUMesh::New("Source mesh Proc0",1);
110           mesh->allocateCells(2);
111           mesh->insertNextCell(INTERP_KERNEL::NORM_SEG2,2,conn);
112           mesh->insertNextCell(INTERP_KERNEL::NORM_SEG2,2,conn+2);
113           mesh->finishInsertingCells();
114           DataArrayDouble *myCoords=DataArrayDouble::New();
115           myCoords->alloc(4,1);
116           std::copy(coords,coords+4,myCoords->getPointer());
117           mesh->setCoords(myCoords);
118           myCoords->decrRef();
119         }
120       if(rank==1)
121         {
122           double coords[2]={0.7,0.9};
123           mcIdType conn[2]={0,1};
124           mesh=MEDCouplingUMesh::New("Source mesh Proc1",1);
125           mesh->allocateCells(1);
126           mesh->insertNextCell(INTERP_KERNEL::NORM_SEG2,2,conn);
127           mesh->finishInsertingCells();
128           DataArrayDouble *myCoords=DataArrayDouble::New();
129           myCoords->alloc(2,1);
130           std::copy(coords,coords+2,myCoords->getPointer());
131           mesh->setCoords(myCoords);
132           myCoords->decrRef();
133         }
134       if(rank==2)
135         {
136           double coords[2]={1.,1.12};
137           mcIdType conn[2]={0,1};
138           mesh=MEDCouplingUMesh::New("Source mesh Proc2",1);
139           mesh->allocateCells(1);
140           mesh->insertNextCell(INTERP_KERNEL::NORM_SEG2,2,conn);
141           mesh->finishInsertingCells();
142           DataArrayDouble *myCoords=DataArrayDouble::New();
143           myCoords->alloc(2,1);
144           std::copy(coords,coords+2,myCoords->getPointer());
145           mesh->setCoords(myCoords);
146           myCoords->decrRef();
147         }
148       paramesh=new ParaMESH(mesh,*source_group,"source mesh");
149       MEDCoupling::ComponentTopology comptopo;
150       parafieldP0 = new ParaFIELD(ON_CELLS,NO_TIME,paramesh, comptopo);
151       double *valueP0=parafieldP0->getField()->getArray()->getPointer();
152       parafieldP0->getField()->setNature(IntensiveMaximum);
153       if(rank==0)
154         {
155           valueP0[0]=7.; valueP0[1]=8.;
156         }
157       if(rank==1)
158         {
159           valueP0[0]=9.;
160         }
161       if(rank==2)
162         {
163           valueP0[0]=10.;
164         }
165     }
166   else
167     {
168       const char targetMeshName[]="target mesh";
169       if(rank==3)
170         {
171           double coords[2]={0.5,0.75};
172           mcIdType conn[2]={0,1};
173           mesh=MEDCouplingUMesh::New("Target mesh Proc3",1);
174           mesh->allocateCells(1);
175           mesh->insertNextCell(INTERP_KERNEL::NORM_SEG2,2,conn);
176           mesh->finishInsertingCells();
177           DataArrayDouble *myCoords=DataArrayDouble::New();
178           myCoords->alloc(2,1);
179           std::copy(coords,coords+2,myCoords->getPointer());
180           mesh->setCoords(myCoords);
181           myCoords->decrRef();
182           paramesh=new ParaMESH(mesh,*target_group,targetMeshName);
183         }
184       if(rank==4)
185         {
186           double coords[2]={0.75,1.2};
187           mcIdType conn[2]={0,1};
188           mesh=MEDCouplingUMesh::New("Target mesh Proc4",1);
189           mesh->allocateCells(1);
190           mesh->insertNextCell(INTERP_KERNEL::NORM_SEG2,2,conn);
191           mesh->finishInsertingCells();
192           DataArrayDouble *myCoords=DataArrayDouble::New();
193           myCoords->alloc(2,1);
194           std::copy(coords,coords+2,myCoords->getPointer());
195           mesh->setCoords(myCoords);
196           myCoords->decrRef();
197           paramesh=new ParaMESH(mesh,*target_group,targetMeshName);
198         }
199       MEDCoupling::ComponentTopology comptopo;
200       parafieldP0 = new ParaFIELD(ON_CELLS,NO_TIME,paramesh, comptopo);
201       parafieldP0->getField()->setNature(IntensiveMaximum);
202     }
203   // test 1
204   MEDCoupling::InterpKernelDEC dec(*source_group,*target_group);
205   if (source_group->containsMyRank())
206     { 
207       dec.setMethod("P0");
208       dec.attachLocalField(parafieldP0);
209       dec.synchronize();
210       dec.setForcedRenormalization(false);
211       dec.sendData();
212       dec.recvData();
213       const double *valueP0=parafieldP0->getField()->getArray()->getPointer();
214       if(rank==0)
215         {
216           CPPUNIT_ASSERT_DOUBLES_EQUAL(7.4,valueP0[0],1e-7);
217           CPPUNIT_ASSERT_DOUBLES_EQUAL(9.0540540540540544,valueP0[1],1e-7);
218         }
219       if(rank==1)
220         {
221           CPPUNIT_ASSERT_DOUBLES_EQUAL(8.64054054054054,valueP0[0],1e-7);
222         }
223       if(rank==2)
224         {
225           CPPUNIT_ASSERT_DOUBLES_EQUAL(9.0540540540540544,valueP0[0],1e-7);
226         }
227     }
228   else
229     {
230       dec.setMethod("P0");
231       dec.attachLocalField(parafieldP0);
232       dec.synchronize();
233       dec.setForcedRenormalization(false);
234       dec.recvData();
235       const double *res=parafieldP0->getField()->getArray()->getConstPointer();
236       if(rank==3)
237         {
238           CPPUNIT_ASSERT_EQUAL(1,(int)parafieldP0->getField()->getNumberOfTuples());
239           CPPUNIT_ASSERT_EQUAL(1,(int)parafieldP0->getField()->getNumberOfComponents());
240           CPPUNIT_ASSERT_DOUBLES_EQUAL(7.4,res[0],1e-12);
241         }
242       if(rank==4)
243         {
244           CPPUNIT_ASSERT_EQUAL(1,(int)parafieldP0->getField()->getNumberOfTuples());
245           CPPUNIT_ASSERT_EQUAL(1,(int)parafieldP0->getField()->getNumberOfComponents());
246           CPPUNIT_ASSERT_DOUBLES_EQUAL(9.0540540540540526,res[0],1e-12);
247         }
248       dec.sendData();
249     }
250   //
251   delete parafieldP0;
252   mesh->decrRef();
253   delete paramesh;
254   delete self_group;
255   delete target_group;
256   delete source_group;
257   //
258   MPI_Barrier(MPI_COMM_WORLD);
259 }
260
261 void ParaMEDMEMTest::testInterpKernelDEC_2DCurve()
262 {
263   int size;
264   int rank;
265   MPI_Comm_size(MPI_COMM_WORLD,&size);
266   MPI_Comm_rank(MPI_COMM_WORLD,&rank);
267   //
268   if(size!=5)
269     return ;
270   int nproc_source = 3;
271   set<int> self_procs;
272   set<int> procs_source;
273   set<int> procs_target;
274   
275   for (int i=0; i<nproc_source; i++)
276     procs_source.insert(i);
277   for (int i=nproc_source; i<size; i++)
278     procs_target.insert(i);
279   self_procs.insert(rank);
280   //
281   MEDCoupling::MEDCouplingUMesh *mesh=0;
282   MEDCoupling::ParaMESH *paramesh=0;
283   MEDCoupling::ParaFIELD *parafieldP0=0;
284   //
285   MEDCoupling::CommInterface interface;
286   //
287   ProcessorGroup* self_group = new MEDCoupling::MPIProcessorGroup(interface,self_procs);
288   ProcessorGroup* target_group = new MEDCoupling::MPIProcessorGroup(interface,procs_target);
289   ProcessorGroup* source_group = new MEDCoupling::MPIProcessorGroup(interface,procs_source);
290   //
291   MPI_Barrier(MPI_COMM_WORLD);
292   if(source_group->containsMyRank())
293     {
294       if(rank==0)
295         {
296           double coords[8]={0.3,0.3,0.7,0.7, 0.9,0.9,1.0,1.0};
297           mcIdType conn[4]={0,1,2,3};
298           mesh=MEDCouplingUMesh::New("Source mesh Proc0",1);
299           mesh->allocateCells(2);
300           mesh->insertNextCell(INTERP_KERNEL::NORM_SEG2,2,conn);
301           mesh->insertNextCell(INTERP_KERNEL::NORM_SEG2,2,conn+2);
302           mesh->finishInsertingCells();
303           DataArrayDouble *myCoords=DataArrayDouble::New();
304           myCoords->alloc(4,2);
305           std::copy(coords,coords+8,myCoords->getPointer());
306           mesh->setCoords(myCoords);
307           myCoords->decrRef();
308         }
309       if(rank==1)
310         {
311           double coords[4]={0.7,0.7,0.9,0.9};
312           mcIdType conn[2]={0,1};
313           mesh=MEDCouplingUMesh::New("Source mesh Proc1",1);
314           mesh->allocateCells(1);
315           mesh->insertNextCell(INTERP_KERNEL::NORM_SEG2,2,conn);
316           mesh->finishInsertingCells();
317           DataArrayDouble *myCoords=DataArrayDouble::New();
318           myCoords->alloc(2,2);
319           std::copy(coords,coords+4,myCoords->getPointer());
320           mesh->setCoords(myCoords);
321           myCoords->decrRef();
322         }
323       if(rank==2)
324         {
325           double coords[4]={1.,1.,1.12,1.12};
326           mcIdType conn[2]={0,1};
327           mesh=MEDCouplingUMesh::New("Source mesh Proc2",1);
328           mesh->allocateCells(1);
329           mesh->insertNextCell(INTERP_KERNEL::NORM_SEG2,2,conn);
330           mesh->finishInsertingCells();
331           DataArrayDouble *myCoords=DataArrayDouble::New();
332           myCoords->alloc(2,2);
333           std::copy(coords,coords+4,myCoords->getPointer());
334           mesh->setCoords(myCoords);
335           myCoords->decrRef();
336         }
337       paramesh=new ParaMESH(mesh,*source_group,"source mesh");
338       MEDCoupling::ComponentTopology comptopo;
339       parafieldP0 = new ParaFIELD(ON_CELLS,NO_TIME,paramesh, comptopo);
340       double *valueP0=parafieldP0->getField()->getArray()->getPointer();
341       parafieldP0->getField()->setNature(IntensiveMaximum);
342       if(rank==0)
343         {
344           valueP0[0]=7.; valueP0[1]=8.;
345         }
346       if(rank==1)
347         {
348           valueP0[0]=9.;
349         }
350       if(rank==2)
351         {
352           valueP0[0]=10.;
353         }
354     }
355   else
356     {
357       const char targetMeshName[]="target mesh";
358       if(rank==3)
359         {
360           double coords[4]={0.5,0.5,0.75,0.75};
361           mcIdType conn[2]={0,1};
362           mesh=MEDCouplingUMesh::New("Target mesh Proc3",1);
363           mesh->allocateCells(1);
364           mesh->insertNextCell(INTERP_KERNEL::NORM_SEG2,2,conn);
365           mesh->finishInsertingCells();
366           DataArrayDouble *myCoords=DataArrayDouble::New();
367           myCoords->alloc(2,2);
368           std::copy(coords,coords+4,myCoords->getPointer());
369           mesh->setCoords(myCoords);
370           myCoords->decrRef();
371           paramesh=new ParaMESH(mesh,*target_group,targetMeshName);
372         }
373       if(rank==4)
374         {
375           double coords[4]={0.75,0.75,1.2,1.2};
376           mcIdType conn[2]={0,1};
377           mesh=MEDCouplingUMesh::New("Target mesh Proc4",1);
378           mesh->allocateCells(1);
379           mesh->insertNextCell(INTERP_KERNEL::NORM_SEG2,2,conn);
380           mesh->finishInsertingCells();
381           DataArrayDouble *myCoords=DataArrayDouble::New();
382           myCoords->alloc(2,2);
383           std::copy(coords,coords+4,myCoords->getPointer());
384           mesh->setCoords(myCoords);
385           myCoords->decrRef();
386           paramesh=new ParaMESH(mesh,*target_group,targetMeshName);
387         }
388       MEDCoupling::ComponentTopology comptopo;
389       parafieldP0 = new ParaFIELD(ON_CELLS,NO_TIME,paramesh, comptopo);
390       parafieldP0->getField()->setNature(IntensiveMaximum);
391     }
392   // test 1
393   MEDCoupling::InterpKernelDEC dec(*source_group,*target_group);
394   if (source_group->containsMyRank())
395     { 
396       dec.setMethod("P0");
397       dec.attachLocalField(parafieldP0);
398       dec.synchronize();
399       dec.setForcedRenormalization(false);
400       dec.sendData();
401       dec.recvData();
402       const double *valueP0=parafieldP0->getField()->getArray()->getPointer();
403       if(rank==0)
404         {
405           CPPUNIT_ASSERT_DOUBLES_EQUAL(7.4,valueP0[0],1e-7);
406           CPPUNIT_ASSERT_DOUBLES_EQUAL(9.0540540540540544,valueP0[1],1e-7);
407         }
408       if(rank==1)
409         {
410           CPPUNIT_ASSERT_DOUBLES_EQUAL(8.64054054054054,valueP0[0],1e-7);
411         }
412       if(rank==2)
413         {
414           CPPUNIT_ASSERT_DOUBLES_EQUAL(9.0540540540540544,valueP0[0],1e-7);
415         }
416     }
417   else
418     {
419       dec.setMethod("P0");
420       dec.attachLocalField(parafieldP0);
421       dec.synchronize();
422       dec.setForcedRenormalization(false);
423       dec.recvData();
424       const double *res=parafieldP0->getField()->getArray()->getConstPointer();
425       if(rank==3)
426         {
427           CPPUNIT_ASSERT_EQUAL(1,(int)parafieldP0->getField()->getNumberOfTuples());
428           CPPUNIT_ASSERT_EQUAL(1,(int)parafieldP0->getField()->getNumberOfComponents());
429           CPPUNIT_ASSERT_DOUBLES_EQUAL(7.4,res[0],1e-12);
430         }
431       if(rank==4)
432         {
433           CPPUNIT_ASSERT_EQUAL(1,(int)parafieldP0->getField()->getNumberOfTuples());
434           CPPUNIT_ASSERT_EQUAL(1,(int)parafieldP0->getField()->getNumberOfComponents());
435           CPPUNIT_ASSERT_DOUBLES_EQUAL(9.0540540540540526,res[0],1e-12);
436         }
437       dec.sendData();
438     }
439   //
440   delete parafieldP0;
441   mesh->decrRef();
442   delete paramesh;
443   delete self_group;
444   delete target_group;
445   delete source_group;
446   //
447   MPI_Barrier(MPI_COMM_WORLD);
448 }
449
450
451 /*
452  * Check methods defined in InterpKernelDEC.hxx
453  *
454  InterpKernelDEC();
455  InterpKernelDEC(ProcessorGroup& local_group, ProcessorGroup& distant_group);
456  virtual ~InterpKernelDEC();
457  void synchronize();
458  void recvData();
459  void sendData();
460 */
461  
462 void ParaMEDMEMTest::testInterpKernelDEC_2D_(const char *srcMeth, const char *targetMeth)
463 {
464   std::string srcM(srcMeth);
465   std::string targetM(targetMeth);
466   int size;
467   int rank;
468   MPI_Comm_size(MPI_COMM_WORLD,&size);
469   MPI_Comm_rank(MPI_COMM_WORLD,&rank);
470
471   //the test is meant to run on five processors
472   if (size !=5) return ;
473    
474   int nproc_source = 3;
475   set<int> self_procs;
476   set<int> procs_source;
477   set<int> procs_target;
478   
479   for (int i=0; i<nproc_source; i++)
480     procs_source.insert(i);
481   for (int i=nproc_source; i<size; i++)
482     procs_target.insert(i);
483   self_procs.insert(rank);
484   
485   MEDCoupling::CommInterface interface;
486     
487   MEDCoupling::ProcessorGroup* self_group = new MEDCoupling::MPIProcessorGroup(interface,self_procs);
488   MEDCoupling::ProcessorGroup* target_group = new MEDCoupling::MPIProcessorGroup(interface,procs_target);
489   MEDCoupling::ProcessorGroup* source_group = new MEDCoupling::MPIProcessorGroup(interface,procs_source);
490   
491   //loading the geometry for the source group
492
493   MEDCoupling::InterpKernelDEC dec (*source_group,*target_group);
494
495   MEDCoupling::MEDCouplingUMesh* mesh = nullptr;
496   MEDCoupling::ParaMESH* paramesh = nullptr;
497   MEDCoupling::ParaFIELD* parafield = nullptr;
498   ICoCo::MEDField* icocofield = nullptr;
499   
500   string filename_xml1              = "square1_split";
501   string filename_xml2              = "square2_split";
502   //string filename_seq_wr            = makeTmpFile("");
503   //string filename_seq_med           = makeTmpFile("myWrField_seq_pointe221.med");
504   
505   // To remove tmp files from disk
506   ParaMEDMEMTest_TmpFilesRemover aRemover;
507   
508   MPI_Barrier(MPI_COMM_WORLD);
509   if (source_group->containsMyRank())
510     {
511       string master = filename_xml1;
512       
513       ostringstream strstream;
514       strstream <<master<<rank+1<<".med";
515       string fName = INTERP_TEST::getResourceFile(strstream.str());
516       ostringstream meshname ;
517       meshname<< "Mesh_2_"<< rank+1;
518       
519       mesh=ReadUMeshFromFile(fName.c_str(),meshname.str().c_str(),0);
520       
521     
522       paramesh=new ParaMESH (mesh,*source_group,"source mesh");
523     
524       //      MEDCoupling::ParaSUPPORT* parasupport=new UnstructuredParaSUPPORT( support,*source_group);
525       MEDCoupling::ComponentTopology comptopo;
526       if(srcM=="P0")
527         {
528           parafield = new ParaFIELD(ON_CELLS,NO_TIME,paramesh, comptopo);
529           parafield->getField()->setNature(IntensiveMaximum);
530         }
531       else
532         parafield = new ParaFIELD(ON_NODES,NO_TIME,paramesh, comptopo);
533       int nb_local;
534       if(srcM=="P0")
535         nb_local=mesh->getNumberOfCells();
536       else
537         nb_local=mesh->getNumberOfNodes();
538       //      double * value= new double[nb_local];
539       double *value=parafield->getField()->getArray()->getPointer();
540       for(int ielem=0; ielem<nb_local;ielem++)
541         value[ielem]=1.0;
542     
543       //      ICoCo::Field* icocofield=new ICoCo::MEDField(paramesh,parafield);
544       icocofield=new ICoCo::MEDField(parafield->getField());
545       dec.setMethod(srcMeth);
546       dec.attachLocalField(icocofield);
547     }
548   
549   //loading the geometry for the target group
550   if (target_group->containsMyRank())
551     {
552       string master= filename_xml2;
553       ostringstream strstream;
554       strstream << master<<(rank-nproc_source+1)<<".med";
555       string fName = INTERP_TEST::getResourceFile(strstream.str());
556       ostringstream meshname ;
557       meshname<< "Mesh_3_"<<rank-nproc_source+1;
558       mesh = ReadUMeshFromFile(fName.c_str(),meshname.str().c_str(),0);
559       
560       paramesh=new ParaMESH (mesh,*target_group,"target mesh");
561       //      MEDCoupling::ParaSUPPORT* parasupport=new UnstructuredParaSUPPORT(support,*target_group);
562       MEDCoupling::ComponentTopology comptopo;
563       if(targetM=="P0")
564         {
565           parafield = new ParaFIELD(ON_CELLS,NO_TIME,paramesh, comptopo);
566           parafield->getField()->setNature(IntensiveMaximum);
567         }
568       else
569         parafield = new ParaFIELD(ON_NODES,NO_TIME,paramesh, comptopo);
570       int nb_local;
571       if(targetM=="P0")
572         nb_local=mesh->getNumberOfCells();
573       else
574         nb_local=mesh->getNumberOfNodes();
575       //      double * value= new double[nb_local];
576       double *value=parafield->getField()->getArray()->getPointer();
577       for(int ielem=0; ielem<nb_local;ielem++)
578         value[ielem]=0.0;
579       //      ICoCo::Field* icocofield=new ICoCo::MEDField(paramesh,parafield);
580       icocofield=new ICoCo::MEDField(parafield->getField());
581       dec.setMethod(targetMeth);
582       dec.attachLocalField(icocofield);
583     }
584     
585   
586   //attaching a DEC to the source group 
587   double field_before_int;
588   double field_after_int;
589   
590   if (source_group->containsMyRank())
591     { 
592       field_before_int = parafield->getVolumeIntegral(0,true);
593       dec.synchronize();
594       cout<<"DEC usage"<<endl;
595       dec.setForcedRenormalization(false);
596
597       dec.sendData();
598       ParaMEDLoader::WriteParaMesh("./sourcesquareb",paramesh);
599       if (source_group->myRank()==0)
600         aRemover.Register("./sourcesquareb");
601       ostringstream filename;
602       filename<<"./sourcesquareb_"<<source_group->myRank()+1;
603       aRemover.Register(filename.str().c_str());
604       //WriteField("./sourcesquareb",parafield->getField());
605    
606       dec.recvData();
607       cout <<"writing"<<endl;
608       ParaMEDLoader::WriteParaMesh("./sourcesquare",paramesh);
609       if (source_group->myRank()==0)
610         aRemover.Register("./sourcesquare");
611       //WriteField("./sourcesquare",parafield->getField());
612       
613      
614       filename<<"./sourcesquare_"<<source_group->myRank()+1;
615       aRemover.Register(filename.str().c_str());
616       field_after_int = parafield->getVolumeIntegral(0,true);
617       
618       
619       //      MPI_Bcast(&field_before_int,1,MPI_DOUBLE,0,MPI_COMM_WORLD);
620       //       MPI_Bcast(&field_after_int,1,MPI_DOUBLE,0,MPI_COMM_WORLD);
621
622       CPPUNIT_ASSERT_DOUBLES_EQUAL(field_before_int, field_after_int, 1e-6);
623     
624     }
625   
626   //attaching a DEC to the target group
627   if (target_group->containsMyRank())
628     {
629       dec.synchronize();
630       dec.setForcedRenormalization(false);
631
632       dec.recvData();
633       ParaMEDLoader::WriteParaMesh("./targetsquareb",paramesh);
634       //WriteField("./targetsquareb",parafield->getField());
635       if (target_group->myRank()==0)
636         aRemover.Register("./targetsquareb");
637       ostringstream filename;
638       filename<<"./targetsquareb_"<<target_group->myRank()+1;
639       aRemover.Register(filename.str().c_str());
640       dec.sendData();
641       ParaMEDLoader::WriteParaMesh("./targetsquare",paramesh);
642       //WriteField("./targetsquare",parafield->getField());
643       
644       if (target_group->myRank()==0)
645         aRemover.Register("./targetsquareb");
646       
647       filename<<"./targetsquareb_"<<target_group->myRank()+1;
648       aRemover.Register(filename.str().c_str());
649       //    double field_before_int, field_after_int;
650       //       MPI_Bcast(&field_before_int,1,MPI_DOUBLE,0,MPI_COMM_WORLD);
651       //       MPI_Bcast(&field_after_int,1,MPI_DOUBLE,0,MPI_COMM_WORLD);
652       
653       //      CPPUNIT_ASSERT_DOUBLES_EQUAL(field_before_int, field_after_int, 1e-6);
654     
655     }
656   
657   delete source_group;
658   delete target_group;
659   delete self_group;
660   delete parafield;
661   delete paramesh;
662   mesh->decrRef();
663
664   delete icocofield;
665
666   MPI_Barrier(MPI_COMM_WORLD);
667   cout << "end of InterpKernelDEC_2D test"<<endl;
668 }
669
670 void ParaMEDMEMTest::testInterpKernelDEC2_2D_(const char *srcMeth, const char *targetMeth)
671 {
672   std::string srcM(srcMeth);
673   std::string targetM(targetMeth);
674   int size;
675   int rank;
676   MPI_Comm_size(MPI_COMM_WORLD,&size);
677   MPI_Comm_rank(MPI_COMM_WORLD,&rank);
678
679   //the test is meant to run on five processors
680   if (size !=5) return ;
681    
682   int nproc_source = 3;
683   set<int> self_procs;
684   set<int> procs_source;
685   set<int> procs_target;
686   
687   for (int i=0; i<nproc_source; i++)
688     procs_source.insert(i);
689   for (int i=nproc_source; i<size; i++)
690     procs_target.insert(i);
691   self_procs.insert(rank);
692   
693   MEDCoupling::CommInterface interface;
694     
695   MEDCoupling::ProcessorGroup* self_group = new MEDCoupling::MPIProcessorGroup(interface,self_procs);
696   MEDCoupling::ProcessorGroup* target_group = new MEDCoupling::MPIProcessorGroup(interface,procs_target);
697   MEDCoupling::ProcessorGroup* source_group = new MEDCoupling::MPIProcessorGroup(interface,procs_source);
698   
699   //loading the geometry for the source group
700
701   MEDCoupling::InterpKernelDEC dec (*source_group,*target_group);
702
703   MEDCoupling::MEDCouplingUMesh* mesh = nullptr;
704   MEDCoupling::MEDCouplingFieldDouble* mcfield = nullptr;
705   
706   string filename_xml1              = "square1_split";
707   string filename_xml2              = "square2_split";
708   
709   // To remove tmp files from disk
710   ParaMEDMEMTest_TmpFilesRemover aRemover;
711   
712   MPI_Barrier(MPI_COMM_WORLD);
713   if (source_group->containsMyRank())
714     {
715       string master = filename_xml1;
716       
717       ostringstream strstream;
718       strstream <<master<<rank+1<<".med";
719       string fName = INTERP_TEST::getResourceFile(strstream.str());
720       ostringstream meshname ;
721       meshname<< "Mesh_2_"<< rank+1;
722       
723       mesh=ReadUMeshFromFile(fName.c_str(),meshname.str().c_str(),0);
724       MEDCoupling::ComponentTopology comptopo;
725       if(srcM=="P0")
726         {
727           mcfield = MEDCouplingFieldDouble::New(ON_CELLS,NO_TIME);
728           mcfield->setMesh(mesh);
729           DataArrayDouble *array=DataArrayDouble::New();
730           array->alloc(mcfield->getNumberOfTuples(),1);
731           mcfield->setArray(array);
732           array->decrRef();
733           mcfield->setNature(IntensiveMaximum);
734         }
735       else
736         {
737           mcfield = MEDCouplingFieldDouble::New(ON_CELLS,NO_TIME);
738           mcfield->setMesh(mesh);
739           DataArrayDouble *array=DataArrayDouble::New();
740           array->alloc(mcfield->getNumberOfTuples(),1);
741           mcfield->setArray(array);
742           array->decrRef();
743         }
744       int nb_local;
745       if(srcM=="P0")
746         nb_local=mesh->getNumberOfCells();
747       else
748         nb_local=mesh->getNumberOfNodes();
749       double *value=mcfield->getArray()->getPointer();
750       for(int ielem=0; ielem<nb_local;ielem++)
751         value[ielem]=1.0;
752       dec.setMethod(srcMeth);
753       dec.attachLocalField(mcfield);
754       dec.attachLocalField(mcfield);
755     }
756   
757   //loading the geometry for the target group
758   if (target_group->containsMyRank())
759     {
760       string master= filename_xml2;
761       ostringstream strstream;
762       strstream << master<<(rank-nproc_source+1)<<".med";
763       string fName = INTERP_TEST::getResourceFile(strstream.str());
764       ostringstream meshname ;
765       meshname<< "Mesh_3_"<<rank-nproc_source+1;
766       mesh = ReadUMeshFromFile(fName.c_str(),meshname.str().c_str(),0);
767       MEDCoupling::ComponentTopology comptopo;
768       if(targetM=="P0")
769         {
770           mcfield = MEDCouplingFieldDouble::New(ON_CELLS,NO_TIME);
771           mcfield->setMesh(mesh);
772           DataArrayDouble *array=DataArrayDouble::New();
773           array->alloc(mcfield->getNumberOfTuples(),1);
774           mcfield->setArray(array);
775           array->decrRef();
776           mcfield->setNature(IntensiveMaximum);
777         }
778       else
779         {
780           mcfield = MEDCouplingFieldDouble::New(ON_NODES,NO_TIME);
781           mcfield->setMesh(mesh);
782           DataArrayDouble *array=DataArrayDouble::New();
783           array->alloc(mcfield->getNumberOfTuples(),1);
784           mcfield->setArray(array);
785           array->decrRef();
786         }
787       int nb_local;
788       if(targetM=="P0")
789         nb_local=mesh->getNumberOfCells();
790       else
791         nb_local=mesh->getNumberOfNodes();
792       double *value=mcfield->getArray()->getPointer();
793       for(int ielem=0; ielem<nb_local;ielem++)
794         value[ielem]=0.0;
795       dec.setMethod(targetMeth);
796       dec.attachLocalField(mcfield);
797       dec.attachLocalField(mcfield);
798     }
799     
800   
801   //attaching a DEC to the source group 
802
803   if (source_group->containsMyRank())
804     { 
805       dec.synchronize();
806       dec.setForcedRenormalization(false);
807       dec.sendData();
808       dec.recvData();
809     }
810   
811   //attaching a DEC to the target group
812   if (target_group->containsMyRank())
813     {
814       dec.synchronize();
815       dec.setForcedRenormalization(false);
816       dec.recvData();
817       dec.sendData();
818     }
819   delete source_group;
820   delete target_group;
821   delete self_group;
822   mcfield->decrRef();
823   mesh->decrRef();
824
825   MPI_Barrier(MPI_COMM_WORLD);
826   cout << "end of InterpKernelDEC2_2D test"<<endl;
827 }
828
829 void ParaMEDMEMTest::testInterpKernelDEC_3D_(const char *srcMeth, const char *targetMeth)
830 {
831   std::string srcM(srcMeth);
832   std::string targetM(targetMeth);
833   int size;
834   int rank;
835   MPI_Comm_size(MPI_COMM_WORLD,&size);
836   MPI_Comm_rank(MPI_COMM_WORLD,&rank);
837
838   //the test is meant to run on five processors
839   if (size !=3) return ;
840    
841   int nproc_source = 2;
842   set<int> self_procs;
843   set<int> procs_source;
844   set<int> procs_target;
845   
846   for (int i=0; i<nproc_source; i++)
847     procs_source.insert(i);
848   for (int i=nproc_source; i<size; i++)
849     procs_target.insert(i);
850   self_procs.insert(rank);
851   
852   MEDCoupling::CommInterface interface;
853     
854   MEDCoupling::ProcessorGroup* self_group = new MEDCoupling::MPIProcessorGroup(interface,self_procs);
855   MEDCoupling::ProcessorGroup* target_group = new MEDCoupling::MPIProcessorGroup(interface,procs_target);
856   MEDCoupling::ProcessorGroup* source_group = new MEDCoupling::MPIProcessorGroup(interface,procs_source);
857   
858   //loading the geometry for the source group
859
860   MEDCoupling::InterpKernelDEC dec (*source_group,*target_group);
861
862   MEDCoupling::MEDCouplingUMesh* mesh = nullptr;
863   MEDCoupling::ParaMESH* paramesh = nullptr;
864   MEDCoupling::ParaFIELD* parafield = nullptr;
865   ICoCo::MEDField* icocofield = nullptr;
866   
867   char * tmp_dir_c                    = getenv("TMP");
868   string tmp_dir;
869   if (tmp_dir_c != NULL)
870     tmp_dir = string(tmp_dir_c);
871   else
872     tmp_dir = "/tmp";
873   string filename_xml1              = "Mesh3D_10_2d";
874   string filename_xml2              = "Mesh3D_11";
875   //string filename_seq_wr            = makeTmpFile("");
876   //string filename_seq_med           = makeTmpFile("myWrField_seq_pointe221.med");
877   
878   // To remove tmp files from disk
879   ParaMEDMEMTest_TmpFilesRemover aRemover;
880   
881   MPI_Barrier(MPI_COMM_WORLD);
882   if (source_group->containsMyRank())
883     {
884       string master = filename_xml1;
885       
886       ostringstream strstream;
887       strstream <<master<<rank+1<<".med";
888       std::string fName = INTERP_TEST::getResourceFile(strstream.str());
889       ostringstream meshname ;
890       meshname<< "Mesh_3_"<< rank+1;
891       
892       mesh=ReadUMeshFromFile(fName.c_str(),meshname.str().c_str(),0);
893       
894     
895       paramesh=new ParaMESH (mesh,*source_group,"source mesh");
896     
897       //      MEDCoupling::ParaSUPPORT* parasupport=new UnstructuredParaSUPPORT( support,*source_group);
898       MEDCoupling::ComponentTopology comptopo;
899       if(srcM=="P0")
900         {
901           parafield = new ParaFIELD(ON_CELLS,NO_TIME,paramesh, comptopo);
902           parafield->getField()->setNature(IntensiveMaximum);
903         }
904       else
905         parafield = new ParaFIELD(ON_NODES,NO_TIME,paramesh, comptopo);
906       int nb_local;
907       if(srcM=="P0")
908         nb_local=mesh->getNumberOfCells();
909       else
910         nb_local=mesh->getNumberOfNodes();
911       //      double * value= new double[nb_local];
912       double *value=parafield->getField()->getArray()->getPointer();
913       for(int ielem=0; ielem<nb_local;ielem++)
914         value[ielem]=1.0;
915     
916       //      ICoCo::Field* icocofield=new ICoCo::MEDField(paramesh,parafield);
917       icocofield=new ICoCo::MEDField(parafield->getField());
918       dec.setMethod(srcMeth);
919       dec.attachLocalField(icocofield);
920     }
921   
922   //loading the geometry for the target group
923   if (target_group->containsMyRank())
924     {
925       string master= filename_xml2;
926       ostringstream strstream;
927       strstream << master << ".med";
928       std::string fName = INTERP_TEST::getResourceFile(strstream.str());
929       ostringstream meshname ;
930       meshname<< "Mesh_6";
931       mesh = ReadUMeshFromFile(fName.c_str(),meshname.str().c_str(),0);
932       
933       paramesh=new ParaMESH (mesh,*target_group,"target mesh");
934       //      MEDCoupling::ParaSUPPORT* parasupport=new UnstructuredParaSUPPORT(support,*target_group);
935       MEDCoupling::ComponentTopology comptopo;
936       if(targetM=="P0")
937         {
938           parafield = new ParaFIELD(ON_CELLS,NO_TIME,paramesh, comptopo);
939           parafield->getField()->setNature(IntensiveMaximum);
940         }
941       else
942         parafield = new ParaFIELD(ON_NODES,NO_TIME,paramesh, comptopo);
943       int nb_local;
944       if(targetM=="P0")
945         nb_local=mesh->getNumberOfCells();
946       else
947         nb_local=mesh->getNumberOfNodes();
948       //      double * value= new double[nb_local];
949       double *value=parafield->getField()->getArray()->getPointer();
950       for(int ielem=0; ielem<nb_local;ielem++)
951         value[ielem]=0.0;
952       //      ICoCo::Field* icocofield=new ICoCo::MEDField(paramesh,parafield);
953       icocofield=new ICoCo::MEDField(parafield->getField());
954       dec.setMethod(targetMeth);
955       dec.attachLocalField(icocofield);
956     }  
957   //attaching a DEC to the source group 
958   double field_before_int;
959   double field_after_int;
960   
961   if (source_group->containsMyRank())
962     { 
963       field_before_int = parafield->getVolumeIntegral(0,true);
964       dec.synchronize();
965       cout<<"DEC usage"<<endl;
966       dec.setForcedRenormalization(false);
967
968       dec.sendData();
969       ParaMEDLoader::WriteParaMesh("./sourcesquareb",paramesh);
970       if (source_group->myRank()==0)
971         aRemover.Register("./sourcesquareb");
972       ostringstream filename;
973       filename<<"./sourcesquareb_"<<source_group->myRank()+1;
974       aRemover.Register(filename.str().c_str());
975       //WriteField("./sourcesquareb",parafield->getField());
976    
977       dec.recvData();
978       cout <<"writing"<<endl;
979       ParaMEDLoader::WriteParaMesh("./sourcesquare",paramesh);
980       if (source_group->myRank()==0)
981         aRemover.Register("./sourcesquare");
982       //WriteField("./sourcesquare",parafield->getField());
983       
984      
985       filename<<"./sourcesquare_"<<source_group->myRank()+1;
986       aRemover.Register(filename.str().c_str());
987       field_after_int = parafield->getVolumeIntegral(0,true);
988
989       CPPUNIT_ASSERT_DOUBLES_EQUAL(field_before_int, field_after_int, 1e-6);
990     
991     }
992   
993   //attaching a DEC to the target group
994   if (target_group->containsMyRank())
995     {
996       dec.synchronize();
997       dec.setForcedRenormalization(false);
998
999       dec.recvData();
1000       ParaMEDLoader::WriteParaMesh("./targetsquareb",paramesh);
1001       //WriteField("./targetsquareb",parafield->getField());
1002       if (target_group->myRank()==0)
1003         aRemover.Register("./targetsquareb");
1004       ostringstream filename;
1005       filename<<"./targetsquareb_"<<target_group->myRank()+1;
1006       aRemover.Register(filename.str().c_str());
1007       dec.sendData();
1008       ParaMEDLoader::WriteParaMesh("./targetsquare",paramesh);
1009       //WriteField("./targetsquare",parafield->getField());
1010       
1011       if (target_group->myRank()==0)
1012         aRemover.Register("./targetsquareb");
1013       
1014       filename<<"./targetsquareb_"<<target_group->myRank()+1;
1015       aRemover.Register(filename.str().c_str());
1016     }
1017   delete source_group;
1018   delete target_group;
1019   delete self_group;
1020   delete parafield;
1021   delete paramesh;
1022   mesh->decrRef();
1023
1024   delete icocofield;
1025
1026   MPI_Barrier(MPI_COMM_WORLD);
1027   cout << "end of InterpKernelDEC_3D test"<<endl;
1028 }
1029
1030 //Synchronous tests without interpolation with native mode (AllToAll(v) from lam/MPI:
1031 void ParaMEDMEMTest::testSynchronousEqualInterpKernelWithoutInterpNativeDEC_2D()
1032 {
1033   testAsynchronousInterpKernelDEC_2D(0.1,1,0.1,1,false,false,false,"P0","P0");
1034 }
1035
1036 //Synchronous tests without interpolation :
1037 void ParaMEDMEMTest::testSynchronousEqualInterpKernelWithoutInterpDEC_2D()
1038 {
1039   testAsynchronousInterpKernelDEC_2D(0.1,1,0.1,1,true,false,false,"P0","P0");
1040 }
1041
1042 //Synchronous tests with interpolation :
1043 void ParaMEDMEMTest::testSynchronousEqualInterpKernelDEC_2D()
1044 {
1045   testAsynchronousInterpKernelDEC_2D(0.1,1,0.1,1,true,false,true,"P0","P0");
1046 }
1047 void ParaMEDMEMTest::testSynchronousFasterSourceInterpKernelDEC_2D()
1048 {
1049   testAsynchronousInterpKernelDEC_2D(0.09,1,0.1,1,true,false,true,"P0","P0");
1050 }
1051 void ParaMEDMEMTest::testSynchronousSlowerSourceInterpKernelDEC_2D()
1052 {
1053   testAsynchronousInterpKernelDEC_2D(0.11,1,0.1,1,true,false,true,"P0","P0");
1054 }
1055 void ParaMEDMEMTest::testSynchronousSlowSourceInterpKernelDEC_2D()
1056 {
1057   testAsynchronousInterpKernelDEC_2D(0.11,1,0.01,1,true,false,true,"P0","P0");
1058 }
1059 void ParaMEDMEMTest::testSynchronousFastSourceInterpKernelDEC_2D()
1060 {
1061   testAsynchronousInterpKernelDEC_2D(0.01,1,0.11,1,true,false,true,"P0","P0");
1062 }
1063
1064 //Asynchronous tests with interpolation :
1065 void ParaMEDMEMTest::testAsynchronousEqualInterpKernelDEC_2D()
1066 {
1067   testAsynchronousInterpKernelDEC_2D(0.1,1,0.1,1,true,true,true,"P0","P0");
1068 }
1069 void ParaMEDMEMTest::testAsynchronousFasterSourceInterpKernelDEC_2D()
1070 {
1071   testAsynchronousInterpKernelDEC_2D(0.09,1,0.1,1,true,true,true,"P0","P0");
1072 }
1073 void ParaMEDMEMTest::testAsynchronousSlowerSourceInterpKernelDEC_2D()
1074 {
1075   testAsynchronousInterpKernelDEC_2D(0.11,1,0.1,1,true,true,true,"P0","P0");
1076 }
1077 void ParaMEDMEMTest::testAsynchronousSlowSourceInterpKernelDEC_2D()
1078 {
1079   testAsynchronousInterpKernelDEC_2D(0.11,1,0.01,1,true,true,true,"P0","P0");
1080 }
1081 void ParaMEDMEMTest::testAsynchronousFastSourceInterpKernelDEC_2D()
1082 {
1083   testAsynchronousInterpKernelDEC_2D(0.01,1,0.11,1,true,true,true,"P0","P0");
1084 }
1085
1086 void ParaMEDMEMTest::testInterpKernelDECNonOverlapp_2D_P0P0()
1087 {
1088   //
1089   const double sourceCoordsAll[2][8]={{0.4,0.5,0.4,1.5,1.6,1.5,1.6,0.5},
1090                                       {0.3,-0.5,1.6,-0.5,1.6,-1.5,0.3,-1.5}};
1091   const double targetCoordsAll[3][16]={{0.7,1.45,0.7,1.65,0.9,1.65,0.9,1.45,  1.1,1.4,1.1,1.6,1.3,1.6,1.3,1.4},
1092                                        {0.7,-0.6,0.7,0.7,0.9,0.7,0.9,-0.6,  1.1,-0.7,1.1,0.6,1.3,0.6,1.3,-0.7},
1093                                        {0.7,-1.55,0.7,-1.35,0.9,-1.35,0.9,-1.55,  1.1,-1.65,1.1,-1.45,1.3,-1.45,1.3,-1.65}};
1094   mcIdType conn4All[8]={0,1,2,3,4,5,6,7};
1095   double targetResults[3][2]={{34.,34.},{38.333333333333336,42.666666666666664},{47.,47.}};
1096   double targetResults2[3][2]={{0.28333333333333344,0.56666666666666687},{1.8564102564102569,2.0128205128205132},{1.0846153846153845,0.36153846153846159}};
1097   double targetResults3[3][2]={{3.7777777777777781,7.5555555555555562},{24.511111111111113,26.355555555555558},{14.1,4.7}};
1098   double targetResults4[3][2]={{8.5,17},{8.8461538461538431, 9.8461538461538449},{35.25,11.75}};
1099   //
1100   int size;
1101   int rank;
1102   MPI_Comm_size(MPI_COMM_WORLD,&size);
1103   MPI_Comm_rank(MPI_COMM_WORLD,&rank);
1104   //
1105   if(size!=5)
1106     return ;
1107   int nproc_source = 2;
1108   set<int> self_procs;
1109   set<int> procs_source;
1110   set<int> procs_target;
1111   
1112   for (int i=0; i<nproc_source; i++)
1113     procs_source.insert(i);
1114   for (int i=nproc_source; i<size; i++)
1115     procs_target.insert(i);
1116   self_procs.insert(rank);
1117   //
1118   MEDCoupling::MEDCouplingUMesh *mesh=0;
1119   MEDCoupling::ParaMESH *paramesh=0;
1120   MEDCoupling::ParaFIELD* parafield=0;
1121   //
1122   MEDCoupling::CommInterface interface;
1123   //
1124   ProcessorGroup* self_group = new MEDCoupling::MPIProcessorGroup(interface,self_procs);
1125   ProcessorGroup* target_group = new MEDCoupling::MPIProcessorGroup(interface,procs_target);
1126   ProcessorGroup* source_group = new MEDCoupling::MPIProcessorGroup(interface,procs_source);
1127   //
1128   MPI_Barrier(MPI_COMM_WORLD);
1129   if(source_group->containsMyRank())
1130     {
1131       std::ostringstream stream; stream << "sourcemesh2D proc " << rank;
1132       mesh=MEDCouplingUMesh::New(stream.str().c_str(),2);
1133       mesh->allocateCells(2);
1134       mesh->insertNextCell(INTERP_KERNEL::NORM_QUAD4,4,conn4All);
1135       mesh->finishInsertingCells();
1136       DataArrayDouble *myCoords=DataArrayDouble::New();
1137       myCoords->alloc(4,2);
1138       const double *sourceCoords=sourceCoordsAll[rank];
1139       std::copy(sourceCoords,sourceCoords+8,myCoords->getPointer());
1140       mesh->setCoords(myCoords);
1141       myCoords->decrRef();
1142       paramesh=new ParaMESH(mesh,*source_group,"source mesh");
1143       MEDCoupling::ComponentTopology comptopo;
1144       parafield = new ParaFIELD(ON_CELLS,NO_TIME,paramesh, comptopo);
1145       double *value=parafield->getField()->getArray()->getPointer();
1146       value[0]=34+13*((double)rank);
1147     }
1148   else
1149     {
1150       std::ostringstream stream; stream << "targetmesh2D proc " << rank-nproc_source;
1151       mesh=MEDCouplingUMesh::New(stream.str().c_str(),2);
1152       mesh->allocateCells(2);
1153       mesh->insertNextCell(INTERP_KERNEL::NORM_QUAD4,4,conn4All);
1154       mesh->insertNextCell(INTERP_KERNEL::NORM_QUAD4,4,conn4All+4);
1155       mesh->finishInsertingCells();
1156       DataArrayDouble *myCoords=DataArrayDouble::New();
1157       myCoords->alloc(8,2);
1158       const double *targetCoords=targetCoordsAll[rank-nproc_source];
1159       std::copy(targetCoords,targetCoords+16,myCoords->getPointer());
1160       mesh->setCoords(myCoords);
1161       myCoords->decrRef();
1162       paramesh=new ParaMESH (mesh,*target_group,"target mesh");
1163       MEDCoupling::ComponentTopology comptopo;
1164       parafield = new ParaFIELD(ON_CELLS,NO_TIME,paramesh, comptopo);
1165     }
1166   //test 1 - Conservative volumic
1167   MEDCoupling::InterpKernelDEC dec(*source_group,*target_group);
1168   parafield->getField()->setNature(IntensiveMaximum);
1169   if (source_group->containsMyRank())
1170     { 
1171       dec.setMethod("P0");
1172       dec.attachLocalField(parafield);
1173       dec.synchronize();
1174       dec.setForcedRenormalization(false);
1175       dec.sendData();
1176     }
1177   else
1178     {
1179       dec.setMethod("P0");
1180       dec.attachLocalField(parafield);
1181       dec.synchronize();
1182       dec.setForcedRenormalization(false);
1183       dec.recvData();
1184       const double *res=parafield->getField()->getArray()->getConstPointer();
1185       const double *expected=targetResults[rank-nproc_source];
1186       CPPUNIT_ASSERT_DOUBLES_EQUAL(expected[0],res[0],1e-13);
1187       CPPUNIT_ASSERT_DOUBLES_EQUAL(expected[1],res[1],1e-13);
1188     }
1189   //test 2 - ExtensiveMaximum
1190   MEDCoupling::InterpKernelDEC dec2(*source_group,*target_group);
1191   parafield->getField()->setNature(ExtensiveMaximum);
1192   if (source_group->containsMyRank())
1193     { 
1194       dec2.setMethod("P0");
1195       dec2.attachLocalField(parafield);
1196       dec2.synchronize();
1197       dec2.setForcedRenormalization(false);
1198       dec2.sendData();
1199     }
1200   else
1201     {
1202       dec2.setMethod("P0");
1203       dec2.attachLocalField(parafield);
1204       dec2.synchronize();
1205       dec2.setForcedRenormalization(false);
1206       dec2.recvData();
1207       const double *res=parafield->getField()->getArray()->getConstPointer();
1208       const double *expected=targetResults2[rank-nproc_source];
1209       CPPUNIT_ASSERT_DOUBLES_EQUAL(expected[0],res[0],1e-13);
1210       CPPUNIT_ASSERT_DOUBLES_EQUAL(expected[1],res[1],1e-13);
1211     }
1212   //test 3 - ExtensiveMaximum with global constraint
1213   MEDCoupling::InterpKernelDEC dec3(*source_group,*target_group);
1214   parafield->getField()->setNature(ExtensiveConservation);
1215   if (source_group->containsMyRank())
1216     { 
1217       dec3.setMethod("P0");
1218       dec3.attachLocalField(parafield);
1219       dec3.synchronize();
1220       dec3.setForcedRenormalization(false);
1221       dec3.sendData();
1222     }
1223   else
1224     {
1225       dec3.setMethod("P0");
1226       dec3.attachLocalField(parafield);
1227       dec3.synchronize();
1228       dec3.setForcedRenormalization(false);
1229       dec3.recvData();
1230       const double *res=parafield->getField()->getArray()->getConstPointer();
1231       const double *expected=targetResults3[rank-nproc_source];
1232       CPPUNIT_ASSERT_DOUBLES_EQUAL(expected[0],res[0],1e-13);
1233       CPPUNIT_ASSERT_DOUBLES_EQUAL(expected[1],res[1],1e-13);
1234     }
1235   //test 4 - IntensiveConservation
1236   MEDCoupling::InterpKernelDEC dec4(*source_group,*target_group);
1237   parafield->getField()->setNature(IntensiveConservation);
1238   if (source_group->containsMyRank())
1239     { 
1240       dec4.setMethod("P0");
1241       dec4.attachLocalField(parafield);
1242       dec4.synchronize();
1243       dec4.setForcedRenormalization(false);
1244       dec4.sendData();
1245     }
1246   else
1247     {
1248       dec4.setMethod("P0");
1249       dec4.attachLocalField(parafield);
1250       dec4.synchronize();
1251       dec4.setForcedRenormalization(false);
1252       dec4.recvData();
1253       const double *res=parafield->getField()->getArray()->getConstPointer();
1254       const double *expected=targetResults4[rank-nproc_source];
1255       CPPUNIT_ASSERT_DOUBLES_EQUAL(expected[0],res[0],1e-13);
1256       CPPUNIT_ASSERT_DOUBLES_EQUAL(expected[1],res[1],1e-13);
1257     }
1258   //test 5 - Conservative volumic reversed
1259   MEDCoupling::InterpKernelDEC dec5(*source_group,*target_group);
1260   parafield->getField()->setNature(IntensiveMaximum);
1261   if (source_group->containsMyRank())
1262     { 
1263       dec5.setMethod("P0");
1264       dec5.attachLocalField(parafield);
1265       dec5.synchronize();
1266       dec5.setForcedRenormalization(false);
1267       dec5.recvData();
1268       const double *res=parafield->getField()->getArray()->getConstPointer();
1269       CPPUNIT_ASSERT_EQUAL(1,(int)parafield->getField()->getNumberOfTuples());
1270       const double expected[]={37.8518518518519,43.5333333333333};
1271       CPPUNIT_ASSERT_DOUBLES_EQUAL(expected[rank],res[0],1e-13);
1272     }
1273   else
1274     {
1275       dec5.setMethod("P0");
1276       dec5.attachLocalField(parafield);
1277       dec5.synchronize();
1278       dec5.setForcedRenormalization(false);
1279       double *res=parafield->getField()->getArray()->getPointer();
1280       const double *toSet=targetResults[rank-nproc_source];
1281       res[0]=toSet[0];
1282       res[1]=toSet[1];
1283       dec5.sendData();
1284     }
1285   //test 6 - ExtensiveMaximum reversed
1286   MEDCoupling::InterpKernelDEC dec6(*source_group,*target_group);
1287   parafield->getField()->setNature(ExtensiveMaximum);
1288   if (source_group->containsMyRank())
1289     { 
1290       dec6.setMethod("P0");
1291       dec6.attachLocalField(parafield);
1292       dec6.synchronize();
1293       dec6.setForcedRenormalization(false);
1294       dec6.recvData();
1295       const double *res=parafield->getField()->getArray()->getConstPointer();
1296       CPPUNIT_ASSERT_EQUAL(1,(int)parafield->getField()->getNumberOfTuples());
1297       const double expected[]={0.794600591715977,1.35631163708087};
1298       CPPUNIT_ASSERT_DOUBLES_EQUAL(expected[rank],res[0],1e-13);
1299     }
1300   else
1301     {
1302       dec6.setMethod("P0");
1303       dec6.attachLocalField(parafield);
1304       dec6.synchronize();
1305       dec6.setForcedRenormalization(false);
1306       double *res=parafield->getField()->getArray()->getPointer();
1307       const double *toSet=targetResults2[rank-nproc_source];
1308       res[0]=toSet[0];
1309       res[1]=toSet[1];
1310       dec6.sendData();
1311     }
1312   //test 7 - ExtensiveMaximum with global constraint reversed
1313   MEDCoupling::InterpKernelDEC dec7(*source_group,*target_group);
1314   parafield->getField()->setNature(ExtensiveConservation);
1315   if (source_group->containsMyRank())
1316     { 
1317       dec7.setMethod("P0");
1318       dec7.attachLocalField(parafield);
1319       dec7.synchronize();
1320       dec7.setForcedRenormalization(false);
1321       dec7.recvData();
1322       const double *res=parafield->getField()->getArray()->getConstPointer();
1323       CPPUNIT_ASSERT_EQUAL(1,(int)parafield->getField()->getNumberOfTuples());
1324       const double expected[]={36.4592592592593,44.5407407407407};
1325       CPPUNIT_ASSERT_DOUBLES_EQUAL(expected[rank],res[0],1e-13);
1326     }
1327   else
1328     {
1329       dec7.setMethod("P0");
1330       dec7.attachLocalField(parafield);
1331       dec7.synchronize();
1332       dec7.setForcedRenormalization(false);
1333       double *res=parafield->getField()->getArray()->getPointer();
1334       const double *toSet=targetResults3[rank-nproc_source];
1335       res[0]=toSet[0];
1336       res[1]=toSet[1];
1337       dec7.sendData();
1338     }
1339   //test 8 - ExtensiveMaximum with IntensiveConservation reversed
1340   MEDCoupling::InterpKernelDEC dec8(*source_group,*target_group);
1341   parafield->getField()->setNature(IntensiveConservation);
1342   if (source_group->containsMyRank())
1343     { 
1344       dec8.setMethod("P0");
1345       dec8.attachLocalField(parafield);
1346       dec8.synchronize();
1347       dec8.setForcedRenormalization(false);
1348       dec8.recvData();
1349       const double *res=parafield->getField()->getArray()->getConstPointer();
1350       CPPUNIT_ASSERT_EQUAL(1,(int)parafield->getField()->getNumberOfTuples());
1351       const double expected[]={0.81314102564102553,1.3428994082840233};
1352       CPPUNIT_ASSERT_DOUBLES_EQUAL(expected[rank],res[0],1e-13);
1353     }
1354   else
1355     {
1356       dec8.setMethod("P0");
1357       dec8.attachLocalField(parafield);
1358       dec8.synchronize();
1359       dec8.setForcedRenormalization(false);
1360       double *res=parafield->getField()->getArray()->getPointer();
1361       const double *toSet=targetResults4[rank-nproc_source];
1362       res[0]=toSet[0];
1363       res[1]=toSet[1];
1364       dec8.sendData();
1365     }
1366   //
1367   delete parafield;
1368   mesh->decrRef();
1369   delete paramesh;
1370   delete self_group;
1371   delete target_group;
1372   delete source_group;
1373   //
1374   MPI_Barrier(MPI_COMM_WORLD);
1375 }
1376
1377 void ParaMEDMEMTest::testInterpKernelDECNonOverlapp_2D_P0P1P1P0()
1378 {
1379   int size;
1380   int rank;
1381   MPI_Comm_size(MPI_COMM_WORLD,&size);
1382   MPI_Comm_rank(MPI_COMM_WORLD,&rank);
1383   //
1384   if(size!=5)
1385     return ;
1386   int nproc_source = 2;
1387   set<int> self_procs;
1388   set<int> procs_source;
1389   set<int> procs_target;
1390   
1391   for (int i=0; i<nproc_source; i++)
1392     procs_source.insert(i);
1393   for (int i=nproc_source; i<size; i++)
1394     procs_target.insert(i);
1395   self_procs.insert(rank);
1396   //
1397   MEDCoupling::MEDCouplingUMesh *mesh=0;
1398   MEDCoupling::ParaMESH *paramesh=0;
1399   MEDCoupling::ParaFIELD *parafieldP0=0,*parafieldP1=0;
1400   //
1401   MEDCoupling::CommInterface interface;
1402   //
1403   ProcessorGroup* self_group = new MEDCoupling::MPIProcessorGroup(interface,self_procs);
1404   ProcessorGroup* target_group = new MEDCoupling::MPIProcessorGroup(interface,procs_target);
1405   ProcessorGroup* source_group = new MEDCoupling::MPIProcessorGroup(interface,procs_source);
1406   //
1407   MPI_Barrier(MPI_COMM_WORLD);
1408   if(source_group->containsMyRank())
1409     {
1410       if(rank==0)
1411         {
1412           double coords[6]={-0.3,-0.3, 0.7,0.7, 0.7,-0.3};
1413           mcIdType conn[3]={0,1,2};
1414           //int globalNode[3]={1,2,0};
1415           mesh=MEDCouplingUMesh::New("Source mesh Proc0",2);
1416           mesh->allocateCells(1);
1417           mesh->insertNextCell(INTERP_KERNEL::NORM_TRI3,3,conn);
1418           mesh->finishInsertingCells();
1419           DataArrayDouble *myCoords=DataArrayDouble::New();
1420           myCoords->alloc(3,2);
1421           std::copy(coords,coords+6,myCoords->getPointer());
1422           mesh->setCoords(myCoords);
1423           myCoords->decrRef();
1424         }
1425       if(rank==1)
1426         {
1427           double coords[6]={-0.3,-0.3, -0.3,0.7, 0.7,0.7};
1428           mcIdType conn[3]={0,1,2};
1429           //int globalNode[3]={1,3,2};
1430           mesh=MEDCouplingUMesh::New("Source mesh Proc1",2);
1431           mesh->allocateCells(1);
1432           mesh->insertNextCell(INTERP_KERNEL::NORM_TRI3,3,conn);
1433           mesh->finishInsertingCells();
1434           DataArrayDouble *myCoords=DataArrayDouble::New();
1435           myCoords->alloc(3,2);
1436           std::copy(coords,coords+6,myCoords->getPointer());
1437           mesh->setCoords(myCoords);
1438           myCoords->decrRef();
1439         }
1440       paramesh=new ParaMESH(mesh,*source_group,"source mesh");
1441       MEDCoupling::ComponentTopology comptopo;
1442       parafieldP0 = new ParaFIELD(ON_CELLS,NO_TIME,paramesh, comptopo);
1443       parafieldP1 = new ParaFIELD(ON_NODES,NO_TIME,paramesh, comptopo);
1444       double *valueP0=parafieldP0->getField()->getArray()->getPointer();
1445       double *valueP1=parafieldP1->getField()->getArray()->getPointer();
1446       parafieldP0->getField()->setNature(IntensiveMaximum);
1447       parafieldP1->getField()->setNature(IntensiveMaximum);
1448       if(rank==0)
1449         {
1450           valueP0[0]=31.;
1451           valueP1[0]=34.; valueP1[1]=77.; valueP1[2]=53.;
1452         }
1453       if(rank==1)
1454         {
1455           valueP0[0]=47.;
1456           valueP1[0]=34.; valueP1[1]=57.; valueP1[2]=77.;
1457         }
1458     }
1459   else
1460     {
1461       const char targetMeshName[]="target mesh";
1462       if(rank==2)
1463         {
1464           double coords[10]={-0.3,-0.3, 0.2,-0.3, 0.7,-0.3, -0.3,0.2, 0.2,0.2 };
1465           mcIdType conn[7]={0,3,4,1, 1,4,2};
1466           //int globalNode[5]={4,3,0,2,1};
1467           mesh=MEDCouplingUMesh::New("Target mesh Proc2",2);
1468           mesh->allocateCells(2);
1469           mesh->insertNextCell(INTERP_KERNEL::NORM_QUAD4,4,conn);
1470           mesh->insertNextCell(INTERP_KERNEL::NORM_TRI3,3,conn+4);
1471           mesh->finishInsertingCells();
1472           DataArrayDouble *myCoords=DataArrayDouble::New();
1473           myCoords->alloc(5,2);
1474           std::copy(coords,coords+10,myCoords->getPointer());
1475           mesh->setCoords(myCoords);
1476           myCoords->decrRef();
1477           paramesh=new ParaMESH(mesh,*target_group,targetMeshName);
1478           std::vector<mcIdType> globalNumberingP2 = {0,1,2,3,4};
1479           MCAuto<DataArrayIdType> da=DataArrayIdType::New(); da->alloc(5,1);
1480           std::copy(globalNumberingP2.begin(), globalNumberingP2.end(), da->rwBegin());
1481           paramesh->setNodeGlobal(da);
1482         }
1483       if(rank==3)
1484         {
1485           double coords[6]={0.2,0.2, 0.7,-0.3, 0.7,0.2};
1486           mcIdType conn[3]={0,2,1};
1487           //int globalNode[3]={1,0,5};
1488           mesh=MEDCouplingUMesh::New("Target mesh Proc3",2);
1489           mesh->allocateCells(1);
1490           mesh->insertNextCell(INTERP_KERNEL::NORM_TRI3,3,conn);
1491           mesh->finishInsertingCells();
1492           DataArrayDouble *myCoords=DataArrayDouble::New();
1493           myCoords->alloc(3,2);
1494           std::copy(coords,coords+6,myCoords->getPointer());
1495           mesh->setCoords(myCoords);
1496           myCoords->decrRef();
1497           paramesh=new ParaMESH(mesh,*target_group,targetMeshName);
1498           std::vector<mcIdType> globalNumberingP3 = {4,2,5};
1499           MCAuto<DataArrayIdType> da=DataArrayIdType::New(); da->alloc(3,1);
1500           std::copy(globalNumberingP3.begin(), globalNumberingP3.end(), da->rwBegin());
1501           paramesh->setNodeGlobal(da);
1502         }
1503       if(rank==4)
1504         {
1505           double coords[12]={-0.3,0.2, -0.3,0.7, 0.2,0.7, 0.2,0.2, 0.7,0.7, 0.7,0.2};
1506           mcIdType conn[8]={0,1,2,3, 3,2,4,5};
1507           //int globalNode[6]={2,6,7,1,8,5};
1508           mesh=MEDCouplingUMesh::New("Target mesh Proc4",2);
1509           mesh->allocateCells(2);
1510           mesh->insertNextCell(INTERP_KERNEL::NORM_QUAD4,4,conn);
1511           mesh->insertNextCell(INTERP_KERNEL::NORM_QUAD4,4,conn+4);
1512           mesh->finishInsertingCells();
1513           DataArrayDouble *myCoords=DataArrayDouble::New();
1514           myCoords->alloc(6,2);
1515           std::copy(coords,coords+12,myCoords->getPointer());
1516           mesh->setCoords(myCoords);
1517           myCoords->decrRef();
1518           paramesh=new ParaMESH(mesh,*target_group,targetMeshName);
1519           std::vector<mcIdType> globalNumberingP4 = {3,6,7,4,8,5};
1520           MCAuto<DataArrayIdType> da=DataArrayIdType::New(); da->alloc(6,1);
1521           std::copy(globalNumberingP4.begin(), globalNumberingP4.end(), da->rwBegin());
1522           paramesh->setNodeGlobal(da);
1523         }
1524       MEDCoupling::ComponentTopology comptopo;
1525       parafieldP0 = new ParaFIELD(ON_CELLS,NO_TIME,paramesh, comptopo);
1526       parafieldP1 = new ParaFIELD(ON_NODES,NO_TIME,paramesh, comptopo);
1527       parafieldP0->getField()->setNature(IntensiveMaximum);
1528       parafieldP1->getField()->setNature(IntensiveMaximum);
1529     }
1530   // test 1 - P0 P1
1531   MEDCoupling::InterpKernelDEC dec(*source_group,*target_group);
1532   if (source_group->containsMyRank())
1533     { 
1534       dec.setMethod("P0");
1535       dec.attachLocalField(parafieldP0);
1536       dec.synchronize();
1537       dec.setForcedRenormalization(false);
1538       dec.sendData();
1539       dec.recvData();
1540       const double *valueP0=parafieldP0->getField()->getArray()->getPointer();
1541       if(rank==0)
1542         {
1543           CPPUNIT_ASSERT_DOUBLES_EQUAL(34.42857143,valueP0[0],1e-7);
1544         }
1545       if(rank==1)
1546         {
1547           CPPUNIT_ASSERT_DOUBLES_EQUAL(44.,valueP0[0],1e-7);
1548         }
1549     }
1550   else
1551     {
1552       dec.setMethod("P1");
1553       dec.attachLocalField(parafieldP1);
1554       dec.synchronize();
1555       dec.setForcedRenormalization(false);
1556       dec.recvData();
1557       const double *res=parafieldP1->getField()->getArray()->getConstPointer();
1558       if(rank==2)
1559         {
1560           const double expectP2[5]={39.0, 31.0, 31.0, 47.0, 39.0};
1561           CPPUNIT_ASSERT_EQUAL(5,(int)parafieldP1->getField()->getNumberOfTuples());
1562           CPPUNIT_ASSERT_EQUAL(1,(int)parafieldP1->getField()->getNumberOfComponents());
1563           for(int kk=0;kk<5;kk++)
1564             CPPUNIT_ASSERT_DOUBLES_EQUAL(expectP2[kk],res[kk],1e-12);
1565         }
1566       if(rank==3)
1567         {
1568           const double expectP3[3]={39.0, 31.0, 31.0};
1569           CPPUNIT_ASSERT_EQUAL(3,(int)parafieldP1->getField()->getNumberOfTuples());
1570           CPPUNIT_ASSERT_EQUAL(1,(int)parafieldP1->getField()->getNumberOfComponents());
1571           for(int kk=0;kk<3;kk++)
1572             CPPUNIT_ASSERT_DOUBLES_EQUAL(expectP3[kk],res[kk],1e-12);
1573         }
1574       if(rank==4)
1575         {
1576           const double expectP4[6]={47.0, 47.0, 47.0, 39.0, 39.0, 31.0};
1577           CPPUNIT_ASSERT_EQUAL(6,(int)parafieldP1->getField()->getNumberOfTuples());
1578           CPPUNIT_ASSERT_EQUAL(1,(int)parafieldP1->getField()->getNumberOfComponents());
1579           for(int kk=0;kk<6;kk++)
1580             CPPUNIT_ASSERT_DOUBLES_EQUAL(expectP4[kk],res[kk],1e-12);
1581         }
1582       dec.sendData();
1583     }
1584   //
1585   delete parafieldP0;
1586   delete parafieldP1;
1587   mesh->decrRef();
1588   delete paramesh;
1589   delete self_group;
1590   delete target_group;
1591   delete source_group;
1592   //
1593   MPI_Barrier(MPI_COMM_WORLD);
1594 }
1595
1596 void ParaMEDMEMTest::testInterpKernelDEC2DM1D_P0P0()
1597 {
1598   int size;
1599   int rank;
1600   MPI_Comm_size(MPI_COMM_WORLD,&size);
1601   MPI_Comm_rank(MPI_COMM_WORLD,&rank);
1602   //
1603   if(size!=3)
1604     return ;
1605   int nproc_source=2;
1606   set<int> procs_source;
1607   set<int> procs_target;
1608   //
1609   for (int i=0; i<nproc_source; i++)
1610     procs_source.insert(i);
1611   for (int i=nproc_source;i<size; i++)
1612     procs_target.insert(i);
1613   //
1614   MEDCoupling::MEDCouplingUMesh *mesh=0;
1615   MEDCoupling::ParaMESH *paramesh=0;
1616   MEDCoupling::ParaFIELD *parafield=0;
1617   //
1618   MEDCoupling::CommInterface interface;
1619   //
1620   ProcessorGroup* target_group = new MEDCoupling::MPIProcessorGroup(interface,procs_target);
1621   ProcessorGroup* source_group = new MEDCoupling::MPIProcessorGroup(interface,procs_source);
1622   //
1623   MPI_Barrier(MPI_COMM_WORLD);
1624   if(source_group->containsMyRank())
1625     {
1626       double targetCoords[18]={-0.3,-0.3, 0.2,-0.3, 0.7,-0.3, -0.3,0.2, 0.2,0.2, 0.7,0.2, -0.3,0.7, 0.2,0.7, 0.7,0.7 };
1627       mesh=MEDCouplingUMesh::New();
1628       mesh->setMeshDimension(2);
1629       DataArrayDouble *myCoords=DataArrayDouble::New();
1630       myCoords->alloc(9,2);
1631       std::copy(targetCoords,targetCoords+18,myCoords->getPointer());
1632       mesh->setCoords(myCoords);
1633       myCoords->decrRef();
1634       if(rank==0)
1635         {
1636           mcIdType targetConn[7]={0,3,4,1, 1,4,2};
1637           mesh->allocateCells(2);
1638           mesh->insertNextCell(INTERP_KERNEL::NORM_QUAD4,4,targetConn);
1639           mesh->insertNextCell(INTERP_KERNEL::NORM_TRI3,3,targetConn+4);
1640           mesh->finishInsertingCells();
1641         }
1642       else
1643         { 
1644           mcIdType targetConn[11]={4,5,2, 6,7,4,3, 7,8,5,4};
1645           mesh->allocateCells(3);
1646           mesh->insertNextCell(INTERP_KERNEL::NORM_TRI3,3,targetConn);
1647           mesh->insertNextCell(INTERP_KERNEL::NORM_QUAD4,4,targetConn+3);
1648           mesh->insertNextCell(INTERP_KERNEL::NORM_QUAD4,4,targetConn+7);
1649           mesh->finishInsertingCells();
1650         }
1651       MEDCoupling::ComponentTopology comptopo;
1652       paramesh=new ParaMESH(mesh,*source_group,"source mesh");
1653       parafield=new ParaFIELD(ON_CELLS,NO_TIME,paramesh, comptopo);
1654       parafield->getField()->setNature(IntensiveMaximum);
1655       double *vals=parafield->getField()->getArray()->getPointer();
1656       if(rank==0)
1657         { vals[0]=7.; vals[1]=8.; }
1658       else
1659         { vals[0]=9.; vals[1]=10.; vals[2]=11.; }
1660     }
1661   else
1662     {
1663       mesh=MEDCouplingUMesh::New("an example of -1 D mesh",-1);
1664       MEDCoupling::ComponentTopology comptopo;
1665       paramesh=new ParaMESH(mesh,*target_group,"target mesh");
1666       parafield=new ParaFIELD(ON_CELLS,NO_TIME,paramesh, comptopo);
1667       parafield->getField()->setNature(IntensiveMaximum);
1668     }
1669   MEDCoupling::InterpKernelDEC dec(*source_group,*target_group);
1670   if(source_group->containsMyRank())
1671     {
1672       dec.setMethod("P0");
1673       dec.attachLocalField(parafield);
1674       dec.synchronize();
1675       dec.setForcedRenormalization(false);
1676       dec.sendData();
1677       dec.recvData();
1678       const double *res=parafield->getField()->getArray()->getConstPointer();
1679       if(rank==0)
1680         {
1681           CPPUNIT_ASSERT_DOUBLES_EQUAL(9.125,res[0],1e-12);
1682           CPPUNIT_ASSERT_DOUBLES_EQUAL(9.125,res[1],1e-12);
1683         }
1684       else
1685         {
1686           CPPUNIT_ASSERT_DOUBLES_EQUAL(9.125,res[0],1e-12);
1687           CPPUNIT_ASSERT_DOUBLES_EQUAL(9.125,res[1],1e-12);
1688           CPPUNIT_ASSERT_DOUBLES_EQUAL(9.125,res[2],1e-12);
1689         }
1690     }
1691   else
1692     {
1693       dec.setMethod("P0");
1694       dec.attachLocalField(parafield);
1695       dec.synchronize();
1696       dec.setForcedRenormalization(false);
1697       dec.recvData();
1698       const double *res=parafield->getField()->getArray()->getConstPointer();
1699       CPPUNIT_ASSERT_DOUBLES_EQUAL(9.125,res[0],1e-12);
1700       dec.sendData();
1701     }
1702   MEDCoupling::InterpKernelDEC dec2(*source_group,*target_group);
1703   dec2.setMethod("P0");
1704   parafield->getField()->setNature(ExtensiveConservation);
1705   if(source_group->containsMyRank())
1706     {
1707       double *vals=parafield->getField()->getArray()->getPointer();
1708       if(rank==0)
1709         { vals[0]=7.; vals[1]=8.; }
1710       else
1711         { vals[0]=9.; vals[1]=10.; vals[2]=11.; }
1712       dec2.attachLocalField(parafield);
1713       dec2.synchronize();
1714       dec2.sendData();
1715       dec2.recvData();
1716       const double *res=parafield->getField()->getArray()->getConstPointer();
1717       if(rank==0)
1718         {
1719           CPPUNIT_ASSERT_DOUBLES_EQUAL(11.25,res[0],1e-12);
1720           CPPUNIT_ASSERT_DOUBLES_EQUAL(5.625,res[1],1e-12);
1721         }
1722       else
1723         {
1724           CPPUNIT_ASSERT_DOUBLES_EQUAL(5.625,res[0],1e-12);
1725           CPPUNIT_ASSERT_DOUBLES_EQUAL(11.25,res[1],1e-12);
1726           CPPUNIT_ASSERT_DOUBLES_EQUAL(11.25,res[2],1e-12);
1727         }
1728     }
1729   else
1730     {
1731       dec2.attachLocalField(parafield);
1732       dec2.synchronize();
1733       dec2.recvData();
1734       const double *res=parafield->getField()->getArray()->getConstPointer();
1735       CPPUNIT_ASSERT_DOUBLES_EQUAL(45.,res[0],1e-12);
1736       dec2.sendData();
1737     }
1738   //
1739   MEDCoupling::InterpKernelDEC dec3(*source_group,*target_group);
1740   dec3.setMethod("P0");
1741   parafield->getField()->setNature(ExtensiveMaximum);
1742   if(source_group->containsMyRank())
1743     {
1744       double *vals=parafield->getField()->getArray()->getPointer();
1745       if(rank==0)
1746         { vals[0]=7.; vals[1]=8.; }
1747       else
1748         { vals[0]=9.; vals[1]=10.; vals[2]=11.; }
1749       dec3.attachLocalField(parafield);
1750       dec3.synchronize();
1751       dec3.sendData();
1752       dec3.recvData();
1753       const double *res=parafield->getField()->getArray()->getConstPointer();
1754       if(rank==0)
1755         {
1756           CPPUNIT_ASSERT_DOUBLES_EQUAL(11.25,res[0],1e-12);
1757           CPPUNIT_ASSERT_DOUBLES_EQUAL(5.625,res[1],1e-12);
1758         }
1759       else
1760         {
1761           CPPUNIT_ASSERT_DOUBLES_EQUAL(5.625,res[0],1e-12);
1762           CPPUNIT_ASSERT_DOUBLES_EQUAL(11.25,res[1],1e-12);
1763           CPPUNIT_ASSERT_DOUBLES_EQUAL(11.25,res[2],1e-12);
1764         }
1765     }
1766   else
1767     {
1768       dec3.attachLocalField(parafield);
1769       dec3.synchronize();
1770       dec3.recvData();
1771       const double *res=parafield->getField()->getArray()->getConstPointer();
1772       CPPUNIT_ASSERT_DOUBLES_EQUAL(45.,res[0],1e-12);
1773       dec3.sendData();
1774     }
1775   //
1776   MEDCoupling::InterpKernelDEC dec4(*source_group,*target_group);
1777   dec4.setMethod("P0");
1778   parafield->getField()->setNature(IntensiveConservation);
1779   if(source_group->containsMyRank())
1780     {
1781       double *vals=parafield->getField()->getArray()->getPointer();
1782       if(rank==0)
1783         { vals[0]=7.; vals[1]=8.; }
1784       else
1785         { vals[0]=9.; vals[1]=10.; vals[2]=11.; }
1786       dec4.attachLocalField(parafield);
1787       dec4.synchronize();
1788       dec4.sendData();
1789       dec4.recvData();
1790       const double *res=parafield->getField()->getArray()->getConstPointer();
1791        if(rank==0)
1792         {
1793           CPPUNIT_ASSERT_DOUBLES_EQUAL(9.125,res[0],1e-12);
1794           CPPUNIT_ASSERT_DOUBLES_EQUAL(9.125,res[1],1e-12);
1795         }
1796       else
1797         {
1798           CPPUNIT_ASSERT_DOUBLES_EQUAL(9.125,res[0],1e-12);
1799           CPPUNIT_ASSERT_DOUBLES_EQUAL(9.125,res[1],1e-12);
1800           CPPUNIT_ASSERT_DOUBLES_EQUAL(9.125,res[2],1e-12);
1801         }
1802     }
1803   else
1804     {
1805       dec4.attachLocalField(parafield);
1806       dec4.synchronize();
1807       dec4.recvData();
1808       const double *res=parafield->getField()->getArray()->getConstPointer();
1809       CPPUNIT_ASSERT_DOUBLES_EQUAL(9.125,res[0],1e-12);
1810       dec4.sendData();
1811     }
1812   delete parafield;
1813   delete paramesh;
1814   mesh->decrRef();
1815   delete target_group;
1816   delete source_group;
1817   //
1818   MPI_Barrier(MPI_COMM_WORLD);
1819 }
1820
1821 void ParaMEDMEMTest::testInterpKernelDECPartialProcs()
1822 {
1823   int size;
1824   int rank;
1825   MPI_Comm_size(MPI_COMM_WORLD,&size);
1826   MPI_Comm_rank(MPI_COMM_WORLD,&rank);
1827   //
1828   if(size!=3)
1829     return ;
1830   set<int> procs_source;
1831   set<int> procs_target;
1832   //
1833   procs_source.insert(0);
1834   procs_target.insert(1);
1835   //
1836   MEDCoupling::MEDCouplingUMesh *mesh=0;
1837   MEDCoupling::ParaMESH *paramesh=0;
1838   MEDCoupling::ParaFIELD *parafield=0;
1839   //
1840   MEDCoupling::CommInterface interface;
1841   //
1842   MPI_Barrier(MPI_COMM_WORLD);
1843   double targetCoords[8]={ 0.,0., 1., 0., 0., 1., 1., 1. };
1844   CommInterface comm;
1845   int grpIds[2]={0,1};
1846   MPI_Group grp,group_world;
1847   comm.commGroup(MPI_COMM_WORLD,&group_world);
1848   comm.groupIncl(group_world,2,grpIds,&grp);
1849   MPI_Comm partialComm;
1850   comm.commCreate(MPI_COMM_WORLD,grp,&partialComm);
1851   //
1852   ProcessorGroup* target_group=0;
1853   ProcessorGroup* source_group=0;
1854   //
1855   MEDCoupling::InterpKernelDEC *dec=0;
1856   if(rank==0 || rank==1)
1857     {
1858       target_group = new MEDCoupling::MPIProcessorGroup(interface,procs_target,partialComm);
1859       source_group = new MEDCoupling::MPIProcessorGroup(interface,procs_source,partialComm);
1860       if(source_group->containsMyRank())
1861         {    
1862           mesh=MEDCouplingUMesh::New();
1863           mesh->setMeshDimension(2);
1864           DataArrayDouble *myCoords=DataArrayDouble::New();
1865           myCoords->alloc(4,2);
1866           std::copy(targetCoords,targetCoords+8,myCoords->getPointer());
1867           mesh->setCoords(myCoords);
1868           myCoords->decrRef();
1869           mcIdType targetConn[4]={0,2,3,1};
1870           mesh->allocateCells(1);
1871           mesh->insertNextCell(INTERP_KERNEL::NORM_QUAD4,4,targetConn);
1872           mesh->finishInsertingCells();
1873           MEDCoupling::ComponentTopology comptopo;
1874           paramesh=new ParaMESH(mesh,*source_group,"source mesh");
1875           parafield=new ParaFIELD(ON_CELLS,NO_TIME,paramesh, comptopo);
1876           parafield->getField()->setNature(IntensiveMaximum);
1877           double *vals=parafield->getField()->getArray()->getPointer();
1878           vals[0]=7.;
1879           dec=new MEDCoupling::InterpKernelDEC(*source_group,*target_group);
1880           dec->attachLocalField(parafield);
1881           dec->synchronize();
1882           dec->sendData();
1883           dec->recvData();
1884         }
1885       else
1886         {
1887           mesh=MEDCouplingUMesh::New();
1888           mesh->setMeshDimension(2);
1889           DataArrayDouble *myCoords=DataArrayDouble::New();
1890           myCoords->alloc(4,2);
1891           std::copy(targetCoords,targetCoords+8,myCoords->getPointer());
1892           mesh->setCoords(myCoords);
1893           myCoords->decrRef();
1894           mcIdType targetConn[6]={0,2,1,2,3,1};
1895           mesh->allocateCells(2);
1896           mesh->insertNextCell(INTERP_KERNEL::NORM_TRI3,3,targetConn);
1897           mesh->insertNextCell(INTERP_KERNEL::NORM_TRI3,3,targetConn+3);
1898           mesh->finishInsertingCells();
1899           MEDCoupling::ComponentTopology comptopo;
1900           paramesh=new ParaMESH(mesh,*target_group,"target mesh");
1901           parafield=new ParaFIELD(ON_CELLS,NO_TIME,paramesh, comptopo);
1902           parafield->getField()->setNature(IntensiveMaximum);
1903           dec=new MEDCoupling::InterpKernelDEC(*source_group,*target_group);
1904           dec->attachLocalField(parafield);
1905           dec->synchronize();
1906           dec->recvData();
1907           dec->sendData();
1908         }
1909     }
1910   delete parafield;
1911   delete paramesh;
1912   if(mesh)
1913     mesh->decrRef();
1914   delete target_group;
1915   delete source_group;
1916   delete dec;
1917   if(partialComm != MPI_COMM_NULL)
1918     comm.commFree(&partialComm);
1919   comm.groupFree(&grp);
1920   comm.groupFree(&group_world);
1921   MPI_Barrier(MPI_COMM_WORLD);
1922 }
1923
1924 /*!
1925  * This test reproduces bug of Gauthier on 13/9/2010 concerning 3DSurf meshes.
1926  * It is possible to lead to dead lock in InterpKernelDEC when 3DSurfMeshes global bounding boxes intersects whereas cell bounding box intersecting only on one side.
1927  */
1928 void ParaMEDMEMTest::testInterpKernelDEC3DSurfEmptyBBox()
1929 {
1930   int size;
1931   int rank;
1932   MPI_Comm_size(MPI_COMM_WORLD,&size);
1933   MPI_Comm_rank(MPI_COMM_WORLD,&rank);
1934   //
1935   if(size!=3)
1936     return ;
1937   int nproc_source = 1;
1938   set<int> self_procs;
1939   set<int> procs_source;
1940   set<int> procs_target;
1941   
1942   for (int i=0; i<nproc_source; i++)
1943     procs_source.insert(i);
1944   for (int i=nproc_source; i<size; i++)
1945     procs_target.insert(i);
1946   self_procs.insert(rank);
1947   //
1948   MEDCoupling::MEDCouplingUMesh *mesh=0;
1949   MEDCoupling::ParaMESH *paramesh=0;
1950   MEDCoupling::ParaFIELD *parafieldP0=0;
1951   //
1952   MEDCoupling::CommInterface interface;
1953   //
1954   ProcessorGroup* self_group = new MEDCoupling::MPIProcessorGroup(interface,self_procs);
1955   ProcessorGroup* target_group = new MEDCoupling::MPIProcessorGroup(interface,procs_target);
1956   ProcessorGroup* source_group = new MEDCoupling::MPIProcessorGroup(interface,procs_source);
1957   //
1958   MPI_Barrier(MPI_COMM_WORLD);
1959   if(source_group->containsMyRank())
1960     {
1961       double coords[15]={1.,0.,0., 2.,0.,0., 2.,2.,0., 0.,2.,0., 0.5,0.5,1.};
1962       mcIdType conn[7]={0,1,2,3,0,3,4};
1963       mesh=MEDCouplingUMesh::New("Source mesh Proc0",2);
1964       mesh->allocateCells(2);
1965       mesh->insertNextCell(INTERP_KERNEL::NORM_QUAD4,4,conn);
1966       mesh->insertNextCell(INTERP_KERNEL::NORM_TRI3,3,conn+4);
1967       mesh->finishInsertingCells();
1968       DataArrayDouble *myCoords=DataArrayDouble::New();
1969       myCoords->alloc(5,3);
1970       std::copy(coords,coords+15,myCoords->getPointer());
1971       mesh->setCoords(myCoords);
1972       myCoords->decrRef();
1973       //
1974       paramesh=new ParaMESH(mesh,*source_group,"source mesh");
1975       MEDCoupling::ComponentTopology comptopo;
1976       parafieldP0 = new ParaFIELD(ON_CELLS,NO_TIME,paramesh, comptopo);
1977       double *valueP0=parafieldP0->getField()->getArray()->getPointer();
1978       parafieldP0->getField()->setNature(IntensiveMaximum);
1979       valueP0[0]=7.; valueP0[1]=8.;
1980     }
1981   else
1982     {
1983       const char targetMeshName[]="target mesh";
1984       if(rank==1)
1985         {
1986           double coords[12]={0.25,0.25,0.5, 0.,0.25,0.5, 0.,0.,0.5, 0.25,0.,0.5};
1987           mcIdType conn[4]={0,1,2,3};
1988           mesh=MEDCouplingUMesh::New("Target mesh Proc1",2);
1989           mesh->allocateCells(1);
1990           mesh->insertNextCell(INTERP_KERNEL::NORM_QUAD4,4,conn);
1991           mesh->finishInsertingCells();
1992           DataArrayDouble *myCoords=DataArrayDouble::New();
1993           myCoords->alloc(4,3);
1994           std::copy(coords,coords+12,myCoords->getPointer());
1995           mesh->setCoords(myCoords);
1996           myCoords->decrRef();
1997           paramesh=new ParaMESH(mesh,*target_group,targetMeshName);
1998         }
1999       if(rank==2)
2000         {
2001           double coords[12]={0.,0.25,0.5, 0.,0.,0.5, -1.,0.,0.5, -1.,0.25,0.5};
2002           mcIdType conn[4]={0,1,2,3};
2003           mesh=MEDCouplingUMesh::New("Target mesh Proc2",2);
2004           mesh->allocateCells(1);
2005           mesh->insertNextCell(INTERP_KERNEL::NORM_QUAD4,4,conn);
2006           mesh->finishInsertingCells();
2007           DataArrayDouble *myCoords=DataArrayDouble::New();
2008           myCoords->alloc(4,3);
2009           std::copy(coords,coords+12,myCoords->getPointer());
2010           mesh->setCoords(myCoords);
2011           myCoords->decrRef();
2012           paramesh=new ParaMESH(mesh,*target_group,targetMeshName);
2013         }
2014       MEDCoupling::ComponentTopology comptopo;
2015       parafieldP0 = new ParaFIELD(ON_CELLS,NO_TIME,paramesh, comptopo);
2016       parafieldP0->getField()->setNature(IntensiveMaximum);
2017     }
2018   // test 1
2019   MEDCoupling::InterpKernelDEC dec(*source_group,*target_group);
2020   if (source_group->containsMyRank())
2021     { 
2022       dec.setMethod("P0");
2023       dec.attachLocalField(parafieldP0);
2024       dec.synchronize();
2025       // dec.setForcedRenormalization(false);
2026       // dec.sendData();
2027       // dec.recvData();
2028       // const double *valueP0=parafieldP0->getField()->getArray()->getPointer();
2029       // if(rank==0)
2030       //   {
2031       //     CPPUNIT_ASSERT_DOUBLES_EQUAL(7.4,valueP0[0],1e-7);
2032       //     CPPUNIT_ASSERT_DOUBLES_EQUAL(9.0540540540540544,valueP0[1],1e-7);
2033       //   }
2034       // if(rank==1)
2035       //   {
2036       //     CPPUNIT_ASSERT_DOUBLES_EQUAL(8.64054054054054,valueP0[0],1e-7);
2037       //   }
2038       // if(rank==2)
2039       //   {
2040       //     CPPUNIT_ASSERT_DOUBLES_EQUAL(9.0540540540540544,valueP0[0],1e-7);
2041       //   }
2042     }
2043   else
2044     {
2045       dec.setMethod("P0");
2046       dec.attachLocalField(parafieldP0);
2047       dec.synchronize();
2048       // dec.setForcedRenormalization(false);
2049       // dec.recvData();
2050       // const double *res=parafieldP0->getField()->getArray()->getConstPointer();
2051       // if(rank==3)
2052       //   {
2053       //     CPPUNIT_ASSERT_EQUAL(1,parafieldP0->getField()->getNumberOfTuples());
2054       //     CPPUNIT_ASSERT_EQUAL(1,parafieldP0->getField()->getNumberOfComponents());
2055       //     CPPUNIT_ASSERT_DOUBLES_EQUAL(7.4,res[0],1e-12);
2056       //   }
2057       // if(rank==4)
2058       //   {
2059       //     CPPUNIT_ASSERT_EQUAL(1,parafieldP0->getField()->getNumberOfTuples());
2060       //     CPPUNIT_ASSERT_EQUAL(1,parafieldP0->getField()->getNumberOfComponents());
2061       //     CPPUNIT_ASSERT_DOUBLES_EQUAL(9.0540540540540526,res[0],1e-12);
2062       //   }
2063       // dec.sendData();
2064     }
2065   //
2066   delete parafieldP0;
2067   mesh->decrRef();
2068   delete paramesh;
2069   delete self_group;
2070   delete target_group;
2071   delete source_group;
2072   //
2073   MPI_Barrier(MPI_COMM_WORLD);
2074 }
2075
2076 /*!
2077  * Tests an asynchronous exchange between two codes
2078  * one sends data with dtA as an interval, the max time being tmaxA
2079  * the other one receives with dtB as an interval, the max time being tmaxB
2080  */
2081 void ParaMEDMEMTest::testAsynchronousInterpKernelDEC_2D(double dtA, double tmaxA, 
2082                                                         double dtB, double tmaxB, bool WithPointToPoint, bool Asynchronous,
2083                                                         bool WithInterp, const char *srcMeth, const char *targetMeth)
2084 {
2085   std::string srcM(srcMeth);
2086   std::string targetM(targetMeth);
2087   int size;
2088   int rank;
2089   MPI_Comm_size(MPI_COMM_WORLD,&size);
2090   MPI_Comm_rank(MPI_COMM_WORLD,&rank);
2091  
2092   //the test is meant to run on five processors
2093   if (size !=5) return ;
2094    
2095   int nproc_source = 3;
2096   set<int> self_procs;
2097   set<int> procs_source;
2098   set<int> procs_target;
2099   
2100   for (int i=0; i<nproc_source; i++)
2101     procs_source.insert(i);
2102   for (int i=nproc_source; i<size; i++)
2103     procs_target.insert(i);
2104   self_procs.insert(rank);
2105   
2106   MEDCoupling::CommInterface interface;
2107     
2108   MEDCoupling::ProcessorGroup* self_group = new MEDCoupling::MPIProcessorGroup(interface,self_procs);
2109   MEDCoupling::ProcessorGroup* target_group = new MEDCoupling::MPIProcessorGroup(interface,procs_target);
2110   MEDCoupling::ProcessorGroup* source_group = new MEDCoupling::MPIProcessorGroup(interface,procs_source);
2111     
2112   //loading the geometry for the source group
2113
2114   MEDCoupling::InterpKernelDEC dec (*source_group,*target_group);
2115   
2116   MEDCoupling::MEDCouplingUMesh* mesh = nullptr;
2117   MEDCoupling::ParaMESH* paramesh = nullptr;
2118   MEDCoupling::ParaFIELD* parafield = nullptr;
2119   ICoCo::MEDField* icocofield = nullptr;
2120
2121   char * tmp_dir_c                    = getenv("TMP");
2122   string tmp_dir;
2123   if (tmp_dir_c != NULL)
2124     tmp_dir = string(tmp_dir_c);
2125   else
2126     tmp_dir = "/tmp";
2127   string filename_xml1              = "square1_split";
2128   string filename_xml2              = "square2_split";
2129   //string filename_seq_wr            = makeTmpFile("");
2130   //string filename_seq_med           = makeTmpFile("myWrField_seq_pointe221.med");
2131   
2132   // To remove tmp files from disk
2133   ParaMEDMEMTest_TmpFilesRemover aRemover;
2134   
2135   MPI_Barrier(MPI_COMM_WORLD);
2136
2137   if (source_group->containsMyRank())
2138     {
2139       string master = filename_xml1;
2140       
2141       ostringstream strstream;
2142       strstream <<master<<rank+1<<".med";
2143       string fName = INTERP_TEST::getResourceFile(strstream.str());
2144       ostringstream meshname ;
2145       meshname<< "Mesh_2_"<< rank+1;
2146       
2147       mesh=ReadUMeshFromFile(fName.c_str(),meshname.str().c_str(),0);
2148
2149       paramesh=new ParaMESH (mesh,*source_group,"source mesh");
2150     
2151       //      MEDCoupling::ParaSUPPORT* parasupport=new UnstructuredParaSUPPORT( support,*source_group);
2152       MEDCoupling::ComponentTopology comptopo;
2153       if(srcM=="P0")
2154         {
2155           parafield = new ParaFIELD(ON_CELLS,NO_TIME,paramesh, comptopo);
2156           parafield->getField()->setNature(IntensiveMaximum);//InvertIntegral);//IntensiveMaximum);
2157         }
2158       else
2159         parafield = new ParaFIELD(ON_NODES,NO_TIME,paramesh, comptopo);
2160
2161       int nb_local;
2162       if(srcM=="P0")
2163         nb_local=mesh->getNumberOfCells();
2164       else
2165         nb_local=mesh->getNumberOfNodes();
2166       //      double * value= new double[nb_local];
2167       double *value=parafield->getField()->getArray()->getPointer();
2168       for(int ielem=0; ielem<nb_local;ielem++)
2169         value[ielem]=0.0;
2170     
2171       //      ICoCo::Field* icocofield=new ICoCo::MEDField(paramesh,parafield);
2172       icocofield=new ICoCo::MEDField(parafield->getField());
2173      
2174       dec.attachLocalField(icocofield);
2175
2176
2177     }
2178   
2179   //loading the geometry for the target group
2180   if (target_group->containsMyRank())
2181     {
2182       string master= filename_xml2;
2183       ostringstream strstream;
2184       strstream << master<<(rank-nproc_source+1)<<".med";
2185       string fName = INTERP_TEST::getResourceFile(strstream.str());
2186       ostringstream meshname ;
2187       meshname<< "Mesh_3_"<<rank-nproc_source+1;
2188       
2189       mesh = ReadUMeshFromFile(fName.c_str(),meshname.str().c_str(),0);
2190
2191       paramesh=new ParaMESH (mesh,*target_group,"target mesh");
2192       //      MEDCoupling::ParaSUPPORT* parasupport=new UnstructuredParaSUPPORT(support,*target_group);
2193       MEDCoupling::ComponentTopology comptopo;
2194       if(targetM=="P0")
2195         {
2196           parafield = new ParaFIELD(ON_CELLS,NO_TIME,paramesh, comptopo);
2197           parafield->getField()->setNature(IntensiveMaximum);//InvertIntegral);//IntensiveMaximum);
2198         }
2199       else
2200         parafield = new ParaFIELD(ON_NODES,NO_TIME,paramesh, comptopo);
2201       
2202       int nb_local;
2203       if(targetM=="P0")
2204         nb_local=mesh->getNumberOfCells();
2205       else
2206         nb_local=mesh->getNumberOfNodes();
2207                         
2208       double *value=parafield->getField()->getArray()->getPointer();
2209       for(int ielem=0; ielem<nb_local;ielem++)
2210         value[ielem]=0.0;
2211       //      ICoCo::Field* icocofield=new ICoCo::MEDField(paramesh,parafield);
2212       icocofield=new ICoCo::MEDField(parafield->getField());
2213       
2214       dec.attachLocalField(icocofield);
2215     }
2216     
2217   
2218   //attaching a DEC to the source group 
2219   
2220   if (source_group->containsMyRank())
2221     { 
2222       cout<<"DEC usage"<<endl;
2223       dec.setAsynchronous(Asynchronous);
2224       if ( WithInterp ) {
2225         dec.setTimeInterpolationMethod(LinearTimeInterp);
2226       }
2227       if ( WithPointToPoint ) {
2228         dec.setAllToAllMethod(PointToPoint);
2229       }
2230       else {
2231         dec.setAllToAllMethod(Native);
2232       }
2233       dec.synchronize();
2234       dec.setForcedRenormalization(false);
2235       for (double time=0; time<tmaxA+1e-10; time+=dtA)
2236         {
2237           cout << "testAsynchronousInterpKernelDEC_2D" << rank << " time " << time
2238                << " dtA " << dtA << " tmaxA " << tmaxA << endl ;
2239           if ( time+dtA < tmaxA+1e-7 ) {
2240             dec.sendData( time , dtA );
2241           }
2242           else {
2243             dec.sendData( time , 0 );
2244           }
2245           double* value = parafield->getField()->getArray()->getPointer();
2246           int nb_local=parafield->getField()->getMesh()->getNumberOfCells();
2247           for (int i=0; i<nb_local;i++)
2248             value[i]= time+dtA;
2249
2250        
2251         }
2252     }
2253   
2254   //attaching a DEC to the target group
2255   if (target_group->containsMyRank())
2256     {
2257       cout<<"DEC usage"<<endl;
2258       dec.setAsynchronous(Asynchronous);
2259       if ( WithInterp ) {
2260         dec.setTimeInterpolationMethod(LinearTimeInterp);
2261       }
2262       if ( WithPointToPoint ) {
2263         dec.setAllToAllMethod(PointToPoint);
2264       }
2265       else {
2266         dec.setAllToAllMethod(Native);
2267       }
2268       dec.synchronize();
2269       dec.setForcedRenormalization(false);
2270       vector<double> times;
2271       for (double time=0; time<tmaxB+1e-10; time+=dtB)
2272         {
2273           cout << "testAsynchronousInterpKernelDEC_2D" << rank << " time " << time
2274                << " dtB " << dtB << " tmaxB " << tmaxB << endl ;
2275           dec.recvData( time );
2276           double vi = parafield->getVolumeIntegral(0,true);
2277           cout << "testAsynchronousInterpKernelDEC_2D" << rank << " time " << time
2278                << " VolumeIntegral " << vi
2279                << " time*10000 " << time*10000 << endl ;
2280           
2281           CPPUNIT_ASSERT_DOUBLES_EQUAL(vi,time*10000,0.001);
2282         }
2283       
2284     }
2285   
2286   delete source_group;
2287   delete target_group;
2288   delete self_group;
2289   delete parafield ;
2290   delete paramesh ;
2291   mesh->decrRef() ;
2292   delete icocofield ;
2293
2294   cout << "testAsynchronousInterpKernelDEC_2D" << rank << " MPI_Barrier " << endl ;
2295  
2296   if (Asynchronous) MPI_Barrier(MPI_COMM_WORLD);
2297   cout << "end of InterpKernelDEC_2D test"<<endl;
2298 }