LifeV
region_marker_id/basic_test.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 /*!
28  @file basic_test.cpp
29  @brief Test the consistency of the mesh data structure
30 
31  @author Alessio Fumagalli <alessio.fumagalli@mail.polimi.it>
32  @contributor
33  @maintainer
34 
35  @date 2012-09-14
36 
37  Colour a mesh with two different colours, count the number
38  of elements equal of one of the two.
39 
40  */
41 
42 // ===================================================
43 //! Includes
44 // ===================================================
45 
46 #include <Epetra_ConfigDefs.h>
47 #ifdef HAVE_MPI
48 #include <mpi.h>
49 #include <Epetra_MpiComm.h>
50 #else
51 #include <Epetra_SerialComm.h>
52 #endif
53 
54 
55 #include <lifev/core/mesh/RegionMesh2DStructured.hpp>
56 #include <lifev/core/mesh/MeshUtility.hpp>
57 
58 #ifdef HAVE_HDF5
59 #include <lifev/core/filter/ExporterHDF5.hpp>
60 #else
61 #include <lifev/core/filter/ExporterEmpty.hpp>
62 #endif
63 
64 using namespace LifeV;
65 
66 UInt colour_fun ( const Vector3D& barycentre )
67 {
68  if ( barycentre[0] < 0.5 && barycentre[1] < 0.5 )
69  {
70  return 2;
71  }
72  return 3;
73 }
74 
75 int main (int argc, char* argv[])
76 {
77 #ifdef HAVE_MPI
78  MPI_Init (&argc, &argv);
79  std::shared_ptr<Epetra_Comm> comm ( new Epetra_MpiComm ( MPI_COMM_WORLD ) );
80 #else
81  std::shared_ptr<Epetra_Comm> comm ( new Epetra_SerialComm );
82 #endif
83 
84  typedef RegionMesh<LinearTriangle> mesh_Type;
85  typedef std::shared_ptr< mesh_Type > meshPtr_Type;
86 
87  // Create the mesh.
88  meshPtr_Type mesh ( new mesh_Type ( comm ) );
89 
90  // Fill the mesh with a structured mesh.
91  regularMesh2D ( *mesh, 0, 8, 11 );
92 
93  // Colour the mesh according to a function.
94  MeshUtility::assignRegionMarkerID ( *mesh, colour_fun );
95 
96  // Count the number of elements with colour 2
97  const UInt colourElements = mesh->elementList().countElementsWithMarkerID ( 2, std::equal_to<markerID_Type>() );
98 
99  // Number of elements with colour 2
100  const UInt exactNumber = 44;
101 
102  {
103  // Needed to correctly destroy the exporterHDF5
104 
105  // Set the exporter for the mesh region.
106 #ifdef HAVE_HDF5
107  ExporterHDF5< mesh_Type > exporter;
108 #else
109  ExporterEmpty< mesh_Type > exporter;
110 #endif
111 
112  // Set the mesh.
113  exporter.setMeshProcId ( mesh, comm->MyPID() );
114 
115  // Export the region marker ID.
116  exporter.exportRegionMarkerID ( mesh, comm );
117 
118  // Do the export.
119  exporter.postProcess ( 0 );
120 
121  } // Needed to correctly destroy the exporterHDF5
122 
123 #ifdef HAVE_MPI
124  MPI_Finalize();
125 #endif
126 
127  if ( colourElements == exactNumber )
128  {
129  return ( EXIT_SUCCESS );
130  }
131  else
132  {
133  return ( EXIT_FAILURE );
134  }
135 }
Real const & operator[](UInt const &i) const
Operator [].
UInt colour_fun(const Vector3D &barycentre)
VectorSmall< 3 > Vector3D
int main(int argc, char **argv)
uint32_type UInt
generic unsigned integer (used mainly for addressing)
Definition: LifeV.hpp:191