OpenFCST: The open-source Fuel Cell Simulation Toolbox
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
List of all members
FuelCellShop::Equation::IncompressibleSingleComponentNSEquations< dim > Class Template Reference

This class deals with Navier-Stokes fluid transport equations. More...

#include <incompressible_single_component_NS_equations.h>

Inheritance diagram for FuelCellShop::Equation::IncompressibleSingleComponentNSEquations< dim >:
Inheritance graph
[legend]
Collaboration diagram for FuelCellShop::Equation::IncompressibleSingleComponentNSEquations< dim >:
Collaboration graph
[legend]

Public Member Functions

Constructors, destructor, and initialization
 IncompressibleSingleComponentNSEquations (FuelCell::SystemManagement &system_management, boost::shared_ptr< FuelCell::ApplicationCore::ApplicationData > data=boost::shared_ptr< FuelCell::ApplicationCore::ApplicationData >())
 Constructor. More...
 
virtual ~IncompressibleSingleComponentNSEquations ()
 Destructor. More...
 
virtual void declare_parameters (ParameterHandler &param) const
 Declare parameters. More...
 
virtual void initialize (ParameterHandler &param)
 Initialize parameters. More...
 
Local CG FEM based assemblers
virtual void assemble_cell_matrix (FuelCell::ApplicationCore::MatrixVector &cell_matrices, const typename FuelCell::ApplicationCore::DoFApplication< dim >::CellInfo &cell_info, FuelCellShop::Layer::BaseLayer< dim > *const layer)
 Assemble local cell matrix. More...
 
virtual void assemble_cell_residual (FuelCell::ApplicationCore::FEVector &cell_residual, const typename FuelCell::ApplicationCore::DoFApplication< dim >::CellInfo &cell_info, FuelCellShop::Layer::BaseLayer< dim > *const layer)
 Assemble local cell residual. More...
 
virtual void assemble_bdry_matrix (FuelCell::ApplicationCore::MatrixVector &bdry_matrices, const typename FuelCell::ApplicationCore::DoFApplication< dim >::FaceInfo &bdry_info, FuelCellShop::Layer::BaseLayer< dim > *const layer)
 Assemble local boundary matrix. More...
 
virtual void assemble_bdry_residual (FuelCell::ApplicationCore::FEVector &bdry_residual, const typename FuelCell::ApplicationCore::DoFApplication< dim >::FaceInfo &bdry_info, FuelCellShop::Layer::BaseLayer< dim > *const layer)
 Assemble local boundary residual. More...
 
Accessors and info
const double & get_inlet_outlet_velocity_max () const
 This function returns inlet_outlet_velocity_max. More...
 
const std::string get_inlet_outlet_equation () const
 User can specify an equation to use for density/velocity profile at boundary. More...
 
const unsigned int get_inlet_outlet_boundary_ID () const
 Outputs the Boundary ID to apply Dirichlet velocity or Dirichlet pressure to. More...
 
const std::string get_press_vel_comp_apply_to () const
 Outputs which equation to apply equation to (Pressure, VelocityX, VelocityY, VelocityZ). More...
 
virtual void print_equation_info () const
 This function prints out the equations info. More...
 
- Public Member Functions inherited from FuelCellShop::Equation::EquationBase< dim >
const couplings_mapget_internal_cell_couplings () const
 This function returns internal_cell_couplings of a derived equation class. More...
 
const couplings_mapget_internal_flux_couplings () const
 This function returns internal_flux_couplings (DG FEM only) of a derived equation class. More...
 
const
component_materialID_value_map
get_component_materialID_value () const
 This function returns component_materialID_value of a derived equation class. More...
 
const
component_boundaryID_value_map
get_component_boundaryID_value () const
 This function returns component_boundaryID_value of a derived equation class. More...
 
const std::vector< BoundaryType > & get_boundary_types () const
 This function returns boundary_types of a derived equation class. More...
 
const std::vector< std::vector
< BoundaryType > > & 
get_multi_boundary_types () const
 This function returns multi_boundary_types of a derived equation class. More...
 
const std::vector< OutputType > & get_output_types () const
 This function returns output_types of a derived equation class. More...
 
