1 // Copyright (C) 2007-2008 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.
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
19 #include "ParaMEDMEMTest.hxx"
20 #include <cppunit/TestAssert.h>
22 #include "MEDMEM_Exception.hxx"
23 #include "CommInterface.hxx"
24 #include "ProcessorGroup.hxx"
25 #include "MPIProcessorGroup.hxx"
26 #include "Topology.hxx"
28 #include "NonCoincidentDEC.hxx"
29 #include "ParaMESH.hxx"
30 #include "ParaFIELD.hxx"
31 #include "UnstructuredParaSUPPORT.hxx"
32 #include "ICoCoMEDField.hxx"
36 // use this define to enable lines, execution of which leads to Segmentation Fault
39 // use this define to enable CPPUNIT asserts and fails, showing bugs
40 #define ENABLE_FORCED_FAILURES
44 using namespace ParaMEDMEM;
45 using namespace MEDMEM;
48 * Check methods defined in IntersectionDEC.hxx
51 IntersectionDEC(ProcessorGroup& local_group, ProcessorGroup& distant_group);
52 virtual ~IntersectionDEC();
58 void ParaMEDMEMTest::testNonCoincidentDEC_2D()
62 MPI_Comm_size(MPI_COMM_WORLD,&size);
64 //the test is meant to run on five processors
65 if (size !=5) return ;
67 testNonCoincidentDEC( "/share/salome/resources/med/square1_split",
69 "/share/salome/resources/med/square2_split",
75 void ParaMEDMEMTest::testNonCoincidentDEC_3D()
78 MPI_Comm_size(MPI_COMM_WORLD,&size);
80 //the test is meant to run on five processors
81 if (size !=4) return ;
83 testNonCoincidentDEC( "/share/salome/resources/med/blade_12000_split2",
85 "/share/salome/resources/med/blade_3000_split2",
91 void ParaMEDMEMTest::testNonCoincidentDEC(const string& filename1,
92 const string& meshname1,
93 const string& filename2,
94 const string& meshname2,
100 MPI_Comm_size(MPI_COMM_WORLD,&size);
101 MPI_Comm_rank(MPI_COMM_WORLD,&rank);
104 set<int> procs_source;
105 set<int> procs_target;
107 for (int i=0; i<nproc_source; i++)
108 procs_source.insert(i);
109 for (int i=nproc_source; i<size; i++)
110 procs_target.insert(i);
111 self_procs.insert(rank);
113 ParaMEDMEM::CommInterface interface;
115 ParaMEDMEM::ProcessorGroup* self_group = new ParaMEDMEM::MPIProcessorGroup(interface,self_procs);
116 ParaMEDMEM::ProcessorGroup* target_group = new ParaMEDMEM::MPIProcessorGroup(interface,procs_target);
117 ParaMEDMEM::ProcessorGroup* source_group = new ParaMEDMEM::MPIProcessorGroup(interface,procs_source);
119 ParaMEDMEM::ParaMESH* source_mesh=0;
120 ParaMEDMEM::ParaMESH* target_mesh=0;
121 ParaMEDMEM::ParaSUPPORT* parasupport=0;
122 //loading the geometry for the source group
124 ParaMEDMEM::NonCoincidentDEC dec (*source_group,*target_group);
127 MEDMEM::SUPPORT* support;
128 MEDMEM::FIELD<double>* field;
129 ParaMEDMEM::ParaMESH* paramesh;
130 ParaMEDMEM::ParaFIELD* parafield;
132 string data_dir = getenv("MED_ROOT_DIR") + "/share/salome/resources/med/";
133 string tmp_dir = getenv("TMP");
136 string filename_xml1 = data_dir + filename1;
137 string filename_xml2 = data_dir + filename2;
138 string filename_seq_wr = tmp_dir + "/";
139 string filename_seq_med = tmp_dir + "/myWrField_seq_pointe221.med";
141 // To remove tmp files from disk
142 ParaMEDMEMTest_TmpFilesRemover aRemover;
143 //aRemover.Register(filename_seq_wr);
144 //aRemover.Register(filename_seq_med);
145 MPI_Barrier(MPI_COMM_WORLD);
146 ICoCo::Field* icocofield;
147 if (source_group->containsMyRank())
149 string master = filename_xml1;
151 ostringstream strstream;
152 strstream <<master<<rank+1<<".med";
153 ostringstream meshname ;
154 meshname<< meshname1<<"_"<< rank+1;
156 CPPUNIT_ASSERT_NO_THROW(mesh = new MESH(MED_DRIVER,strstream.str(),meshname.str()));
157 support=new MEDMEM::SUPPORT(mesh,"all elements",MED_EN::MED_CELL);
159 paramesh=new ParaMESH (*mesh,*source_group,"source mesh");
161 parasupport=new UnstructuredParaSUPPORT( support,*source_group);
162 ParaMEDMEM::ComponentTopology comptopo;
163 parafield = new ParaFIELD(parasupport, comptopo);
166 int nb_local=support->getNumberOfElements(MED_EN::MED_ALL_ELEMENTS);
167 double * value= new double[nb_local];
168 for(int ielem=0; ielem<nb_local;ielem++)
170 parafield->getField()->setValue(value);
172 icocofield=new ICoCo::MEDField(paramesh,parafield);
174 dec.attachLocalField(icocofield);
178 //loading the geometry for the target group
179 if (target_group->containsMyRank())
181 string master= filename_xml2;
182 ostringstream strstream;
183 strstream << master<<(rank-nproc_source+1)<<".med";
184 ostringstream meshname ;
185 meshname<< meshname2<<"_"<<rank-nproc_source+1;
187 CPPUNIT_ASSERT_NO_THROW(mesh = new MESH(MED_DRIVER,strstream.str(),meshname.str()));
188 support=new MEDMEM::SUPPORT(mesh,"all elements",MED_EN::MED_CELL);
190 paramesh=new ParaMESH (*mesh,*target_group,"target mesh");
191 parasupport=new UnstructuredParaSUPPORT(support,*target_group);
192 ParaMEDMEM::ComponentTopology comptopo;
193 parafield = new ParaFIELD(parasupport, comptopo);
196 int nb_local=support->getNumberOfElements(MED_EN::MED_ALL_ELEMENTS);
197 double * value= new double[nb_local];
198 for(int ielem=0; ielem<nb_local;ielem++)
200 parafield->getField()->setValue(value);
201 icocofield=new ICoCo::MEDField(paramesh,parafield);
203 dec.attachLocalField(icocofield);
208 //attaching a DEC to the source group
209 double field_before_int;
210 double field_after_int;
212 if (source_group->containsMyRank())
214 field_before_int = parafield->getVolumeIntegral(1);
215 MPI_Bcast(&field_before_int, 1,MPI_DOUBLE, 0,MPI_COMM_WORLD);
217 cout<<"DEC usage"<<endl;
218 dec.setOption("ForcedRenormalization",false);
221 // paramesh->write(MED_DRIVER,"./sourcesquarenc");
222 //parafield->write(MED_DRIVER,"./sourcesquarenc","boundary");
227 //attaching a DEC to the target group
228 if (target_group->containsMyRank())
230 MPI_Bcast(&field_before_int, 1,MPI_DOUBLE, 0,MPI_COMM_WORLD);
233 dec.setOption("ForcedRenormalization",false);
235 //paramesh->write(MED_DRIVER, "./targetsquarenc");
236 //parafield->write(MED_DRIVER, "./targetsquarenc", "boundary");
237 field_after_int = parafield->getVolumeIntegral(1);
240 MPI_Bcast(&field_before_int,1,MPI_DOUBLE,0,MPI_COMM_WORLD);
241 MPI_Bcast(&field_after_int, 1,MPI_DOUBLE, size-1,MPI_COMM_WORLD);
243 CPPUNIT_ASSERT_DOUBLES_EQUAL(field_before_int, field_after_int, epsilon);
254 MPI_Barrier(MPI_COMM_WORLD);