LifeV
core/testsuite/offline_partition_io/main_write.cpp
Go to the documentation of this file.
1 //@HEADER
2 /*
3 *******************************************************************************
4 
5 Copyright (C) 2004, 2005, 2007 EPFL, Politecnico di Milano, INRIA
6 Copyright (C) 2010, 2011, 2012 EPFL, Politecnico di Milano, Emory University
7 
8 This file is part of LifeV.
9 
10 LifeV is free software; you can redistribute it and/or modify
11 it under the terms of the GNU Lesser General Public License as published by
12 the Free Software Foundation, either version 3 of the License, or
13 (at your option) any later version.
14 
15 LifeV is distributed in the hope that it will be useful,
16 but WITHOUT ANY WARRANTY; without even the implied warranty of
17 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
18 Lesser General Public License for more details.
19 
20 You should have received a copy of the GNU Lesser General Public License
21 along with LifeV. If not, see <http://www.gnu.org/licenses/>.
22 
23 *******************************************************************************
24 */
25 //@HEADER
26 
27 /*!
28  @file
29  @brief Test for PartitionIO class - cut and write
30 
31  @author Radu Popescu <radu.popescu@epfl.ch>
32  @maintainer Radu Popescu <radu.popescu@epfl.ch>
33  @date 10-05-2012
34 
35  Partition a mesh using a single (MPI) process and save mesh parts
36  to an HDF5 file.
37  */
38 
39 #include <lifev/core/LifeV.hpp>
40 
41 #ifdef LIFEV_HAS_HDF5
42 
43 
44 #include "Epetra_config.h"
45 
46 #ifdef HAVE_MPI
47 
48 #include <Epetra_MpiComm.h>
49 
50 
51 #include <lifev/core/filter/GetPot.hpp>
52 #include <lifev/core/mesh/MeshPartitioner.hpp>
53 #include <lifev/core/mesh/MeshPartitionTool.hpp>
54 #include <lifev/core/filter/PartitionIO.hpp>
55 #include <lifev/core/mesh/RegionMesh3DStructured.hpp>
56 
57 using namespace LifeV;
58 
59 #endif /* HAVE_MPI */
60 #endif /* LIFEV_HAS_HDF5 */
61 
62 typedef MeshPartitionTool < RegionMesh<LinearTetra> > meshCutter_Type;
63 
64 int main (int argc, char** argv)
65 {
66 #ifdef LIFEV_HAS_HDF5
67 #ifdef HAVE_MPI
68 
69  typedef RegionMesh<LinearTetra> mesh_Type;
70  typedef std::shared_ptr<mesh_Type> meshPtr_Type;
71  typedef std::vector<meshPtr_Type> meshParts_Type;
72  typedef std::shared_ptr<meshParts_Type> meshPartsPtr_Type;
73 
74  MPI_Init (&argc, &argv);
75  std::shared_ptr<Epetra_Comm> comm (new Epetra_MpiComm (MPI_COMM_WORLD) );
76 
77  if (comm->NumProc() != 1)
78  {
79  std::cout << "This test needs to be run "
80  << "with a single process. Aborting."
81  << std::endl;
82  return EXIT_FAILURE;
83  }
84 
85  GetPot cl (argc, argv);
86  const UInt numElements = cl.follow (9, "--num-elem");
87  const Int numParts = cl.follow (3, "--num-parts");
88  // partitionerType should be MeshPartitioner, MeshPartitionTool_ParMETIS or
89  // MeshPartitionTool_Zoltan
90  const std::string partitionerType = cl.follow ("MeshPartitioner",
91  "--partitioner-type");
92  std::string partsFile;
93  partsFile.reserve (50);
94  partsFile += "cube_";
95  partsFile += partitionerType;
96  partsFile += ".h5";
97 
98  std::cout << "Number of elements in mesh: " << numElements << std::endl;
99  std::cout << "Number of parts: " << numParts << std::endl;
100  std::cout << "Mesh partitioner type: " << partitionerType << std::endl;
101  std::cout << "Name of HDF5 container: " << partsFile << std::endl;
102 
103  meshPtr_Type fullMeshPtr (new mesh_Type ( comm ) );
104  meshPartsPtr_Type meshPartPtr;
105  regularMesh3D (*fullMeshPtr, 1, numElements, numElements, numElements,
106  false, 2.0, 2.0, 2.0, -1.0, -1.0, -1.0);
107 
108  if (partitionerType == "MeshPartitioner")
109  {
110  // Using old MeshPartitioner class
111  MeshPartitioner<mesh_Type> meshCutter;
112  meshCutter.setup (numParts, comm);
113  meshCutter.attachUnpartitionedMesh (fullMeshPtr);
114  meshCutter.doPartitionGraph();
115  meshCutter.fillEntityPID();
116  meshCutter.doPartitionMesh();
117  meshCutter.releaseUnpartitionedMesh();
118  meshPartPtr = meshCutter.meshPartitions();
119  }
120  else if (partitionerType == "MeshPartitionTool_ParMETIS")
121  {
122  // Using new MeshPartitionTool class with ParMETIS
123  Teuchos::ParameterList meshParameters;
124  meshParameters.set ("num-parts", numParts, "");
125  meshParameters.set ("offline-mode", true, "");
126  meshParameters.set ("graph-lib", "parmetis", "");
127  meshCutter_Type meshCutter (fullMeshPtr, comm, meshParameters);
128  if (! meshCutter.success() )
129  {
130  std::cout << "Mesh partition failed.";
131  return EXIT_FAILURE;
132  }
133  meshPartPtr = meshCutter.allMeshParts();
134  }
135  else if (partitionerType == "MeshPartitionTool_Zoltan")
136  {
137  // Using new MeshPartitionTool class with Zoltan
138  Teuchos::ParameterList meshParameters;
139  meshParameters.set ("num-parts", numParts, "");
140  meshParameters.set ("offline-mode", true, "");
141  meshParameters.set ("graph-lib", "zoltan", "");
142  meshCutter_Type meshCutter (fullMeshPtr, comm, meshParameters);
143  if (! meshCutter.success() )
144  {
145  std::cout << "Mesh partition failed.";
146  return EXIT_FAILURE;
147  }
148  meshPartPtr = meshCutter.allMeshParts();
149  }
150 
151  // delete the RegionMesh object
152  fullMeshPtr.reset();
153 
154  // Write mesh parts to HDF5 container
155  {
156  std::shared_ptr<Epetra_MpiComm> mpiComm =
157  std::dynamic_pointer_cast<Epetra_MpiComm> (comm);
158  PartitionIO<mesh_Type> partitionIO (partsFile, mpiComm);
159  partitionIO.write (meshPartPtr);
160  }
161 
162  MPI_Finalize();
163 
164 #else
165  std::cout << "This test needs MPI to run. Aborting." << std::endl;
166  return (EXIT_FAILURE);
167 #endif /* HAVE_MPI */
168 #else
169  std::cout << "This test needs HDF5 to run. Aborting." << std::endl;
170  return (EXIT_FAILURE);
171 #endif /* LIFEV_HAS_HDF5 */
172 
173  return EXIT_SUCCESS;
174 }
int main(int argc, char **argv)
Definition: dummy.cpp:5