const std::vector< std::vector
< OutputType > > & 
get_multi_output_types () const
 This function returns multi_output_types of a derived equation class. More...
 
const std::string & get_equation_name () const
 This function returns equation_name of a derived equation class. More...
 
const std::vector< unsigned int > & get_matrix_block_indices () const
 This function returns matrix_block_indices of a derived equation class. More...
 
const std::vector< unsigned int > & get_residual_indices () const
 This function returns residual_indices of a derived equation class. More...
 

Protected Member Functions

Other - make_ functions
virtual void make_internal_cell_couplings ()
 This function fills out internal_cell_couplings. More...
 
virtual void make_matrix_block_indices ()
 This function fills out matrix_block_indices. More...
 
virtual void make_residual_indices ()
 This function fills out residual_indices. More...
 
- Protected Member Functions inherited from FuelCellShop::Equation::EquationBase< dim >
 EquationBase (FuelCell::SystemManagement &sys_management, boost::shared_ptr< FuelCell::ApplicationCore::ApplicationData > data=boost::shared_ptr< FuelCell::ApplicationCore::ApplicationData >())
 Constructor. More...
 
virtual ~EquationBase ()
 Destructor. More...
 
virtual void set_parameters (const std::vector< std::string > &name_dvar, const std::vector< double > &value_dvar, ParameterHandler &param)
 Set parameters using the parameter file, in order to run parametric/optimization studies. More...
 
virtual void make_assemblers_generic_constant_data ()
 Function used to initialize variable information that will be needed to assemble matrix and residual and that will remain constant throughout the simulation. More...
 
virtual void make_assemblers_cell_constant_data (const typename FuelCell::ApplicationCore::DoFApplication< dim >::CellInfo &cell_info)
 Function used to initialize cell speciific information that remains constant regardless of the cell being computed. More...
 
virtual void make_assemblers_bdry_constant_data (const typename FuelCell::ApplicationCore::DoFApplication< dim >::FaceInfo &bdry_info)
 
virtual void make_assemblers_cell_variable_data (const typename FuelCell::ApplicationCore::DoFApplication< dim >::CellInfo &cell_info, FuelCellShop::Layer::BaseLayer< dim > *const layer)
 Function used to compute cell specific information such as shape functions, shape function gradients, solution at each quadrature point. More...
 
virtual void make_assemblers_bdry_variable_data (const typename FuelCell::ApplicationCore::DoFApplication< dim >::FaceInfo &bdry_info, FuelCellShop::Layer::BaseLayer< dim > *const layer)
 
void select_cell_assemblers (const typename FuelCell::ApplicationCore::DoFApplication< dim >::CellInfo &cell_info, FuelCellShop::Layer::BaseLayer< dim > *const layer)
 This routine is used to select the make_assembly routines that need to be called inside assemble_cell_matrix to compute. More...
 
virtual void make_internal_flux_couplings ()
 This function fills out internal_flux_couplings (DG FEM only) of a derived equation class. More...
 
virtual void make_component_materialID_value ()
 This function fills out component_materialID_value of a derived equation class. More...
 
virtual void make_component_boundaryID_value ()
 This function fills out component_boundaryID_value of a derived equation class. More...
 
virtual void make_boundary_types ()
 This function fills out boundary_types of a derived equation class. More...
 
virtual void make_multi_boundary_types ()
 This function fills out multi_boundary_types of a derived equation class. More...
 
virtual void make_output_types ()
 This function fills out output_types of a derived equation class. More...
 
virtual void make_multi_output_types ()
 This function fills out multi_output_types of a derived equation class. More...
 
virtual void assemble_cell_Jacobian_matrix (FuelCell::ApplicationCore::MatrixVector &cell_matrices, const typename FuelCell::ApplicationCore::DoFApplication< dim >::CellInfo &cell_info, FuelCellShop::Layer::BaseLayer< dim > *const layer)
 Assemble the local Jacobian Matrix for Non-Linear problems. More...
 
