LifeV
TimeIterationPolicyNonlinearIncremental.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 TimeIterationNonlinear 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 TIMEITERATIONPOLICYNONLINEARINCREMENTAL_HPP
37 #define TIMEITERATIONPOLICYNONLINEARINCREMENTAL_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 TimeIterationPolicyNonlinearIncremental : 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
99 
104 
105  virtual Displayer displayer() = 0;
106  virtual fespacePtr_Type uFESpace() const = 0;
107  virtual fespacePtr_Type pFESpace() const = 0;
108 };
109 
110 template< class mesh_Type, class AssemblyPolicy, class SolverPolicy >
111 void
112 TimeIterationPolicyNonlinearIncremental<mesh_Type, AssemblyPolicy, SolverPolicy>::
113 initTimeIteration ( Teuchos::ParameterList& list )
114 {
115  // Parameters
116  M_computeResidual = list.get ( "Compute exact residual", false );
117  M_nonLinearTolerance = list.get ( "Nonlinear tolerance", 1e-6 );
118 
119  // Initialization
120  M_solutionMap.reset ( new map_Type ( uFESpace()->map() + pFESpace()->map() ) );
121 
122  // Init the assembler
123  Teuchos::ParameterList assemblyList = list.sublist ( "Assembly: Parameter list" );
124  AssemblyPolicy::initAssembly ( assemblyList );
125 
126  // Init the solver
127  Teuchos::ParameterList solverList = list.sublist ( "Solver: Parameter list" );
128  SolverPolicy::initSolver ( solverList );
129  M_rhs.reset ( new vector_Type ( *M_solutionMap, Unique ) );
130  M_deltaSolution.reset ( new vector_Type ( *M_solutionMap, Unique ) );
131 }
132 
133 template< class mesh_Type, class AssemblyPolicy, class SolverPolicy >
134 void
135 TimeIterationPolicyNonlinearIncremental<mesh_Type, AssemblyPolicy, SolverPolicy>::
137  bcContainerPtr_Type bchandler,
138  const Real& currentTime )
139 {
140  int subiter = 0;
141 
142  Real normRhs ( 0.0 );
143  Real nonLinearResidual ( 0.0 );
144  Real rhsIterNorm ( 0.0 );
145 
146  do
147  {
148  //
149  // STEP 1: Updating the system
150  //
151  displayer().leaderPrint ( "Updating the system... " );
152  *M_rhs = 0.0;
153  M_systemMatrix.reset ( new matrix_Type ( *M_solutionMap ) );
154  AssemblyPolicy::assembleSystem ( M_systemMatrix, M_rhs, solution, SolverPolicy::preconditioner() );
155  displayer().leaderPrint ( "done\n" );
156 
157  //
158  // STEP 2: Applying the boundary conditions
159  //
160  displayer().leaderPrint ( "Applying BC... " );
161  bcManage ( *M_systemMatrix, *M_rhs, *uFESpace()->mesh(), uFESpace()->dof(), *bchandler, uFESpace()->feBd(), 1.0, currentTime );
162  M_systemMatrix->globalAssemble();
163  displayer().leaderPrint ( "done\n" );
164 
165  // Norm of the rhs needed for the nonlinear convergence test
166  if ( subiter == 0 )
167  {
168  normRhs = M_rhs->norm2();
169  }
170 
171  //
172  // STEP 3: Computing the residual
173  //
174 
175  // Computing the RHS as RHS=b-Ax_k
176  vector_Type Ax ( solution->map() );
177  M_systemMatrix->matrixPtr()->Apply ( solution->epetraVector(), Ax.epetraVector() );
178 
179  M_rhs->epetraVector().Update ( -1, Ax.epetraVector(), 1 );
180  nonLinearResidual = M_rhs->norm2();
181 
182  displayer().leaderPrint ( "Nonlinear residual : ", nonLinearResidual, "\n" );
183  displayer().leaderPrint ( "Nonlinear residual (scaled) : ", nonLinearResidual / normRhs, "\n" );
184 
185  if ( nonLinearResidual > M_nonLinearTolerance * normRhs )
186  {
187  displayer().leaderPrint ( "---\nSubiteration [", ++subiter, "]\n" );
188 
189  // Extra information if we want to know the exact residual
190  if ( M_computeResidual )
191  {
192  rhsIterNorm = M_rhs->norm2();
193  }
194 
195  //
196  // Solving the system
197  //
198  displayer().leaderPrint ( "Solving the system... \n" );
199  *M_deltaSolution = 0.0;
200  SolverPolicy::solve ( M_systemMatrix, M_rhs, M_deltaSolution );
201  // int numIter = SolverPolicy::solve( M_systemMatrix, M_rhs, M_deltaSolution );
202  // numIterSum += numIter; //
203 
204  if ( M_computeResidual )
205  {
206  vector_Type Ax ( M_deltaSolution->map() );
207  vector_Type res ( *M_rhs );
208  M_systemMatrix->matrixPtr()->Apply ( M_deltaSolution->epetraVector(), Ax.epetraVector() );
209  res.epetraVector().Update ( -1, Ax.epetraVector(), 1 );
210  Real residual;
211  res.norm2 ( &residual );
212  residual /= rhsIterNorm;
213  displayer().leaderPrint ( "Scaled residual: ", residual, "\n" );
214  }
215 
216  // Update the solution
217  *solution += *M_deltaSolution;
218  }
219  }
220  while ( nonLinearResidual > M_nonLinearTolerance * normRhs );
221 
222  displayer().leaderPrint ( "Nonlinear iterations : ", subiter, "\n" );
223 }
224 
225 } // namespace LifeV
226 
227 #endif /* TIMEITERATIONPOLICYNONLINEARINCREMENTAL_HPP */
VectorEpetra - The Epetra Vector format Wrapper.
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
void updateInverseJacobian(const UInt &iQuadPt)
Epetra_Import const & importer()
Getter for the Epetra_Import.
Definition: MapEpetra.cpp:394
virtual fespacePtr_Type pFESpace() const =0
void norm2(Real *result) const
Compute and store the norm 2 in the given pointed variable.
double Real
Generic real data.
Definition: LifeV.hpp:175
virtual fespacePtr_Type uFESpace() const =0
Displayer - This class is used to display messages in parallel simulations.
Definition: Displayer.hpp:62