LifeV
test_gmsh_parser.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 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 #include <Epetra_ConfigDefs.h>
28 #ifdef HAVE_MPI
29 #include <mpi.h>
30 #include <Epetra_MpiComm.h>
31 #else
32 #include <Epetra_SerialComm.h>
33 #endif
34 
35 #include <boost/shared_ptr.hpp>
36 
37 #include <lifev/core/LifeV.cpp>
38 #include <lifev/core/mesh/MeshPartitioner.hpp>
39 #include <lifev/core/filter/ParserGmsh.hpp>
40 #include <lifev/core/mesh/ConvertBareMesh.hpp>
41 
42 #include <cstdlib>
43 #include <iostream>
44 
46 
47 template <int dim>
49 {
50  template<typename M>
51  static void dopartition (std::shared_ptr<M>& mesh_p, comm_t comm)
52  {
53  mesh_p->updateElementFacets (true, true);
54  mesh_p->updateElementRidges (true, true);
55 
56  LifeV::MeshPartitioner<M> part (mesh_p, comm);
57  }
58 };
59 
60 template <>
61 struct partitioner<2>
62 {
63  template<typename M>
64  static void dopartition (std::shared_ptr<M>&, comm_t)
65  {
66  // TODO: not working yet, problem with boundary ridges (aka points).
67  /*
68  mesh_p->updateElementFacets(true, true);
69 
70  LifeV::MeshPartitioner<M> part(mesh_p, comm);
71  */
72  }
73 };
74 
75 template <>
76 struct partitioner<1>
77 {
78  template<typename M>
79  static void dopartition (std::shared_ptr<M>&, comm_t)
80  {}
81 };
82 
83 
84 template <typename S>
85 struct tester
86 {
87  typedef typename LifeV::BareMesh<S> baremesh_t;
88  typedef typename LifeV::RegionMesh<S> mesh_t;
89 
90  static bool test (std::string& filename, comm_t comm, bool part = true)
91  {
92  baremesh_t baremesh;
93 
94  if (!LifeV::MeshIO::ReadGmshFile (filename, baremesh, 1, true) )
95  {
96  return false;
97  }
98 
99  // Convert the mesh
100  std::shared_ptr<mesh_t> mesh (new mesh_t (comm) );
101  LifeV::convertBareMesh (baremesh, *mesh, true);
102 
103  // Partitioning (if possible)
104  if (part)
105  {
106  partitioner<S::S_nDimensions>::dopartition (mesh, comm);
107  }
108 
109  return true;
110  }
111 };
112 
113 int main (int argc, char* argv[])
114 {
115 
116 #ifdef HAVE_MPI
117  MPI_Init (&argc, &argv);
118  std::shared_ptr<Epetra_Comm> comm (new Epetra_MpiComm (MPI_COMM_WORLD) );
119 #else
120  std::shared_ptr<Epetra_Comm> comm (new Epetra_SerialComm);
121 #endif
122 
123  bool ilead = (comm->MyPID() == 0);
124 
125  // Parse command-line
126  if (argc < 5)
127  {
128  std::cout << "Usage: " << argv[0] << " [mesh_1D_P1] [mesh_2D_Q2] [mesh_3D_P1] [mesh_3D_Q2]" << std::endl;
129  return EXIT_SUCCESS;
130  }
131 
132  // Read the file
133  std::vector<std::string> files;
134  std::copy (argv + 1, argv + 5, std::back_inserter (files) );
135 
136  if (ilead)
137  {
138  std::cout << "\n [[ LINEAR LINE ]] \n\n";
139  }
140  tester<LifeV::LinearLine>::test (files[0], comm);
141 
142  if (ilead)
143  {
144  std::cout << "\n [[ QUADRATIC QUAD ]] \n\n";
145  }
146  tester<LifeV::QuadraticQuad>::test (files[1], comm);
147 
148  if (ilead)
149  {
150  std::cout << "\n [[ LINEAR TETRA ]] \n\n";
151  }
152  tester<LifeV::LinearTetra>::test (files[2], comm);
153 
154  if (ilead)
155  {
156  std::cout << "\n [[ QUADRATIC HEXA ]] \n\n";
157  }
158  // Online partitioning of Q2 Hexa not working yet.
159  tester<LifeV::QuadraticHexa>::test (files[3], comm, false);
160 
161 #ifdef HAVE_MPI
162  MPI_Finalize();
163 #endif
164 
165  return EXIT_SUCCESS;
166 }
static void dopartition(std::shared_ptr< M > &mesh_p, comm_t comm)
std::shared_ptr< Epetra_Comm > comm_t
static bool test(std::string &filename, comm_t comm, bool part=true)
static void dopartition(std::shared_ptr< M > &, comm_t)
int main(int argc, char **argv)
static void dopartition(std::shared_ptr< M > &, comm_t)
elms_t::const_iterator elms_citer