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