LifeV
InitPolicySolver.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 InitPolicySolver 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 INITPOLICYSOLVER_HPP
36 #define INITPOLICYSOLVER_HPP
37 
38 #include <iostream>
39 #include <boost/shared_ptr.hpp>
40 
41 
42 #ifdef HAVE_MPI
43 #include <Epetra_MpiComm.h>
44 #else
45 #include <Epetra_SerialComm.h>
46 #endif
47 
48 #include <Teuchos_ParameterList.hpp>
49 #include <Teuchos_XMLParameterListHelpers.hpp>
50 #include <Teuchos_RCP.hpp>
51 
52 
53 #include <lifev/core/LifeV.hpp>
54 #include <lifev/core/array/VectorEpetra.hpp>
55 #include <lifev/core/util/Displayer.hpp>
56 #include <lifev/core/mesh/RegionMesh.hpp>
57 #include <lifev/core/fem/TimeAdvanceBDF.hpp>
58 #include <lifev/core/fem/BCHandler.hpp>
59 #include <lifev/navier_stokes/solver/NavierStokesSolver/NavierStokesProblem.hpp>
60 
61 
62 namespace LifeV
63 {
64 
65 template< class mesh_Type, class TimeIterationPolicy >
66 struct InitPolicySolver : public virtual TimeIterationPolicy
67 {
70  typedef MeshPartitioner< mesh_Type > meshPartitioner_Type;
73  typedef TimeAdvanceBDF<vector_Type> bdf_Type;
78 
80  virtual ~InitPolicySolver() {}
81 
82  void setupInit ( Teuchos::ParameterList& list );
83 
84  void initSimulation ( bcContainerPtr_Type bchandler,
85  vectorPtr_Type solution );
86 
87  virtual Displayer displayer() = 0;
88  virtual Real timestep() const = 0;
89  virtual Real initialTime() const = 0;
90  virtual bdfPtr_Type bdf() const = 0;
91 };
92 
93 template< class mesh_Type, class TimeIterationPolicy >
94 void
95 InitPolicySolver< mesh_Type, TimeIterationPolicy >::
96 setupInit ( Teuchos::ParameterList& list )
97 {
98  TimeIterationPolicy::initTimeIteration ( list );
99 }
100 
101 template< class mesh_Type, class TimeIterationPolicy >
102 void
103 InitPolicySolver< mesh_Type, TimeIterationPolicy >::
105  vectorPtr_Type solution )
106 {
107  displayer().leaderPrint ( "\n[Initializing the problem]\n" );
108 
109  LifeChrono nsTimeLoopChrono;
110  nsTimeLoopChrono.start();
111 
112  Real currentTime = initialTime() - timestep() * ( bdf()->order() - 1 );
113 
114  *solution = 0.0;
115  bdf()->setInitialCondition ( *solution );
116 
117  while ( currentTime <= initialTime() + timestep() / 2.0 )
118  {
119  LifeChrono iterChrono;
120  iterChrono.start();
121 
122  displayer().leaderPrint ( "\n[t = ", currentTime, " s.]\n" );
123 
124  // if( !M_usePreviousSolutionAsGuess ) *solution = 0;
125 
126  TimeIterationPolicy::iterate ( solution,
127  bchandler,
128  currentTime );
129 
130  // Updating the BDF scheme
131  bdf()->shiftRight ( *solution );
132 
133  iterChrono.stop();
134  displayer().leaderPrintMax ( "Iteration time: ", iterChrono.diff(), " s.\n" );
135 
136  currentTime += timestep();
137 
138 #ifdef HAVE_MPI
139  MPI_Barrier ( MPI_COMM_WORLD );
140 #endif
141  }
142 
143  nsTimeLoopChrono.stop();
144  displayer().leaderPrintMax ("Time of the temporal loop: ", nsTimeLoopChrono.diff(), " s.\n");
145 }
146 
147 } // namespace LifeV
148 
149 #endif /* INITPOLICYSTOKES_HPP */
VectorEpetra - The Epetra Vector format Wrapper.
void start()
Start the timer.
Definition: LifeChrono.hpp:93
std::shared_ptr< map_Type > mapPtr_Type
virtual bdfPtr_Type bdf() const =0
BCHandler - class for handling boundary conditions.
Definition: BCHandler.hpp:100
void initSimulation(bcContainerPtr_Type bchandler, vectorPtr_Type solution)
void setupInit(Teuchos::ParameterList &list)
std::shared_ptr< bdf_Type > bdfPtr_Type
void updateInverseJacobian(const UInt &iQuadPt)
TimeAdvanceBDF< vector_Type > bdf_Type
Epetra_Import const & importer()
Getter for the Epetra_Import.
Definition: MapEpetra.cpp:394
virtual Real timestep() const =0
std::shared_ptr< bcContainer_Type > bcContainerPtr_Type
Real diff()
Compute the difference in time between start and stop.
Definition: LifeChrono.hpp:111
MeshPartitioner< mesh_Type > meshPartitioner_Type
virtual Real initialTime() const =0
double Real
Generic real data.
Definition: LifeV.hpp:175
virtual Displayer displayer()=0
std::shared_ptr< VectorEpetra > vectorPtr_Type
void stop()
Stop the timer.
Definition: LifeChrono.hpp:100
Displayer - This class is used to display messages in parallel simulations.
Definition: Displayer.hpp:62
std::shared_ptr< NavierStokesProblem< mesh_Type > > NSProblemPtr_Type