LifeV
cavity_stokes.cpp
Go to the documentation of this file.
1 /* -*- mode: c++ -*-
2  This program is part of the LifeV library
3  Copyright (C) 2001,2002,2003,2004 EPFL, INRIA, Politechnico di Milano
4 
5  This program is free software; you can redistribute it and/or
6  modify it under the terms of the GNU General Public License
7  as published by the Free Software Foundation; either version 2
8  of the License, or (at your option) any later version.
9 
10  This program is distributed in the hope that it will be useful,
11  but WITHOUT ANY WARRANTY; without even the implied warranty of
12  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13  GNU General Public License for more details.
14 
15  You should have received a copy of the GNU General Public License
16  along with this program; if not, write to the Free Software
17  Foundation, Inc., 59 Temple Place - Suite 330, Boston, MA 02111-1307, USA.
18 */
19 /**
20  \file cavity.cpp
21  \author Gwenol Grandperrin <gwenol.grandperrin@epfl.ch>
22  \date 2010-10-25
23  */
24 
25 // LifeV definition files
26 #include <lifev/core/LifeV.hpp>
27 #include <lifev/core/util/LifeChrono.hpp>
28 #include <lifev/core/array/MapEpetra.hpp>
29 #include <lifev/core/mesh/MeshData.hpp>
30 #include <lifev/core/mesh/MeshPartitioner.hpp>
31 #include <lifev/core/fem/FESpace.hpp>
32 #include <lifev/navier_stokes/solver/OseenData.hpp>
33 #include <lifev/navier_stokes/solver/OseenSolver.hpp>
34 #include <lifev/core/filter/ExporterEnsight.hpp>
35 #include <lifev/core/filter/ExporterHDF5.hpp>
36 #include <lifev/core/filter/ExporterEmpty.hpp>
37 
38 
39 
40 //#include <fstream> // To create an output for the flux
41 
42 // Trilinos-MPI communication definitions
43 #include <Epetra_ConfigDefs.h>
44 #ifdef HAVE_MPI
45 #include "Epetra_MpiComm.h"
46 #else
47 #include "Epetra_SerialComm.h"
48 #endif
49 
50 
51 using namespace LifeV;
52 
53 
54 // Object type definitions
61 
62 // +-----------------------------------------------+
63 // | Data and functions for the boundary conditions|
64 // +-----------------------------------------------+
65 const int BACK = 1;
66 const int FRONT = 2;
67 const int LEFT = 3;
68 const int RIGHT = 4;
69 const int BOTTOM = 5;
70 const int TOP = 6;
71 
72 Real lidBC (const Real& /*t*/,
73  const Real& /*x*/,
74  const Real& /*y*/,
75  const Real& /*z*/,
76  const ID& i)
77 {
78  switch (i)
79  {
80  case 1:
81  return 1.0;
82  default:
83  return 0.0;
84  }
85 }
86 
87 Real fZero ( const Real& /* t */,
88  const Real& /* x */,
89  const Real& /* y */,
90  const Real& /* z */,
91  const ID& /* i */ )
92 {
93  return 0.0;
94 }
95 
96 int main (int argc, char** argv)
97 {
98 
99  // +-----------------------------------------------+
100  // | Initialization of MPI |
101  // +-----------------------------------------------+
102 #ifdef HAVE_MPI
103  MPI_Init (&argc, &argv);
104 #endif
105 
106  std::shared_ptr<Epetra_Comm> comm;
107 #ifdef EPETRA_MPI
108  comm.reset ( new Epetra_MpiComm ( MPI_COMM_WORLD ) );
109  int nproc;
110  MPI_Comm_size (MPI_COMM_WORLD, &nproc);
111 #else
112  comm.reset ( new Epetra_SerialComm() );
113 #endif
114 
115  bool verbose (false);
116  if (comm->MyPID() == 0)
117  {
118  verbose = true;
119  std::cout
120  << " +-----------------------------------------------+" << std::endl
121  << " | Cavity example for LifeV |" << std::endl
122  << " +-----------------------------------------------+" << std::endl
123  << std::endl
124  << " +-----------------------------------------------+" << std::endl
125  << " | Author: Gwenol Grandperrin |" << std::endl
126  << " | Date: October 25, 2010 |" << std::endl
127  << " +-----------------------------------------------+" << std::endl
128  << std::endl;
129 
130  std::cout << "[Initilization of MPI]" << std::endl;
131 #ifdef HAVE_MPI
132  std::cout << "Using MPI (" << nproc << " proc.)" << std::endl;
133 #else
134  std::cout << "Using serial version" << std::endl;
135 #endif
136  }
137 
138  // Chronometer
139  LifeChrono globalChrono;
140  globalChrono.start();
141 
142  // +-----------------------------------------------+
143  // | Loading the data |
144  // +-----------------------------------------------+
145  if (verbose)
146  {
147  std::cout << std::endl << "[Loading the data]" << std::endl;
148  }
149  GetPot command_line (argc, argv);
150  const std::string dataFileName = command_line.follow ("data", 2, "-f", "--file");
151  GetPot dataFile (dataFileName);
152 
153  // +-----------------------------------------------+
154  // | Loading the mesh |
155  // +-----------------------------------------------+
156  if (verbose)
157  {
158  std::cout << std::endl << "[Loading the mesh]" << std::endl;
159  }
160  MeshData meshData;
161  meshData.setup (dataFile, "fluid/space_discretization");
162  if (verbose)
163  {
164  std::cout << "Mesh file: " << meshData.meshDir() << meshData.meshFile() << std::endl;
165  }
166  std::shared_ptr< mesh_Type > fullMeshPtr (new mesh_Type);
167  readMesh (*fullMeshPtr, meshData);
168  // Split the mesh between processors
169  MeshPartitioner< mesh_Type > meshPart (fullMeshPtr, comm);
170 
171  // +-----------------------------------------------+
172  // | Creating the FE spaces |
173  // +-----------------------------------------------+
174  if (verbose)
175  {
176  std::cout << std::endl << "[Creating the FE spaces]" << std::endl;
177  }
178  std::string uOrder = dataFile ( "fluid/space_discretization/vel_order", "P2");
179  std::string pOrder = dataFile ( "fluid/space_discretization/press_order", "P1");
180  if (verbose) std::cout << "FE for the velocity: " << uOrder << std::endl
181  << "FE for the pressure: " << pOrder << std::endl;
182 
183  if (verbose)
184  {
185  std::cout << "Building the velocity FE space... " << std::flush;
186  }
187  feSpacePtr_Type uFESpacePtr ( new feSpace_Type (meshPart.meshPartition(), uOrder, 3, comm) );
188  if (verbose)
189  {
190  std::cout << "ok." << std::endl;
191  }
192 
193  if (verbose)
194  {
195  std::cout << "Building the pressure FE space... " << std::flush;
196  }
197  feSpacePtr_Type pFESpacePtr ( new feSpace_Type (meshPart.meshPartition(), pOrder, 1, comm) );
198  if (verbose)
199  {
200  std::cout << "ok." << std::endl;
201  }
202 
203  // Total degrees of freedom (elements of matrix)
204  UInt totalVelDof = uFESpacePtr->map().map (Unique)->NumGlobalElements();
205  UInt totalPressDof = pFESpacePtr->map().map (Unique)->NumGlobalElements();
206 
207  if (verbose)
208  {
209  std::cout << "Total Velocity Dof: " << totalVelDof << std::endl;
210  }
211  if (verbose)
212  {
213  std::cout << "Total Pressure Dof: " << totalPressDof << std::endl;
214  }
215 
216  // +-----------------------------------------------+
217  // | Boundary conditions |
218  // +-----------------------------------------------+
219  if (verbose)
220  {
221  std::cout << std::endl << "[Boundary conditions]" << std::endl;
222  }
223 
224  BCFunctionBase uZero (fZero);
225  BCFunctionBase uLid (lidBC);
226 
227  std::vector<ID> xComp (1);
228  xComp[0] = 1;
229 
230  BCHandler bcH;
231  // A boundary condition in every face
232  bcH.addBC ( "Top" , TOP , Essential, Full , uLid , 3 );
233  bcH.addBC ( "Left" , LEFT , Essential, Full , uZero, 3 );
234  bcH.addBC ( "Front" , FRONT , Essential, Component, uZero, xComp );
235  bcH.addBC ( "Right" , RIGHT , Essential, Full , uZero, 3 );
236  bcH.addBC ( "Back" , BACK , Essential, Component, uZero, xComp );
237  bcH.addBC ( "Bottom", BOTTOM, Essential, Full , uZero, 3 );
238 
239  // Get the number of Lagrange Multiplyers (LM) and set the offsets
240  std::vector<bcName_Type> fluxVector = bcH.findAllBCWithType ( Flux );
241  UInt numLM = static_cast<UInt> ( fluxVector.size() );
242 
243  UInt offset = uFESpacePtr->map().map (Unique)->NumGlobalElements()
244  + pFESpacePtr->map().map (Unique)->NumGlobalElements();
245 
246  for ( UInt i = 0; i < numLM; ++i )
247  {
248  bcH.setOffset ( fluxVector[i], offset + i );
249  }
250 
251  // +-----------------------------------------------+
252  // | Creating the problem |
253  // +-----------------------------------------------+
254  if (verbose)
255  {
256  std::cout << std::endl << "[Creating the problem]" << std::endl;
257  }
258  std::shared_ptr< OseenData > oseenData (new OseenData);
259  oseenData->setup ( dataFile );
260 
261  // The problem (matrix and rhs) is packed in an object called fluid
262  OseenSolver< mesh_Type > fluid (oseenData,
263  *uFESpacePtr,
264  *pFESpacePtr,
265  comm,
266  numLM);
267 
268  // Gets inputs from the data file
269  fluid.setUp (dataFile);
270 
271  // Assemble the matrices
272  fluid.buildSystem();
273 
274  // Communication map
275  MapEpetra fullMap (fluid.getMap() );
276 
277  vector_Type beta ( fullMap );
278  vector_Type rhs ( fullMap );
279  double alpha = 0.;
280  beta *= 0.;
281  rhs *= 0.;
282  fluid.updateSystem (alpha, beta, rhs);
283  fluid.iterate (bcH);
284 
285  std::shared_ptr< ExporterHDF5<mesh_Type> > exporter;
286 
287  std::string const exporterType = dataFile ( "exporter/type", "ensight");
288 
289  exporter.reset ( new ExporterHDF5<mesh_Type> ( dataFile, "cavity_stokes_example" ) );
290  exporter->setPostDir ( "./" ); // This is a test to see if M_post_dir is working
291  exporter->setMeshProcId ( meshPart.meshPartition(), comm->MyPID() );
292 
293  vectorPtr_Type velAndPressure;
294  velAndPressure.reset ( new vector_Type (*fluid.solution(), exporter->mapType() ) );
295 
296  exporter->addVariable ( ExporterData<mesh_Type>::VectorField, "velocity", uFESpacePtr,
297  velAndPressure, UInt (0) );
298 
299  exporter->addVariable ( ExporterData<mesh_Type>::ScalarField, "pressure", pFESpacePtr,
300  velAndPressure, UInt (3 * uFESpacePtr->dof().numTotalDof() ) );
301  exporter->postProcess ( 0 );
302 
303  globalChrono.stop();
304  if (verbose)
305  {
306  std::cout << "Total simulation time: " << globalChrono.diff() << " s." << std::endl;
307  }
308 
309  exporter->closeFile();
310 
311 #ifdef HAVE_MPI
312  MPI_Finalize();
313 #endif
314 
315  return 0;
316 }
const int RIGHT
Real fZero(const Real &, const Real &, const Real &, const Real &, const ID &)
std::shared_ptr< vector_Type > vectorPtr_Type
void start()
Start the timer.
Definition: LifeChrono.hpp:93
FESpace - Short description here please!
Definition: FESpace.hpp:78
RegionMesh< LinearTetra > mesh_Type
const int FRONT
GetPot(const int argc_, char **argv_, const char *FieldSeparator=0x0)
Definition: GetPot.hpp:507
fluid_Type::vector_Type vector_Type
const int BACK
Real lidBC(const Real &, const Real &, const Real &, const Real &, const ID &i)
const int LEFT
OseenSolver< mesh_Type > fluid_Type
const int TOP
Epetra_Import const & importer()
Getter for the Epetra_Import.
Definition: MapEpetra.cpp:394
uint32_type ID
IDs.
Definition: LifeV.hpp:194
const int BOTTOM
FESpace< mesh_Type, MapEpetra > feSpace_Type
std::shared_ptr< std::vector< Int > > M_isOnProc
MeshData - class for handling spatial discretization.
Definition: MeshData.hpp:72
double Real
Generic real data.
Definition: LifeV.hpp:175
void stop()
Stop the timer.
Definition: LifeChrono.hpp:100
std::shared_ptr< feSpace_Type > feSpacePtr_Type
int main(int argc, char **argv)
uint32_type UInt
generic unsigned integer (used mainly for addressing)
Definition: LifeV.hpp:191