pi-DoMUS: Parallel Deal.II MUltiphysics Solver
NavierStokes< dim, spacedim, LAC > Class Template Reference

This interface solves a Navier Stokes Equation:

\[ \begin{cases} \partial_t u + (u\cdot\nabla)u - \nu\textrm{div} \epsilon(u) + \frac{1}{\rho}\nabla p = f \\ \textrm{div}u=0 \end{cases} \]

where \( \epsilon(u) = \frac{\nabla u + [\nabla u]^t}{2}. \). More...

#include <navier_stokes.h>

Inheritance diagram for NavierStokes< dim, spacedim, LAC >:
PDESystemInterface< dim, spacedim, NavierStokes< dim, spacedim, LAC >, LAC >

Public Member Functions

 ~NavierStokes ()
 
 NavierStokes (bool dynamic, bool convection)
 
void declare_parameters (ParameterHandler &prm)
 
void parse_parameters_call_back ()
 
template<typename EnergyType , typename ResidualType >
void energies_and_residuals (const typename DoFHandler< dim, spacedim >::active_cell_iterator &cell, FEValuesCache< dim, spacedim > &scratch, std::vector< EnergyType > &energies, std::vector< std::vector< ResidualType >> &residuals, bool compute_only_system_terms) const
 
void compute_system_operators (const DoFHandler< dim, spacedim > &, const std::vector< shared_ptr< LATrilinos::BlockMatrix >>, LinearOperator< LATrilinos::VectorType > &, LinearOperator< LATrilinos::VectorType > &) const
 
void set_matrix_couplings (std::vector< std::string > &couplings) const
 
template<typename EnergyType , typename ResidualType >
void energies_and_residuals (const typename DoFHandler< dim, spacedim >::active_cell_iterator &cell, FEValuesCache< dim, spacedim > &fe_cache, std::vector< EnergyType > &, std::vector< std::vector< ResidualType > > &residual, bool compute_only_system_terms) const
 
- Public Member Functions inherited from PDESystemInterface< dim, spacedim, NavierStokes< dim, spacedim, LAC >, LAC >
virtual ~PDESystemInterface ()
 
 PDESystemInterface (const std::string &name="", const unsigned int &n_components=0, const unsigned int &n_matrices=0, const std::string &default_fe="FE_Q(1)", const std::string &default_component_names="u", const std::string &default_differential_components="")
 Pass initializers to the base class constructor. More...
 
virtual void assemble_energies_and_residuals (const typename DoFHandler< dim, spacedim >::active_cell_iterator &cell, FEValuesCache< dim, spacedim > &scratch, std::vector< Sdouble > &energies, std::vector< std::vector< double > > &local_residuals, bool compute_only_system_terms) const
 
virtual void assemble_energies_and_residuals (const typename DoFHandler< dim, spacedim >::active_cell_iterator &cell, FEValuesCache< dim, spacedim > &scratch, std::vector< SSdouble > &energies, std::vector< std::vector< Sdouble > > &local_residuals, bool compute_only_system_terms) const
 

Detailed Description

template<int dim, int spacedim = dim, typename LAC = LATrilinos>
class NavierStokes< dim, spacedim, LAC >

This interface solves a Navier Stokes Equation:

\[ \begin{cases} \partial_t u + (u\cdot\nabla)u - \nu\textrm{div} \epsilon(u) + \frac{1}{\rho}\nabla p = f \\ \textrm{div}u=0 \end{cases} \]

where \( \epsilon(u) = \frac{\nabla u + [\nabla u]^t}{2}. \).

Non time-depending Navier Stokes Equation:

\[ \begin{cases} (u\cdot\nabla)u - \nu\textrm{div} \epsilon(u) + \frac{1}{\rho}\nabla p = f \\ \textrm{div}u=0 \end{cases} \]

can be recoverd setting dynamic = false

Dynamic Stokes Equation:

\[ \begin{cases} \partial_t u - \nu\textrm{div} \epsilon(u) + \frac{1}{\rho}\nabla p = f \\ \textrm{div}u=0 \end{cases} \]

can be recoverd setting convection = false

Stokes Equation:

\[ \begin{cases} - \nu\textrm{div} \epsilon(u) + \frac{1}{\rho}\nabla p = f \\ \textrm{div}u=0 \end{cases} \]

can be recoverd setting dynamic = false and convection = false

In the code we adopt the following notations:

  • Mp := block resulting from \( ( \partial_t p, q ) \)
  • Ap := block resulting from \( \nu ( \nabla p,\nabla q ) \)
  • Np := block resulting from \( ( u \cdot \nabla p, q) \)
  • Fp := Mp + Ap + Np

where:

  • p = pressure
  • q = test function for the pressure
  • u = velocity
  • v = test function for the velocity

Notes on preconditioners:

  • default: This preconditioner uses the mass matrix of pressure block as inverse for the Schur block. This is a preconditioner suitable for problems wher the viscosity is higher than the density.

    \[ S^-1 = \nu M_p \]

  • identity: Identity matrix preconditioner
  • low-nu: This preconditioner uses the stifness matrix of pressure block as inverse for the Schur block.

    \[ S^-1 = \rho \frac{1}{\Delta t} A_p \]

  • cah-cha: Preconditioner suggested by J. Cahouet and J.-P. Chabard.

    \[ S^-1 = \nu M_p + \rho \frac{1}{\Delta t} A_p. \]

Constructor & Destructor Documentation

template<int dim, int spacedim = dim, typename LAC = LATrilinos>
NavierStokes< dim, spacedim, LAC >::~NavierStokes ( )
inline

Member Function Documentation

template<int dim, int spacedim = dim, typename LAC = LATrilinos>
template<typename EnergyType , typename ResidualType >
void NavierStokes< dim, spacedim, LAC >::energies_and_residuals ( const typename DoFHandler< dim, spacedim >::active_cell_iterator cell,
FEValuesCache< dim, spacedim > &  scratch,
std::vector< EnergyType > &  energies,
std::vector< std::vector< ResidualType >> &  residuals,
bool  compute_only_system_terms 
) const

The documentation for this class was generated from the following file: