LifeV
StabilizationSUPG_semi_implicit_ale.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  @file
28  @brief SUPG stabilization, semi-implicit case.
29  @author Davide Forti <davide.forti@epfl.ch>
30  @contributor Luca Dede <luca.dede@epfl.ch>
31  @maintainer Davide Forti <davide.forti@epfl.ch>
32  @date 15-04-2015
33 
34  This file contains an ETA implementation of SUPG, for semi-implicit treatment of the convective term.
35 
36  For reference, see Paper
37 
38  * D. Forti, L. Dede, Semi–implicit BDF time discretization of the Navier–Stokes equations with VMS–LES modeling in
39  a High Performance Computing framework. Submitted, available as MATHICSE report, 2015.
40 
41  * Y. Bazilevs, V.M. Calo, J.A. Cottrell, T.J.R. Hughes, A. Reali, and G. Scovazzi. <i>Variational multiscale
42  residual-based turbulence modeling for large eddy simulation of incompressible flows</i>. Comput. Methods Appl. Mech. Engr. 197(1):173–201, 2007.
43 
44 */
45 
46 #ifndef _STABILIZATIONSUPG_SEMI_IMPLICIT_ALE_HPP_
47 #define _STABILIZATIONSUPG_SEMI_IMPLICIT_ALE_HPP_ 1
48 
49 // Tell the compiler to ignore specific kind of warnings:
50 #pragma GCC diagnostic ignored "-Wunused-variable"
51 #pragma GCC diagnostic ignored "-Wunused-parameter"
52 
53 #include <Epetra_ConfigDefs.h>
54 #ifdef EPETRA_MPI
55  #include <mpi.h>
56  #include <Epetra_MpiComm.h>
57 #else
58  #include <Epetra_SerialComm.h>
59 #endif
60 
61 #include <Epetra_FECrsMatrix.h>
62 
63 //Tell the compiler to restore the warning previously silented
64 #pragma GCC diagnostic warning "-Wunused-variable"
65 #pragma GCC diagnostic warning "-Wunused-parameter"
66 
67 #include <lifev/navier_stokes_blocks/solver/Stabilization.hpp>
68 
69 #include <lifev/core/fem/FESpace.hpp>
70 #include <lifev/core/fem/ReferenceFE.hpp>
71 
72 #include <lifev/eta/fem/ETFESpace.hpp>
73 
74 // includes for building the matrix graph
75 #include <Epetra_FECrsGraph.h>
76 #include <lifev/eta/expression/Integrate.hpp>
77 #include <lifev/eta/expression/BuildGraph.hpp>
78 
79 #include <lifev/core/fem/TimeAndExtrapolationHandlerQuadPts.hpp>
80 
81 namespace LifeV
82 {
83 
85 {
86 public:
87  typedef Real return_Type;
88 
89  inline return_Type operator() (const Real& a)
90  {
91  return std::sqrt(a);
92  }
93 
97 };
98 
100 {
101 public:
102 
103  //@name Public Types
104  //@{
107 
108  typedef VectorEpetra vector_Type;
110 
113 
116 
119 
122 
125 
126  //@}
127 
128  //! @name Constructor and Destructor
129  //@{
130  //! Default Constructor
132 
133  //! ~Destructor
135 
136  //@}
137 
138  //! Build the graphs of each single block
139  void buildGraphs();
140 
141  //! Updates the system matrix in Navier-Stokes simulations in fixed coordinates
142  // with semi-implicit treatment of the convective term.
143  /*!
144  @param velocityExtrapolated extrapolation of the fluid velocity
145  @param aleVelocity velocity of the domain
146  */
147  void apply_matrix( const vector_Type& velocityExtrapolated, const vector_Type& velocityALE );
148 
149  //! Adds to the right hand side the contribution coming from the SUPG stabilization
150  // in Navier-Stokes simulations in fixed coordinates. Used for NS semi-implicit
151  /*!
152  @param rhs_velocity velocity component of the right hand side
153  @param rhs_pressure pressure component of the right hand side
154  @param velocity_extrapolated velocity extrapolated
155  @param velocity_rhs velocity term from approximation time derivative
156  */
157  void apply_vector( vectorPtr_Type& rhs_velocity,
158  vectorPtr_Type& rhs_pressure,
159  const vector_Type& velocityExtrapolated,
160  const vector_Type& velocityALE,
161  const vector_Type& velocity_rhs);
162 
163  void setVelocitySpace(fespacePtr_Type velocityFESpace){ M_uFESpace = velocityFESpace;}
164 
165  void setPressureSpace(fespacePtr_Type pressureFESpace){ M_pFESpace = pressureFESpace;}
166 
167  //! Set the constant C_I for the supg
168  void setConstant (const int & value);
169 
170  //! Set the fluid density
171  void setDensity (const Real & density) { M_density = density;}
172 
173  //! Set the bdf order
174  void setBDForder (const Real & bdfOrder) { M_bdfOrder = bdfOrder;}
175 
176  //! Set the bdf order
177  void setAlpha (const Real & alpha) { M_alpha = alpha;}
178 
179  //! Set the fluid dynamic viscosity
180  void setViscosity (const Real & viscosity) { M_viscosity = viscosity;}
181 
182  //! Set the Epetra communicator
183  void setCommunicator (std::shared_ptr<Epetra_Comm> comm) { M_comm = comm;}
184 
185  //! Set the time step size
186  void setTimeStep (const Real & timestep) { M_timestep = timestep;}
187 
188  void setETvelocitySpace(const ETFESpacePtr_velocity & velocityEta_fespace){ M_fespaceUETA = velocityEta_fespace;}
189 
190  void setETpressureSpace(const ETFESpacePtr_pressure & pressureEta_fespace){ M_fespacePETA = pressureEta_fespace;}
191 
192  void setUseGraph (const bool& useGraph) { M_useGraph = useGraph; }
193 
194  //! @name Getters
195  //@{
196 
197  matrixPtr_Type const& block_00() const
198  {
199  return M_block_00;
200  }
201 
202  matrixPtr_Type const& block_01() const
203  {
204  return M_block_01;
205  }
206 
207  matrixPtr_Type const& block_10() const
208  {
209  return M_block_10;
210  }
211 
212  matrixPtr_Type const& block_11() const
213  {
214  return M_block_11;
215  }
216 
217  std::string label () { return M_label; }
218 
219  void setUseODEfineScale ( const bool& useODEfineScale );
220 
221  void updateODEfineScale ( const vectorPtr_Type& velocity, const vectorPtr_Type& pressure );
222 
223  void setExportFineScaleVelocity ( ExporterHDF5<mesh_Type> & exporter, const int& numElementsTotal);
224 
225  //@}
226 
227 private:
228 
229  void setupODEfineScale();
230 
231  void computeFineScales ( const vectorPtr_Type& velocity, const vectorPtr_Type& pressure );
232 
233  void computeFineScalesForVisualization ( const vectorPtr_Type& velocity, const vectorPtr_Type& pressure );
234 
235  //! @name Private Attributes
236  //@{
237 
238  //! finite element spaces for velocity and pressure
241 
244 
245  //! Epetra communicator
247 
248  //! fluid dynamic viscosity @f$\nu@f$
250 
251  //! fluid density @f$\nu@f$
253 
254  //! stabilization parameters for the momentum and continuity equations
259 
261 
262  // Graphs
267 
268  // Matrices
273 
275 
277 
279 
281 
285 
289 
290  //@}
291 }; // class StabilizationSUPG_semi_implicit
292 
293 //! Factory create function
295 {
297 }
298 namespace
299 {
302 }
303 
304 } // namespace LifeV
305 
306 #endif /* _STABILIZATIONSUPG_SEMI_IMPLICIT_HPP_ */
void setPressureSpace(fespacePtr_Type pressureFESpace)
Set Expression Template FE space for velocity.
void computeFineScales(const vectorPtr_Type &velocity, const vectorPtr_Type &pressure)
SquareRoot_supg_semi_implicit_ale(const SquareRoot_supg_semi_implicit_ale &)
void apply_matrix(const vector_Type &velocityExtrapolated, const vector_Type &velocityALE)
Updates the system matrix in Navier-Stokes simulations in fixed coordinates.
void setConstant(const int &value)
Set the constant C_I for the supg.
void setupODEfineScale()
Setup of the fine scale component.
void setBDForder(const Real &bdfOrder)
Set the bdf order.
void setExportFineScaleVelocity(ExporterHDF5< mesh_Type > &exporter, const int &numElementsTotal)
Set if the user wants to export the fine scale component.
void setCommunicator(std::shared_ptr< Epetra_Comm > comm)
Set the Epetra communicator.
FESpace - Short description here please!
Definition: FESpace.hpp:78
Real M_timestep
stabilization parameters for the momentum and continuity equations
void importFromHDF5(std::string const &fileName, std::string const &matrixName="matrix")
Read a matrix from a HDF5 (.h5) file.
std::vector< std::vector< VectorSmall< 1 > > > M_fineScalePressure
matrixPtr_Type const & block_10() const
Get block10 of the stabilization matrix.
StabilizationSUPG_semi_implicit_ale * createStabilizationSUPG_semi_implicit_ale()
Factory create function.
void setTimeStep(const Real &timestep)
Set the time step size.
void setETvelocitySpace(const ETFESpacePtr_velocity &velocityEta_fespace)
Set Expression Template FE space for velocity.
void setUseODEfineScale(const bool &useODEfineScale)
Set if using dynamic fine scale model.
void apply_vector(vectorPtr_Type &rhs_velocity, vectorPtr_Type &rhs_pressure, const vector_Type &velocityExtrapolated, const vector_Type &velocityALE, const vector_Type &velocity_rhs)
Adds to the right hand side the contribution coming from the SUPG stabilization.
fespacePtr_Type M_uFESpace
finite element spaces for velocity and pressure
std::shared_ptr< Epetra_Comm > M_comm
Epetra communicator.
Epetra_Import const & importer()
Getter for the Epetra_Import.
Definition: MapEpetra.cpp:394
matrixPtr_Type const & block_11() const
Get block11 of the stabilization matrix.
matrixPtr_Type const & block_00() const
Get block00 of the stabilization matrix.
void computeFineScalesForVisualization(const vectorPtr_Type &velocity, const vectorPtr_Type &pressure)
void setETpressureSpace(const ETFESpacePtr_pressure &pressureEta_fespace)
Set Expression Template FE space for velocity.
void buildGraphs()
Build the graphs of each single block.
double Real
Generic real data.
Definition: LifeV.hpp:175
void setDensity(const Real &density)
Set the fluid density.
std::string label()
Get name of stabilization being used.
void setUseGraph(const bool &useGraph)
Set if using the graph.
void setViscosity(const Real &viscosity)
Set the fluid dynamic viscosity.
class ETFESpace A light, templated version of the FESpace
Definition: ETFESpace.hpp:93
matrixPtr_Type const & block_01() const
Get block01 of the stabilization matrix.
void setVelocitySpace(fespacePtr_Type velocityFESpace)
Set FE space for velocity.
void updateODEfineScale(const vectorPtr_Type &velocity, const vectorPtr_Type &pressure)
Updates the fine scale component.