virtual void assemble_bdry_Jacobian_matrix (FuelCell::ApplicationCore::MatrixVector &bdry_matrices, const typename FuelCell::ApplicationCore::DoFApplication< dim >::FaceInfo &bdry_info, FuelCellShop::Layer::BaseLayer< dim > *const layer)
 Assemble local Jacobian boundary matrix for Non-Linear problems. More...
 
virtual void assemble_cell_linear_matrix (FuelCell::ApplicationCore::MatrixVector &cell_matrices, const typename FuelCell::ApplicationCore::DoFApplication< dim >::CellInfo &cell_info, FuelCellShop::Layer::BaseLayer< dim > *const layer)
 Assemble the local cell matrix for Linear problems. More...
 
virtual void assemble_cell_residual_rhs (FuelCell::ApplicationCore::FEVector &cell_residual, const typename FuelCell::ApplicationCore::DoFApplication< dim >::CellInfo &cell_info, FuelCellShop::Layer::BaseLayer< dim > *const layer)
 Assemble local cell RHS for nonlinear problems. More...
 
virtual void assemble_cell_linear_rhs (FuelCell::ApplicationCore::FEVector &cell_residual, const typename FuelCell::ApplicationCore::DoFApplication< dim >::CellInfo &cell_info, FuelCellShop::Layer::BaseLayer< dim > *const layer)
 Assemble local cell RHS for Linear problems. More...
 
virtual void assemble_bdry_linear_matrix (FuelCell::ApplicationCore::MatrixVector &bdry_matrices, const typename FuelCell::ApplicationCore::DoFApplication< dim >::FaceInfo &bdry_info, FuelCellShop::Layer::BaseLayer< dim > *const layer)
 Assemble local boundary matrix for linear problems. More...
 
virtual void assemble_bdry_linear_rhs (FuelCell::ApplicationCore::FEVector &bdry_residual, const typename FuelCell::ApplicationCore::DoFApplication< dim >::FaceInfo &bdry_info, FuelCellShop::Layer::BaseLayer< dim > *const layer)
 Assemble local boundary RHS for linear problems. More...
 
void standard_to_block_wise (FullMatrix< double > &target) const
 This function changes the order of dealii::FullMatrix<double> target from standard to block-wise. More...
 
void standard_to_block_wise (Vector< double > &target) const
 This function changes the order of dealii::Vector<double> target from standard to block-wise. More...
 
void dealII_to_appframe (FuelCell::ApplicationCore::MatrixVector &dst, const FullMatrix< double > &src, const std::vector< unsigned int > &matrix_block_indices) const
 This function converts the standard ordered structure dealii::FullMatrix<double> src into the block-wise ordered structure FuelCell::ApplicationCore::MatrixVector dst. More...
 
void dealII_to_appframe (FuelCell::ApplicationCore::FEVector &dst, const Vector< double > &src, const std::vector< unsigned int > &residual_indices) const
 This function converts the standard ordered structure dealii::Vector<double> src into the block-wise ordered structure FuelCell::ApplicationCore::FEVector dst. More...
 
bool belongs_to_boundary (const unsigned int &tria_boundary_id, const unsigned int &param_boundary_id) const
 This function returns true if a boundary indicator of an external face on the triangulation coincides with a boundary indicator defined in the parameters file of a derived equation class. More...
 
void print_caller_name (const std::string &caller_name) const
 This function is used to print out the name of another function that has been declared in the scope of this class, but not yet been implemented. More...
 

Protected Attributes

Boolean constants and form of the drag force in porous media
bool inertia_in_channels
 This object indicates whether the inertia term is ON or OFF in channels. More...
 
bool shear_stress_in_channels
 This object indicates whether the shear stress term is ON or OFF in channels. More...
 
bool gravity_in_channels
 This object indicates whether the gravity term is ON or OFF in channels. More...
 
bool inertia_in_porous_media
 This object indicates whether the inertia term is ON or OFF in porous media. More...
 
bool shear_stress_in_porous_media
 This object indicates whether the shear stress term is ON or OFF in porous media. More...
 
bool gravity_in_porous_media
 This object indicates whether the gravity term is ON or OFF in porous media. More...
 
std::string drag_in_porous_media
 This object indicates which form of the drag term is supposed to be chosen. More...
 
bool normal_velocity_is_suppressed_weakly
 Sometimes we implement different types of slip boundary conditions on the curved impermeable walls of the domain or the walls that are not perfectly aligned with the coordinate axes. More...
 
Computational constants
Tensor< 1, dimgravity_acceleration
 Gravitational acceleration, $ \mathbf{g} = \{ g_{\alpha} \}_{\alpha = 1}^d \quad \text{such that} \quad \forall \alpha \neq d : \quad g_{\alpha} = 0 \quad \left[ \frac{\text{m}}{\text{sec}^2} \right] \quad \text{and} \quad g_d = -9.81 \quad \left[ \frac{\text{m}}{\text{sec}^2} \right] $. More...
 
SymmetricTensor< 2, dimunit
 Unit tensor, $ \hat{ \mathbf{I} } = \{ \delta_{\alpha \beta} \}_{\alpha,\beta = 1}^d $. More...
 
double theta
 Navier slip coefficient, $ \theta $. More...
 
double eta
 Normal velocity suppression coefficient, $ \eta $. More...
 
double inlet_outlet_velocity_max
 Maximum inlet-outlet velocity, $ \mathbf{u}_{\text{in,out}}^{\text{max}} $ $ \quad \left[ \frac{\text{m}}{\text{sec}} \right] $. More...
 
std::string inlet_outlet_equation
 
unsigned int inlet_outlet_boundary_ID
 
std::string press_vel_comp_apply_to
 
- Protected Attributes inherited from FuelCellShop::Equation::EquationBase< dim >
unsigned int dofs_per_cell
 Number of degrees of freedom per cell. More...
 
unsigned int n_q_points_cell
 Number of quadrature points per cell. More...
 
unsigned int n_q_points_bdry
 Number of quadrature points per boundary. More...
 
DoFHandler< dim >
::active_cell_iterator 
cell
 Currently active DoFHandler<dim> active cell iterator. More...
 
DoFHandler< dim >
::active_face_iterator 
bdry
 Currently active DoFHandler<dim> active boundary iterator. More...
 
std::vector< double > JxW_cell
 Jacobian of mapping by Weight in the quadrature points of a cell. More...
 
std::vector< double > JxW_bdry
 Jacobian of mapping by Weight in the quadrature points of a boundary. More...
 
std::vector< Point< dim > > normal_vectors
 Normal vectors in the quadrature points of a boundary. More...
 
std::vector< std::vector
< Point< dim > > > 
tangential_vectors
 Tangential vectors in the quadrature points of a boundary. More...
 
FuelCell::SystemManagementsystem_management
 Pointer to the external YourApplication<dim>::system_management object. More...
 
couplings_map internal_cell_couplings
 This object contains the info on how the equations and solution variables of a derived equation class are coupled. More...
 
couplings_map internal_flux_couplings
 This object contains the info on how the "X" and "Y" of a derived equation class are coupled (DG FEM only). More...
 
component_materialID_value_map component_materialID_value
 This object reflects the following structure (see FuelCell::InitialAndBoundaryData namespace docs): More...
 
component_boundaryID_value_map component_boundaryID_value
 This object reflects the following structure (see FuelCell::InitialAndBoundaryData namespace docs): More...
 
std::vector< BoundaryTypeboundary_types
 The list of boundary types of a derived equation class. More...
 
std::vector< std::vector
< BoundaryType > > 
multi_boundary_types
 The list of multiple boundary types of a derived equation class. More...
 
std::vector< OutputTypeoutput_types
 The list of output types of a derived equation class. More...
 
std::vector< std::vector
< OutputType > > 
multi_output_types
 The list of multiple output types of a derived equation class. More...
 
std::string equation_name
 The name of a derived equation class. More...
 
std::string name_base_variable
 Const std::string member storing name of the base solution variable corresponding to the equation represented by this class. More...
 
std::vector< unsigned int > matrix_block_indices
 The system matrix block indices (a derived equation class) drawn from the global structure (a derived equation class + other active equation classes included into the computation). More...
 
std::vector< unsigned int > residual_indices
 The residual indices (a derived equation class) drawn from the global structure (a derived equation class + other active equation classes included into the computation). More...
 
