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>
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 |
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:
where:
Notes on preconditioners:
\[ S^-1 = \nu M_p \]
\[ S^-1 = \rho \frac{1}{\Delta t} A_p \]
\[ S^-1 = \nu M_p + \rho \frac{1}{\Delta t} A_p. \]
|
inline |
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 |