Salome HOME
Merge from V6_main_20120808 08Aug12
[modules/med.git] / src / ParaMEDMEMTest / test_perf.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 #include <time.h>
21 #include <sys/times.h>
22 #include <sys/time.h>
23 #include "ParaMEDMEMTest.hxx"
24 #include <cppunit/TestAssert.h>
25
26 #include "CommInterface.hxx"
27 #include "ProcessorGroup.hxx"
28 #include "MPIProcessorGroup.hxx"
29 #include "Topology.hxx"
30 #include "DEC.hxx"
31 #include "MxN_Mapping.hxx"
32 #include "InterpKernelDEC.hxx"
33 #include "ParaMESH.hxx"
34 #include "ParaFIELD.hxx"
35 #include "ComponentTopology.hxx"
36 #include "ICoCoMEDField.hxx"
37 #include "MEDLoader.hxx"
38  
39 #include <string>
40 #include <cstring>
41
42 // use this define to enable lines, execution of which leads to Segmentation Fault
43 #define ENABLE_FAULTS
44
45 // use this define to enable CPPUNIT asserts and fails, showing bugs
46 #define ENABLE_FORCED_FAILURES
47
48 #ifndef CLK_TCK 
49 #include <unistd.h>
50 #define CLK_TCK sysconf(_SC_CLK_TCK);
51 #endif 
52
53 using namespace std;
54 using namespace ParaMEDMEM;
55  
56 void testInterpKernelDEC_2D(const string& filename1, const string& meshname1,
57                             const string& filename2, const string& meshname2,
58                             int nproc_source, double epsilon, bool tri, bool all);
59 void get_time( float *telps, float *tuser, float *tsys, float *tcpu );
60
61 int main(int argc, char *argv[])
62 {
63   string filename1, filename2;
64   string meshname1, meshname2;
65   int nproc_source=1, rank;
66   double epsilon=1.e-6;
67   int count=0;
68   bool tri=false;
69   bool all=false;
70
71   MPI_Init(&argc,&argv);
72
73   for(int i=1;i<argc;i++){
74     if( strcmp(argv[i],"-f1") == 0 ){
75       filename1 = argv[++i];
76       count++;
77     }
78     else if( strcmp(argv[i],"-f2") == 0 ){
79       filename2 = argv[++i];
80       count++;
81     }
82     else if( strcmp(argv[i],"-m1") == 0 ){
83       meshname1 = argv[++i];
84       count++;
85     }
86     else if( strcmp(argv[i],"-m2") == 0 ){
87       meshname2 = argv[++i];
88       count++;
89     }
90     else if( strcmp(argv[i],"-ns") == 0 ){
91       nproc_source = atoi(argv[++i]);
92     }
93     else if( strcmp(argv[i],"-eps") == 0 ){
94       epsilon = atof(argv[++i]);
95     }
96     else if( strcmp(argv[i],"-tri") == 0 ){
97       tri = true;
98     }
99     else if( strcmp(argv[i],"-all") == 0 ){
100       all = true;
101     }
102   }
103
104   if( count != 4 ){
105     cout << "usage test_perf -f1 filename1 -m1 meshname1 -f2 filename2 -m2 meshname2 (-ns nproc_source -eps epsilon -tri -all)" << endl;
106     exit(0);
107   }
108
109   MPI_Comm_rank(MPI_COMM_WORLD,&rank);
110   testInterpKernelDEC_2D(filename1,meshname1,filename2,meshname2,nproc_source,epsilon,tri,all);
111
112   MPI_Finalize();
113 }
114
115 void testInterpKernelDEC_2D(const string& filename_xml1, const string& meshname1,
116                             const string& filename_xml2, const string& meshname2,
117                             int nproc_source, double epsilon, bool tri, bool all)
118 {
119   float tcpu, tcpu_u, tcpu_s, telps;
120   int size;
121   int rank;
122   MPI_Comm_size(MPI_COMM_WORLD,&size);
123   MPI_Comm_rank(MPI_COMM_WORLD,&rank);
124  
125   set<int> self_procs;
126   set<int> procs_source;
127   set<int> procs_target;
128   
129   for (int i=0; i<nproc_source; i++)
130     procs_source.insert(i);
131   for (int i=nproc_source; i<size; i++)
132     procs_target.insert(i);
133   self_procs.insert(rank);
134   
135   ParaMEDMEM::CommInterface interface;
136     
137   ParaMEDMEM::ProcessorGroup* self_group = new ParaMEDMEM::MPIProcessorGroup(interface,self_procs);
138   ParaMEDMEM::ProcessorGroup* target_group = new ParaMEDMEM::MPIProcessorGroup(interface,procs_target);
139   ParaMEDMEM::ProcessorGroup* source_group = new ParaMEDMEM::MPIProcessorGroup(interface,procs_source);
140   
141   //loading the geometry for the source group
142
143   ParaMEDMEM::InterpKernelDEC dec (*source_group,*target_group);
144   if(tri)
145     dec.setIntersectionType(INTERP_KERNEL::Triangulation);
146   else
147     dec.setIntersectionType(INTERP_KERNEL::Convex);
148
149   ParaMEDMEM::MEDCouplingUMesh* mesh;
150   ParaMEDMEM::ParaMESH* paramesh;
151   ParaMEDMEM::ParaFIELD* parafield;
152   ICoCo::Field* icocofield ;
153   
154   // To remove tmp files from disk
155   ParaMEDMEMTest_TmpFilesRemover aRemover;
156   
157   MPI_Barrier(MPI_COMM_WORLD);
158   if (source_group->containsMyRank()){
159     string master = filename_xml1;
160       
161     ostringstream strstream;
162     if( nproc_source == 1 )
163       strstream <<master<<".med";
164     else
165       strstream <<master<<rank+1<<".med";
166
167     ostringstream meshname ;
168     if( nproc_source == 1 )
169       meshname<< meshname1;
170     else
171       meshname<< meshname1<<"_"<< rank+1;
172       
173     get_time( &telps, &tcpu_u, &tcpu_s, &tcpu );
174     mesh=MEDLoader::ReadUMeshFromFile(strstream.str().c_str(),meshname.str().c_str(),0);
175     get_time( &telps, &tcpu_u, &tcpu_s, &tcpu );
176     if( rank == 0 )
177       cout << "IO : Telapse = " << telps << " TuserCPU = " << tcpu_u << " TsysCPU = " << tcpu_s << " TCPU = " << tcpu << endl;
178     mesh->incrRef();
179     
180     paramesh=new ParaMESH (mesh,*source_group,"source mesh");
181     
182     ParaMEDMEM::ComponentTopology comptopo;
183     parafield = new ParaFIELD(ON_CELLS, NO_TIME, paramesh, comptopo);
184
185     int nb_local=mesh->getNumberOfCells();
186     double *value=parafield->getField()->getArray()->getPointer();
187     for(int ielem=0; ielem<nb_local;ielem++)
188       value[ielem]=1.0;
189     
190     icocofield=new ICoCo::MEDField((MEDCouplingUMesh *)paramesh->getCellMesh(),parafield->getField());
191      
192     dec.attachLocalField(icocofield);
193   }
194   
195   //loading the geometry for the target group
196   if (target_group->containsMyRank()){
197     string master= filename_xml2;
198     ostringstream strstream;
199     if( (size-nproc_source) == 1 )
200       strstream << master<<".med";
201     else
202       strstream << master<<(rank-nproc_source+1)<<".med";
203     ostringstream meshname ;
204     if( (size-nproc_source) == 1 )
205       meshname<< meshname2;
206     else
207       meshname<< meshname2<<"_"<<rank-nproc_source+1;
208       
209     get_time( &telps, &tcpu_u, &tcpu_s, &tcpu );
210     mesh = MEDLoader::ReadUMeshFromFile(strstream.str().c_str(),meshname.str().c_str(),0);
211     get_time( &telps, &tcpu_u, &tcpu_s, &tcpu );
212     mesh->incrRef();
213
214     paramesh=new ParaMESH (mesh,*target_group,"target mesh");
215     ParaMEDMEM::ComponentTopology comptopo;
216     parafield = new ParaFIELD(ON_CELLS,NO_TIME,paramesh, comptopo);
217
218     int nb_local=mesh->getNumberOfCells();
219     double *value=parafield->getField()->getArray()->getPointer();
220     for(int ielem=0; ielem<nb_local;ielem++)
221       value[ielem]=0.0;
222     icocofield=new ICoCo::MEDField((MEDCouplingUMesh *)paramesh->getCellMesh(),parafield->getField());
223       
224     dec.attachLocalField(icocofield);
225   }
226     
227   
228   //attaching a DEC to the source group 
229   double field_before_int;
230   double field_after_int;
231   
232   if (source_group->containsMyRank()){ 
233     field_before_int = parafield->getVolumeIntegral(0,true);
234     get_time( &telps, &tcpu_u, &tcpu_s, &tcpu );
235     dec.synchronize();
236     get_time( &telps, &tcpu_u, &tcpu_s, &tcpu );
237     if( rank == 0 )
238       cout << "SYNCHRONIZE : Telapse = " << telps << " TuserCPU = " << tcpu_u << " TsysCPU = " << tcpu_s << " TCPU = " << tcpu << endl;
239     cout<<"DEC usage"<<endl;
240     dec.setForcedRenormalization(false);
241     if(all)
242       dec.setAllToAllMethod(PointToPoint);
243
244     get_time( &telps, &tcpu_u, &tcpu_s, &tcpu );
245     dec.sendData();
246     
247     get_time( &telps, &tcpu_u, &tcpu_s, &tcpu );
248     if( rank == 0 )
249       cout << "SEND DATA : Telapse = " << telps << " TuserCPU = " << tcpu_u << " TsysCPU = " << tcpu_s << " TCPU = " << tcpu << endl;
250     dec.recvData();
251      
252     field_after_int = parafield->getVolumeIntegral(0,true);
253 //    CPPUNIT_ASSERT_DOUBLES_EQUAL(field_before_int, field_after_int, epsilon);
254       
255   }
256   
257   //attaching a DEC to the target group
258   if (target_group->containsMyRank()){
259     get_time( &telps, &tcpu_u, &tcpu_s, &tcpu );
260     dec.synchronize();
261     get_time( &telps, &tcpu_u, &tcpu_s, &tcpu );
262     dec.setForcedRenormalization(false);
263     if(all)
264       dec.setAllToAllMethod(PointToPoint);
265
266     get_time( &telps, &tcpu_u, &tcpu_s, &tcpu );
267     dec.recvData();
268     get_time( &telps, &tcpu_u, &tcpu_s, &tcpu );
269     dec.sendData();
270   }
271   
272   get_time( &telps, &tcpu_u, &tcpu_s, &tcpu );
273   if( rank == 0 )
274     cout << "RECV DATA : Telapse = " << telps << " TuserCPU = " << tcpu_u << " TsysCPU = " << tcpu_s << " TCPU = " << tcpu << endl;
275
276   delete source_group;
277   delete target_group;
278   delete self_group;
279   delete paramesh;
280   delete parafield;
281   mesh->decrRef() ;
282   delete icocofield;
283
284   MPI_Barrier(MPI_COMM_WORLD);
285   cout << "end of InterpKernelDEC_2D test"<<endl;
286 }
287
288 void get_time( float *telps, float *tuser, float *tsys, float *tcpu )
289 {
290
291   /* Variables declaration */
292   static time_t zsec = 0;
293   static long zusec = 0;
294   time_t nsec;
295   long nusec;
296   static clock_t zclock = 0;
297   clock_t nclock;
298   static clock_t zuser = 0;
299   static clock_t zsys = 0;
300   clock_t nuser, nsys;
301
302   struct timeval tp;
303   struct timezone tzp;
304   struct tms local;
305
306   MPI_Barrier(MPI_COMM_WORLD);
307
308   /* Elapsed time reading */
309
310   gettimeofday(&tp,&tzp);
311   nsec = tp.tv_sec;
312   nusec = tp.tv_usec;
313   *telps = (float)(nsec-zsec) + (float)(nusec-zusec)/(float)CLOCKS_PER_SEC;
314   
315   zsec = nsec;
316   zusec = nusec;
317
318   /* User and system CPU time reading */
319
320   times(&local);
321   nuser = local.tms_utime;
322   nsys = local.tms_stime;
323   *tuser = (float)(nuser-zuser) / (float)CLK_TCK;
324   *tsys = (float)(nsys-zsys) / (float)CLK_TCK;
325
326   zuser = nuser;
327   zsys = nsys;
328
329   /* CPU time reading */
330
331   nclock = clock();
332   *tcpu = (float)(nclock-zclock) / (float)CLOCKS_PER_SEC;
333   zclock = nclock;
334
335 }
336
337