std::vector< bool > counter
 This vector contains the collection of internal "counters" used by the derived equation classes. More...
 
EquationFlags assemble_flags
 This vector contains a collection of internal flags to tell derived equation classes what needs to be re-computed. More...
 
boost::shared_ptr
< FuelCell::ApplicationCore::ApplicationData
data
 Data object for the application data to be passed to the equation classes. More...
 
std::string solution_vector_name
 The name of the solution vector in FEVectors. More...
 
std::string residual_vector_name
 The name of the residual vector name in FEVectors. More...
 

Additional Inherited Members

- Public Attributes inherited from FuelCellShop::Equation::EquationBase< dim >
bool variable_initial_data
 true, if variable initial data is prescribed on a part of the domain. More...
 
bool variable_boundary_data
 true, if variable Dirichlet boundary conditions are prescribed on a part of the boundary. More...
 

Detailed Description

template<int dim>
class FuelCellShop::Equation::IncompressibleSingleComponentNSEquations< dim >

This class deals with Navier-Stokes fluid transport equations.

According to classification, these equations are

These equations can be written as

$ \nabla \cdot \mathbf{u} = 0 \quad \text{in} \quad \Omega $ $ \quad - \quad $ mass conservation equation

$ \rho \nabla \cdot \left( \frac{1}{\epsilon} \mathbf{u} \otimes \mathbf{u} \right) = \nabla \cdot \left( -p \hat{\mathbf{I}} + \hat{\boldsymbol\sigma} \right) + \left( \mathbf{F} + \epsilon \rho \mathbf{g} \right) \quad \text{in} \quad \Omega $ $ \quad - \quad $ momentum conservation equation

where

$ \hat{\boldsymbol\sigma} = 2 \mu \nabla_s \mathbf{u} $ $ \quad - \quad $ shear stress

$ \mathbf{F} = \begin{cases} \mathbf{0} \quad \text{in} \quad \Omega_c \\ - \mu \hat{ \mathbf{K} }^{-1} \mathbf{u} - \hat{\boldsymbol\beta} \rho \big| \mathbf{u} \big| \mathbf{u} \quad \text{in} \quad \Omega_p \end{cases} $ $ \quad - \quad $ drag force

To be well-posed, these equations are equipped with appropriate boundary conditions. Below we consider several types of boundaries and several types of boundary conditions assigned to those boundaries (here $ k $ denotes some portion of a boundary of a certain type).

* Impermeable walls = [boundary_id]: [type]
*

where

* Symmetry line or plane = [boundary_id]: [type]
*

where

* Inlet-Outlet = [boundary_id]: [type]
*

where

These equations utilize

as layer classes and

as material classes and have the following output:

The whole problem is solved by linearizing the governing equations at each Newton iteration with subsequent CG FEM discretization in space.

The implementation is performed by means of using only one weak formulation for all $ d + 1 $ scalar governing equations. This approach guarantees the dim-independent programming of the equations at hand.

Author
Valentin N. Zingan, 2012

Constructor & Destructor Documentation

Constructor.

Destructor.

Member Function Documentation

template<int dim>
virtual void FuelCellShop::Equation::IncompressibleSingleComponentNSEquations< dim >::assemble_bdry_matrix ( FuelCell::ApplicationCore::MatrixVector bdry_matrices,
const typename FuelCell::ApplicationCore::DoFApplication< dim >::FaceInfo &  bdry_info,
FuelCellShop::Layer::BaseLayer< dim > *const  layer 
)
virtual

Assemble local boundary matrix.

Reimplemented from FuelCellShop::Equation::EquationBase< dim >.

template<int dim>
virtual void FuelCellShop::Equation::IncompressibleSingleComponentNSEquations< dim >::assemble_bdry_residual ( FuelCell::ApplicationCore::FEVector bdry_residual,
const typename FuelCell::ApplicationCore::DoFApplication< dim >::FaceInfo &  bdry_info,
FuelCellShop::Layer::BaseLayer< dim > *const  layer 
)
virtual

Assemble local boundary residual.

Reimplemented from FuelCellShop::Equation::EquationBase< dim >.

