Salome HOME
0f9521ef9c210c27fbc8ff35a84bad40cfb46755
[tools/medcoupling.git] / src / MEDPartitioner / medpartitioner_para.cxx
1 // Copyright (C) 2007-2019  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 /*
21   examples of launch
22   rm ttmp* tttmp*
23   export verb=11
24   mpirun -np 2 medpartitioner_para --input-file=blade.med --output-file=ttmp1_ --ndomains=2 --dump-cpu-memory --verbose=$verb
25   mpirun -np 5 medpartitioner_para --input-file=blade.med --output-file=ttmp1_ --ndomains=2 --dump-cpu-memory --verbose=$verb
26   mpirun -np 4 medpartitioner_para --input-file=ttmp1_.xml --output-file=tttmp1_ --ndomains=4 --verbose=$verb
27
28   mpirun -np 2 medpartitioner_para --ndomains=2 --input-file=tmp_testMesh_20x30x50.med --output-file=ttmp2_ --verbose=$verb
29   mpirun -np 2 medpartitioner_para --ndomains=2 --input-file=ttmp2_.xml --output-file=ttmp3_ --verbose=$verb
30
31   mpirun -np 2 medpartitioner_para --ndomains=2 --input-file=tmp_testMeshWithFaces_20x30x50.med --output-file=ttmp2_ --verbose=$verb
32   mpirun -np 2 medpartitioner_para --ndomains=2 --input-file=ttmp2_.xml --output-file=ttmp3_ --verbose=$verb
33
34   mpirun -np 2 medpartitioner_para --ndomains=2 --input-file=tmp_testMesh_20x30x50_WithVecFieldOnCells.med --output-file=ttmp2_ --verbose=$verb
35   mpirun -np 2 medpartitioner_para --ndomains=2 --input-file=tmp_testMesh_20x30x50_WithVecFieldOnNodes.med --output-file=ttmp2_ --verbose=$verb
36   mpirun -np 2 medpartitioner_para --ndomains=2 --input-file=ttmp2_.xml --output-file=ttmp3_ --verbose=$verb
37
38   mpirun -np 4 medpartitioner_para --ndomains=4 --input-file=tmp_testMeshHuge_20x30x50_6.xml --output-file=ttmp3_ --verbose=1
39 */
40
41
42 #include "MEDPARTITIONER_MeshCollection.hxx"
43 #include "MEDPARTITIONER_ParallelTopology.hxx"
44 #include "MEDPARTITIONER_ParaDomainSelector.hxx"
45 #include "MEDPARTITIONER_Utils.hxx"
46
47 #include "MEDLoader.hxx"
48
49 #include <fstream>
50 #include <iostream>
51 #include <iomanip>
52 #include <sstream>
53 #include <string>
54 #include <cstring>
55
56 #ifdef HAVE_MPI
57 #include <mpi.h>
58 #endif
59
60 using namespace std;
61 using namespace MEDPARTITIONER;
62
63 int main(int argc, char** argv)
64 {
65 #if !defined(MED_ENABLE_PARMETIS)
66   cout << "Sorry, no one split method is available. Please, compile with ParMETIS."<<endl;
67   return 1;
68 #else
69
70   // Defining options
71   // by parsing the command line
72   
73   bool split_family=false;
74   bool empty_groups=false;
75   bool mesure_memory=false;
76   bool filter_face=true;
77
78   string input;
79   string output;
80   string meshname;
81   string library="metis";  //default
82   int ndomains;
83   int help=0;
84   int test=0;
85   MPI_Init(&argc,&argv);
86   MPI_Comm_size(MPI_COMM_WORLD, &MyGlobals::_World_Size);
87   MPI_Comm_rank(MPI_COMM_WORLD, &MyGlobals::_Rank);
88   
89   //cout<<"proc "<<MyGlobals::_Rank<<" of "<<MyGlobals::_World_Size<<endl; //for debug
90   //TestVectorOfStringMpi(); //cvw
91   //TestRandomize();
92   
93   //Primitive parsing of command-line options
94   string desc ("Available options of medpartitioner_para V1.0.3:\n"
95                "\t--help                   : produces this help message\n"
96                "\t--verbose                : echoes arguments\n"
97                "\t--input-file=<string>    : name of the input .med file or .xml master file\n"
98                "\t--output-file=<string>   : name of the resulting file (without extension)\n"
99                "\t--ndomains=<number>      : number of subdomains in the output file, default is 1\n"
100 #if defined(MED_ENABLE_PARMETIS) 
101 // || defined(MED_ENABLE_PTSCOTCH)
102 //user can choose! (not yet)
103                "\t--split-method=<string>  : name of the splitting library (metis/scotch), default is metis\n"
104 #endif
105                "\t--create-boundary-faces : creates boundary faces mesh in the output files\n"
106                "\t--dump-cpu-memory        : dumps passed CPU time and maximal increase of used memory\n"
107                //"\t--randomize=<number>     : random seed for other partitionning (only on one proc)\n"
108                //"\t--atomize                : do the opposite of a good partitionner (only on one proc)\n"
109                );
110
111   if (argc<=1) help=1;
112   string value;
113   for (int i = 1; i < argc; i++)
114     {
115       if (strlen(argv[i]) < 3) 
116         {
117           if (MyGlobals::_Rank==0) cerr << "bad argument : "<< argv[i] << endl;
118           MPI_Finalize(); return 1;
119         }
120     
121       if (TestArg(argv[i],"--verbose",value)) 
122         {
123           MyGlobals::_Verbose=1;
124           if (value!="") MyGlobals::_Verbose = atoi(value.c_str());
125         }
126       else if (TestArg(argv[i],"--help",value)) help=1;
127       else if (TestArg(argv[i],"--test",value)) test=1;
128       else if (TestArg(argv[i],"--input-file",value)) input=value;
129       else if (TestArg(argv[i],"--output-file",value)) output=value;
130       else if (TestArg(argv[i],"--split-method",value)) library=value;
131       else if (TestArg(argv[i],"--ndomains",value)) ndomains=atoi(value.c_str());
132       else if (TestArg(argv[i],"--randomize",value)) MyGlobals::_Randomize=atoi(value.c_str());
133       else if (TestArg(argv[i],"--atomize",value)) MyGlobals::_Atomize=atoi(value.c_str());
134       else if (TestArg(argv[i],"--create-boundary-faces",value)) MyGlobals::_Create_Boundary_Faces=1;
135       else if (TestArg(argv[i],"--dump-cpu-memory",value)) mesure_memory=true;
136       else 
137         {
138           if (MyGlobals::_Rank==0) cerr << "unknown argument : "<< argv[i] << endl;
139           MPI_Finalize(); return 1;
140         }
141     }
142
143
144   if (MyGlobals::_Randomize!=0 && MyGlobals::_World_Size!=1)
145     {
146       if (MyGlobals::_Rank==0) cerr << "randomize only available on 1 proc (mpirun -np 1)" << endl;
147       MyGlobals::_Randomize=0;
148     }
149
150 //no choice
151 #if defined(MED_ENABLE_PARMETIS) && !defined(MED_ENABLE_SCOTCH)
152   library = "metis";
153 #endif
154 #if !defined(MED_ENABLE_PARMETIS) && defined(MED_ENABLE_SCOTCH)
155   library = "scotch";
156 #endif
157 //user choice
158 #if defined(MED_ENABLE_PARMETIS) && defined(MED_ENABLE_SCOTCH)
159   if ((library!="metis") && (library!="scotch"))
160     {
161       if (MyGlobals::_Rank==0) cerr << "split-method only available : metis, scotch" << endl;
162       MPI_Finalize(); return 1;
163     }
164 #endif
165  
166   if (help==1)
167     {
168       if (MyGlobals::_Rank==0) cout<<desc<<"\n";
169       MPI_Finalize(); return 0;
170     }
171   
172   MyGlobals::_Is0verbose=0;
173   if (MyGlobals::_Rank==0) MyGlobals::_Is0verbose=MyGlobals::_Verbose;
174   //MyGlobals::_Is0verbose=((MyGlobals::_Rank==0) && MyGlobals::_Verbose);
175   if (test==1) //only for debugger
176     {
177       if (MyGlobals::_Is0verbose>0) cout<<"tests on "<<MyGlobals::_Atomize<<" "<<ndomains<<endl;
178       //TestPersistantMpi0To1(MyGlobals::_Atomize, ndomains);
179       //TestPersistantMpiRing(MyGlobals::_Atomize, ndomains);
180       TestPersistantMpiRingOnCommSplit(MyGlobals::_Atomize, ndomains);
181       //MPI_Barrier(MPI_COMM_WORLD);
182       MPI_Finalize();
183       return 0;
184       TestVectorOfStringMpi();
185       TestMapOfStringIntMpi();
186       TestMapOfStringVectorOfStringMpi();
187       TestDataArrayMpi();
188       MPI_Finalize();
189       return 1;
190     }
191   
192   if (MyGlobals::_Is0verbose)
193     {
194       cout << "medpartitioner_para V1.0 :" << endl;
195       cout << "  input-file = " << input << endl;
196       cout << "  output-file = " << output << endl;
197       cout << "  split-method = " << library << endl;
198       cout << "  ndomains = " << ndomains << endl;
199       cout << "  create_boundary_faces = " << MyGlobals::_Create_Boundary_Faces << endl;
200       cout << "  dump-cpu-memory = " << mesure_memory<< endl;
201       cout << "  verbose = " << MyGlobals::_Verbose << endl;
202     }
203   //testing whether it is possible to write a file at the specified location
204   if (MyGlobals::_Rank==0)
205     {
206       string outputtest = output + ".testioms.";
207       ofstream testfile (outputtest.c_str());
208       if (testfile.fail())
209         { 
210           cerr << "output-file directory does not exist or is in read-only access" << endl;
211           MPI_Finalize(); return 1;
212         }
213       //deletes test file
214       remove(outputtest.c_str());
215     }
216   
217   // Beginning of the computation
218
219   // Loading the mesh collection
220   if (MyGlobals::_Is0verbose) cout << "Reading input files "<<endl;
221   
222   try
223     {
224       MEDPARTITIONER::ParaDomainSelector parallelizer(mesure_memory);
225       MEDPARTITIONER::MeshCollection collection(input,parallelizer);
226       MEDPARTITIONER::ParallelTopology* aPT = (MEDPARTITIONER::ParallelTopology*) collection.getTopology();
227       aPT->setGlobalNumerotationDefault(collection.getParaDomainSelector());
228       //int nbfiles=MyGlobals::_fileMedNames->size(); //nb domains
229       //to have unique valid fields names/pointers/descriptions for partitionning
230       collection.prepareFieldDescriptions();
231       //int nbfields=collection.getFieldDescriptions().size(); //on all domains
232       //cout<<ReprVectorOfString(collection.getFieldDescriptions());
233     
234       if (MyGlobals::_Is0verbose)
235         {
236           cout<<"fileNames :"<<endl
237               <<ReprVectorOfString(MyGlobals::_File_Names);
238           cout<<"fieldDescriptions :"<<endl
239               <<ReprFieldDescriptions(collection.getFieldDescriptions()," ");
240           cout<<"familyInfo :\n"
241               <<ReprMapOfStringInt(collection.getFamilyInfo())<<endl;
242           cout<<"groupInfo :\n"
243               <<ReprMapOfStringVectorOfString(collection.getGroupInfo())<<endl;
244         }
245     
246       //Creating the graph and partitioning it
247       if (MyGlobals::_Is0verbose) cout << "Computing partition "<< endl;
248   
249       auto_ptr< MEDPARTITIONER::Topology > new_topo;
250       if (library == "metis")
251         new_topo.reset( collection.createPartition(ndomains,MEDPARTITIONER::Graph::METIS));
252       else
253         new_topo.reset( collection.createPartition(ndomains,MEDPARTITIONER::Graph::SCOTCH));
254       parallelizer.evaluateMemory();
255   
256       //Creating a new mesh collection from the partitioning
257       if (MyGlobals::_Is0verbose)  cout << "Creating new meshes"<< endl;
258       MEDPARTITIONER::MeshCollection new_collection(collection,new_topo.get(),split_family,empty_groups);
259       //cout<<"proc "<<MyGlobals::_Rank<<" : new_collection done"<<endl;
260       parallelizer.evaluateMemory();
261   
262       //if (!xml_output_master) new_collection.setDriverType(MEDPARTITIONER::MedAscii);
263       if (filter_face) new_collection.filterFaceOnCell();
264     
265       //to get infos on all procs
266     
267       //see meshName
268       vector<string> finalInformations;
269       vector<string> r1,r2;
270       r1=AllgathervVectorOfString(MyGlobals::_General_Informations);
271       //if (MyGlobals::_Is0verbose>1000) cout << "generalInformations : \n"<<ReprVectorOfString(r1);
272       r2=SelectTagsInVectorOfString(r1,"ioldDomain=");
273       r2=SelectTagsInVectorOfString(r2,"meshName=");
274       if (r2.size()==(collection.getMesh()).size())
275         {
276           for (std::size_t i=0; i<r2.size(); i++)
277             r2[i]=EraseTagSerialized(r2[i],"ioldDomain=");
278           r2=DeleteDuplicatesInVectorOfString(r2);
279           if (r2.size()==1)
280             {
281               string finalMesh="finalMeshName="+ExtractFromDescription(r2[0], "meshName=");
282               finalInformations.push_back(SerializeFromString(finalMesh));
283             }
284         }
285       if (finalInformations.size()==0)
286         {
287           if (MyGlobals::_Rank==0)
288             cerr<<"Problem on final meshName : set at 'Merge'"<<endl;
289           finalInformations.push_back(SerializeFromString("finalMeshName=Merge"));
290         }
291     
292       //see field info nbComponents & componentInfo (if fields present)
293       r2=SelectTagsInVectorOfString(r1,"fieldName=");
294       r2=SelectTagsInVectorOfString(r2,"nbComponents=");
295       //may be yes? or not?
296       for (std::size_t i=0; i<r2.size(); i++)
297         r2[i]=EraseTagSerialized(r2[i],"ioldFieldDouble=");
298       r2=DeleteDuplicatesInVectorOfString(r2);
299       for (std::size_t i=0; i<r2.size(); i++)
300         finalInformations.push_back(r2[i]);
301     
302       MyGlobals::_General_Informations=finalInformations;
303       if (MyGlobals::_Is0verbose) 
304         cout << "generalInformations : \n"<<ReprVectorOfString(finalInformations);
305     
306       //new_collection.setSubdomainBoundaryCreates(create_boundary_faces);
307       if (MyGlobals::_Is0verbose) cout << "Writing "<<ndomains<<" output files "<<output<<"xx.med"<<" and "<<output<<".xml"<<endl;
308       new_collection.write(output);
309   
310       if ( mesure_memory )
311         if ( parallelizer.isOnDifferentHosts() || MyGlobals::_Rank==0 )
312           {
313             cout << "Elapsed time = " << parallelizer.getPassedTime()
314                  << ", max memory usage = " << parallelizer.evaluateMemory() << " KB"
315                  << endl;
316           }
317       if(MyGlobals::_World_Size>1)
318         MPI_Barrier(MPI_COMM_WORLD);
319       if (MyGlobals::_Is0verbose>0) cout<<"OK END"<< endl;
320       MPI_Finalize();
321       return 0;
322     }
323   catch(const char *mess)
324     {
325       cerr<<"proc "<<MyGlobals::_Rank<<" : "<<mess<<endl;
326       fflush(stderr);
327       MPI_Finalize();
328       return 1;
329     }
330   catch(INTERP_KERNEL::Exception& e)
331     {
332       cerr<<"proc "<<MyGlobals::_Rank<<" : INTERP_KERNEL_Exception : "<<e.what()<<endl;
333       fflush(stderr);
334       MPI_Finalize();
335       return 1;
336     }
337   catch(std::exception& e)
338     {
339       cerr<<"proc "<<MyGlobals::_Rank<<" : std_Exception : "<<e.what()<<endl;
340       fflush(stderr);
341       MPI_Finalize();
342       return 1;
343     }
344   catch(...)
345     {
346       cerr<<"proc "<<MyGlobals::_Rank<<" : an unknown type exception error has occurred"<<endl;
347       fflush(stderr);
348       MPI_Finalize();
349       return 1;
350     }
351 #endif
352 }
353