1 // Copyright (C) 2007-2015 CEA/DEN, EDF R&D
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.
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.
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
17 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
20 #include "ParaMEDMEMTest.hxx"
21 #include <cppunit/TestAssert.h>
23 #include "CommInterface.hxx"
24 #include "ProcessorGroup.hxx"
25 #include "MPIProcessorGroup.hxx"
26 #include "Topology.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"
41 // use this define to enable lines, execution of which leads to Segmentation Fault
44 // use this define to enable CPPUNIT asserts and fails, showing bugs
45 #define ENABLE_FORCED_FAILURES
49 using namespace ParaMEDMEM;
51 void ParaMEDMEMTest::testInterpKernelDEC_2D()
53 testInterpKernelDEC_2D_("P0","P0");
56 void ParaMEDMEMTest::testInterpKernelDEC2_2D()
58 testInterpKernelDEC2_2D_("P0","P0");
61 void ParaMEDMEMTest::testInterpKernelDEC_3D()
63 testInterpKernelDEC_3D_("P0","P0");
66 void ParaMEDMEMTest::testInterpKernelDEC_2DP0P1()
68 //testInterpKernelDEC_2D_("P0","P1");
71 void ParaMEDMEMTest::testInterpKernelDEC_1D()
75 MPI_Comm_size(MPI_COMM_WORLD,&size);
76 MPI_Comm_rank(MPI_COMM_WORLD,&rank);
82 set<int> procs_source;
83 set<int> procs_target;
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);
91 ParaMEDMEM::MEDCouplingUMesh *mesh=0;
92 ParaMEDMEM::ParaMESH *paramesh=0;
93 ParaMEDMEM::ParaFIELD *parafieldP0=0;
95 ParaMEDMEM::CommInterface interface;
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);
101 MPI_Barrier(MPI_COMM_WORLD);
102 if(source_group->containsMyRank())
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);
121 double coords[2]={0.7,0.9};
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);
135 double coords[2]={1.,1.12};
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);
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);
154 valueP0[0]=7.; valueP0[1]=8.;
167 const char targetMeshName[]="target mesh";
170 double coords[2]={0.5,0.75};
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);
181 paramesh=new ParaMESH(mesh,*target_group,targetMeshName);
185 double coords[2]={0.75,1.2};
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);
196 paramesh=new ParaMESH(mesh,*target_group,targetMeshName);
198 ParaMEDMEM::ComponentTopology comptopo;
199 parafieldP0 = new ParaFIELD(ON_CELLS,NO_TIME,paramesh, comptopo);
200 parafieldP0->getField()->setNature(ConservativeVolumic);
203 ParaMEDMEM::InterpKernelDEC dec(*source_group,*target_group);
204 if (source_group->containsMyRank())
207 dec.attachLocalField(parafieldP0);
209 dec.setForcedRenormalization(false);
212 const double *valueP0=parafieldP0->getField()->getArray()->getPointer();
215 CPPUNIT_ASSERT_DOUBLES_EQUAL(7.4,valueP0[0],1e-7);
216 CPPUNIT_ASSERT_DOUBLES_EQUAL(9.0540540540540544,valueP0[1],1e-7);
220 CPPUNIT_ASSERT_DOUBLES_EQUAL(8.64054054054054,valueP0[0],1e-7);
224 CPPUNIT_ASSERT_DOUBLES_EQUAL(9.0540540540540544,valueP0[0],1e-7);
230 dec.attachLocalField(parafieldP0);
232 dec.setForcedRenormalization(false);
234 const double *res=parafieldP0->getField()->getArray()->getConstPointer();
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);
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);
257 MPI_Barrier(MPI_COMM_WORLD);
260 void ParaMEDMEMTest::testInterpKernelDEC_2DCurve()
264 MPI_Comm_size(MPI_COMM_WORLD,&size);
265 MPI_Comm_rank(MPI_COMM_WORLD,&rank);
269 int nproc_source = 3;
271 set<int> procs_source;
272 set<int> procs_target;
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);
280 ParaMEDMEM::MEDCouplingUMesh *mesh=0;
281 ParaMEDMEM::ParaMESH *paramesh=0;
282 ParaMEDMEM::ParaFIELD *parafieldP0=0;
284 ParaMEDMEM::CommInterface interface;
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);
290 MPI_Barrier(MPI_COMM_WORLD);
291 if(source_group->containsMyRank())
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);
310 double coords[4]={0.7,0.7,0.9,0.9};
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);
324 double coords[4]={1.,1.,1.12,1.12};
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);
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);
343 valueP0[0]=7.; valueP0[1]=8.;
356 const char targetMeshName[]="target mesh";
359 double coords[4]={0.5,0.5,0.75,0.75};
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);
370 paramesh=new ParaMESH(mesh,*target_group,targetMeshName);
374 double coords[4]={0.75,0.75,1.2,1.2};
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);
385 paramesh=new ParaMESH(mesh,*target_group,targetMeshName);
387 ParaMEDMEM::ComponentTopology comptopo;
388 parafieldP0 = new ParaFIELD(ON_CELLS,NO_TIME,paramesh, comptopo);
389 parafieldP0->getField()->setNature(ConservativeVolumic);
392 ParaMEDMEM::InterpKernelDEC dec(*source_group,*target_group);
393 if (source_group->containsMyRank())
396 dec.attachLocalField(parafieldP0);
398 dec.setForcedRenormalization(false);
401 const double *valueP0=parafieldP0->getField()->getArray()->getPointer();
404 CPPUNIT_ASSERT_DOUBLES_EQUAL(7.4,valueP0[0],1e-7);
405 CPPUNIT_ASSERT_DOUBLES_EQUAL(9.0540540540540544,valueP0[1],1e-7);
409 CPPUNIT_ASSERT_DOUBLES_EQUAL(8.64054054054054,valueP0[0],1e-7);
413 CPPUNIT_ASSERT_DOUBLES_EQUAL(9.0540540540540544,valueP0[0],1e-7);
419 dec.attachLocalField(parafieldP0);
421 dec.setForcedRenormalization(false);
423 const double *res=parafieldP0->getField()->getArray()->getConstPointer();
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);
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);
446 MPI_Barrier(MPI_COMM_WORLD);
451 * Check methods defined in InterpKernelDEC.hxx
454 InterpKernelDEC(ProcessorGroup& local_group, ProcessorGroup& distant_group);
455 virtual ~InterpKernelDEC();
461 void ParaMEDMEMTest::testInterpKernelDEC_2D_(const char *srcMeth, const char *targetMeth)
463 std::string srcM(srcMeth);
464 std::string targetM(targetMeth);
467 MPI_Comm_size(MPI_COMM_WORLD,&size);
468 MPI_Comm_rank(MPI_COMM_WORLD,&rank);
470 //the test is meant to run on five processors
471 if (size !=5) return ;
473 int nproc_source = 3;
475 set<int> procs_source;
476 set<int> procs_target;
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);
484 ParaMEDMEM::CommInterface interface;
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);
490 //loading the geometry for the source group
492 ParaMEDMEM::InterpKernelDEC dec (*source_group,*target_group);
494 ParaMEDMEM::MEDCouplingUMesh* mesh;
495 ParaMEDMEM::ParaMESH* paramesh;
496 ParaMEDMEM::ParaFIELD* parafield;
497 ICoCo::MEDField* icocofield ;
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");
504 // To remove tmp files from disk
505 ParaMEDMEMTest_TmpFilesRemover aRemover;
507 MPI_Barrier(MPI_COMM_WORLD);
508 if (source_group->containsMyRank())
510 string master = filename_xml1;
512 ostringstream strstream;
513 strstream <<master<<rank+1<<".med";
514 ostringstream meshname ;
515 meshname<< "Mesh_2_"<< rank+1;
517 mesh=MEDLoader::ReadUMeshFromFile(strstream.str().c_str(),meshname.str().c_str(),0);
520 paramesh=new ParaMESH (mesh,*source_group,"source mesh");
522 // ParaMEDMEM::ParaSUPPORT* parasupport=new UnstructuredParaSUPPORT( support,*source_group);
523 ParaMEDMEM::ComponentTopology comptopo;
526 parafield = new ParaFIELD(ON_CELLS,NO_TIME,paramesh, comptopo);
527 parafield->getField()->setNature(ConservativeVolumic);
530 parafield = new ParaFIELD(ON_NODES,NO_TIME,paramesh, comptopo);
533 nb_local=mesh->getNumberOfCells();
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++)
541 // ICoCo::Field* icocofield=new ICoCo::MEDField(paramesh,parafield);
542 icocofield=new ICoCo::MEDField(parafield->getField());
543 dec.setMethod(srcMeth);
544 dec.attachLocalField(icocofield);
547 //loading the geometry for the target group
548 if (target_group->containsMyRank())
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);
557 paramesh=new ParaMESH (mesh,*target_group,"target mesh");
558 // ParaMEDMEM::ParaSUPPORT* parasupport=new UnstructuredParaSUPPORT(support,*target_group);
559 ParaMEDMEM::ComponentTopology comptopo;
562 parafield = new ParaFIELD(ON_CELLS,NO_TIME,paramesh, comptopo);
563 parafield->getField()->setNature(ConservativeVolumic);
566 parafield = new ParaFIELD(ON_NODES,NO_TIME,paramesh, comptopo);
569 nb_local=mesh->getNumberOfCells();
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++)
576 // ICoCo::Field* icocofield=new ICoCo::MEDField(paramesh,parafield);
577 icocofield=new ICoCo::MEDField(parafield->getField());
578 dec.setMethod(targetMeth);
579 dec.attachLocalField(icocofield);
583 //attaching a DEC to the source group
584 double field_before_int;
585 double field_after_int;
587 if (source_group->containsMyRank())
589 field_before_int = parafield->getVolumeIntegral(0,true);
591 cout<<"DEC usage"<<endl;
592 dec.setForcedRenormalization(false);
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());
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());
611 filename<<"./sourcesquare_"<<source_group->myRank()+1;
612 aRemover.Register(filename.str().c_str());
613 field_after_int = parafield->getVolumeIntegral(0,true);
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);
619 CPPUNIT_ASSERT_DOUBLES_EQUAL(field_before_int, field_after_int, 1e-6);
623 //attaching a DEC to the target group
624 if (target_group->containsMyRank())
627 dec.setForcedRenormalization(false);
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());
638 ParaMEDLoader::WriteParaMesh("./targetsquare",paramesh);
639 //MEDLoader::WriteField("./targetsquare",parafield->getField());
641 if (target_group->myRank()==0)
642 aRemover.Register("./targetsquareb");
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);
650 // CPPUNIT_ASSERT_DOUBLES_EQUAL(field_before_int, field_after_int, 1e-6);
663 MPI_Barrier(MPI_COMM_WORLD);
664 cout << "end of InterpKernelDEC_2D test"<<endl;
667 void ParaMEDMEMTest::testInterpKernelDEC2_2D_(const char *srcMeth, const char *targetMeth)
669 std::string srcM(srcMeth);
670 std::string targetM(targetMeth);
673 MPI_Comm_size(MPI_COMM_WORLD,&size);
674 MPI_Comm_rank(MPI_COMM_WORLD,&rank);
676 //the test is meant to run on five processors
677 if (size !=5) return ;
679 int nproc_source = 3;
681 set<int> procs_source;
682 set<int> procs_target;
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);
690 ParaMEDMEM::CommInterface interface;
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);
696 //loading the geometry for the source group
698 ParaMEDMEM::InterpKernelDEC dec (*source_group,*target_group);
700 ParaMEDMEM::MEDCouplingUMesh* mesh;
701 ParaMEDMEM::MEDCouplingFieldDouble* mcfield;
703 string filename_xml1 = getResourceFile("square1_split");
704 string filename_xml2 = getResourceFile("square2_split");
706 // To remove tmp files from disk
707 ParaMEDMEMTest_TmpFilesRemover aRemover;
709 MPI_Barrier(MPI_COMM_WORLD);
710 if (source_group->containsMyRank())
712 string master = filename_xml1;
714 ostringstream strstream;
715 strstream <<master<<rank+1<<".med";
716 ostringstream meshname ;
717 meshname<< "Mesh_2_"<< rank+1;
719 mesh=MEDLoader::ReadUMeshFromFile(strstream.str().c_str(),meshname.str().c_str(),0);
720 ParaMEDMEM::ComponentTopology comptopo;
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);
729 mcfield->setNature(ConservativeVolumic);
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);
742 nb_local=mesh->getNumberOfCells();
744 nb_local=mesh->getNumberOfNodes();
745 double *value=mcfield->getArray()->getPointer();
746 for(int ielem=0; ielem<nb_local;ielem++)
748 dec.setMethod(srcMeth);
749 dec.attachLocalField(mcfield);
750 dec.attachLocalField(mcfield);
753 //loading the geometry for the target group
754 if (target_group->containsMyRank())
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;
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);
771 mcfield->setNature(ConservativeVolumic);
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);
784 nb_local=mesh->getNumberOfCells();
786 nb_local=mesh->getNumberOfNodes();
787 double *value=mcfield->getArray()->getPointer();
788 for(int ielem=0; ielem<nb_local;ielem++)
790 dec.setMethod(targetMeth);
791 dec.attachLocalField(mcfield);
792 dec.attachLocalField(mcfield);
796 //attaching a DEC to the source group
798 if (source_group->containsMyRank())
801 dec.setForcedRenormalization(false);
806 //attaching a DEC to the target group
807 if (target_group->containsMyRank())
810 dec.setForcedRenormalization(false);
820 MPI_Barrier(MPI_COMM_WORLD);
821 cout << "end of InterpKernelDEC2_2D test"<<endl;
824 void ParaMEDMEMTest::testInterpKernelDEC_3D_(const char *srcMeth, const char *targetMeth)
826 std::string srcM(srcMeth);
827 std::string targetM(targetMeth);
830 MPI_Comm_size(MPI_COMM_WORLD,&size);
831 MPI_Comm_rank(MPI_COMM_WORLD,&rank);
833 //the test is meant to run on five processors
834 if (size !=3) return ;
836 int nproc_source = 2;
838 set<int> procs_source;
839 set<int> procs_target;
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);
847 ParaMEDMEM::CommInterface interface;
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);
853 //loading the geometry for the source group
855 ParaMEDMEM::InterpKernelDEC dec (*source_group,*target_group);
857 ParaMEDMEM::MEDCouplingUMesh* mesh;
858 ParaMEDMEM::ParaMESH* paramesh;
859 ParaMEDMEM::ParaFIELD* parafield;
860 ICoCo::MEDField* icocofield ;
862 char * tmp_dir_c = getenv("TMP");
864 if (tmp_dir_c != NULL)
865 tmp_dir = string(tmp_dir_c);
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");
873 // To remove tmp files from disk
874 ParaMEDMEMTest_TmpFilesRemover aRemover;
876 MPI_Barrier(MPI_COMM_WORLD);
877 if (source_group->containsMyRank())
879 string master = filename_xml1;
881 ostringstream strstream;
882 strstream <<master<<rank+1<<".med";
883 ostringstream meshname ;
884 meshname<< "Mesh_3_"<< rank+1;
886 mesh=MEDLoader::ReadUMeshFromFile(strstream.str().c_str(),meshname.str().c_str(),0);
889 paramesh=new ParaMESH (mesh,*source_group,"source mesh");
891 // ParaMEDMEM::ParaSUPPORT* parasupport=new UnstructuredParaSUPPORT( support,*source_group);
892 ParaMEDMEM::ComponentTopology comptopo;
895 parafield = new ParaFIELD(ON_CELLS,NO_TIME,paramesh, comptopo);
896 parafield->getField()->setNature(ConservativeVolumic);
899 parafield = new ParaFIELD(ON_NODES,NO_TIME,paramesh, comptopo);
902 nb_local=mesh->getNumberOfCells();
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++)
910 // ICoCo::Field* icocofield=new ICoCo::MEDField(paramesh,parafield);
911 icocofield=new ICoCo::MEDField(parafield->getField());
912 dec.setMethod(srcMeth);
913 dec.attachLocalField(icocofield);
916 //loading the geometry for the target group
917 if (target_group->containsMyRank())
919 string master= filename_xml2;
920 ostringstream strstream;
921 strstream << master << ".med";
922 ostringstream meshname ;
924 mesh = MEDLoader::ReadUMeshFromFile(strstream.str().c_str(),meshname.str().c_str(),0);
926 paramesh=new ParaMESH (mesh,*target_group,"target mesh");
927 // ParaMEDMEM::ParaSUPPORT* parasupport=new UnstructuredParaSUPPORT(support,*target_group);
928 ParaMEDMEM::ComponentTopology comptopo;
931 parafield = new ParaFIELD(ON_CELLS,NO_TIME,paramesh, comptopo);
932 parafield->getField()->setNature(ConservativeVolumic);
935 parafield = new ParaFIELD(ON_NODES,NO_TIME,paramesh, comptopo);
938 nb_local=mesh->getNumberOfCells();
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++)
945 // ICoCo::Field* icocofield=new ICoCo::MEDField(paramesh,parafield);
946 icocofield=new ICoCo::MEDField(parafield->getField());
947 dec.setMethod(targetMeth);
948 dec.attachLocalField(icocofield);
950 //attaching a DEC to the source group
951 double field_before_int;
952 double field_after_int;
954 if (source_group->containsMyRank())
956 field_before_int = parafield->getVolumeIntegral(0,true);
958 cout<<"DEC usage"<<endl;
959 dec.setForcedRenormalization(false);
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());
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());
978 filename<<"./sourcesquare_"<<source_group->myRank()+1;
979 aRemover.Register(filename.str().c_str());
980 field_after_int = parafield->getVolumeIntegral(0,true);
982 CPPUNIT_ASSERT_DOUBLES_EQUAL(field_before_int, field_after_int, 1e-6);
986 //attaching a DEC to the target group
987 if (target_group->containsMyRank())
990 dec.setForcedRenormalization(false);
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());
1001 ParaMEDLoader::WriteParaMesh("./targetsquare",paramesh);
1002 //MEDLoader::WriteField("./targetsquare",parafield->getField());
1004 if (target_group->myRank()==0)
1005 aRemover.Register("./targetsquareb");
1007 filename<<"./targetsquareb_"<<target_group->myRank()+1;
1008 aRemover.Register(filename.str().c_str());
1010 delete source_group;
1011 delete target_group;
1019 MPI_Barrier(MPI_COMM_WORLD);
1020 cout << "end of InterpKernelDEC_3D test"<<endl;
1023 //Synchronous tests without interpolation with native mode (AllToAll(v) from lam/MPI:
1024 void ParaMEDMEMTest::testSynchronousEqualInterpKernelWithoutInterpNativeDEC_2D()
1026 testAsynchronousInterpKernelDEC_2D(0.1,1,0.1,1,false,false,false,"P0","P0");
1029 //Synchronous tests without interpolation :
1030 void ParaMEDMEMTest::testSynchronousEqualInterpKernelWithoutInterpDEC_2D()
1032 testAsynchronousInterpKernelDEC_2D(0.1,1,0.1,1,true,false,false,"P0","P0");
1035 //Synchronous tests with interpolation :
1036 void ParaMEDMEMTest::testSynchronousEqualInterpKernelDEC_2D()
1038 testAsynchronousInterpKernelDEC_2D(0.1,1,0.1,1,true,false,true,"P0","P0");
1040 void ParaMEDMEMTest::testSynchronousFasterSourceInterpKernelDEC_2D()
1042 testAsynchronousInterpKernelDEC_2D(0.09,1,0.1,1,true,false,true,"P0","P0");
1044 void ParaMEDMEMTest::testSynchronousSlowerSourceInterpKernelDEC_2D()
1046 testAsynchronousInterpKernelDEC_2D(0.11,1,0.1,1,true,false,true,"P0","P0");
1048 void ParaMEDMEMTest::testSynchronousSlowSourceInterpKernelDEC_2D()
1050 testAsynchronousInterpKernelDEC_2D(0.11,1,0.01,1,true,false,true,"P0","P0");
1052 void ParaMEDMEMTest::testSynchronousFastSourceInterpKernelDEC_2D()
1054 testAsynchronousInterpKernelDEC_2D(0.01,1,0.11,1,true,false,true,"P0","P0");
1057 //Asynchronous tests with interpolation :
1058 void ParaMEDMEMTest::testAsynchronousEqualInterpKernelDEC_2D()
1060 testAsynchronousInterpKernelDEC_2D(0.1,1,0.1,1,true,true,true,"P0","P0");
1062 void ParaMEDMEMTest::testAsynchronousFasterSourceInterpKernelDEC_2D()
1064 testAsynchronousInterpKernelDEC_2D(0.09,1,0.1,1,true,true,true,"P0","P0");
1066 void ParaMEDMEMTest::testAsynchronousSlowerSourceInterpKernelDEC_2D()
1068 testAsynchronousInterpKernelDEC_2D(0.11,1,0.1,1,true,true,true,"P0","P0");
1070 void ParaMEDMEMTest::testAsynchronousSlowSourceInterpKernelDEC_2D()
1072 testAsynchronousInterpKernelDEC_2D(0.11,1,0.01,1,true,true,true,"P0","P0");
1074 void ParaMEDMEMTest::testAsynchronousFastSourceInterpKernelDEC_2D()
1076 testAsynchronousInterpKernelDEC_2D(0.01,1,0.11,1,true,true,true,"P0","P0");
1079 void ParaMEDMEMTest::testInterpKernelDECNonOverlapp_2D_P0P0()
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}};
1095 MPI_Comm_size(MPI_COMM_WORLD,&size);
1096 MPI_Comm_rank(MPI_COMM_WORLD,&rank);
1100 int nproc_source = 2;
1101 set<int> self_procs;
1102 set<int> procs_source;
1103 set<int> procs_target;
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);
1111 ParaMEDMEM::MEDCouplingUMesh *mesh=0;
1112 ParaMEDMEM::ParaMESH *paramesh=0;
1113 ParaMEDMEM::ParaFIELD* parafield=0;
1115 ParaMEDMEM::CommInterface interface;
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);
1121 MPI_Barrier(MPI_COMM_WORLD);
1122 if(source_group->containsMyRank())
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);
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);
1159 //test 1 - Conservative volumic
1160 ParaMEDMEM::InterpKernelDEC dec(*source_group,*target_group);
1161 parafield->getField()->setNature(ConservativeVolumic);
1162 if (source_group->containsMyRank())
1164 dec.setMethod("P0");
1165 dec.attachLocalField(parafield);
1167 dec.setForcedRenormalization(false);
1172 dec.setMethod("P0");
1173 dec.attachLocalField(parafield);
1175 dec.setForcedRenormalization(false);
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);
1183 ParaMEDMEM::InterpKernelDEC dec2(*source_group,*target_group);
1184 parafield->getField()->setNature(Integral);
1185 if (source_group->containsMyRank())
1187 dec2.setMethod("P0");
1188 dec2.attachLocalField(parafield);
1190 dec2.setForcedRenormalization(false);
1195 dec2.setMethod("P0");
1196 dec2.attachLocalField(parafield);
1198 dec2.setForcedRenormalization(false);
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);
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())
1210 dec3.setMethod("P0");
1211 dec3.attachLocalField(parafield);
1213 dec3.setForcedRenormalization(false);
1218 dec3.setMethod("P0");
1219 dec3.attachLocalField(parafield);
1221 dec3.setForcedRenormalization(false);
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);
1228 //test 4 - RevIntegral
1229 ParaMEDMEM::InterpKernelDEC dec4(*source_group,*target_group);
1230 parafield->getField()->setNature(RevIntegral);
1231 if (source_group->containsMyRank())
1233 dec4.setMethod("P0");
1234 dec4.attachLocalField(parafield);
1236 dec4.setForcedRenormalization(false);
1241 dec4.setMethod("P0");
1242 dec4.attachLocalField(parafield);
1244 dec4.setForcedRenormalization(false);
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);
1251 //test 5 - Conservative volumic reversed
1252 ParaMEDMEM::InterpKernelDEC dec5(*source_group,*target_group);
1253 parafield->getField()->setNature(ConservativeVolumic);
1254 if (source_group->containsMyRank())
1256 dec5.setMethod("P0");
1257 dec5.attachLocalField(parafield);
1259 dec5.setForcedRenormalization(false);
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);
1268 dec5.setMethod("P0");
1269 dec5.attachLocalField(parafield);
1271 dec5.setForcedRenormalization(false);
1272 double *res=parafield->getField()->getArray()->getPointer();
1273 const double *toSet=targetResults[rank-nproc_source];
1278 //test 6 - Integral reversed
1279 ParaMEDMEM::InterpKernelDEC dec6(*source_group,*target_group);
1280 parafield->getField()->setNature(Integral);
1281 if (source_group->containsMyRank())
1283 dec6.setMethod("P0");
1284 dec6.attachLocalField(parafield);
1286 dec6.setForcedRenormalization(false);
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);
1295 dec6.setMethod("P0");
1296 dec6.attachLocalField(parafield);
1298 dec6.setForcedRenormalization(false);
1299 double *res=parafield->getField()->getArray()->getPointer();
1300 const double *toSet=targetResults2[rank-nproc_source];
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())
1310 dec7.setMethod("P0");
1311 dec7.attachLocalField(parafield);
1313 dec7.setForcedRenormalization(false);
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);
1322 dec7.setMethod("P0");
1323 dec7.attachLocalField(parafield);
1325 dec7.setForcedRenormalization(false);
1326 double *res=parafield->getField()->getArray()->getPointer();
1327 const double *toSet=targetResults3[rank-nproc_source];
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())
1337 dec8.setMethod("P0");
1338 dec8.attachLocalField(parafield);
1340 dec8.setForcedRenormalization(false);
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);
1349 dec8.setMethod("P0");
1350 dec8.attachLocalField(parafield);
1352 dec8.setForcedRenormalization(false);
1353 double *res=parafield->getField()->getArray()->getPointer();
1354 const double *toSet=targetResults4[rank-nproc_source];
1364 delete target_group;
1365 delete source_group;
1367 MPI_Barrier(MPI_COMM_WORLD);
1370 void ParaMEDMEMTest::testInterpKernelDECNonOverlapp_2D_P0P1P1P0()
1374 MPI_Comm_size(MPI_COMM_WORLD,&size);
1375 MPI_Comm_rank(MPI_COMM_WORLD,&rank);
1379 int nproc_source = 2;
1380 set<int> self_procs;
1381 set<int> procs_source;
1382 set<int> procs_target;
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);
1390 ParaMEDMEM::MEDCouplingUMesh *mesh=0;
1391 ParaMEDMEM::ParaMESH *paramesh=0;
1392 ParaMEDMEM::ParaFIELD *parafieldP0=0,*parafieldP1=0;
1394 ParaMEDMEM::CommInterface interface;
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);
1400 MPI_Barrier(MPI_COMM_WORLD);
1401 if(source_group->containsMyRank())
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();
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();
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);
1444 valueP1[0]=34.; valueP1[1]=77.; valueP1[2]=53.;
1449 valueP1[0]=34.; valueP1[1]=57.; valueP1[2]=77.;
1454 const char targetMeshName[]="target mesh";
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);
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);
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);
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);
1527 ParaMEDMEM::InterpKernelDEC dec(*source_group,*target_group);
1528 if (source_group->containsMyRank())
1530 dec.setMethod("P0");
1531 dec.attachLocalField(parafieldP0);
1533 dec.setForcedRenormalization(false);
1536 const double *valueP0=parafieldP0->getField()->getArray()->getPointer();
1539 CPPUNIT_ASSERT_DOUBLES_EQUAL(34.42857143,valueP0[0],1e-7);
1543 CPPUNIT_ASSERT_DOUBLES_EQUAL(44.,valueP0[0],1e-7);
1548 dec.setMethod("P1");
1549 dec.attachLocalField(parafieldP1);
1551 dec.setForcedRenormalization(false);
1553 const double *res=parafieldP1->getField()->getArray()->getConstPointer();
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);
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);
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);
1586 delete target_group;
1587 delete source_group;
1589 MPI_Barrier(MPI_COMM_WORLD);
1592 void ParaMEDMEMTest::testInterpKernelDEC2DM1D_P0P0()
1596 MPI_Comm_size(MPI_COMM_WORLD,&size);
1597 MPI_Comm_rank(MPI_COMM_WORLD,&rank);
1602 set<int> procs_source;
1603 set<int> procs_target;
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);
1610 ParaMEDMEM::MEDCouplingUMesh *mesh=0;
1611 ParaMEDMEM::ParaMESH *paramesh=0;
1612 ParaMEDMEM::ParaFIELD *parafield=0;
1614 ParaMEDMEM::CommInterface interface;
1616 ProcessorGroup* target_group = new ParaMEDMEM::MPIProcessorGroup(interface,procs_target);
1617 ProcessorGroup* source_group = new ParaMEDMEM::MPIProcessorGroup(interface,procs_source);
1619 MPI_Barrier(MPI_COMM_WORLD);
1620 if(source_group->containsMyRank())
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();
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();
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();
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();
1653 { vals[0]=7.; vals[1]=8.; }
1655 { vals[0]=9.; vals[1]=10.; vals[2]=11.; }
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);
1665 ParaMEDMEM::InterpKernelDEC dec(*source_group,*target_group);
1666 if(source_group->containsMyRank())
1668 dec.setMethod("P0");
1669 dec.attachLocalField(parafield);
1671 dec.setForcedRenormalization(false);
1674 const double *res=parafield->getField()->getArray()->getConstPointer();
1677 CPPUNIT_ASSERT_DOUBLES_EQUAL(9.125,res[0],1e-12);
1678 CPPUNIT_ASSERT_DOUBLES_EQUAL(9.125,res[1],1e-12);
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);
1689 dec.setMethod("P0");
1690 dec.attachLocalField(parafield);
1692 dec.setForcedRenormalization(false);
1694 const double *res=parafield->getField()->getArray()->getConstPointer();
1695 CPPUNIT_ASSERT_DOUBLES_EQUAL(9.125,res[0],1e-12);
1698 ParaMEDMEM::InterpKernelDEC dec2(*source_group,*target_group);
1699 dec2.setMethod("P0");
1700 parafield->getField()->setNature(IntegralGlobConstraint);
1701 if(source_group->containsMyRank())
1703 double *vals=parafield->getField()->getArray()->getPointer();
1705 { vals[0]=7.; vals[1]=8.; }
1707 { vals[0]=9.; vals[1]=10.; vals[2]=11.; }
1708 dec2.attachLocalField(parafield);
1712 const double *res=parafield->getField()->getArray()->getConstPointer();
1715 CPPUNIT_ASSERT_DOUBLES_EQUAL(11.25,res[0],1e-12);
1716 CPPUNIT_ASSERT_DOUBLES_EQUAL(5.625,res[1],1e-12);
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);
1727 dec2.attachLocalField(parafield);
1730 const double *res=parafield->getField()->getArray()->getConstPointer();
1731 CPPUNIT_ASSERT_DOUBLES_EQUAL(45.,res[0],1e-12);
1735 ParaMEDMEM::InterpKernelDEC dec3(*source_group,*target_group);
1736 dec3.setMethod("P0");
1737 parafield->getField()->setNature(Integral);
1738 if(source_group->containsMyRank())
1740 double *vals=parafield->getField()->getArray()->getPointer();
1742 { vals[0]=7.; vals[1]=8.; }
1744 { vals[0]=9.; vals[1]=10.; vals[2]=11.; }
1745 dec3.attachLocalField(parafield);
1749 const double *res=parafield->getField()->getArray()->getConstPointer();
1752 CPPUNIT_ASSERT_DOUBLES_EQUAL(11.25,res[0],1e-12);
1753 CPPUNIT_ASSERT_DOUBLES_EQUAL(5.625,res[1],1e-12);
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);
1764 dec3.attachLocalField(parafield);
1767 const double *res=parafield->getField()->getArray()->getConstPointer();
1768 CPPUNIT_ASSERT_DOUBLES_EQUAL(45.,res[0],1e-12);
1772 ParaMEDMEM::InterpKernelDEC dec4(*source_group,*target_group);
1773 dec4.setMethod("P0");
1774 parafield->getField()->setNature(RevIntegral);
1775 if(source_group->containsMyRank())
1777 double *vals=parafield->getField()->getArray()->getPointer();
1779 { vals[0]=7.; vals[1]=8.; }
1781 { vals[0]=9.; vals[1]=10.; vals[2]=11.; }
1782 dec4.attachLocalField(parafield);
1786 const double *res=parafield->getField()->getArray()->getConstPointer();
1789 CPPUNIT_ASSERT_DOUBLES_EQUAL(9.125,res[0],1e-12);
1790 CPPUNIT_ASSERT_DOUBLES_EQUAL(9.125,res[1],1e-12);
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);
1801 dec4.attachLocalField(parafield);
1804 const double *res=parafield->getField()->getArray()->getConstPointer();
1805 CPPUNIT_ASSERT_DOUBLES_EQUAL(9.125,res[0],1e-12);
1811 delete target_group;
1812 delete source_group;
1814 MPI_Barrier(MPI_COMM_WORLD);
1817 void ParaMEDMEMTest::testInterpKernelDECPartialProcs()
1821 MPI_Comm_size(MPI_COMM_WORLD,&size);
1822 MPI_Comm_rank(MPI_COMM_WORLD,&rank);
1826 set<int> procs_source;
1827 set<int> procs_target;
1829 procs_source.insert(0);
1830 procs_target.insert(1);
1832 ParaMEDMEM::MEDCouplingUMesh *mesh=0;
1833 ParaMEDMEM::ParaMESH *paramesh=0;
1834 ParaMEDMEM::ParaFIELD *parafield=0;
1836 ParaMEDMEM::CommInterface interface;
1838 MPI_Barrier(MPI_COMM_WORLD);
1839 double targetCoords[8]={ 0.,0., 1., 0., 0., 1., 1., 1. };
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);
1848 ProcessorGroup* target_group=0;
1849 ProcessorGroup* source_group=0;
1851 ParaMEDMEM::InterpKernelDEC *dec=0;
1852 if(rank==0 || rank==1)
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())
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();
1875 dec=new ParaMEDMEM::InterpKernelDEC(*source_group,*target_group);
1876 dec->attachLocalField(parafield);
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);
1910 delete target_group;
1911 delete source_group;
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);
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.
1924 void ParaMEDMEMTest::testInterpKernelDEC3DSurfEmptyBBox()
1928 MPI_Comm_size(MPI_COMM_WORLD,&size);
1929 MPI_Comm_rank(MPI_COMM_WORLD,&rank);
1933 int nproc_source = 1;
1934 set<int> self_procs;
1935 set<int> procs_source;
1936 set<int> procs_target;
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);
1944 ParaMEDMEM::MEDCouplingUMesh *mesh=0;
1945 ParaMEDMEM::ParaMESH *paramesh=0;
1946 ParaMEDMEM::ParaFIELD *parafieldP0=0;
1948 ParaMEDMEM::CommInterface interface;
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);
1954 MPI_Barrier(MPI_COMM_WORLD);
1955 if(source_group->containsMyRank())
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();
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.;
1979 const char targetMeshName[]="target mesh";
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);
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);
2010 ParaMEDMEM::ComponentTopology comptopo;
2011 parafieldP0 = new ParaFIELD(ON_CELLS,NO_TIME,paramesh, comptopo);
2012 parafieldP0->getField()->setNature(ConservativeVolumic);
2015 ParaMEDMEM::InterpKernelDEC dec(*source_group,*target_group);
2016 if (source_group->containsMyRank())
2018 dec.setMethod("P0");
2019 dec.attachLocalField(parafieldP0);
2021 // dec.setForcedRenormalization(false);
2024 // const double *valueP0=parafieldP0->getField()->getArray()->getPointer();
2027 // CPPUNIT_ASSERT_DOUBLES_EQUAL(7.4,valueP0[0],1e-7);
2028 // CPPUNIT_ASSERT_DOUBLES_EQUAL(9.0540540540540544,valueP0[1],1e-7);
2032 // CPPUNIT_ASSERT_DOUBLES_EQUAL(8.64054054054054,valueP0[0],1e-7);
2036 // CPPUNIT_ASSERT_DOUBLES_EQUAL(9.0540540540540544,valueP0[0],1e-7);
2041 dec.setMethod("P0");
2042 dec.attachLocalField(parafieldP0);
2044 // dec.setForcedRenormalization(false);
2046 // const double *res=parafieldP0->getField()->getArray()->getConstPointer();
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);
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);
2066 delete target_group;
2067 delete source_group;
2069 MPI_Barrier(MPI_COMM_WORLD);
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
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)
2081 std::string srcM(srcMeth);
2082 std::string targetM(targetMeth);
2085 MPI_Comm_size(MPI_COMM_WORLD,&size);
2086 MPI_Comm_rank(MPI_COMM_WORLD,&rank);
2088 //the test is meant to run on five processors
2089 if (size !=5) return ;
2091 int nproc_source = 3;
2092 set<int> self_procs;
2093 set<int> procs_source;
2094 set<int> procs_target;
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);
2102 ParaMEDMEM::CommInterface interface;
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);
2108 //loading the geometry for the source group
2110 ParaMEDMEM::InterpKernelDEC dec (*source_group,*target_group);
2112 ParaMEDMEM::MEDCouplingUMesh* mesh;
2113 ParaMEDMEM::ParaMESH* paramesh;
2114 ParaMEDMEM::ParaFIELD* parafield;
2116 ICoCo::MEDField* icocofield ;
2118 char * tmp_dir_c = getenv("TMP");
2120 if (tmp_dir_c != NULL)
2121 tmp_dir = string(tmp_dir_c);
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");
2129 // To remove tmp files from disk
2130 ParaMEDMEMTest_TmpFilesRemover aRemover;
2132 MPI_Barrier(MPI_COMM_WORLD);
2134 if (source_group->containsMyRank())
2136 string master = filename_xml1;
2138 ostringstream strstream;
2139 strstream <<master<<rank+1<<".med";
2140 ostringstream meshname ;
2141 meshname<< "Mesh_2_"<< rank+1;
2143 mesh=MEDLoader::ReadUMeshFromFile(strstream.str().c_str(),meshname.str().c_str(),0);
2145 paramesh=new ParaMESH (mesh,*source_group,"source mesh");
2147 // ParaMEDMEM::ParaSUPPORT* parasupport=new UnstructuredParaSUPPORT( support,*source_group);
2148 ParaMEDMEM::ComponentTopology comptopo;
2151 parafield = new ParaFIELD(ON_CELLS,NO_TIME,paramesh, comptopo);
2152 parafield->getField()->setNature(ConservativeVolumic);//InvertIntegral);//ConservativeVolumic);
2155 parafield = new ParaFIELD(ON_NODES,NO_TIME,paramesh, comptopo);
2159 nb_local=mesh->getNumberOfCells();
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++)
2167 // ICoCo::Field* icocofield=new ICoCo::MEDField(paramesh,parafield);
2168 icocofield=new ICoCo::MEDField(parafield->getField());
2170 dec.attachLocalField(icocofield);
2175 //loading the geometry for the target group
2176 if (target_group->containsMyRank())
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;
2184 mesh = MEDLoader::ReadUMeshFromFile(strstream.str().c_str(),meshname.str().c_str(),0);
2186 paramesh=new ParaMESH (mesh,*target_group,"target mesh");
2187 // ParaMEDMEM::ParaSUPPORT* parasupport=new UnstructuredParaSUPPORT(support,*target_group);
2188 ParaMEDMEM::ComponentTopology comptopo;
2191 parafield = new ParaFIELD(ON_CELLS,NO_TIME,paramesh, comptopo);
2192 parafield->getField()->setNature(ConservativeVolumic);//InvertIntegral);//ConservativeVolumic);
2195 parafield = new ParaFIELD(ON_NODES,NO_TIME,paramesh, comptopo);
2199 nb_local=mesh->getNumberOfCells();
2201 nb_local=mesh->getNumberOfNodes();
2203 double *value=parafield->getField()->getArray()->getPointer();
2204 for(int ielem=0; ielem<nb_local;ielem++)
2206 // ICoCo::Field* icocofield=new ICoCo::MEDField(paramesh,parafield);
2207 icocofield=new ICoCo::MEDField(parafield->getField());
2209 dec.attachLocalField(icocofield);
2213 //attaching a DEC to the source group
2215 if (source_group->containsMyRank())
2217 cout<<"DEC usage"<<endl;
2218 dec.setAsynchronous(Asynchronous);
2220 dec.setTimeInterpolationMethod(LinearTimeInterp);
2222 if ( WithPointToPoint ) {
2223 dec.setAllToAllMethod(PointToPoint);
2226 dec.setAllToAllMethod(Native);
2229 dec.setForcedRenormalization(false);
2230 for (double time=0; time<tmaxA+1e-10; time+=dtA)
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 );
2238 dec.sendData( time , 0 );
2240 double* value = parafield->getField()->getArray()->getPointer();
2241 int nb_local=parafield->getField()->getMesh()->getNumberOfCells();
2242 for (int i=0; i<nb_local;i++)
2249 //attaching a DEC to the target group
2250 if (target_group->containsMyRank())
2252 cout<<"DEC usage"<<endl;
2253 dec.setAsynchronous(Asynchronous);
2255 dec.setTimeInterpolationMethod(LinearTimeInterp);
2257 if ( WithPointToPoint ) {
2258 dec.setAllToAllMethod(PointToPoint);
2261 dec.setAllToAllMethod(Native);
2264 dec.setForcedRenormalization(false);
2265 vector<double> times;
2266 for (double time=0; time<tmaxB+1e-10; time+=dtB)
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 ;
2276 CPPUNIT_ASSERT_DOUBLES_EQUAL(vi,time*10000,0.001);
2281 delete source_group;
2282 delete target_group;
2289 cout << "testAsynchronousInterpKernelDEC_2D" << rank << " MPI_Barrier " << endl ;
2291 if (Asynchronous) MPI_Barrier(MPI_COMM_WORLD);
2292 cout << "end of InterpKernelDEC_2D test"<<endl;