template<int dim>
virtual void FuelCellShop::Equation::IncompressibleSingleComponentNSEquations< dim >::assemble_cell_matrix ( FuelCell::ApplicationCore::MatrixVector cell_matrices,
const typename FuelCell::ApplicationCore::DoFApplication< dim >::CellInfo &  cell_info,
FuelCellShop::Layer::BaseLayer< dim > *const  layer 
)
virtual

Assemble local cell matrix.

Reimplemented from FuelCellShop::Equation::EquationBase< dim >.

template<int dim>
virtual void FuelCellShop::Equation::IncompressibleSingleComponentNSEquations< dim >::assemble_cell_residual ( FuelCell::ApplicationCore::FEVector cell_residual,
const typename FuelCell::ApplicationCore::DoFApplication< dim >::CellInfo &  cell_info,
FuelCellShop::Layer::BaseLayer< dim > *const  layer 
)
virtual

Assemble local cell residual.

Reimplemented from FuelCellShop::Equation::EquationBase< dim >.

template<int dim>
virtual void FuelCellShop::Equation::IncompressibleSingleComponentNSEquations< dim >::declare_parameters ( ParameterHandler &  param) const
virtual

Declare parameters.

Reimplemented from FuelCellShop::Equation::EquationBase< dim >.

template<int dim>
const unsigned int FuelCellShop::Equation::IncompressibleSingleComponentNSEquations< dim >::get_inlet_outlet_boundary_ID ( ) const
inline

Outputs the Boundary ID to apply Dirichlet velocity or Dirichlet pressure to.

References FuelCellShop::Equation::IncompressibleSingleComponentNSEquations< dim >::inlet_outlet_boundary_ID.

template<int dim>
const std::string FuelCellShop::Equation::IncompressibleSingleComponentNSEquations< dim >::get_inlet_outlet_equation ( ) const
inline

User can specify an equation to use for density/velocity profile at boundary.

E. inlet-outlet equation = -625*y*y + 66.25*y - 0.755625

References FuelCellShop::Equation::IncompressibleSingleComponentNSEquations< dim >::inlet_outlet_equation.

template<int dim>
const double& FuelCellShop::Equation::IncompressibleSingleComponentNSEquations< dim >::get_inlet_outlet_velocity_max ( ) const
inline
template<int dim>
const std::string FuelCellShop::Equation::IncompressibleSingleComponentNSEquations< dim >::get_press_vel_comp_apply_to ( ) const
inline

Outputs which equation to apply equation to (Pressure, VelocityX, VelocityY, VelocityZ).

References FuelCellShop::Equation::IncompressibleSingleComponentNSEquations< dim >::press_vel_comp_apply_to.

template<int dim>
virtual void FuelCellShop::Equation::IncompressibleSingleComponentNSEquations< dim >::initialize ( ParameterHandler &  param)
virtual

Initialize parameters.

Reimplemented from FuelCellShop::Equation::EquationBase< dim >.

template<int dim>
virtual void FuelCellShop::Equation::IncompressibleSingleComponentNSEquations< dim >::make_internal_cell_couplings ( )
protectedvirtual

This function fills out internal_cell_couplings.

Reimplemented from FuelCellShop::Equation::EquationBase< dim >.

template<int dim>
virtual void FuelCellShop::Equation::IncompressibleSingleComponentNSEquations< dim >::make_matrix_block_indices ( )
protectedvirtual

This function fills out matrix_block_indices.

Reimplemented from FuelCellShop::Equation::EquationBase< dim >.

template<int dim>
virtual void FuelCellShop::Equation::IncompressibleSingleComponentNSEquations< dim >::make_residual_indices ( )
protectedvirtual

This function fills out residual_indices.

Reimplemented from FuelCellShop::Equation::EquationBase< dim >.

template<int dim>
virtual void FuelCellShop::Equation::IncompressibleSingleComponentNSEquations< dim >::print_equation_info ( ) const
virtual

This function prints out the equations info.

Reimplemented from FuelCellShop::Equation::EquationBase< dim >.

Member Data Documentation

template<int dim>
std::string FuelCellShop::Equation::IncompressibleSingleComponentNSEquations< dim >::drag_in_porous_media
protected

