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 "StructuredCoincidentDEC.hxx"
29 #include "ParaMESH.hxx"
30 #include "ParaFIELD.hxx"
31 #include "ComponentTopology.hxx"
32 #include "ICoCoMEDField.hxx"
33 #include "MEDLoader.hxx"
34 #include "TestInterpKernelUtils.hxx"
38 // use this define to enable lines, execution of which leads to Segmentation Fault
41 // use this define to enable CPPUNIT asserts and fails, showing bugs
42 #define ENABLE_FORCED_FAILURES
45 using namespace ParaMEDMEM;
48 * Check methods defined in StructuredCoincidentDEC.hxx
50 StructuredCoincidentDEC();
51 StructuredCoincidentDEC(ProcessorGroup& local_group, ProcessorGroup& distant_group);
52 virtual ~StructuredCoincidentDEC();
58 void ParaMEDMEMTest::testStructuredCoincidentDEC() {
59 string testname="ParaMEDMEM - testStructured CoincidentDEC";
60 // MPI_Init(&argc, &argv);
63 MPI_Comm_size(MPI_COMM_WORLD, &size);
64 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
69 ParaMEDMEM::CommInterface interface;
71 ParaMEDMEM::MPIProcessorGroup self_group (interface,rank,rank);
72 ParaMEDMEM::MPIProcessorGroup target_group(interface,3,size-1);
73 ParaMEDMEM::MPIProcessorGroup source_group (interface,0,2);
75 ParaMEDMEM::MEDCouplingUMesh* mesh;
76 ParaMEDMEM::ParaMESH* paramesh;
77 ParaMEDMEM::ParaFIELD* parafield;
79 string filename_xml1 = INTERP_TEST::getResourceFile("square1_split");
80 string filename_2 = INTERP_TEST::getResourceFile("square1.med");
81 //string filename_seq_wr = makeTmpFile("");
82 //string filename_seq_med = makeTmpFile("myWrField_seq_pointe221.med");
84 // To remove tmp files from disk
85 ParaMEDMEMTest_TmpFilesRemover aRemover;
87 //loading the geometry for the source group
89 ParaMEDMEM::StructuredCoincidentDEC dec(source_group, target_group);
91 MPI_Barrier(MPI_COMM_WORLD);
92 if (source_group.containsMyRank()) {
93 string master = filename_xml1;
95 ostringstream strstream;
96 strstream <<master<<rank+1<<".med";
97 ostringstream meshname;
98 meshname<< "Mesh_2_"<< rank+1;
100 mesh=MEDLoader::ReadUMeshFromFile(strstream.str().c_str(),meshname.str().c_str(),0);
103 paramesh=new ParaMESH (mesh,source_group,"source mesh");
105 ParaMEDMEM::ComponentTopology comptopo(6);
106 parafield = new ParaFIELD(ON_CELLS,NO_TIME,paramesh, comptopo);
108 int nb_local=mesh->getNumberOfCells();
109 const int* global_numbering = paramesh->getGlobalNumberingCell();
111 double *value=parafield->getField()->getArray()->getPointer();
112 for(int ielem=0; ielem<nb_local;ielem++)
113 for (int icomp=0; icomp<6; icomp++)
114 value[ielem*6+icomp]=global_numbering[ielem]*6+icomp;
116 //ICoCo::Field* icocofield=new ICoCo::MEDField((MEDCouplingUMesh *)paramesh->getCellMesh(),parafield->getField());
118 dec.attachLocalField(parafield);
124 //loading the geometry for the target group
125 if (target_group.containsMyRank()) {
127 string meshname2("Mesh_2");
128 mesh = MEDLoader::ReadUMeshFromFile(filename_2.c_str(),meshname2.c_str(),0);
130 paramesh=new ParaMESH (mesh,self_group,"target mesh");
131 ParaMEDMEM::ComponentTopology comptopo(6, &target_group);
133 parafield = new ParaFIELD(ON_CELLS,NO_TIME,paramesh, comptopo);
135 int nb_local=mesh->getNumberOfCells();
136 double *value=parafield->getField()->getArray()->getPointer();
137 for (int ielem=0; ielem<nb_local; ielem++)
138 for (int icomp=0; icomp<comptopo.nbLocalComponents(); icomp++)
139 value[ielem*comptopo.nbLocalComponents()+icomp]=0.0;
140 //ICoCo::Field* icocofield=new ICoCo::MEDField((MEDCouplingUMesh *)paramesh->getCellMesh(),parafield->getField());
142 dec.attachLocalField(parafield);
146 //checking validity of field
147 const double* recv_value = parafield->getField()->getArray()->getPointer();
148 for (int i=0; i< nb_local; i++) {
149 int first = comptopo.firstLocalComponent();
150 for (int icomp = 0; icomp < comptopo.nbLocalComponents(); icomp++)
151 CPPUNIT_ASSERT_DOUBLES_EQUAL(recv_value[i*comptopo.nbLocalComponents()+icomp],(double)(i*6+icomp+first),1e-12);
159 // MPI_Barrier(MPI_COMM_WORLD);