Salome HOME
Merge from BR_V5_DEV 16Feb09
[modules/med.git] / src / ParaMEDMEM / Test / ParaMEDMEMTest_NonCoincidentDEC.cxx
1 //  Copyright (C) 2007-2008  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.
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 #include "ParaMEDMEMTest.hxx"
20 #include <cppunit/TestAssert.h>
21
22 #include "MEDMEM_Exception.hxx"
23 #include "CommInterface.hxx"
24 #include "ProcessorGroup.hxx"
25 #include "MPIProcessorGroup.hxx"
26 #include "Topology.hxx"
27 #include "DEC.hxx"
28 #include "NonCoincidentDEC.hxx"
29 #include "ParaMESH.hxx"
30 #include "ParaFIELD.hxx"
31 #include "UnstructuredParaSUPPORT.hxx"
32 #include "ICoCoMEDField.hxx"
33  
34 #include <string>
35
36 // use this define to enable lines, execution of which leads to Segmentation Fault
37 #define ENABLE_FAULTS
38
39 // use this define to enable CPPUNIT asserts and fails, showing bugs
40 #define ENABLE_FORCED_FAILURES
41
42
43 using namespace std;
44 using namespace ParaMEDMEM;
45 using namespace MEDMEM;
46  
47 /*
48  * Check methods defined in IntersectionDEC.hxx
49  *
50  IntersectionDEC();
51  IntersectionDEC(ProcessorGroup& local_group, ProcessorGroup& distant_group);
52  virtual ~IntersectionDEC();
53  void synchronize();
54  void recvData();
55  void sendData();
56 */
57
58 void ParaMEDMEMTest::testNonCoincidentDEC_2D()
59 {
60
61   int size;
62   MPI_Comm_size(MPI_COMM_WORLD,&size);
63
64   //the test is meant to run on five processors
65   if (size !=5) return ;
66   
67   testNonCoincidentDEC( "/share/salome/resources/med/square1_split",
68                         "Mesh_2",
69                         "/share/salome/resources/med/square2_split",
70                         "Mesh_3",
71                         3,
72                         1e-6);
73
74
75 void ParaMEDMEMTest::testNonCoincidentDEC_3D()
76 {
77   int size;
78   MPI_Comm_size(MPI_COMM_WORLD,&size);
79   
80   //the test is meant to run on five processors
81   if (size !=4) return ;
82   
83   testNonCoincidentDEC( "/share/salome/resources/med/blade_12000_split2",
84                         "Mesh_1",
85                         "/share/salome/resources/med/blade_3000_split2",
86                         "Mesh_1",
87                         2,
88                         1e4);
89
90
91 void ParaMEDMEMTest::testNonCoincidentDEC(const string& filename1,
92                                           const string& meshname1,
93                                           const string& filename2,
94                                           const string& meshname2,
95                                           int nproc_source,
96                                           double epsilon)
97 {
98   int size;
99   int rank;
100   MPI_Comm_size(MPI_COMM_WORLD,&size);
101   MPI_Comm_rank(MPI_COMM_WORLD,&rank);
102    
103   set<int> self_procs;
104   set<int> procs_source;
105   set<int> procs_target;
106   
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);
112   
113   ParaMEDMEM::CommInterface interface;
114     
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);
118   
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
123
124   ParaMEDMEM::NonCoincidentDEC dec (*source_group,*target_group);
125
126   MEDMEM::MESH* mesh;
127   MEDMEM::SUPPORT* support;
128   MEDMEM::FIELD<double>* field;
129   ParaMEDMEM::ParaMESH* paramesh;
130   ParaMEDMEM::ParaFIELD* parafield;
131   
132   string data_dir                   = getenv("MED_ROOT_DIR") + "/share/salome/resources/med/";
133   string tmp_dir                    = getenv("TMP");
134   if (tmp_dir == "")
135     tmp_dir = "/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";
140   
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())
148     {
149       string master = filename_xml1;
150       
151       ostringstream strstream;
152       strstream <<master<<rank+1<<".med";
153       ostringstream meshname ;
154       meshname<< meshname1<<"_"<< rank+1;
155       
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);
158     
159       paramesh=new ParaMESH (*mesh,*source_group,"source mesh");
160     
161       parasupport=new UnstructuredParaSUPPORT( support,*source_group);
162       ParaMEDMEM::ComponentTopology comptopo;
163       parafield = new ParaFIELD(parasupport, comptopo);
164
165       
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++)
169         value[ielem]=1.0;
170       parafield->getField()->setValue(value);
171
172       icocofield=new ICoCo::MEDField(paramesh,parafield);
173      
174       dec.attachLocalField(icocofield);
175       delete [] value;
176     }
177   
178   //loading the geometry for the target group
179   if (target_group->containsMyRank())
180     {
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;
186       
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);
189       
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);
194
195       
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++)
199         value[ielem]=0.0;
200       parafield->getField()->setValue(value);
201       icocofield=new ICoCo::MEDField(paramesh,parafield);
202       
203       dec.attachLocalField(icocofield);
204       delete [] value;
205     }
206     
207   
208   //attaching a DEC to the source group 
209   double field_before_int;
210   double field_after_int;
211   
212   if (source_group->containsMyRank())
213     { 
214       field_before_int = parafield->getVolumeIntegral(1);
215       MPI_Bcast(&field_before_int, 1,MPI_DOUBLE, 0,MPI_COMM_WORLD);
216       dec.synchronize();
217       cout<<"DEC usage"<<endl;
218       dec.setOption("ForcedRenormalization",false);
219
220       dec.sendData();
221       //      paramesh->write(MED_DRIVER,"./sourcesquarenc");
222       //parafield->write(MED_DRIVER,"./sourcesquarenc","boundary");  
223    
224       
225     }
226   
227   //attaching a DEC to the target group
228   if (target_group->containsMyRank())
229     {
230       MPI_Bcast(&field_before_int, 1,MPI_DOUBLE, 0,MPI_COMM_WORLD);
231      
232       dec.synchronize();
233       dec.setOption("ForcedRenormalization",false);
234       dec.recvData();
235       //paramesh->write(MED_DRIVER, "./targetsquarenc");
236       //parafield->write(MED_DRIVER, "./targetsquarenc", "boundary");
237       field_after_int = parafield->getVolumeIntegral(1);
238       
239     }
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);
242      
243   CPPUNIT_ASSERT_DOUBLES_EQUAL(field_before_int, field_after_int, epsilon);
244     
245   delete source_group;
246   delete target_group;
247   delete self_group;
248   delete icocofield;
249   delete paramesh;
250   delete parafield;
251   delete support;
252   delete parasupport;
253   delete mesh;
254   MPI_Barrier(MPI_COMM_WORLD);
255   
256 }