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