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