LifeV
TimeIterationPolicyLinear.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 TimeIterationLinear class
29  @brief This class contains all the informations necessary
30  to perform a time iteration
31 
32  @author Gwenol Grandperrin <gwenol.grandperrin@epfl.ch>
33  @date 06-12-2012
34  */
35 
36 #ifndef TIMEITERATIONPOLICYLINEAR_HPP
37 #define TIMEITERATIONPOLICYLINEAR_HPP
38 
39 #include <iostream>
40 #include <string>
41 #include <boost/shared_ptr.hpp>
42 
43 
44 #ifdef HAVE_MPI
45 #include <Epetra_MpiComm.h>
46 #else
47 #include <Epetra_SerialComm.h>
48 #endif
49 
50 #include <Teuchos_ParameterList.hpp>
51 #include <Teuchos_XMLParameterListHelpers.hpp>
52 #include <Teuchos_RCP.hpp>
53 
54 
55 #include <lifev/core/LifeV.hpp>
56 #include <lifev/core/array/MatrixEpetra.hpp>
57 #include <lifev/core/array/VectorEpetra.hpp>
58 #include <lifev/core/util/Displayer.hpp>
59 #include <lifev/core/util/LifeChrono.hpp>
60 #include <lifev/core/mesh/RegionMesh.hpp>
61 #include <lifev/core/fem/FESpace.hpp>
62 #include <lifev/core/fem/TimeAdvanceBDF.hpp>
63 #include <lifev/core/fem/BCHandler.hpp>
64 #include <lifev/navier_stokes/solver/NavierStokesSolver/SolverPolicyLinearSolver.hpp>
65 
66 
67 namespace LifeV
68 {
69 
70 template< class mesh_Type, class AssemblyPolicy, class SolverPolicy = SolverPolicyLinearSolver >
71 struct TimeIterationPolicyLinear : private AssemblyPolicy, public virtual SolverPolicy
72 {
73 protected:
78  typedef MeshPartitioner< mesh_Type > meshPartitioner_Type;
81  typedef FESpace< mesh_Type, map_Type > fespace_Type;
83  typedef TimeAdvanceBDF<vector_Type> bdf_Type;
87 
88  enum { BDFOrder = AssemblyPolicy::BDFOrder };
89 
90  void initTimeIteration ( Teuchos::ParameterList& list );
91 
92  void iterate ( vectorPtr_Type solution,
93  bcContainerPtr_Type bchandler,
94  const Real& currentTime );
95 
96  // Parameters
98 
102 
103  virtual Displayer displayer() = 0;
104  virtual fespacePtr_Type uFESpace() const = 0;
105  virtual fespacePtr_Type pFESpace() const = 0;
106 };
107 
108 template< class mesh_Type, class AssemblyPolicy, class SolverPolicy >
109 void
110 TimeIterationPolicyLinear<mesh_Type, AssemblyPolicy, SolverPolicy>::
111 initTimeIteration ( Teuchos::ParameterList& list )
112 {
113  // Parameters
114  M_computeResidual = list.get ( "Compute exact residual", false );
115 
116  // Initialization
117  M_solutionMap.reset ( new map_Type ( uFESpace()->map() + pFESpace()->map() ) );
118 
119  // Init the assembler
120  Teuchos::ParameterList assemblyList = list.sublist ( "Assembly: Parameter list" );
121  AssemblyPolicy::initAssembly ( assemblyList );
122 
123  // Init the solver
124  Teuchos::ParameterList solverList = list.sublist ( "Solver: Parameter list" );
125  SolverPolicy::initSolver ( solverList );
126  M_rhs.reset ( new vector_Type ( *M_solutionMap, Unique ) );
127 }
128 
129 template< class mesh_Type, class AssemblyPolicy, class SolverPolicy >
130 void
131 TimeIterationPolicyLinear<mesh_Type, AssemblyPolicy, SolverPolicy>::
133  bcContainerPtr_Type bchandler,
134  const Real& currentTime )
135 {
136  Real rhsIterNorm ( 0.0 );
137 
138  //
139  // STEP 1: Updating the system
140  //
141  displayer().leaderPrint ( "Updating the system... " );
142  *M_rhs = 0.0;
143  M_systemMatrix.reset ( new matrix_Type ( *M_solutionMap ) );
144  AssemblyPolicy::assembleSystem ( M_systemMatrix, M_rhs, solution, SolverPolicy::preconditioner() );
145  displayer().leaderPrint ( "done\n" );
146 
147  //
148  // STEP 2: Applying the boundary conditions
149  //
150  displayer().leaderPrint ( "Applying BC... " );
151  bcManage ( *M_systemMatrix, *M_rhs, *uFESpace()->mesh(), uFESpace()->dof(), *bchandler, uFESpace()->feBd(), 1.0, currentTime );
152  M_systemMatrix->globalAssemble();
153  displayer().leaderPrint ( "done\n" );
154 
155  // Extra information if we want to know the exact residual
156  if ( M_computeResidual )
157  {
158  rhsIterNorm = M_rhs->norm2();
159  }
160 
161  //
162  // STEP 3: Solving the system
163  //
164  displayer().leaderPrint ( "Solving the system... \n" );
165  SolverPolicy::solve ( M_systemMatrix, M_rhs, solution );
166 
167  if ( M_computeResidual )
168  {
169  vector_Type Ax ( solution->map() );
170  vector_Type res ( *M_rhs );
171  M_systemMatrix->matrixPtr()->Apply ( solution->epetraVector(), Ax.epetraVector() );
172  res.epetraVector().Update ( -1, Ax.epetraVector(), 1 );
173  Real residual;
174  res.norm2 ( &residual );
175  residual /= rhsIterNorm;
176  displayer().leaderPrint ( "Scaled residual: ", residual, "\n" );
177  }
178 }
179 
180 } // namespace LifeV
181 
182 #endif /* TIMEITERATIONPOLICYLINEAR_HPP */
VectorEpetra - The Epetra Vector format Wrapper.
std::shared_ptr< matrix_Type > matrixPtr_Type
std::shared_ptr< VectorEpetra > vectorPtr_Type
virtual Displayer displayer()=0
void iterate(vectorPtr_Type solution, bcContainerPtr_Type bchandler, const Real &currentTime)
void importFromHDF5(std::string const &fileName, std::string const &matrixName="matrix")
Read a matrix from a HDF5 (.h5) file.
BCHandler - class for handling boundary conditions.
Definition: BCHandler.hpp:100
virtual fespacePtr_Type uFESpace() const =0
void updateInverseJacobian(const UInt &iQuadPt)
Epetra_Import const & importer()
Getter for the Epetra_Import.
Definition: MapEpetra.cpp:394
void norm2(Real *result) const
Compute and store the norm 2 in the given pointed variable.
std::shared_ptr< fespace_Type > fespacePtr_Type
void initTimeIteration(Teuchos::ParameterList &list)
std::shared_ptr< bcContainer_Type > bcContainerPtr_Type
double Real
Generic real data.
Definition: LifeV.hpp:175
FESpace< mesh_Type, map_Type > fespace_Type
TimeAdvanceBDF< vector_Type > bdf_Type
MeshPartitioner< mesh_Type > meshPartitioner_Type
virtual fespacePtr_Type pFESpace() const =0
Displayer - This class is used to display messages in parallel simulations.
Definition: Displayer.hpp:62