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
22 #include "ParaMEDMEMTest.hxx"
23 #include <cppunit/TestAssert.h>
25 #include "MEDMEM_Exception.hxx"
26 #include "CommInterface.hxx"
27 #include "ProcessorGroup.hxx"
28 #include "MPIProcessorGroup.hxx"
29 #include "Topology.hxx"
31 #include "NonCoincidentDEC.hxx"
32 #include "ParaMESH.hxx"
33 #include "ParaFIELD.hxx"
34 #include "UnstructuredParaSUPPORT.hxx"
35 #include "ICoCoMEDField.hxx"
39 // use this define to enable lines, execution of which leads to Segmentation Fault
42 // use this define to enable CPPUNIT asserts and fails, showing bugs
43 #define ENABLE_FORCED_FAILURES
47 using namespace ParaMEDMEM;
48 using namespace MEDMEM;
51 * Check methods defined in InterpKernelDEC.hxx
54 InterpKernelDEC(ProcessorGroup& local_group, ProcessorGroup& distant_group);
55 virtual ~InterpKernelDEC();
61 void ParaMEDMEMTest::testNonCoincidentDEC_2D()
65 MPI_Comm_size(MPI_COMM_WORLD,&size);
67 //the test is meant to run on five processors
68 if (size !=5) return ;
70 testNonCoincidentDEC( "/share/salome/resources/med/square1_split",
72 "/share/salome/resources/med/square2_split",
78 void ParaMEDMEMTest::testNonCoincidentDEC_3D()
81 MPI_Comm_size(MPI_COMM_WORLD,&size);
83 //the test is meant to run on five processors
84 if (size !=4) return ;
86 testNonCoincidentDEC( "/share/salome/resources/med/blade_12000_split2",
88 "/share/salome/resources/med/blade_3000_split2",
94 void ParaMEDMEMTest::testNonCoincidentDEC(const string& filename1,
95 const string& meshname1,
96 const string& filename2,
97 const string& meshname2,
103 MPI_Comm_size(MPI_COMM_WORLD,&size);
104 MPI_Comm_rank(MPI_COMM_WORLD,&rank);
107 set<int> procs_source;
108 set<int> procs_target;
110 for (int i=0; i<nproc_source; i++)
111 procs_source.insert(i);
112 for (int i=nproc_source; i<size; i++)
113 procs_target.insert(i);
114 self_procs.insert(rank);
116 ParaMEDMEM::CommInterface interface;
118 ParaMEDMEM::ProcessorGroup* self_group = new ParaMEDMEM::MPIProcessorGroup(interface,self_procs);
119 ParaMEDMEM::ProcessorGroup* target_group = new ParaMEDMEM::MPIProcessorGroup(interface,procs_target);
120 ParaMEDMEM::ProcessorGroup* source_group = new ParaMEDMEM::MPIProcessorGroup(interface,procs_source);
122 ParaMEDMEM::ParaMESH* source_mesh=0;
123 ParaMEDMEM::ParaMESH* target_mesh=0;
124 ParaMEDMEM::ParaSUPPORT* parasupport=0;
125 //loading the geometry for the source group
127 ParaMEDMEM::NonCoincidentDEC dec (*source_group,*target_group);
130 MEDMEM::SUPPORT* support;
131 MEDMEM::FIELD<double>* field;
132 ParaMEDMEM::ParaMESH* paramesh;
133 ParaMEDMEM::ParaFIELD* parafield;
135 string filename_xml1 = getResourceFile(filename1);
136 string filename_xml2 = getResourceFile(filename2);
137 //string filename_seq_wr = makeTmpFile("");
138 //string filename_seq_med = makeTmpFile("myWrField_seq_pointe221.med");
140 // To remove tmp files from disk
141 ParaMEDMEMTest_TmpFilesRemover aRemover;
142 //aRemover.Register(filename_seq_wr);
143 //aRemover.Register(filename_seq_med);
144 MPI_Barrier(MPI_COMM_WORLD);
145 ICoCo::Field* icocofield;
146 if (source_group->containsMyRank())
148 string master = filename_xml1;
150 ostringstream strstream;
151 strstream <<master<<rank+1<<".med";
152 ostringstream meshname ;
153 meshname<< meshname1<<"_"<< rank+1;
155 CPPUNIT_ASSERT_NO_THROW(mesh = new MESH(MED_DRIVER,strstream.str(),meshname.str()));
156 support=new MEDMEM::SUPPORT(mesh,"all elements",MED_EN::MED_CELL);
158 paramesh=new ParaMESH (*mesh,*source_group,"source mesh");
160 parasupport=new UnstructuredParaSUPPORT( support,*source_group);
161 ParaMEDMEM::ComponentTopology comptopo;
162 parafield = new ParaFIELD(parasupport, comptopo);
165 int nb_local=support->getNumberOfElements(MED_EN::MED_ALL_ELEMENTS);
166 double * value= new double[nb_local];
167 for(int ielem=0; ielem<nb_local;ielem++)
169 parafield->getField()->setValue(value);
171 icocofield=new ICoCo::MEDField(paramesh,parafield);
173 dec.attachLocalField(icocofield);
177 //loading the geometry for the target group
178 if (target_group->containsMyRank())
180 string master= filename_xml2;
181 ostringstream strstream;
182 strstream << master<<(rank-nproc_source+1)<<".med";
183 ostringstream meshname ;
184 meshname<< meshname2<<"_"<<rank-nproc_source+1;
186 CPPUNIT_ASSERT_NO_THROW(mesh = new MESH(MED_DRIVER,strstream.str(),meshname.str()));
187 support=new MEDMEM::SUPPORT(mesh,"all elements",MED_EN::MED_CELL);
189 paramesh=new ParaMESH (*mesh,*target_group,"target mesh");
190 parasupport=new UnstructuredParaSUPPORT(support,*target_group);
191 ParaMEDMEM::ComponentTopology comptopo;
192 parafield = new ParaFIELD(parasupport, comptopo);
195 int nb_local=support->getNumberOfElements(MED_EN::MED_ALL_ELEMENTS);
196 double * value= new double[nb_local];
197 for(int ielem=0; ielem<nb_local;ielem++)
199 parafield->getField()->setValue(value);
200 icocofield=new ICoCo::MEDField(paramesh,parafield);
202 dec.attachLocalField(icocofield);
207 //attaching a DEC to the source group
208 double field_before_int;
209 double field_after_int;
211 if (source_group->containsMyRank())
213 field_before_int = parafield->getVolumeIntegral(1);
214 MPI_Bcast(&field_before_int, 1,MPI_DOUBLE, 0,MPI_COMM_WORLD);
216 cout<<"DEC usage"<<endl;
217 dec.setOption("ForcedRenormalization",false);
220 // paramesh->write(MED_DRIVER,"./sourcesquarenc");
221 //parafield->write(MED_DRIVER,"./sourcesquarenc","boundary");
226 //attaching a DEC to the target group
227 if (target_group->containsMyRank())
229 MPI_Bcast(&field_before_int, 1,MPI_DOUBLE, 0,MPI_COMM_WORLD);
232 dec.setOption("ForcedRenormalization",false);
234 //paramesh->write(MED_DRIVER, "./targetsquarenc");
235 //parafield->write(MED_DRIVER, "./targetsquarenc", "boundary");
236 field_after_int = parafield->getVolumeIntegral(1);
239 MPI_Bcast(&field_before_int,1,MPI_DOUBLE,0,MPI_COMM_WORLD);
240 MPI_Bcast(&field_after_int, 1,MPI_DOUBLE, size-1,MPI_COMM_WORLD);
242 CPPUNIT_ASSERT_DOUBLES_EQUAL(field_before_int, field_after_int, epsilon);
253 MPI_Barrier(MPI_COMM_WORLD);