]> SALOME platform Git repositories - tools/medcoupling.git/blob - src/ParaMEDMEMTest/ParaMEDMEMTest_InterpKernelDEC.cxx
Salome HOME
9ca93d3d69541f0d7b56f59ea608583b3e7ff280
[tools/medcoupling.git] / src / ParaMEDMEMTest / ParaMEDMEMTest_InterpKernelDEC.cxx
1 // Copyright (C) 2007-2015  CEA/DEN, EDF R&D
2 //
3 // This library is free software; you can redistribute it and/or
4 // modify it under the terms of the GNU Lesser General Public
5 // License as published by the Free Software Foundation; either
6 // version 2.1 of the License, or (at your option) any later version.
7 //
8 // This library is distributed in the hope that it will be useful,
9 // but WITHOUT ANY WARRANTY; without even the implied warranty of
10 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
11 // Lesser General Public License for more details.
12 //
13 // You should have received a copy of the GNU Lesser General Public
14 // License along with this library; if not, write to the Free Software
15 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA
16 //
17 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
18 //
19
20 #include "ParaMEDMEMTest.hxx"
21 #include <cppunit/TestAssert.h>
22
23 #include "CommInterface.hxx"
24 #include "ProcessorGroup.hxx"
25 #include "MPIProcessorGroup.hxx"
26 #include "Topology.hxx"
27 #include "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::MEDField* 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(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(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::MEDField* icocofield ;
861   
862   char * tmp_dir_c                    = getenv("TMP");
863   string tmp_dir;
864   if (tmp_dir_c != NULL)
865     tmp_dir = string(tmp_dir_c);
866   else
867     tmp_dir = "/tmp";
868   string filename_xml1              = getResourceFile("Mesh3D_10_2d");
869   string filename_xml2              = getResourceFile("Mesh3D_11");
870   //string filename_seq_wr            = makeTmpFile("");
871   //string filename_seq_med           = makeTmpFile("myWrField_seq_pointe221.med");
872   
873   // To remove tmp files from disk
874   ParaMEDMEMTest_TmpFilesRemover aRemover;
875   
876   MPI_Barrier(MPI_COMM_WORLD);
877   if (source_group->containsMyRank())
878     {
879       string master = filename_xml1;
880       
881       ostringstream strstream;
882       strstream <<master<<rank+1<<".med";
883       ostringstream meshname ;
884       meshname<< "Mesh_3_"<< rank+1;
885       
886       mesh=MEDLoader::ReadUMeshFromFile(strstream.str().c_str(),meshname.str().c_str(),0);
887       
888     
889       paramesh=new ParaMESH (mesh,*source_group,"source mesh");
890     
891       //      ParaMEDMEM::ParaSUPPORT* parasupport=new UnstructuredParaSUPPORT( support,*source_group);
892       ParaMEDMEM::ComponentTopology comptopo;
893       if(srcM=="P0")
894         {
895           parafield = new ParaFIELD(ON_CELLS,NO_TIME,paramesh, comptopo);
896           parafield->getField()->setNature(ConservativeVolumic);
897         }
898       else
899         parafield = new ParaFIELD(ON_NODES,NO_TIME,paramesh, comptopo);
900       int nb_local;
901       if(srcM=="P0")
902         nb_local=mesh->getNumberOfCells();
903       else
904         nb_local=mesh->getNumberOfNodes();
905       //      double * value= new double[nb_local];
906       double *value=parafield->getField()->getArray()->getPointer();
907       for(int ielem=0; ielem<nb_local;ielem++)
908         value[ielem]=1.0;
909     
910       //      ICoCo::Field* icocofield=new ICoCo::MEDField(paramesh,parafield);
911       icocofield=new ICoCo::MEDField(parafield->getField());
912       dec.setMethod(srcMeth);
913       dec.attachLocalField(icocofield);
914     }
915   
916   //loading the geometry for the target group
917   if (target_group->containsMyRank())
918     {
919       string master= filename_xml2;
920       ostringstream strstream;
921       strstream << master << ".med";
922       ostringstream meshname ;
923       meshname<< "Mesh_6";
924       mesh = MEDLoader::ReadUMeshFromFile(strstream.str().c_str(),meshname.str().c_str(),0);
925       
926       paramesh=new ParaMESH (mesh,*target_group,"target mesh");
927       //      ParaMEDMEM::ParaSUPPORT* parasupport=new UnstructuredParaSUPPORT(support,*target_group);
928       ParaMEDMEM::ComponentTopology comptopo;
929       if(targetM=="P0")
930         {
931           parafield = new ParaFIELD(ON_CELLS,NO_TIME,paramesh, comptopo);
932           parafield->getField()->setNature(ConservativeVolumic);
933         }
934       else
935         parafield = new ParaFIELD(ON_NODES,NO_TIME,paramesh, comptopo);
936       int nb_local;
937       if(targetM=="P0")
938         nb_local=mesh->getNumberOfCells();
939       else
940         nb_local=mesh->getNumberOfNodes();
941       //      double * value= new double[nb_local];
942       double *value=parafield->getField()->getArray()->getPointer();
943       for(int ielem=0; ielem<nb_local;ielem++)
944         value[ielem]=0.0;
945       //      ICoCo::Field* icocofield=new ICoCo::MEDField(paramesh,parafield);
946       icocofield=new ICoCo::MEDField(parafield->getField());
947       dec.setMethod(targetMeth);
948       dec.attachLocalField(icocofield);
949     }  
950   //attaching a DEC to the source group 
951   double field_before_int;
952   double field_after_int;
953   
954   if (source_group->containsMyRank())
955     { 
956       field_before_int = parafield->getVolumeIntegral(0,true);
957       dec.synchronize();
958       cout<<"DEC usage"<<endl;
959       dec.setForcedRenormalization(false);
960
961       dec.sendData();
962       ParaMEDLoader::WriteParaMesh("./sourcesquareb",paramesh);
963       if (source_group->myRank()==0)
964         aRemover.Register("./sourcesquareb");
965       ostringstream filename;
966       filename<<"./sourcesquareb_"<<source_group->myRank()+1;
967       aRemover.Register(filename.str().c_str());
968       //MEDLoader::WriteField("./sourcesquareb",parafield->getField());
969    
970       dec.recvData();
971       cout <<"writing"<<endl;
972       ParaMEDLoader::WriteParaMesh("./sourcesquare",paramesh);
973       if (source_group->myRank()==0)
974         aRemover.Register("./sourcesquare");
975       //MEDLoader::WriteField("./sourcesquare",parafield->getField());
976       
977      
978       filename<<"./sourcesquare_"<<source_group->myRank()+1;
979       aRemover.Register(filename.str().c_str());
980       field_after_int = parafield->getVolumeIntegral(0,true);
981
982       CPPUNIT_ASSERT_DOUBLES_EQUAL(field_before_int, field_after_int, 1e-6);
983     
984     }
985   
986   //attaching a DEC to the target group
987   if (target_group->containsMyRank())
988     {
989       dec.synchronize();
990       dec.setForcedRenormalization(false);
991
992       dec.recvData();
993       ParaMEDLoader::WriteParaMesh("./targetsquareb",paramesh);
994       //MEDLoader::WriteField("./targetsquareb",parafield->getField());
995       if (target_group->myRank()==0)
996         aRemover.Register("./targetsquareb");
997       ostringstream filename;
998       filename<<"./targetsquareb_"<<target_group->myRank()+1;
999       aRemover.Register(filename.str().c_str());
1000       dec.sendData();
1001       ParaMEDLoader::WriteParaMesh("./targetsquare",paramesh);
1002       //MEDLoader::WriteField("./targetsquare",parafield->getField());
1003       
1004       if (target_group->myRank()==0)
1005         aRemover.Register("./targetsquareb");
1006       
1007       filename<<"./targetsquareb_"<<target_group->myRank()+1;
1008       aRemover.Register(filename.str().c_str());
1009     }
1010   delete source_group;
1011   delete target_group;
1012   delete self_group;
1013   delete parafield;
1014   delete paramesh;
1015   mesh->decrRef();
1016
1017   delete icocofield;
1018
1019   MPI_Barrier(MPI_COMM_WORLD);
1020   cout << "end of InterpKernelDEC_3D test"<<endl;
1021 }
1022
1023 //Synchronous tests without interpolation with native mode (AllToAll(v) from lam/MPI:
1024 void ParaMEDMEMTest::testSynchronousEqualInterpKernelWithoutInterpNativeDEC_2D()
1025 {
1026   testAsynchronousInterpKernelDEC_2D(0.1,1,0.1,1,false,false,false,"P0","P0");
1027 }
1028
1029 //Synchronous tests without interpolation :
1030 void ParaMEDMEMTest::testSynchronousEqualInterpKernelWithoutInterpDEC_2D()
1031 {
1032   testAsynchronousInterpKernelDEC_2D(0.1,1,0.1,1,true,false,false,"P0","P0");
1033 }
1034
1035 //Synchronous tests with interpolation :
1036 void ParaMEDMEMTest::testSynchronousEqualInterpKernelDEC_2D()
1037 {
1038   testAsynchronousInterpKernelDEC_2D(0.1,1,0.1,1,true,false,true,"P0","P0");
1039 }
1040 void ParaMEDMEMTest::testSynchronousFasterSourceInterpKernelDEC_2D()
1041 {
1042   testAsynchronousInterpKernelDEC_2D(0.09,1,0.1,1,true,false,true,"P0","P0");
1043 }
1044 void ParaMEDMEMTest::testSynchronousSlowerSourceInterpKernelDEC_2D()
1045 {
1046   testAsynchronousInterpKernelDEC_2D(0.11,1,0.1,1,true,false,true,"P0","P0");
1047 }
1048 void ParaMEDMEMTest::testSynchronousSlowSourceInterpKernelDEC_2D()
1049 {
1050   testAsynchronousInterpKernelDEC_2D(0.11,1,0.01,1,true,false,true,"P0","P0");
1051 }
1052 void ParaMEDMEMTest::testSynchronousFastSourceInterpKernelDEC_2D()
1053 {
1054   testAsynchronousInterpKernelDEC_2D(0.01,1,0.11,1,true,false,true,"P0","P0");
1055 }
1056
1057 //Asynchronous tests with interpolation :
1058 void ParaMEDMEMTest::testAsynchronousEqualInterpKernelDEC_2D()
1059 {
1060   testAsynchronousInterpKernelDEC_2D(0.1,1,0.1,1,true,true,true,"P0","P0");
1061 }
1062 void ParaMEDMEMTest::testAsynchronousFasterSourceInterpKernelDEC_2D()
1063 {
1064   testAsynchronousInterpKernelDEC_2D(0.09,1,0.1,1,true,true,true,"P0","P0");
1065 }
1066 void ParaMEDMEMTest::testAsynchronousSlowerSourceInterpKernelDEC_2D()
1067 {
1068   testAsynchronousInterpKernelDEC_2D(0.11,1,0.1,1,true,true,true,"P0","P0");
1069 }
1070 void ParaMEDMEMTest::testAsynchronousSlowSourceInterpKernelDEC_2D()
1071 {
1072   testAsynchronousInterpKernelDEC_2D(0.11,1,0.01,1,true,true,true,"P0","P0");
1073 }
1074 void ParaMEDMEMTest::testAsynchronousFastSourceInterpKernelDEC_2D()
1075 {
1076   testAsynchronousInterpKernelDEC_2D(0.01,1,0.11,1,true,true,true,"P0","P0");
1077 }
1078
1079 void ParaMEDMEMTest::testInterpKernelDECNonOverlapp_2D_P0P0()
1080 {
1081   //
1082   const double sourceCoordsAll[2][8]={{0.4,0.5,0.4,1.5,1.6,1.5,1.6,0.5},
1083                                       {0.3,-0.5,1.6,-0.5,1.6,-1.5,0.3,-1.5}};
1084   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},
1085                                        {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},
1086                                        {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}};
1087   int conn4All[8]={0,1,2,3,4,5,6,7};
1088   double targetResults[3][2]={{34.,34.},{38.333333333333336,42.666666666666664},{47.,47.}};
1089   double targetResults2[3][2]={{0.28333333333333344,0.56666666666666687},{1.8564102564102569,2.0128205128205132},{1.0846153846153845,0.36153846153846159}};
1090   double targetResults3[3][2]={{3.7777777777777781,7.5555555555555562},{24.511111111111113,26.355555555555558},{14.1,4.7}};
1091   double targetResults4[3][2]={{8.5,17},{8.8461538461538431, 9.8461538461538449},{35.25,11.75}};
1092   //
1093   int size;
1094   int rank;
1095   MPI_Comm_size(MPI_COMM_WORLD,&size);
1096   MPI_Comm_rank(MPI_COMM_WORLD,&rank);
1097   //
1098   if(size!=5)
1099     return ;
1100   int nproc_source = 2;
1101   set<int> self_procs;
1102   set<int> procs_source;
1103   set<int> procs_target;
1104   
1105   for (int i=0; i<nproc_source; i++)
1106     procs_source.insert(i);
1107   for (int i=nproc_source; i<size; i++)
1108     procs_target.insert(i);
1109   self_procs.insert(rank);
1110   //
1111   ParaMEDMEM::MEDCouplingUMesh *mesh=0;
1112   ParaMEDMEM::ParaMESH *paramesh=0;
1113   ParaMEDMEM::ParaFIELD* parafield=0;
1114   //
1115   ParaMEDMEM::CommInterface interface;
1116   //
1117   ProcessorGroup* self_group = new ParaMEDMEM::MPIProcessorGroup(interface,self_procs);
1118   ProcessorGroup* target_group = new ParaMEDMEM::MPIProcessorGroup(interface,procs_target);
1119   ProcessorGroup* source_group = new ParaMEDMEM::MPIProcessorGroup(interface,procs_source);
1120   //
1121   MPI_Barrier(MPI_COMM_WORLD);
1122   if(source_group->containsMyRank())
1123     {
1124       std::ostringstream stream; stream << "sourcemesh2D proc " << rank;
1125       mesh=MEDCouplingUMesh::New(stream.str().c_str(),2);
1126       mesh->allocateCells(2);
1127       mesh->insertNextCell(INTERP_KERNEL::NORM_QUAD4,4,conn4All);
1128       mesh->finishInsertingCells();
1129       DataArrayDouble *myCoords=DataArrayDouble::New();
1130       myCoords->alloc(4,2);
1131       const double *sourceCoords=sourceCoordsAll[rank];
1132       std::copy(sourceCoords,sourceCoords+8,myCoords->getPointer());
1133       mesh->setCoords(myCoords);
1134       myCoords->decrRef();
1135       paramesh=new ParaMESH(mesh,*source_group,"source mesh");
1136       ParaMEDMEM::ComponentTopology comptopo;
1137       parafield = new ParaFIELD(ON_CELLS,NO_TIME,paramesh, comptopo);
1138       double *value=parafield->getField()->getArray()->getPointer();
1139       value[0]=34+13*((double)rank);
1140     }
1141   else
1142     {
1143       std::ostringstream stream; stream << "targetmesh2D proc " << rank-nproc_source;
1144       mesh=MEDCouplingUMesh::New(stream.str().c_str(),2);
1145       mesh->allocateCells(2);
1146       mesh->insertNextCell(INTERP_KERNEL::NORM_QUAD4,4,conn4All);
1147       mesh->insertNextCell(INTERP_KERNEL::NORM_QUAD4,4,conn4All+4);
1148       mesh->finishInsertingCells();
1149       DataArrayDouble *myCoords=DataArrayDouble::New();
1150       myCoords->alloc(8,2);
1151       const double *targetCoords=targetCoordsAll[rank-nproc_source];
1152       std::copy(targetCoords,targetCoords+16,myCoords->getPointer());
1153       mesh->setCoords(myCoords);
1154       myCoords->decrRef();
1155       paramesh=new ParaMESH (mesh,*target_group,"target mesh");
1156       ParaMEDMEM::ComponentTopology comptopo;
1157       parafield = new ParaFIELD(ON_CELLS,NO_TIME,paramesh, comptopo);
1158     }
1159   //test 1 - Conservative volumic
1160   ParaMEDMEM::InterpKernelDEC dec(*source_group,*target_group);
1161   parafield->getField()->setNature(ConservativeVolumic);
1162   if (source_group->containsMyRank())
1163     { 
1164       dec.setMethod("P0");
1165       dec.attachLocalField(parafield);
1166       dec.synchronize();
1167       dec.setForcedRenormalization(false);
1168       dec.sendData();
1169     }
1170   else
1171     {
1172       dec.setMethod("P0");
1173       dec.attachLocalField(parafield);
1174       dec.synchronize();
1175       dec.setForcedRenormalization(false);
1176       dec.recvData();
1177       const double *res=parafield->getField()->getArray()->getConstPointer();
1178       const double *expected=targetResults[rank-nproc_source];
1179       CPPUNIT_ASSERT_DOUBLES_EQUAL(expected[0],res[0],1e-13);
1180       CPPUNIT_ASSERT_DOUBLES_EQUAL(expected[1],res[1],1e-13);
1181     }
1182   //test 2 - Integral
1183   ParaMEDMEM::InterpKernelDEC dec2(*source_group,*target_group);
1184   parafield->getField()->setNature(Integral);
1185   if (source_group->containsMyRank())
1186     { 
1187       dec2.setMethod("P0");
1188       dec2.attachLocalField(parafield);
1189       dec2.synchronize();
1190       dec2.setForcedRenormalization(false);
1191       dec2.sendData();
1192     }
1193   else
1194     {
1195       dec2.setMethod("P0");
1196       dec2.attachLocalField(parafield);
1197       dec2.synchronize();
1198       dec2.setForcedRenormalization(false);
1199       dec2.recvData();
1200       const double *res=parafield->getField()->getArray()->getConstPointer();
1201       const double *expected=targetResults2[rank-nproc_source];
1202       CPPUNIT_ASSERT_DOUBLES_EQUAL(expected[0],res[0],1e-13);
1203       CPPUNIT_ASSERT_DOUBLES_EQUAL(expected[1],res[1],1e-13);
1204     }
1205   //test 3 - Integral with global constraint
1206   ParaMEDMEM::InterpKernelDEC dec3(*source_group,*target_group);
1207   parafield->getField()->setNature(IntegralGlobConstraint);
1208   if (source_group->containsMyRank())
1209     { 
1210       dec3.setMethod("P0");
1211       dec3.attachLocalField(parafield);
1212       dec3.synchronize();
1213       dec3.setForcedRenormalization(false);
1214       dec3.sendData();
1215     }
1216   else
1217     {
1218       dec3.setMethod("P0");
1219       dec3.attachLocalField(parafield);
1220       dec3.synchronize();
1221       dec3.setForcedRenormalization(false);
1222       dec3.recvData();
1223       const double *res=parafield->getField()->getArray()->getConstPointer();
1224       const double *expected=targetResults3[rank-nproc_source];
1225       CPPUNIT_ASSERT_DOUBLES_EQUAL(expected[0],res[0],1e-13);
1226       CPPUNIT_ASSERT_DOUBLES_EQUAL(expected[1],res[1],1e-13);
1227     }
1228   //test 4 - RevIntegral
1229   ParaMEDMEM::InterpKernelDEC dec4(*source_group,*target_group);
1230   parafield->getField()->setNature(RevIntegral);
1231   if (source_group->containsMyRank())
1232     { 
1233       dec4.setMethod("P0");
1234       dec4.attachLocalField(parafield);
1235       dec4.synchronize();
1236       dec4.setForcedRenormalization(false);
1237       dec4.sendData();
1238     }
1239   else
1240     {
1241       dec4.setMethod("P0");
1242       dec4.attachLocalField(parafield);
1243       dec4.synchronize();
1244       dec4.setForcedRenormalization(false);
1245       dec4.recvData();
1246       const double *res=parafield->getField()->getArray()->getConstPointer();
1247       const double *expected=targetResults4[rank-nproc_source];
1248       CPPUNIT_ASSERT_DOUBLES_EQUAL(expected[0],res[0],1e-13);
1249       CPPUNIT_ASSERT_DOUBLES_EQUAL(expected[1],res[1],1e-13);
1250     }
1251   //test 5 - Conservative volumic reversed
1252   ParaMEDMEM::InterpKernelDEC dec5(*source_group,*target_group);
1253   parafield->getField()->setNature(ConservativeVolumic);
1254   if (source_group->containsMyRank())
1255     { 
1256       dec5.setMethod("P0");
1257       dec5.attachLocalField(parafield);
1258       dec5.synchronize();
1259       dec5.setForcedRenormalization(false);
1260       dec5.recvData();
1261       const double *res=parafield->getField()->getArray()->getConstPointer();
1262       CPPUNIT_ASSERT_EQUAL(1,parafield->getField()->getNumberOfTuples());
1263       const double expected[]={37.8518518518519,43.5333333333333};
1264       CPPUNIT_ASSERT_DOUBLES_EQUAL(expected[rank],res[0],1e-13);
1265     }
1266   else
1267     {
1268       dec5.setMethod("P0");
1269       dec5.attachLocalField(parafield);
1270       dec5.synchronize();
1271       dec5.setForcedRenormalization(false);
1272       double *res=parafield->getField()->getArray()->getPointer();
1273       const double *toSet=targetResults[rank-nproc_source];
1274       res[0]=toSet[0];
1275       res[1]=toSet[1];
1276       dec5.sendData();
1277     }
1278   //test 6 - Integral reversed
1279   ParaMEDMEM::InterpKernelDEC dec6(*source_group,*target_group);
1280   parafield->getField()->setNature(Integral);
1281   if (source_group->containsMyRank())
1282     { 
1283       dec6.setMethod("P0");
1284       dec6.attachLocalField(parafield);
1285       dec6.synchronize();
1286       dec6.setForcedRenormalization(false);
1287       dec6.recvData();
1288       const double *res=parafield->getField()->getArray()->getConstPointer();
1289       CPPUNIT_ASSERT_EQUAL(1,parafield->getField()->getNumberOfTuples());
1290       const double expected[]={0.794600591715977,1.35631163708087};
1291       CPPUNIT_ASSERT_DOUBLES_EQUAL(expected[rank],res[0],1e-13);
1292     }
1293   else
1294     {
1295       dec6.setMethod("P0");
1296       dec6.attachLocalField(parafield);
1297       dec6.synchronize();
1298       dec6.setForcedRenormalization(false);
1299       double *res=parafield->getField()->getArray()->getPointer();
1300       const double *toSet=targetResults2[rank-nproc_source];
1301       res[0]=toSet[0];
1302       res[1]=toSet[1];
1303       dec6.sendData();
1304     }
1305   //test 7 - Integral with global constraint reversed
1306   ParaMEDMEM::InterpKernelDEC dec7(*source_group,*target_group);
1307   parafield->getField()->setNature(IntegralGlobConstraint);
1308   if (source_group->containsMyRank())
1309     { 
1310       dec7.setMethod("P0");
1311       dec7.attachLocalField(parafield);
1312       dec7.synchronize();
1313       dec7.setForcedRenormalization(false);
1314       dec7.recvData();
1315       const double *res=parafield->getField()->getArray()->getConstPointer();
1316       CPPUNIT_ASSERT_EQUAL(1,parafield->getField()->getNumberOfTuples());
1317       const double expected[]={36.4592592592593,44.5407407407407};
1318       CPPUNIT_ASSERT_DOUBLES_EQUAL(expected[rank],res[0],1e-13);
1319     }
1320   else
1321     {
1322       dec7.setMethod("P0");
1323       dec7.attachLocalField(parafield);
1324       dec7.synchronize();
1325       dec7.setForcedRenormalization(false);
1326       double *res=parafield->getField()->getArray()->getPointer();
1327       const double *toSet=targetResults3[rank-nproc_source];
1328       res[0]=toSet[0];
1329       res[1]=toSet[1];
1330       dec7.sendData();
1331     }
1332   //test 8 - Integral with RevIntegral reversed
1333   ParaMEDMEM::InterpKernelDEC dec8(*source_group,*target_group);
1334   parafield->getField()->setNature(RevIntegral);
1335   if (source_group->containsMyRank())
1336     { 
1337       dec8.setMethod("P0");
1338       dec8.attachLocalField(parafield);
1339       dec8.synchronize();
1340       dec8.setForcedRenormalization(false);
1341       dec8.recvData();
1342       const double *res=parafield->getField()->getArray()->getConstPointer();
1343       CPPUNIT_ASSERT_EQUAL(1,parafield->getField()->getNumberOfTuples());
1344       const double expected[]={0.81314102564102553,1.3428994082840233};
1345       CPPUNIT_ASSERT_DOUBLES_EQUAL(expected[rank],res[0],1e-13);
1346     }
1347   else
1348     {
1349       dec8.setMethod("P0");
1350       dec8.attachLocalField(parafield);
1351       dec8.synchronize();
1352       dec8.setForcedRenormalization(false);
1353       double *res=parafield->getField()->getArray()->getPointer();
1354       const double *toSet=targetResults4[rank-nproc_source];
1355       res[0]=toSet[0];
1356       res[1]=toSet[1];
1357       dec8.sendData();
1358     }
1359   //
1360   delete parafield;
1361   mesh->decrRef();
1362   delete paramesh;
1363   delete self_group;
1364   delete target_group;
1365   delete source_group;
1366   //
1367   MPI_Barrier(MPI_COMM_WORLD);
1368 }
1369
1370 void ParaMEDMEMTest::testInterpKernelDECNonOverlapp_2D_P0P1P1P0()
1371 {
1372   int size;
1373   int rank;
1374   MPI_Comm_size(MPI_COMM_WORLD,&size);
1375   MPI_Comm_rank(MPI_COMM_WORLD,&rank);
1376   //
1377   if(size!=5)
1378     return ;
1379   int nproc_source = 2;
1380   set<int> self_procs;
1381   set<int> procs_source;
1382   set<int> procs_target;
1383   
1384   for (int i=0; i<nproc_source; i++)
1385     procs_source.insert(i);
1386   for (int i=nproc_source; i<size; i++)
1387     procs_target.insert(i);
1388   self_procs.insert(rank);
1389   //
1390   ParaMEDMEM::MEDCouplingUMesh *mesh=0;
1391   ParaMEDMEM::ParaMESH *paramesh=0;
1392   ParaMEDMEM::ParaFIELD *parafieldP0=0,*parafieldP1=0;
1393   //
1394   ParaMEDMEM::CommInterface interface;
1395   //
1396   ProcessorGroup* self_group = new ParaMEDMEM::MPIProcessorGroup(interface,self_procs);
1397   ProcessorGroup* target_group = new ParaMEDMEM::MPIProcessorGroup(interface,procs_target);
1398   ProcessorGroup* source_group = new ParaMEDMEM::MPIProcessorGroup(interface,procs_source);
1399   //
1400   MPI_Barrier(MPI_COMM_WORLD);
1401   if(source_group->containsMyRank())
1402     {
1403       if(rank==0)
1404         {
1405           double coords[6]={-0.3,-0.3, 0.7,0.7, 0.7,-0.3};
1406           int conn[3]={0,1,2};
1407           //int globalNode[3]={1,2,0};
1408           mesh=MEDCouplingUMesh::New("Source mesh Proc0",2);
1409           mesh->allocateCells(1);
1410           mesh->insertNextCell(INTERP_KERNEL::NORM_TRI3,3,conn);
1411           mesh->finishInsertingCells();
1412           DataArrayDouble *myCoords=DataArrayDouble::New();
1413           myCoords->alloc(3,2);
1414           std::copy(coords,coords+6,myCoords->getPointer());
1415           mesh->setCoords(myCoords);
1416           myCoords->decrRef();
1417         }
1418       if(rank==1)
1419         {
1420           double coords[6]={-0.3,-0.3, -0.3,0.7, 0.7,0.7};
1421           int conn[3]={0,1,2};
1422           //int globalNode[3]={1,3,2};
1423           mesh=MEDCouplingUMesh::New("Source mesh Proc1",2);
1424           mesh->allocateCells(1);
1425           mesh->insertNextCell(INTERP_KERNEL::NORM_TRI3,3,conn);
1426           mesh->finishInsertingCells();
1427           DataArrayDouble *myCoords=DataArrayDouble::New();
1428           myCoords->alloc(3,2);
1429           std::copy(coords,coords+6,myCoords->getPointer());
1430           mesh->setCoords(myCoords);
1431           myCoords->decrRef();
1432         }
1433       paramesh=new ParaMESH(mesh,*source_group,"source mesh");
1434       ParaMEDMEM::ComponentTopology comptopo;
1435       parafieldP0 = new ParaFIELD(ON_CELLS,NO_TIME,paramesh, comptopo);
1436       parafieldP1 = new ParaFIELD(ON_NODES,NO_TIME,paramesh, comptopo);
1437       double *valueP0=parafieldP0->getField()->getArray()->getPointer();
1438       double *valueP1=parafieldP1->getField()->getArray()->getPointer();
1439       parafieldP0->getField()->setNature(ConservativeVolumic);
1440       parafieldP1->getField()->setNature(ConservativeVolumic);
1441       if(rank==0)
1442         {
1443           valueP0[0]=31.;
1444           valueP1[0]=34.; valueP1[1]=77.; valueP1[2]=53.;
1445         }
1446       if(rank==1)
1447         {
1448           valueP0[0]=47.;
1449           valueP1[0]=34.; valueP1[1]=57.; valueP1[2]=77.;
1450         }
1451     }
1452   else
1453     {
1454       const char targetMeshName[]="target mesh";
1455       if(rank==2)
1456         {
1457           double coords[10]={-0.3,-0.3, 0.2,-0.3, 0.7,-0.3, -0.3,0.2, 0.2,0.2 };
1458           int conn[7]={0,3,4,1, 1,4,2};
1459           //int globalNode[5]={4,3,0,2,1};
1460           mesh=MEDCouplingUMesh::New("Target mesh Proc2",2);
1461           mesh->allocateCells(2);
1462           mesh->insertNextCell(INTERP_KERNEL::NORM_QUAD4,4,conn);
1463           mesh->insertNextCell(INTERP_KERNEL::NORM_TRI3,3,conn+4);
1464           mesh->finishInsertingCells();
1465           DataArrayDouble *myCoords=DataArrayDouble::New();
1466           myCoords->alloc(5,2);
1467           std::copy(coords,coords+10,myCoords->getPointer());
1468           mesh->setCoords(myCoords);
1469           myCoords->decrRef();
1470           paramesh=new ParaMESH(mesh,*target_group,targetMeshName);
1471           DataArrayInt *da=DataArrayInt::New();
1472           const int globalNumberingP2[5]={0,1,2,3,4};
1473           da->useArray(globalNumberingP2,false,CPP_DEALLOC,5,1);
1474           paramesh->setNodeGlobal(da);
1475           da->decrRef();
1476         }
1477       if(rank==3)
1478         {
1479           double coords[6]={0.2,0.2, 0.7,-0.3, 0.7,0.2};
1480           int conn[3]={0,2,1};
1481           //int globalNode[3]={1,0,5};
1482           mesh=MEDCouplingUMesh::New("Target mesh Proc3",2);
1483           mesh->allocateCells(1);
1484           mesh->insertNextCell(INTERP_KERNEL::NORM_TRI3,3,conn);
1485           mesh->finishInsertingCells();
1486           DataArrayDouble *myCoords=DataArrayDouble::New();
1487           myCoords->alloc(3,2);
1488           std::copy(coords,coords+6,myCoords->getPointer());
1489           mesh->setCoords(myCoords);
1490           myCoords->decrRef();
1491           paramesh=new ParaMESH(mesh,*target_group,targetMeshName);
1492           DataArrayInt *da=DataArrayInt::New();
1493           const int globalNumberingP3[3]={4,2,5};
1494           da->useArray(globalNumberingP3,false,CPP_DEALLOC,3,1);
1495           paramesh->setNodeGlobal(da);
1496           da->decrRef();
1497         }
1498       if(rank==4)
1499         {
1500           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};
1501           int conn[8]={0,1,2,3, 3,2,4,5};
1502           //int globalNode[6]={2,6,7,1,8,5};
1503           mesh=MEDCouplingUMesh::New("Target mesh Proc4",2);
1504           mesh->allocateCells(2);
1505           mesh->insertNextCell(INTERP_KERNEL::NORM_QUAD4,4,conn);
1506           mesh->insertNextCell(INTERP_KERNEL::NORM_QUAD4,4,conn+4);
1507           mesh->finishInsertingCells();
1508           DataArrayDouble *myCoords=DataArrayDouble::New();
1509           myCoords->alloc(6,2);
1510           std::copy(coords,coords+12,myCoords->getPointer());
1511           mesh->setCoords(myCoords);
1512           myCoords->decrRef();
1513           paramesh=new ParaMESH(mesh,*target_group,targetMeshName);
1514           DataArrayInt *da=DataArrayInt::New();
1515           const int globalNumberingP4[6]={3,6,7,4,8,5};
1516           da->useArray(globalNumberingP4,false,CPP_DEALLOC,6,1);
1517           paramesh->setNodeGlobal(da);
1518           da->decrRef();
1519         }
1520       ParaMEDMEM::ComponentTopology comptopo;
1521       parafieldP0 = new ParaFIELD(ON_CELLS,NO_TIME,paramesh, comptopo);
1522       parafieldP1 = new ParaFIELD(ON_NODES,NO_TIME,paramesh, comptopo);
1523       parafieldP0->getField()->setNature(ConservativeVolumic);
1524       parafieldP1->getField()->setNature(ConservativeVolumic);
1525     }
1526   // test 1 - P0 P1
1527   ParaMEDMEM::InterpKernelDEC dec(*source_group,*target_group);
1528   if (source_group->containsMyRank())
1529     { 
1530       dec.setMethod("P0");
1531       dec.attachLocalField(parafieldP0);
1532       dec.synchronize();
1533       dec.setForcedRenormalization(false);
1534       dec.sendData();
1535       dec.recvData();
1536       const double *valueP0=parafieldP0->getField()->getArray()->getPointer();
1537       if(rank==0)
1538         {
1539           CPPUNIT_ASSERT_DOUBLES_EQUAL(34.42857143,valueP0[0],1e-7);
1540         }
1541       if(rank==1)
1542         {
1543           CPPUNIT_ASSERT_DOUBLES_EQUAL(44.,valueP0[0],1e-7);
1544         }
1545     }
1546   else
1547     {
1548       dec.setMethod("P1");
1549       dec.attachLocalField(parafieldP1);
1550       dec.synchronize();
1551       dec.setForcedRenormalization(false);
1552       dec.recvData();
1553       const double *res=parafieldP1->getField()->getArray()->getConstPointer();
1554       if(rank==2)
1555         {
1556           const double expectP2[5]={39.0, 31.0, 31.0, 47.0, 39.0};
1557           CPPUNIT_ASSERT_EQUAL(5,parafieldP1->getField()->getNumberOfTuples());
1558           CPPUNIT_ASSERT_EQUAL(1,parafieldP1->getField()->getNumberOfComponents());
1559           for(int kk=0;kk<5;kk++)
1560             CPPUNIT_ASSERT_DOUBLES_EQUAL(expectP2[kk],res[kk],1e-12);
1561         }
1562       if(rank==3)
1563         {
1564           const double expectP3[3]={39.0, 31.0, 31.0};
1565           CPPUNIT_ASSERT_EQUAL(3,parafieldP1->getField()->getNumberOfTuples());
1566           CPPUNIT_ASSERT_EQUAL(1,parafieldP1->getField()->getNumberOfComponents());
1567           for(int kk=0;kk<3;kk++)
1568             CPPUNIT_ASSERT_DOUBLES_EQUAL(expectP3[kk],res[kk],1e-12);
1569         }
1570       if(rank==4)
1571         {
1572           const double expectP4[6]={47.0, 47.0, 47.0, 39.0, 39.0, 31.0};
1573           CPPUNIT_ASSERT_EQUAL(6,parafieldP1->getField()->getNumberOfTuples());
1574           CPPUNIT_ASSERT_EQUAL(1,parafieldP1->getField()->getNumberOfComponents());
1575           for(int kk=0;kk<6;kk++)
1576             CPPUNIT_ASSERT_DOUBLES_EQUAL(expectP4[kk],res[kk],1e-12);
1577         }
1578       dec.sendData();
1579     }
1580   //
1581   delete parafieldP0;
1582   delete parafieldP1;
1583   mesh->decrRef();
1584   delete paramesh;
1585   delete self_group;
1586   delete target_group;
1587   delete source_group;
1588   //
1589   MPI_Barrier(MPI_COMM_WORLD);
1590 }
1591
1592 void ParaMEDMEMTest::testInterpKernelDEC2DM1D_P0P0()
1593 {
1594   int size;
1595   int rank;
1596   MPI_Comm_size(MPI_COMM_WORLD,&size);
1597   MPI_Comm_rank(MPI_COMM_WORLD,&rank);
1598   //
1599   if(size!=3)
1600     return ;
1601   int nproc_source=2;
1602   set<int> procs_source;
1603   set<int> procs_target;
1604   //
1605   for (int i=0; i<nproc_source; i++)
1606     procs_source.insert(i);
1607   for (int i=nproc_source;i<size; i++)
1608     procs_target.insert(i);
1609   //
1610   ParaMEDMEM::MEDCouplingUMesh *mesh=0;
1611   ParaMEDMEM::ParaMESH *paramesh=0;
1612   ParaMEDMEM::ParaFIELD *parafield=0;
1613   //
1614   ParaMEDMEM::CommInterface interface;
1615   //
1616   ProcessorGroup* target_group = new ParaMEDMEM::MPIProcessorGroup(interface,procs_target);
1617   ProcessorGroup* source_group = new ParaMEDMEM::MPIProcessorGroup(interface,procs_source);
1618   //
1619   MPI_Barrier(MPI_COMM_WORLD);
1620   if(source_group->containsMyRank())
1621     {
1622       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 };
1623       mesh=MEDCouplingUMesh::New();
1624       mesh->setMeshDimension(2);
1625       DataArrayDouble *myCoords=DataArrayDouble::New();
1626       myCoords->alloc(9,2);
1627       std::copy(targetCoords,targetCoords+18,myCoords->getPointer());
1628       mesh->setCoords(myCoords);
1629       myCoords->decrRef();
1630       if(rank==0)
1631         {
1632           int targetConn[7]={0,3,4,1, 1,4,2};
1633           mesh->allocateCells(2);
1634           mesh->insertNextCell(INTERP_KERNEL::NORM_QUAD4,4,targetConn);
1635           mesh->insertNextCell(INTERP_KERNEL::NORM_TRI3,3,targetConn+4);
1636           mesh->finishInsertingCells();
1637         }
1638       else
1639         { 
1640           int targetConn[11]={4,5,2, 6,7,4,3, 7,8,5,4};
1641           mesh->allocateCells(3);
1642           mesh->insertNextCell(INTERP_KERNEL::NORM_TRI3,3,targetConn);
1643           mesh->insertNextCell(INTERP_KERNEL::NORM_QUAD4,4,targetConn+3);
1644           mesh->insertNextCell(INTERP_KERNEL::NORM_QUAD4,4,targetConn+7);
1645           mesh->finishInsertingCells();
1646         }
1647       ParaMEDMEM::ComponentTopology comptopo;
1648       paramesh=new ParaMESH(mesh,*source_group,"source mesh");
1649       parafield=new ParaFIELD(ON_CELLS,NO_TIME,paramesh, comptopo);
1650       parafield->getField()->setNature(ConservativeVolumic);
1651       double *vals=parafield->getField()->getArray()->getPointer();
1652       if(rank==0)
1653         { vals[0]=7.; vals[1]=8.; }
1654       else
1655         { vals[0]=9.; vals[1]=10.; vals[2]=11.; }
1656     }
1657   else
1658     {
1659       mesh=MEDCouplingUMesh::New("an example of -1 D mesh",-1);
1660       ParaMEDMEM::ComponentTopology comptopo;
1661       paramesh=new ParaMESH(mesh,*target_group,"target mesh");
1662       parafield=new ParaFIELD(ON_CELLS,NO_TIME,paramesh, comptopo);
1663       parafield->getField()->setNature(ConservativeVolumic);
1664     }
1665   ParaMEDMEM::InterpKernelDEC dec(*source_group,*target_group);
1666   if(source_group->containsMyRank())
1667     {
1668       dec.setMethod("P0");
1669       dec.attachLocalField(parafield);
1670       dec.synchronize();
1671       dec.setForcedRenormalization(false);
1672       dec.sendData();
1673       dec.recvData();
1674       const double *res=parafield->getField()->getArray()->getConstPointer();
1675       if(rank==0)
1676         {
1677           CPPUNIT_ASSERT_DOUBLES_EQUAL(9.125,res[0],1e-12);
1678           CPPUNIT_ASSERT_DOUBLES_EQUAL(9.125,res[1],1e-12);
1679         }
1680       else
1681         {
1682           CPPUNIT_ASSERT_DOUBLES_EQUAL(9.125,res[0],1e-12);
1683           CPPUNIT_ASSERT_DOUBLES_EQUAL(9.125,res[1],1e-12);
1684           CPPUNIT_ASSERT_DOUBLES_EQUAL(9.125,res[2],1e-12);
1685         }
1686     }
1687   else
1688     {
1689       dec.setMethod("P0");
1690       dec.attachLocalField(parafield);
1691       dec.synchronize();
1692       dec.setForcedRenormalization(false);
1693       dec.recvData();
1694       const double *res=parafield->getField()->getArray()->getConstPointer();
1695       CPPUNIT_ASSERT_DOUBLES_EQUAL(9.125,res[0],1e-12);
1696       dec.sendData();
1697     }
1698   ParaMEDMEM::InterpKernelDEC dec2(*source_group,*target_group);
1699   dec2.setMethod("P0");
1700   parafield->getField()->setNature(IntegralGlobConstraint);
1701   if(source_group->containsMyRank())
1702     {
1703       double *vals=parafield->getField()->getArray()->getPointer();
1704       if(rank==0)
1705         { vals[0]=7.; vals[1]=8.; }
1706       else
1707         { vals[0]=9.; vals[1]=10.; vals[2]=11.; }
1708       dec2.attachLocalField(parafield);
1709       dec2.synchronize();
1710       dec2.sendData();
1711       dec2.recvData();
1712       const double *res=parafield->getField()->getArray()->getConstPointer();
1713       if(rank==0)
1714         {
1715           CPPUNIT_ASSERT_DOUBLES_EQUAL(11.25,res[0],1e-12);
1716           CPPUNIT_ASSERT_DOUBLES_EQUAL(5.625,res[1],1e-12);
1717         }
1718       else
1719         {
1720           CPPUNIT_ASSERT_DOUBLES_EQUAL(5.625,res[0],1e-12);
1721           CPPUNIT_ASSERT_DOUBLES_EQUAL(11.25,res[1],1e-12);
1722           CPPUNIT_ASSERT_DOUBLES_EQUAL(11.25,res[2],1e-12);
1723         }
1724     }
1725   else
1726     {
1727       dec2.attachLocalField(parafield);
1728       dec2.synchronize();
1729       dec2.recvData();
1730       const double *res=parafield->getField()->getArray()->getConstPointer();
1731       CPPUNIT_ASSERT_DOUBLES_EQUAL(45.,res[0],1e-12);
1732       dec2.sendData();
1733     }
1734   //
1735   ParaMEDMEM::InterpKernelDEC dec3(*source_group,*target_group);
1736   dec3.setMethod("P0");
1737   parafield->getField()->setNature(Integral);
1738   if(source_group->containsMyRank())
1739     {
1740       double *vals=parafield->getField()->getArray()->getPointer();
1741       if(rank==0)
1742         { vals[0]=7.; vals[1]=8.; }
1743       else
1744         { vals[0]=9.; vals[1]=10.; vals[2]=11.; }
1745       dec3.attachLocalField(parafield);
1746       dec3.synchronize();
1747       dec3.sendData();
1748       dec3.recvData();
1749       const double *res=parafield->getField()->getArray()->getConstPointer();
1750       if(rank==0)
1751         {
1752           CPPUNIT_ASSERT_DOUBLES_EQUAL(11.25,res[0],1e-12);
1753           CPPUNIT_ASSERT_DOUBLES_EQUAL(5.625,res[1],1e-12);
1754         }
1755       else
1756         {
1757           CPPUNIT_ASSERT_DOUBLES_EQUAL(5.625,res[0],1e-12);
1758           CPPUNIT_ASSERT_DOUBLES_EQUAL(11.25,res[1],1e-12);
1759           CPPUNIT_ASSERT_DOUBLES_EQUAL(11.25,res[2],1e-12);
1760         }
1761     }
1762   else
1763     {
1764       dec3.attachLocalField(parafield);
1765       dec3.synchronize();
1766       dec3.recvData();
1767       const double *res=parafield->getField()->getArray()->getConstPointer();
1768       CPPUNIT_ASSERT_DOUBLES_EQUAL(45.,res[0],1e-12);
1769       dec3.sendData();
1770     }
1771   //
1772   ParaMEDMEM::InterpKernelDEC dec4(*source_group,*target_group);
1773   dec4.setMethod("P0");
1774   parafield->getField()->setNature(RevIntegral);
1775   if(source_group->containsMyRank())
1776     {
1777       double *vals=parafield->getField()->getArray()->getPointer();
1778       if(rank==0)
1779         { vals[0]=7.; vals[1]=8.; }
1780       else
1781         { vals[0]=9.; vals[1]=10.; vals[2]=11.; }
1782       dec4.attachLocalField(parafield);
1783       dec4.synchronize();
1784       dec4.sendData();
1785       dec4.recvData();
1786       const double *res=parafield->getField()->getArray()->getConstPointer();
1787        if(rank==0)
1788         {
1789           CPPUNIT_ASSERT_DOUBLES_EQUAL(9.125,res[0],1e-12);
1790           CPPUNIT_ASSERT_DOUBLES_EQUAL(9.125,res[1],1e-12);
1791         }
1792       else
1793         {
1794           CPPUNIT_ASSERT_DOUBLES_EQUAL(9.125,res[0],1e-12);
1795           CPPUNIT_ASSERT_DOUBLES_EQUAL(9.125,res[1],1e-12);
1796           CPPUNIT_ASSERT_DOUBLES_EQUAL(9.125,res[2],1e-12);
1797         }
1798     }
1799   else
1800     {
1801       dec4.attachLocalField(parafield);
1802       dec4.synchronize();
1803       dec4.recvData();
1804       const double *res=parafield->getField()->getArray()->getConstPointer();
1805       CPPUNIT_ASSERT_DOUBLES_EQUAL(9.125,res[0],1e-12);
1806       dec4.sendData();
1807     }
1808   delete parafield;
1809   delete paramesh;
1810   mesh->decrRef();
1811   delete target_group;
1812   delete source_group;
1813   //
1814   MPI_Barrier(MPI_COMM_WORLD);
1815 }
1816
1817 void ParaMEDMEMTest::testInterpKernelDECPartialProcs()
1818 {
1819   int size;
1820   int rank;
1821   MPI_Comm_size(MPI_COMM_WORLD,&size);
1822   MPI_Comm_rank(MPI_COMM_WORLD,&rank);
1823   //
1824   if(size!=3)
1825     return ;
1826   set<int> procs_source;
1827   set<int> procs_target;
1828   //
1829   procs_source.insert(0);
1830   procs_target.insert(1);
1831   //
1832   ParaMEDMEM::MEDCouplingUMesh *mesh=0;
1833   ParaMEDMEM::ParaMESH *paramesh=0;
1834   ParaMEDMEM::ParaFIELD *parafield=0;
1835   //
1836   ParaMEDMEM::CommInterface interface;
1837   //
1838   MPI_Barrier(MPI_COMM_WORLD);
1839   double targetCoords[8]={ 0.,0., 1., 0., 0., 1., 1., 1. };
1840   CommInterface comm;
1841   int grpIds[2]={0,1};
1842   MPI_Group grp,group_world;
1843   comm.commGroup(MPI_COMM_WORLD,&group_world);
1844   comm.groupIncl(group_world,2,grpIds,&grp);
1845   MPI_Comm partialComm;
1846   comm.commCreate(MPI_COMM_WORLD,grp,&partialComm);
1847   //
1848   ProcessorGroup* target_group=0;
1849   ProcessorGroup* source_group=0;
1850   //
1851   ParaMEDMEM::InterpKernelDEC *dec=0;
1852   if(rank==0 || rank==1)
1853     {
1854       target_group = new ParaMEDMEM::MPIProcessorGroup(interface,procs_target,partialComm);
1855       source_group = new ParaMEDMEM::MPIProcessorGroup(interface,procs_source,partialComm);
1856       if(source_group->containsMyRank())
1857         {    
1858           mesh=MEDCouplingUMesh::New();
1859           mesh->setMeshDimension(2);
1860           DataArrayDouble *myCoords=DataArrayDouble::New();
1861           myCoords->alloc(4,2);
1862           std::copy(targetCoords,targetCoords+8,myCoords->getPointer());
1863           mesh->setCoords(myCoords);
1864           myCoords->decrRef();
1865           int targetConn[4]={0,2,3,1};
1866           mesh->allocateCells(1);
1867           mesh->insertNextCell(INTERP_KERNEL::NORM_QUAD4,4,targetConn);
1868           mesh->finishInsertingCells();
1869           ParaMEDMEM::ComponentTopology comptopo;
1870           paramesh=new ParaMESH(mesh,*source_group,"source mesh");
1871           parafield=new ParaFIELD(ON_CELLS,NO_TIME,paramesh, comptopo);
1872           parafield->getField()->setNature(ConservativeVolumic);
1873           double *vals=parafield->getField()->getArray()->getPointer();
1874           vals[0]=7.;
1875           dec=new ParaMEDMEM::InterpKernelDEC(*source_group,*target_group);
1876           dec->attachLocalField(parafield);
1877           dec->synchronize();
1878           dec->sendData();
1879           dec->recvData();
1880         }
1881       else
1882         {
1883           mesh=MEDCouplingUMesh::New();
1884           mesh->setMeshDimension(2);
1885           DataArrayDouble *myCoords=DataArrayDouble::New();
1886           myCoords->alloc(4,2);
1887           std::copy(targetCoords,targetCoords+8,myCoords->getPointer());
1888           mesh->setCoords(myCoords);
1889           myCoords->decrRef();
1890           int targetConn[6]={0,2,1,2,3,1};
1891           mesh->allocateCells(2);
1892           mesh->insertNextCell(INTERP_KERNEL::NORM_TRI3,3,targetConn);
1893           mesh->insertNextCell(INTERP_KERNEL::NORM_TRI3,3,targetConn+3);
1894           mesh->finishInsertingCells();
1895           ParaMEDMEM::ComponentTopology comptopo;
1896           paramesh=new ParaMESH(mesh,*target_group,"target mesh");
1897           parafield=new ParaFIELD(ON_CELLS,NO_TIME,paramesh, comptopo);
1898           parafield->getField()->setNature(ConservativeVolumic);
1899           dec=new ParaMEDMEM::InterpKernelDEC(*source_group,*target_group);
1900           dec->attachLocalField(parafield);
1901           dec->synchronize();
1902           dec->recvData();
1903           dec->sendData();
1904         }
1905     }
1906   delete parafield;
1907   delete paramesh;
1908   if(mesh)
1909     mesh->decrRef();
1910   delete target_group;
1911   delete source_group;
1912   delete dec;
1913   if(partialComm != MPI_COMM_NULL)
1914     comm.commFree(&partialComm);
1915   comm.groupFree(&grp);
1916   comm.groupFree(&group_world);
1917   MPI_Barrier(MPI_COMM_WORLD);
1918 }
1919
1920 /*!
1921  * This test reproduces bug of Gauthier on 13/9/2010 concerning 3DSurf meshes.
1922  * 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.
1923  */
1924 void ParaMEDMEMTest::testInterpKernelDEC3DSurfEmptyBBox()
1925 {
1926   int size;
1927   int rank;
1928   MPI_Comm_size(MPI_COMM_WORLD,&size);
1929   MPI_Comm_rank(MPI_COMM_WORLD,&rank);
1930   //
1931   if(size!=3)
1932     return ;
1933   int nproc_source = 1;
1934   set<int> self_procs;
1935   set<int> procs_source;
1936   set<int> procs_target;
1937   
1938   for (int i=0; i<nproc_source; i++)
1939     procs_source.insert(i);
1940   for (int i=nproc_source; i<size; i++)
1941     procs_target.insert(i);
1942   self_procs.insert(rank);
1943   //
1944   ParaMEDMEM::MEDCouplingUMesh *mesh=0;
1945   ParaMEDMEM::ParaMESH *paramesh=0;
1946   ParaMEDMEM::ParaFIELD *parafieldP0=0;
1947   //
1948   ParaMEDMEM::CommInterface interface;
1949   //
1950   ProcessorGroup* self_group = new ParaMEDMEM::MPIProcessorGroup(interface,self_procs);
1951   ProcessorGroup* target_group = new ParaMEDMEM::MPIProcessorGroup(interface,procs_target);
1952   ProcessorGroup* source_group = new ParaMEDMEM::MPIProcessorGroup(interface,procs_source);
1953   //
1954   MPI_Barrier(MPI_COMM_WORLD);
1955   if(source_group->containsMyRank())
1956     {
1957       double coords[15]={1.,0.,0., 2.,0.,0., 2.,2.,0., 0.,2.,0., 0.5,0.5,1.};
1958       int conn[7]={0,1,2,3,0,3,4};
1959       mesh=MEDCouplingUMesh::New("Source mesh Proc0",2);
1960       mesh->allocateCells(2);
1961       mesh->insertNextCell(INTERP_KERNEL::NORM_QUAD4,4,conn);
1962       mesh->insertNextCell(INTERP_KERNEL::NORM_TRI3,3,conn+4);
1963       mesh->finishInsertingCells();
1964       DataArrayDouble *myCoords=DataArrayDouble::New();
1965       myCoords->alloc(5,3);
1966       std::copy(coords,coords+15,myCoords->getPointer());
1967       mesh->setCoords(myCoords);
1968       myCoords->decrRef();
1969       //
1970       paramesh=new ParaMESH(mesh,*source_group,"source mesh");
1971       ParaMEDMEM::ComponentTopology comptopo;
1972       parafieldP0 = new ParaFIELD(ON_CELLS,NO_TIME,paramesh, comptopo);
1973       double *valueP0=parafieldP0->getField()->getArray()->getPointer();
1974       parafieldP0->getField()->setNature(ConservativeVolumic);
1975       valueP0[0]=7.; valueP0[1]=8.;
1976     }
1977   else
1978     {
1979       const char targetMeshName[]="target mesh";
1980       if(rank==1)
1981         {
1982           double coords[12]={0.25,0.25,0.5, 0.,0.25,0.5, 0.,0.,0.5, 0.25,0.,0.5};
1983           int conn[4]={0,1,2,3};
1984           mesh=MEDCouplingUMesh::New("Target mesh Proc1",2);
1985           mesh->allocateCells(1);
1986           mesh->insertNextCell(INTERP_KERNEL::NORM_QUAD4,4,conn);
1987           mesh->finishInsertingCells();
1988           DataArrayDouble *myCoords=DataArrayDouble::New();
1989           myCoords->alloc(4,3);
1990           std::copy(coords,coords+12,myCoords->getPointer());
1991           mesh->setCoords(myCoords);
1992           myCoords->decrRef();
1993           paramesh=new ParaMESH(mesh,*target_group,targetMeshName);
1994         }
1995       if(rank==2)
1996         {
1997           double coords[12]={0.,0.25,0.5, 0.,0.,0.5, -1.,0.,0.5, -1.,0.25,0.5};
1998           int conn[4]={0,1,2,3};
1999           mesh=MEDCouplingUMesh::New("Target mesh Proc2",2);
2000           mesh->allocateCells(1);
2001           mesh->insertNextCell(INTERP_KERNEL::NORM_QUAD4,4,conn);
2002           mesh->finishInsertingCells();
2003           DataArrayDouble *myCoords=DataArrayDouble::New();
2004           myCoords->alloc(4,3);
2005           std::copy(coords,coords+12,myCoords->getPointer());
2006           mesh->setCoords(myCoords);
2007           myCoords->decrRef();
2008           paramesh=new ParaMESH(mesh,*target_group,targetMeshName);
2009         }
2010       ParaMEDMEM::ComponentTopology comptopo;
2011       parafieldP0 = new ParaFIELD(ON_CELLS,NO_TIME,paramesh, comptopo);
2012       parafieldP0->getField()->setNature(ConservativeVolumic);
2013     }
2014   // test 1
2015   ParaMEDMEM::InterpKernelDEC dec(*source_group,*target_group);
2016   if (source_group->containsMyRank())
2017     { 
2018       dec.setMethod("P0");
2019       dec.attachLocalField(parafieldP0);
2020       dec.synchronize();
2021       // dec.setForcedRenormalization(false);
2022       // dec.sendData();
2023       // dec.recvData();
2024       // const double *valueP0=parafieldP0->getField()->getArray()->getPointer();
2025       // if(rank==0)
2026       //   {
2027       //     CPPUNIT_ASSERT_DOUBLES_EQUAL(7.4,valueP0[0],1e-7);
2028       //     CPPUNIT_ASSERT_DOUBLES_EQUAL(9.0540540540540544,valueP0[1],1e-7);
2029       //   }
2030       // if(rank==1)
2031       //   {
2032       //     CPPUNIT_ASSERT_DOUBLES_EQUAL(8.64054054054054,valueP0[0],1e-7);
2033       //   }
2034       // if(rank==2)
2035       //   {
2036       //     CPPUNIT_ASSERT_DOUBLES_EQUAL(9.0540540540540544,valueP0[0],1e-7);
2037       //   }
2038     }
2039   else
2040     {
2041       dec.setMethod("P0");
2042       dec.attachLocalField(parafieldP0);
2043       dec.synchronize();
2044       // dec.setForcedRenormalization(false);
2045       // dec.recvData();
2046       // const double *res=parafieldP0->getField()->getArray()->getConstPointer();
2047       // if(rank==3)
2048       //   {
2049       //     CPPUNIT_ASSERT_EQUAL(1,parafieldP0->getField()->getNumberOfTuples());
2050       //     CPPUNIT_ASSERT_EQUAL(1,parafieldP0->getField()->getNumberOfComponents());
2051       //     CPPUNIT_ASSERT_DOUBLES_EQUAL(7.4,res[0],1e-12);
2052       //   }
2053       // if(rank==4)
2054       //   {
2055       //     CPPUNIT_ASSERT_EQUAL(1,parafieldP0->getField()->getNumberOfTuples());
2056       //     CPPUNIT_ASSERT_EQUAL(1,parafieldP0->getField()->getNumberOfComponents());
2057       //     CPPUNIT_ASSERT_DOUBLES_EQUAL(9.0540540540540526,res[0],1e-12);
2058       //   }
2059       // dec.sendData();
2060     }
2061   //
2062   delete parafieldP0;
2063   mesh->decrRef();
2064   delete paramesh;
2065   delete self_group;
2066   delete target_group;
2067   delete source_group;
2068   //
2069   MPI_Barrier(MPI_COMM_WORLD);
2070 }
2071
2072 /*!
2073  * Tests an asynchronous exchange between two codes
2074  * one sends data with dtA as an interval, the max time being tmaxA
2075  * the other one receives with dtB as an interval, the max time being tmaxB
2076  */
2077 void ParaMEDMEMTest::testAsynchronousInterpKernelDEC_2D(double dtA, double tmaxA, 
2078                                                         double dtB, double tmaxB, bool WithPointToPoint, bool Asynchronous,
2079                                                         bool WithInterp, const char *srcMeth, const char *targetMeth)
2080 {
2081   std::string srcM(srcMeth);
2082   std::string targetM(targetMeth);
2083   int size;
2084   int rank;
2085   MPI_Comm_size(MPI_COMM_WORLD,&size);
2086   MPI_Comm_rank(MPI_COMM_WORLD,&rank);
2087  
2088   //the test is meant to run on five processors
2089   if (size !=5) return ;
2090    
2091   int nproc_source = 3;
2092   set<int> self_procs;
2093   set<int> procs_source;
2094   set<int> procs_target;
2095   
2096   for (int i=0; i<nproc_source; i++)
2097     procs_source.insert(i);
2098   for (int i=nproc_source; i<size; i++)
2099     procs_target.insert(i);
2100   self_procs.insert(rank);
2101   
2102   ParaMEDMEM::CommInterface interface;
2103     
2104   ParaMEDMEM::ProcessorGroup* self_group = new ParaMEDMEM::MPIProcessorGroup(interface,self_procs);
2105   ParaMEDMEM::ProcessorGroup* target_group = new ParaMEDMEM::MPIProcessorGroup(interface,procs_target);
2106   ParaMEDMEM::ProcessorGroup* source_group = new ParaMEDMEM::MPIProcessorGroup(interface,procs_source);
2107     
2108   //loading the geometry for the source group
2109
2110   ParaMEDMEM::InterpKernelDEC dec (*source_group,*target_group);
2111   
2112   ParaMEDMEM::MEDCouplingUMesh* mesh;
2113   ParaMEDMEM::ParaMESH* paramesh;
2114   ParaMEDMEM::ParaFIELD* parafield;
2115   
2116   ICoCo::MEDField* icocofield ;
2117
2118   char * tmp_dir_c                    = getenv("TMP");
2119   string tmp_dir;
2120   if (tmp_dir_c != NULL)
2121     tmp_dir = string(tmp_dir_c);
2122   else
2123     tmp_dir = "/tmp";
2124   string filename_xml1              = getResourceFile("square1_split");
2125   string filename_xml2              = getResourceFile("square2_split"); 
2126   //string filename_seq_wr            = makeTmpFile("");
2127   //string filename_seq_med           = makeTmpFile("myWrField_seq_pointe221.med");
2128   
2129   // To remove tmp files from disk
2130   ParaMEDMEMTest_TmpFilesRemover aRemover;
2131   
2132   MPI_Barrier(MPI_COMM_WORLD);
2133
2134   if (source_group->containsMyRank())
2135     {
2136       string master = filename_xml1;
2137       
2138       ostringstream strstream;
2139       strstream <<master<<rank+1<<".med";
2140       ostringstream meshname ;
2141       meshname<< "Mesh_2_"<< rank+1;
2142       
2143       mesh=MEDLoader::ReadUMeshFromFile(strstream.str().c_str(),meshname.str().c_str(),0);
2144
2145       paramesh=new ParaMESH (mesh,*source_group,"source mesh");
2146     
2147       //      ParaMEDMEM::ParaSUPPORT* parasupport=new UnstructuredParaSUPPORT( support,*source_group);
2148       ParaMEDMEM::ComponentTopology comptopo;
2149       if(srcM=="P0")
2150         {
2151           parafield = new ParaFIELD(ON_CELLS,NO_TIME,paramesh, comptopo);
2152           parafield->getField()->setNature(ConservativeVolumic);//InvertIntegral);//ConservativeVolumic);
2153         }
2154       else
2155         parafield = new ParaFIELD(ON_NODES,NO_TIME,paramesh, comptopo);
2156
2157       int nb_local;
2158       if(srcM=="P0")
2159         nb_local=mesh->getNumberOfCells();
2160       else
2161         nb_local=mesh->getNumberOfNodes();
2162       //      double * value= new double[nb_local];
2163       double *value=parafield->getField()->getArray()->getPointer();
2164       for(int ielem=0; ielem<nb_local;ielem++)
2165         value[ielem]=0.0;
2166     
2167       //      ICoCo::Field* icocofield=new ICoCo::MEDField(paramesh,parafield);
2168       icocofield=new ICoCo::MEDField(parafield->getField());
2169      
2170       dec.attachLocalField(icocofield);
2171
2172
2173     }
2174   
2175   //loading the geometry for the target group
2176   if (target_group->containsMyRank())
2177     {
2178       string master= filename_xml2;
2179       ostringstream strstream;
2180       strstream << master<<(rank-nproc_source+1)<<".med";
2181       ostringstream meshname ;
2182       meshname<< "Mesh_3_"<<rank-nproc_source+1;
2183       
2184       mesh = MEDLoader::ReadUMeshFromFile(strstream.str().c_str(),meshname.str().c_str(),0);
2185
2186       paramesh=new ParaMESH (mesh,*target_group,"target mesh");
2187       //      ParaMEDMEM::ParaSUPPORT* parasupport=new UnstructuredParaSUPPORT(support,*target_group);
2188       ParaMEDMEM::ComponentTopology comptopo;
2189       if(targetM=="P0")
2190         {
2191           parafield = new ParaFIELD(ON_CELLS,NO_TIME,paramesh, comptopo);
2192           parafield->getField()->setNature(ConservativeVolumic);//InvertIntegral);//ConservativeVolumic);
2193         }
2194       else
2195         parafield = new ParaFIELD(ON_NODES,NO_TIME,paramesh, comptopo);
2196       
2197       int nb_local;
2198       if(targetM=="P0")
2199         nb_local=mesh->getNumberOfCells();
2200       else
2201         nb_local=mesh->getNumberOfNodes();
2202                         
2203       double *value=parafield->getField()->getArray()->getPointer();
2204       for(int ielem=0; ielem<nb_local;ielem++)
2205         value[ielem]=0.0;
2206       //      ICoCo::Field* icocofield=new ICoCo::MEDField(paramesh,parafield);
2207       icocofield=new ICoCo::MEDField(parafield->getField());
2208       
2209       dec.attachLocalField(icocofield);
2210     }
2211     
2212   
2213   //attaching a DEC to the source group 
2214   
2215   if (source_group->containsMyRank())
2216     { 
2217       cout<<"DEC usage"<<endl;
2218       dec.setAsynchronous(Asynchronous);
2219       if ( WithInterp ) {
2220         dec.setTimeInterpolationMethod(LinearTimeInterp);
2221       }
2222       if ( WithPointToPoint ) {
2223         dec.setAllToAllMethod(PointToPoint);
2224       }
2225       else {
2226         dec.setAllToAllMethod(Native);
2227       }
2228       dec.synchronize();
2229       dec.setForcedRenormalization(false);
2230       for (double time=0; time<tmaxA+1e-10; time+=dtA)
2231         {
2232           cout << "testAsynchronousInterpKernelDEC_2D" << rank << " time " << time
2233                << " dtA " << dtA << " tmaxA " << tmaxA << endl ;
2234           if ( time+dtA < tmaxA+1e-7 ) {
2235             dec.sendData( time , dtA );
2236           }
2237           else {
2238             dec.sendData( time , 0 );
2239           }
2240           double* value = parafield->getField()->getArray()->getPointer();
2241           int nb_local=parafield->getField()->getMesh()->getNumberOfCells();
2242           for (int i=0; i<nb_local;i++)
2243             value[i]= time+dtA;
2244
2245        
2246         }
2247     }
2248   
2249   //attaching a DEC to the target group
2250   if (target_group->containsMyRank())
2251     {
2252       cout<<"DEC usage"<<endl;
2253       dec.setAsynchronous(Asynchronous);
2254       if ( WithInterp ) {
2255         dec.setTimeInterpolationMethod(LinearTimeInterp);
2256       }
2257       if ( WithPointToPoint ) {
2258         dec.setAllToAllMethod(PointToPoint);
2259       }
2260       else {
2261         dec.setAllToAllMethod(Native);
2262       }
2263       dec.synchronize();
2264       dec.setForcedRenormalization(false);
2265       vector<double> times;
2266       for (double time=0; time<tmaxB+1e-10; time+=dtB)
2267         {
2268           cout << "testAsynchronousInterpKernelDEC_2D" << rank << " time " << time
2269                << " dtB " << dtB << " tmaxB " << tmaxB << endl ;
2270           dec.recvData( time );
2271           double vi = parafield->getVolumeIntegral(0,true);
2272           cout << "testAsynchronousInterpKernelDEC_2D" << rank << " time " << time
2273                << " VolumeIntegral " << vi
2274                << " time*10000 " << time*10000 << endl ;
2275           
2276           CPPUNIT_ASSERT_DOUBLES_EQUAL(vi,time*10000,0.001);
2277         }
2278       
2279     }
2280   
2281   delete source_group;
2282   delete target_group;
2283   delete self_group;
2284   delete parafield ;
2285   delete paramesh ;
2286   mesh->decrRef() ;
2287   delete icocofield ;
2288
2289   cout << "testAsynchronousInterpKernelDEC_2D" << rank << " MPI_Barrier " << endl ;
2290  
2291   if (Asynchronous) MPI_Barrier(MPI_COMM_WORLD);
2292   cout << "end of InterpKernelDEC_2D test"<<endl;
2293 }