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