LifeV
cavity.old.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 <life/lifecore/LifeV.hpp>
27 #include <life/lifecore/LifeChrono.hpp>
28 #include <life/lifesolver/OseenSolver.hpp>
29 #include <life/lifemesh/MeshData.hpp>
30 #include <life/lifesolver/OseenData.hpp>
31 #include <life/lifearray/MapEpetra.hpp>
32 #include <life/lifefem/FESpace.hpp>
33 #include <life/lifefem/TimeAdvanceBDFNavierStokes.hpp>
34 #include <life/lifemesh/MeshPartitioner.hpp>
35 #include <life/lifefilters/ExporterEnsight.hpp>
36 #include <life/lifefilters/ExporterHDF5.hpp>
37 #include <life/lifefilters/ExporterEmpty.hpp>
38 //#include <fstream> // To create an output for the flux
39 
40 // Trilinos-MPI communication definitions
41 #include "Epetra_config.h"
42 #ifdef HAVE_MPI
43 #include "Epetra_MpiComm.h"
44 #else
45 #include "Epetra_SerialComm.h"
46 #endif
47 
48 // Object type definitions
53 
54 // +-----------------------------------------------+
55 // | Data and functions for the boundary conditions|
56 // +-----------------------------------------------+
57 const int BACK = 1;
58 const int FRONT = 2;
59 const int LEFT = 3;
60 const int RIGHT = 4;
61 const int BOTTOM = 5;
62 const int TOP = 6;
63 
64 LifeV::Real lidBC (const LifeV::Real& t, const LifeV::Real& /*x*/, const LifeV::Real& /*y*/, const LifeV::Real& /*z*/, const LifeV::ID& i)
65 {
66  switch (i)
67  {
68  case 2:
69  return 1.0;
70  default:
71  return 0.0;
72  }
73 }
74 
75 LifeV::Real fZero ( const LifeV::Real& /* t */,
76  const LifeV::Real& /* x */,
77  const LifeV::Real& /* y */,
78  const LifeV::Real& /* z */,
79  const LifeV::ID& /* i */ )
80 {
81  return 0.0;
82 }
83 
84 int main (int argc, char** argv)
85 {
86 
87  // +-----------------------------------------------+
88  // | Initialization of MPI |
89  // +-----------------------------------------------+
90 #ifdef HAVE_MPI
91  MPI_Init (&argc, &argv);
92 #endif
93 
94  std::shared_ptr<Epetra_Comm> comm;
95 #ifdef EPETRA_MPI
96  comm.reset ( new Epetra_MpiComm ( MPI_COMM_WORLD ) );
97  int nproc;
98  MPI_Comm_size (MPI_COMM_WORLD, &nproc);
99 #else
100  comm.reset ( new Epetra_SerialComm() );
101 #endif
102 
103  bool verbose (false);
104  if (comm->MyPID() == 0)
105  {
106  verbose = true;
107  std::cout
108  << " +-----------------------------------------------+" << std::endl
109  << " | Cavity example for LifeV |" << std::endl
110  << " +-----------------------------------------------+" << std::endl
111  << std::endl
112  << " +-----------------------------------------------+" << std::endl
113  << " | Author: Gwenol Grandperrin |" << std::endl
114  << " | Date: October 25, 2010 |" << std::endl
115  << " +-----------------------------------------------+" << std::endl
116  << std::endl;
117 
118  std::cout << "[Initilization of MPI]" << std::endl;
119 #ifdef HAVE_MPI
120  std::cout << "Using MPI (" << nproc << " proc.)" << std::endl;
121 #else
122  std::cout << "Using serial version" << std::endl;
123 #endif
124  }
125 
126  // Chronometer
127  LifeV::Chrono globalChrono;
128  LifeV::Chrono initChrono;
129  LifeV::Chrono iterChrono;
130  globalChrono.start();
131  initChrono.start();
132 
133  // +-----------------------------------------------+
134  // | Loading the data |
135  // +-----------------------------------------------+
136  if (verbose)
137  {
138  std::cout << std::endl << "[Loading the data]" << std::endl;
139  }
140  GetPot command_line (argc, argv);
141  const std::string dataFileName = command_line.follow ("data", 2, "-f", "--file");
142  GetPot dataFile (dataFileName);
143 
144  // +-----------------------------------------------+
145  // | Loading the mesh |
146  // +-----------------------------------------------+
147  if (verbose)
148  {
149  std::cout << std::endl << "[Loading the mesh]" << std::endl;
150  }
151  LifeV::MeshData meshData;
152  meshData.setup (dataFile, "fluid/space_discretization");
153  if (verbose)
154  {
155  std::cout << "Mesh file: " << meshData.meshDir() << meshData.meshFile() << std::endl;
156  }
157  std::shared_ptr< LifeV::RegionMesh<LifeV::LinearTetra> > fullMeshPtr (new LifeV::RegionMesh<LifeV::LinearTetra>);
158  LifeV::readMesh (*fullMeshPtr, meshData);
159  // Split the mesh between processors
160  LifeV::partitionMesh< LifeV::RegionMesh<LifeV::LinearTetra> > meshPart (fullMeshPtr, comm);
161 
162  // +-----------------------------------------------+
163  // | Creating the FE spaces |
164  // +-----------------------------------------------+
165  if (verbose)
166  {
167  std::cout << std::endl << "[Creating the FE spaces]" << std::endl;
168  }
169  std::string uOrder = dataFile ( "fluid/space_discretization/vel_order", "P2");
170  std::string pOrder = dataFile ( "fluid/space_discretization/press_order", "P1");
171  if (verbose) std::cout << "FE for the velocity: " << uOrder << std::endl
172  << "FE for the pressure: " << pOrder << std::endl;
173 
174  if (verbose)
175  {
176  std::cout << "Building the velocity FE space... " << std::flush;
177  }
178  LifeV::FESpace< LifeV::RegionMesh<LifeV::LinearTetra>, LifeV::MapEpetra > uFESpace (meshPart, uOrder, 3, comm);
179  if (verbose)
180  {
181  std::cout << "ok." << std::endl;
182  }
183 
184  if (verbose)
185  {
186  std::cout << "Building the pressure FE space... " << std::flush;
187  }
188  LifeV::FESpace< LifeV::RegionMesh<LifeV::LinearTetra>, LifeV::MapEpetra > pFESpace (meshPart, pOrder, 1, comm);
189  if (verbose)
190  {
191  std::cout << "ok." << std::endl;
192  }
193 
194  // Total degrees of freedom (elements of matrix)
195  LifeV::UInt totalVelDof = uFESpace.map().getMap (LifeV::Unique)->NumGlobalElements();
196  LifeV::UInt totalPressDof = pFESpace.map().getMap (LifeV::Unique)->NumGlobalElements();
197 
198  if (verbose)
199  {
200  std::cout << "Total Velocity Dof: " << totalVelDof << std::endl;
201  }
202  if (verbose)
203  {
204  std::cout << "Total Pressure Dof: " << totalPressDof << std::endl;
205  }
206 
207  // +-----------------------------------------------+
208  // | Boundary conditions |
209  // +-----------------------------------------------+
210  if (verbose)
211  {
212  std::cout << std::endl << "[Boundary conditions]" << std::endl;
213  }
214 
215  LifeV::BCFunctionBase uZero (fZero);
216  LifeV::BCFunctionBase uLid (lidBC);
217 
218  LifeV::BCHandler::BCHints hint;
219  hint = LifeV::BCHandler::HINT_BC_NONE;
220 
221  std::vector<LifeV::ID> xComp (1);
222  xComp[0] = 1;
223 
224  LifeV::BCHandler bcH ( 0, hint );
225  // A boundary condition in every face
226  bcH.addBC ( "Top" , TOP , LifeV::Essential, LifeV::Full , uLid , 3 );
227  bcH.addBC ( "Left" , LEFT , LifeV::Essential, LifeV::Full , uZero, 3 );
228  bcH.addBC ( "Front" , FRONT , LifeV::Essential, LifeV::Component, uZero, xComp );
229  bcH.addBC ( "Right" , RIGHT , LifeV::Essential, LifeV::Full , uZero, 3 );
230  bcH.addBC ( "Back" , BACK , LifeV::Essential, LifeV::Component, uZero, xComp );
231  bcH.addBC ( "Bottom", BOTTOM, LifeV::Essential, LifeV::Full , uZero, 3 );
232 
233  // Get the number of Lagrange Multiplyers (LM) and set the offsets
234  std::vector<LifeV::bcName_Type> fluxVector = bcH.getBCWithType ( LifeV::Flux );
235  LifeV::UInt numLM = static_cast<LifeV::UInt> ( fluxVector.size() );
236 
237  LifeV::UInt offset = uFESpace.map().getMap (LifeV::Unique)->NumGlobalElements()
238  + pFESpace.map().getMap (LifeV::Unique)->NumGlobalElements();
239 
240  for ( LifeV::UInt i = 0; i < numLM; ++i )
241  {
242  bcH.setOffset ( fluxVector[i], offset + i );
243  }
244 
245  // +-----------------------------------------------+
246  // | Creating the problem |
247  // +-----------------------------------------------+
248  if (verbose)
249  {
250  std::cout << std::endl << "[Creating the problem]" << std::endl;
251  }
252  std::shared_ptr<LifeV::OseenData> oseenData (new LifeV::OseenData() );
253  oseenData->setup ( dataFile );
254 
255  if (verbose)
256  {
257  std::cout << "Time discretization order " << oseenData->dataTime()->orderBDF() << std::endl;
258  }
259 
260  // The problem (matrix and rhs) is packed in an object called fluid
261  LifeV::OseenSolver< LifeV::RegionMesh<LifeV::LinearTetra> > fluid (oseenData,
262  uFESpace,
263  pFESpace,
264  comm,
265  numLM);
266  // Gets inputs from the data file
267  fluid.setUp (dataFile);
268 
269  // Assemble the matrices
270  fluid.buildSystem();
271 
272  // Communication map
273  LifeV::MapEpetra fullMap (fluid.getMap() );
274 
275  // Synchronization
276  MPI_Barrier (MPI_COMM_WORLD);
277 
278  // +-----------------------------------------------+
279  // | Initialization of the simulation |
280  // +-----------------------------------------------+
281  if (verbose)
282  {
283  std::cout << std::endl << "[Initialization of the simulation]" << std::endl;
284  }
285  LifeV::Real dt = oseenData->dataTime()->timeStep();
286  LifeV::Real t0 = oseenData->dataTime()->initialTime();
287  LifeV::Real tFinal = oseenData->dataTime()->endTime();
288 
289  // bdf object to store the previous solutions
290  LifeV::TimeAdvanceBDFNavierStokes<vector_type> bdf (oseenData->dataTime()->orderBDF() );
291 
292  // Initialization with exact solution: either interpolation or "L2-NS"-projection
293  t0 -= dt * bdf.bdfVelocity().order();
294 
295  vector_type beta ( fullMap );
296  vector_type rhs ( fullMap );
297 
298  MPI_Barrier (MPI_COMM_WORLD);
299 
300  // We get the initial solution using a steady Stokes problem
301  oseenData->dataTime()->setTime (t0);
302 
303  beta *= 0.;
304  rhs *= 0.;
305  fluid.updateSystem (0.0, beta, rhs);
306  fluid.iterate (bcH);
307  bdf.bdfVelocity().setInitialCondition ( *fluid.solution() );
308 
309  LifeV::Real time = t0 + dt;
310  for ( ; time <= oseenData->dataTime()->initialTime() + dt / 2.; time += dt)
311  {
312  oseenData->dataTime()->setTime (time);
313 
314  fluid.updateSystem (0.0, beta, rhs);
315  fluid.iterate (bcH);
316  bdf.bdfVelocity().shiftRight ( *fluid.solution() );
317  }
318 
319  // We erase the preconditioner build for Stokes
320  // (A new one should be built for Navier-Stokes)
321  fluid.resetPrec();
322 
323  std::shared_ptr< LifeV::Exporter<LifeV::RegionMesh<LifeV::LinearTetra> > > exporter;
324 
325  vector_ptrtype velAndPressure;
326 
327  std::string const exporterType = dataFile ( "exporter/type", "ensight");
328 
329 #ifdef HAVE_HDF5
330  if (exporterType.compare ("hdf5") == 0)
331  {
332  exporter.reset ( new LifeV::Hdf5exporter<LifeV::RegionMesh<LifeV::LinearTetra> > ( dataFile, "cavity_example" ) );
333  exporter->setDirectory ( "./" ); // This is a test to see if M_post_dir is working
334  exporter->setMeshProcId ( meshPart.meshPartition(), comm->MyPID() );
335  }
336  else
337 #endif
338  {
339  if (exporterType.compare ("none") == 0)
340  {
341  exporter.reset ( new LifeV::NoExport<LifeV::RegionMesh<LifeV::LinearTetra> > ( dataFile, meshPart.meshPartition(), "cavity_example", comm->MyPID() ) );
342  }
343  else
344  {
345  exporter.reset ( new LifeV::Ensight<LifeV::RegionMesh<LifeV::LinearTetra> > ( dataFile, meshPart.meshPartition(), "cavity_example", comm->MyPID() ) );
346  }
347  }
348 
349  velAndPressure.reset ( new vector_type (*fluid.solution(), exporter->mapType() ) );
350 
351  exporter->addVariable ( LifeV::ExporterData::Vector, "velocity", velAndPressure,
352  LifeV::UInt (0), uFESpace.dof().numTotalDof() );
353 
354  exporter->addVariable ( LifeV::ExporterData::Scalar, "pressure", velAndPressure,
355  LifeV::UInt (3 * uFESpace.dof().numTotalDof() ),
356  LifeV::UInt (pFESpace.dof().numTotalDof() ) );
357  exporter->postProcess ( 0 );
358 
359  initChrono.stop();
360  if (verbose)
361  {
362  std::cout << "Initialization time: " << initChrono.diff() << " s." << std::endl;
363  }
364 
365  // +-----------------------------------------------+
366  // | Solving the problem |
367  // +-----------------------------------------------+
368  if (verbose)
369  {
370  std::cout << std::endl << "[Solving the problem]" << std::endl;
371  }
372  int iter = 1;
373 
374  for ( ; time <= tFinal + dt / 2.; time += dt, iter++)
375  {
376 
377  oseenData->dataTime()->setTime (time);
378 
379  if (verbose)
380  {
381  std::cout << "[t = " << oseenData->dataTime()->time() << " s.]" << std::endl;
382  }
383 
384  iterChrono.start();
385 
386  double alpha = bdf.bdfVelocity().coefficientFirstDerivative ( 0 ) / oseenData->dataTime()->timeStep();
387 
388  beta = bdf.bdfVelocity().extrapolation(); // Extrapolation for the convective term
389 
390  bdf.bdfVelocity().updateRHSContribution ( oseenData->dataTime()->timeStep() );
391  rhs = fluid.matrMass() * bdf.bdfVelocity().rhsContributionFirstDerivative();
392  fluid.getDisplayer().leaderPrint ("alpha ", alpha);
393  fluid.getDisplayer().leaderPrint ("\n");
394  fluid.getDisplayer().leaderPrint ("norm beta ", beta.Norm2() );
395  fluid.getDisplayer().leaderPrint ("\n");
396  fluid.getDisplayer().leaderPrint ("norm rhs ", rhs.Norm2() );
397  fluid.getDisplayer().leaderPrint ("\n");
398 
399  fluid.updateSystem ( alpha, beta, rhs );
400  fluid.iterate ( bcH );
401 
402  bdf.bdfVelocity().shiftRight ( *fluid.solution() );
403 
404  // Computation of the error
405  vector_type vel (uFESpace.map(), LifeV::Repeated);
406  vector_type press (pFESpace.map(), LifeV::Repeated);
407  vector_type velpressure ( *fluid.solution(), LifeV::Repeated );
408 
409  velpressure = *fluid.solution();
410  vel.subset (velpressure);
411  press.subset (velpressure, uFESpace.dim() *uFESpace.fieldDim() );
412 
413 
414  bool verbose = (comm->MyPID() == 0);
415 
416 
417  // Exporting the solution
418  *velAndPressure = *fluid.solution();
419  exporter->postProcess ( time );
420 
421 
422  MPI_Barrier (MPI_COMM_WORLD);
423 
424  iterChrono.stop();
425  if (verbose)
426  {
427  std::cout << "Iteration time: " << iterChrono.diff() << " s." << std::endl << std::endl;
428  }
429  }
430 
431  globalChrono.stop();
432  if (verbose)
433  {
434  std::cout << "Total simulation time: " << globalChrono.diff() << " s." << std::endl;
435  }
436 
437  //exporter->CloseFile();
438 
439 #ifdef HAVE_MPI
440  MPI_Finalize();
441 #endif
442 
443  return 0;
444 }
const int BOTTOM
Definition: cavity.old.cpp:61
const int LEFT
Definition: cavity.old.cpp:59
LifeV::Real fZero(const LifeV::Real &, const LifeV::Real &, const LifeV::Real &, const LifeV::Real &, const LifeV::ID &)
Definition: cavity.old.cpp:75
const int FRONT
Definition: cavity.old.cpp:58
fluid_type::vector_type vector_type
Definition: cavity.old.cpp:51
int main(int argc, char **argv)
Definition: dummy.cpp:5
const int RIGHT
Definition: cavity.old.cpp:60
LifeV::Real lidBC(const LifeV::Real &t, const LifeV::Real &, const LifeV::Real &, const LifeV::Real &, const LifeV::ID &i)
Definition: cavity.old.cpp:64
LifeV::OseenSolver< mesh_type > fluid_type
Definition: cavity.old.cpp:50
LifeV::RegionMesh< LifeV::LinearTetra > mesh_type
Definition: cavity.old.cpp:49
std::shared_ptr< vector_type > vector_ptrtype
Definition: cavity.old.cpp:52
const int TOP
Definition: cavity.old.cpp:62
const int BACK
Definition: cavity.old.cpp:57