LifeV
InitPolicyInterpolation.hpp
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 InitPolicyInterpolation class
29  @brief This class is a strategy to initialize a Navier-Stokes problem
30 
31  @author Gwenol Grandperrin <gwenol.grandperrin@epfl.ch>
32  @date 11-12-2012
33  */
34 
35 #ifndef INITPOLICYINTERPOLATION_HPP
36 #define INITPOLICYINTERPOLATION_HPP
37 
38 #include <iostream>
39 #include <string>
40 #include <boost/shared_ptr.hpp>
41 
42 
43 #ifdef HAVE_MPI
44 #include <Epetra_MpiComm.h>
45 #else
46 #include <Epetra_SerialComm.h>
47 #endif
48 
49 #include <Teuchos_ParameterList.hpp>
50 #include <Teuchos_XMLParameterListHelpers.hpp>
51 #include <Teuchos_RCP.hpp>
52 
53 
54 #include <lifev/core/LifeV.hpp>
55 #include <lifev/core/util/LifeChrono.hpp>
56 #include <lifev/core/array/VectorEpetra.hpp>
57 #include <lifev/core/mesh/RegionMesh.hpp>
58 #include <lifev/core/fem/FESpace.hpp>
59 #include <lifev/core/fem/TimeAdvanceBDF.hpp>
60 #include <lifev/core/fem/BCHandler.hpp>
61 #include <lifev/navier_stokes/solver/NavierStokesSolver/NavierStokesProblem.hpp>
62 
63 namespace LifeV
64 {
65 
66 template< class mesh_Type >
68 {
73  typedef FESpace< mesh_Type, map_Type > fespace_Type;
75  typedef TimeAdvanceBDF<vector_Type> bdf_Type;
80 
81  void setupInit ( Teuchos::ParameterList& list );
82 
83  void initSimulation ( bcContainerPtr_Type bchandler,
84  vectorPtr_Type solution );
85 
86  virtual fespacePtr_Type uFESpace() const = 0;
87  virtual fespacePtr_Type pFESpace() const = 0;
88  virtual Real timestep() const = 0;
89  virtual Real initialTime() const = 0;
90  virtual bdfPtr_Type bdf() const = 0;
91  virtual NSProblemPtr_Type problem() const = 0;
92 };
93 
94 template< class mesh_Type >
95 void
96 InitPolicyInterpolation< mesh_Type >::
98 {
99 
100 }
101 
102 template< class mesh_Type >
103 void
104 InitPolicyInterpolation< mesh_Type >::
106  vectorPtr_Type solution )
107 {
108  ASSERT ( problem()->hasExactSolution(), "Error: You cannot use the interpolation method if the problem has not an analytical solution." );
109 
110  Real currentTime = initialTime() - timestep() * bdf()->order();
111  UInt pressureOffset = uFESpace()->fieldDim() * uFESpace()->dof().numTotalDof();
112 
113  vectorPtr_Type velocity;
114  velocity.reset ( new vector_Type ( uFESpace()->map(), Unique ) );
115 
116  vectorPtr_Type pressure;
117  pressure.reset ( new vector_Type ( pFESpace()->map(), Unique ) );
118 
119  *solution = 0.0;
120  uFESpace()->interpolate ( problem()->uexact(), *velocity, currentTime );
121  pFESpace()->interpolate ( problem()->pexact(), *pressure, currentTime );
122  solution->add ( *velocity );
123  solution->add ( *pressure, pressureOffset );
124  bdf()->setInitialCondition ( *solution );
125 
126  currentTime += timestep();
127  for ( ; currentTime <= initialTime() + timestep() / 2.0; currentTime += timestep() )
128  {
129  *solution = 0.0;
130  *velocity = 0.0;
131  *pressure = 0.0;
132 
133  uFESpace()->interpolate ( problem()->uexact(), *velocity, currentTime );
134  pFESpace()->interpolate ( problem()->pexact(), *pressure, currentTime );
135  solution->add ( *velocity );
136  solution->add ( *pressure, pressureOffset );
137 
138  // Updating bdf
139  bdf()->shiftRight ( *solution );
140  }
141 }
142 
143 } // namespace LifeV
144 
145 #endif /* INITPOLICYINTERPOLATION_HPP */
void initSimulation(bcContainerPtr_Type bchandler, vectorPtr_Type solution)
VectorEpetra - The Epetra Vector format Wrapper.
virtual fespacePtr_Type uFESpace() const =0
std::shared_ptr< fespace_Type > fespacePtr_Type
std::shared_ptr< map_Type > mapPtr_Type
std::shared_ptr< VectorEpetra > vectorPtr_Type
BCHandler - class for handling boundary conditions.
Definition: BCHandler.hpp:100
std::shared_ptr< NavierStokesProblem< mesh_Type > > NSProblemPtr_Type
std::shared_ptr< bdf_Type > bdfPtr_Type
virtual Real initialTime() const =0
virtual fespacePtr_Type pFESpace() const =0
virtual Real timestep() const =0
void updateInverseJacobian(const UInt &iQuadPt)
FESpace< mesh_Type, map_Type > fespace_Type
virtual NSProblemPtr_Type problem() const =0
Epetra_Import const & importer()
Getter for the Epetra_Import.
Definition: MapEpetra.cpp:394
#define ASSERT(X, A)
Definition: LifeAssert.hpp:90
double Real
Generic real data.
Definition: LifeV.hpp:175
TimeAdvanceBDF< vector_Type > bdf_Type
void setupInit(Teuchos::ParameterList &list)
virtual bdfPtr_Type bdf() const =0
std::shared_ptr< bcContainer_Type > bcContainerPtr_Type
uint32_type UInt
generic unsigned integer (used mainly for addressing)
Definition: LifeV.hpp:191