This object indicates which form of the drag term is supposed to be chosen.

There are four options currently implemented in FCST. These options are

  • none,
  • Darcy,
  • Forchheimer,
  • Forchheimer modified.
template<int dim>
double FuelCellShop::Equation::IncompressibleSingleComponentNSEquations< dim >::eta
protected

Normal velocity suppression coefficient, $ \eta $.

template<int dim>
Tensor<1,dim> FuelCellShop::Equation::IncompressibleSingleComponentNSEquations< dim >::gravity_acceleration
protected

Gravitational acceleration, $ \mathbf{g} = \{ g_{\alpha} \}_{\alpha = 1}^d \quad \text{such that} \quad \forall \alpha \neq d : \quad g_{\alpha} = 0 \quad \left[ \frac{\text{m}}{\text{sec}^2} \right] \quad \text{and} \quad g_d = -9.81 \quad \left[ \frac{\text{m}}{\text{sec}^2} \right] $.

template<int dim>
bool FuelCellShop::Equation::IncompressibleSingleComponentNSEquations< dim >::gravity_in_channels
protected

This object indicates whether the gravity term is ON or OFF in channels.

template<int dim>
bool FuelCellShop::Equation::IncompressibleSingleComponentNSEquations< dim >::gravity_in_porous_media
protected

This object indicates whether the gravity term is ON or OFF in porous media.

template<int dim>
bool FuelCellShop::Equation::IncompressibleSingleComponentNSEquations< dim >::inertia_in_channels
protected

This object indicates whether the inertia term is ON or OFF in channels.

template<int dim>
bool FuelCellShop::Equation::IncompressibleSingleComponentNSEquations< dim >::inertia_in_porous_media
protected

This object indicates whether the inertia term is ON or OFF in porous media.

template<int dim>
unsigned int FuelCellShop::Equation::IncompressibleSingleComponentNSEquations< dim >::inlet_outlet_boundary_ID
protected
template<int dim>
std::string FuelCellShop::Equation::IncompressibleSingleComponentNSEquations< dim >::inlet_outlet_equation
protected
template<int dim>
double FuelCellShop::Equation::IncompressibleSingleComponentNSEquations< dim >::inlet_outlet_velocity_max
protected
template<int dim>
bool FuelCellShop::Equation::IncompressibleSingleComponentNSEquations< dim >::normal_velocity_is_suppressed_weakly
protected

Sometimes we implement different types of slip boundary conditions on the curved impermeable walls of the domain or the walls that are not perfectly aligned with the coordinate axes.

In this case, it might be quite problematic (but still possible) to strictly enforce the normal component of velocity be equal to 0.

The alternative approach is doing that by introducing a penalization term $ \displaystyle \int_{ \Gamma^{\text{walls},k} } \frac{1}{\eta} \left( \boldsymbol\omega \cdot \mathbf{n} \right) \left( \mathbf{u} \cdot \mathbf{n} \right) dS $ with $ \eta \sim 10^{-10} - 10^{-12} $ into the weak formulation of the equations at hand.

This object indicates whether the penalty method is used.

template<int dim>
std::string FuelCellShop::Equation::IncompressibleSingleComponentNSEquations< dim >::press_vel_comp_apply_to
protected
template<int dim>
bool FuelCellShop::Equation::IncompressibleSingleComponentNSEquations< dim >::shear_stress_in_channels
protected

This object indicates whether the shear stress term is ON or OFF in channels.

template<int dim>
bool FuelCellShop::Equation::IncompressibleSingleComponentNSEquations< dim >::shear_stress_in_porous_media
protected

This object indicates whether the shear stress term is ON or OFF in porous media.

template<int dim>
double FuelCellShop::Equation::IncompressibleSingleComponentNSEquations< dim >::theta
protected

Navier slip coefficient, $ \theta $.

template<int dim>
SymmetricTensor<2,dim> FuelCellShop::Equation::IncompressibleSingleComponentNSEquations< dim >::unit
protected

Unit tensor, $ \hat{ \mathbf{I} } = \{ \delta_{\alpha \beta} \}_{\alpha,\beta = 1}^d $.


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