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

This is the base class used for all Equation classes. More...

#include <equation_base.h>

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

Public Member Functions

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 for the RHS. 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 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...
 
virtual void print_equation_info () const
 This function prints out the equations info of a derived equation class. More...
 

Public Attributes

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...
 

Protected Member Functions

Constructors, destructor, and initialization
 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 declare_parameters (ParameterHandler &param) const
 Declare parameters. More...
 
virtual void initialize (ParameterHandler &param)
 Initialize parameters. 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...
 
Local CG FEM based assemblers - make_ functions
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...
 
Other - make_ functions
virtual void make_internal_cell_couplings ()
 This function fills out internal_cell_couplings of a derived equation class. 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 make_matrix_block_indices ()
 This function is only needed to provide the last argument to dealII_to_appframe. More...
 
virtual void make_residual_indices ()
 This function is only needed to provide the last argument to dealII_to_appframe. More...
 
FEM assemblers
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...
 
Converters
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...
 
Other functions
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...
 
Minor functions
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

Local CG FEM based assemblers - constant data (generic)
unsigned int dofs_per_cell
 Number of degrees of freedom per cell. More...
 
Local CG FEM based assemblers - constant data (cell)
unsigned int n_q_points_cell
 Number of quadrature points per cell. More...
 
Local CG FEM based assemblers - constant data (boundary)
unsigned int n_q_points_bdry
 Number of quadrature points per boundary. More...
 
Local CG FEM based assemblers - variable data (active mesh iterators)
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...
 
Local CG FEM based assemblers - variable data (other - cell)

Implementation is in the derived equation classes.

std::vector< double > JxW_cell
 Jacobian of mapping by Weight in the quadrature points of a cell. More...
 
Local CG FEM based assemblers - variable data (other - boundary)

Implementation is in the derived equation classes.

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...
 
Generic Data
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...
 

Detailed Description

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

This is the base class used for all Equation classes.

Equation classes are those that implement the weak for of a given governing equation such as Fick's law, Ohm's law or a source term.

This class is used to lock the interface of all derived equation classes under the FuelCellShop::Equation namespace.

This class also contains generic data and methods heavily used by all derived equation classes. The functionality of this class can be extended if needed.

Author
Valentin N. Zingan, 2012
Marc Secanell, 2012-16
Mayank Sabharwal, 2015
Aslan Kosakian, 2015

Constructor & Destructor Documentation

template<int dim>
FuelCellShop::Equation::EquationBase< dim >::EquationBase ( FuelCell::SystemManagement sys_management,
boost::shared_ptr< FuelCell::ApplicationCore::ApplicationData data = boost::shared_ptr< FuelCell::ApplicationCore::ApplicationData >() 
)
protected

Constructor.

Warning
The constructor of the child classes needs to define the following to variables:
  • name_base_variable
  • equation_name for example, for lambda_transport_equation
    * this->equation_name = "Membrane Water Content Transport Equation";
    * this->name_base_variable = "membrane_water_content";
    *

These two variables are used to define the section where the information is stored in the parameter file for the equations.

template<int dim>
virtual FuelCellShop::Equation::EquationBase< dim >::~EquationBase ( )
protectedvirtual

Destructor.

Member Function Documentation

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

Assemble local Jacobian boundary matrix for Non-Linear problems.

Note
This member function is called by assemble_cell_matrix depending on the type of problem
template<int dim>
virtual void FuelCellShop::Equation::EquationBase< dim >::assemble_bdry_linear_matrix ( FuelCell::ApplicationCore::MatrixVector bdry_matrices,
const typename FuelCell::ApplicationCore::DoFApplication< dim >::FaceInfo &  bdry_info,
FuelCellShop::Layer::BaseLayer< dim > *const  layer 
)
inlineprotectedvirtual

Assemble local boundary matrix for linear problems.

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

Assemble local boundary RHS for linear problems.

template<int dim>
virtual void FuelCellShop::Equation::EquationBase< dim >::assemble_bdry_matrix ( FuelCell::ApplicationCore::MatrixVector bdry_matrices,
const typename FuelCell::ApplicationCore::DoFApplication< dim >::FaceInfo &  bdry_info,
FuelCellShop::Layer::BaseLayer< dim > *const  layer 
)
inlinevirtual
template<int dim>
virtual void FuelCellShop::Equation::EquationBase< dim >::assemble_bdry_residual ( FuelCell::ApplicationCore::FEVector bdry_residual,
const typename FuelCell::ApplicationCore::DoFApplication< dim >::FaceInfo &  bdry_info,
FuelCellShop::Layer::BaseLayer< dim > *const  layer 
)
inlinevirtual
template<int dim>
virtual void FuelCellShop::Equation::EquationBase< dim >::assemble_cell_Jacobian_matrix ( FuelCell::ApplicationCore::MatrixVector cell_matrices,
const typename FuelCell::ApplicationCore::DoFApplication< dim >::CellInfo &  cell_info,
FuelCellShop::Layer::BaseLayer< dim > *const  layer 
)
inlineprotectedvirtual

Assemble the local Jacobian Matrix for Non-Linear problems.

Note
This member function is called by assemble_cell_matrix depending on the type of problem
template<int dim>
virtual void FuelCellShop::Equation::EquationBase< dim >::assemble_cell_linear_matrix ( FuelCell::ApplicationCore::MatrixVector cell_matrices,
const typename FuelCell::ApplicationCore::DoFApplication< dim >::CellInfo &  cell_info,
FuelCellShop::Layer::BaseLayer< dim > *const  layer 
)
inlineprotectedvirtual

Assemble the local cell matrix for Linear problems.

Note
This member function is called by assemble_cell_matrix depending on the type of problem
template<int dim>
virtual void FuelCellShop::Equation::EquationBase< dim >::assemble_cell_linear_rhs ( FuelCell::ApplicationCore::FEVector cell_residual,
const typename FuelCell::ApplicationCore::DoFApplication< dim >::CellInfo &  cell_info,
FuelCellShop::Layer::BaseLayer< dim > *const  layer 
)
inlineprotectedvirtual

Assemble local cell RHS for Linear problems.

Note
This member function is called by assemble_cell_residual depending on the type of problem
template<int dim>
virtual void FuelCellShop::Equation::EquationBase< 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
template<int dim>
virtual void FuelCellShop::Equation::EquationBase< 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
template<int dim>
virtual void FuelCellShop::Equation::EquationBase< dim >::assemble_cell_residual_rhs ( FuelCell::ApplicationCore::FEVector cell_residual,
const typename FuelCell::ApplicationCore::DoFApplication< dim >::CellInfo &  cell_info,
FuelCellShop::Layer::BaseLayer< dim > *const  layer 
)
inlineprotectedvirtual

Assemble local cell RHS for nonlinear problems.

Note
This member function is called by assemble_cell_residual depending on the type of problem
template<int dim>
bool FuelCellShop::Equation::EquationBase< dim >::belongs_to_boundary ( const unsigned int &  tria_boundary_id,
const unsigned int &  param_boundary_id 
) const
inlineprotected

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.

template<int dim>
void FuelCellShop::Equation::EquationBase< dim >::dealII_to_appframe ( FuelCell::ApplicationCore::MatrixVector dst,
const FullMatrix< double > &  src,
const std::vector< unsigned int > &  matrix_block_indices 
) const
protected

This function converts the standard ordered structure dealii::FullMatrix<double> src into the block-wise ordered structure FuelCell::ApplicationCore::MatrixVector dst.

matrix_block_indices is 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).

template<int dim>
void FuelCellShop::Equation::EquationBase< dim >::dealII_to_appframe ( FuelCell::ApplicationCore::FEVector dst,
const Vector< double > &  src,
const std::vector< unsigned int > &  residual_indices 
) const
protected

This function converts the standard ordered structure dealii::Vector<double> src into the block-wise ordered structure FuelCell::ApplicationCore::FEVector dst.

residual_indices is the residual indices (a derived equation class) drawn from the global structure (a derived equation class + other active equation classes included into the computation).

template<int dim>
virtual void FuelCellShop::Equation::EquationBase< dim >::declare_parameters ( ParameterHandler &  param) const
protectedvirtual
template<int dim>
const std::vector< BoundaryType >& FuelCellShop::Equation::EquationBase< dim >::get_boundary_types ( ) const
inline

This function returns boundary_types of a derived equation class.

template<int dim>
const component_boundaryID_value_map& FuelCellShop::Equation::EquationBase< dim >::get_component_boundaryID_value ( ) const
inline

This function returns component_boundaryID_value of a derived equation class.

template<int dim>
const component_materialID_value_map& FuelCellShop::Equation::EquationBase< dim >::get_component_materialID_value ( ) const
inline

This function returns component_materialID_value of a derived equation class.

template<int dim>
const std::string& FuelCellShop::Equation::EquationBase< dim >::get_equation_name ( ) const
inline

This function returns equation_name of a derived equation class.

template<int dim>
const couplings_map& FuelCellShop::Equation::EquationBase< dim >::get_internal_cell_couplings ( ) const
inline

This function returns internal_cell_couplings of a derived equation class.

template<int dim>
const couplings_map& FuelCellShop::Equation::EquationBase< dim >::get_internal_flux_couplings ( ) const
inline

This function returns internal_flux_couplings (DG FEM only) of a derived equation class.

template<int dim>
const std::vector<unsigned int>& FuelCellShop::Equation::EquationBase< dim >::get_matrix_block_indices ( ) const
inline

This function returns matrix_block_indices of a derived equation class.

template<int dim>
const std::vector< std::vector< BoundaryType > >& FuelCellShop::Equation::EquationBase< dim >::get_multi_boundary_types ( ) const
inline

This function returns multi_boundary_types of a derived equation class.

template<int dim>
const std::vector< std::vector< OutputType > >& FuelCellShop::Equation::EquationBase< dim >::get_multi_output_types ( ) const
inline

This function returns multi_output_types of a derived equation class.

template<int dim>
const std::vector< OutputType >& FuelCellShop::Equation::EquationBase< dim >::get_output_types ( ) const
inline

This function returns output_types of a derived equation class.

template<int dim>
const std::vector<unsigned int>& FuelCellShop::Equation::EquationBase< dim >::get_residual_indices ( ) const
inline

This function returns residual_indices of a derived equation class.

template<int dim>
virtual void FuelCellShop::Equation::EquationBase< dim >::initialize ( ParameterHandler &  param)
protectedvirtual
template<int dim>
virtual void FuelCellShop::Equation::EquationBase< dim >::make_assemblers_bdry_constant_data ( const typename FuelCell::ApplicationCore::DoFApplication< dim >::FaceInfo &  bdry_info)
inlineprotectedvirtual
template<int dim>
virtual void FuelCellShop::Equation::EquationBase< dim >::make_assemblers_bdry_variable_data ( const typename FuelCell::ApplicationCore::DoFApplication< dim >::FaceInfo &  bdry_info,
FuelCellShop::Layer::BaseLayer< dim > *const  layer 
)
inlineprotectedvirtual
template<int dim>
virtual void FuelCellShop::Equation::EquationBase< dim >::make_assemblers_cell_constant_data ( const typename FuelCell::ApplicationCore::DoFApplication< dim >::CellInfo &  cell_info)
inlineprotectedvirtual
template<int dim>
virtual void FuelCellShop::Equation::EquationBase< dim >::make_assemblers_cell_variable_data ( const typename FuelCell::ApplicationCore::DoFApplication< dim >::CellInfo &  cell_info,
FuelCellShop::Layer::BaseLayer< dim > *const  layer 
)
inlineprotectedvirtual
template<int dim>
virtual void FuelCellShop::Equation::EquationBase< dim >::make_assemblers_generic_constant_data ( )
inlineprotectedvirtual

Function used to initialize variable information that will be needed to assemble matrix and residual and that will remain constant throughout the simulation.

The data should not depend on the cell. If it does, then initialize this informatin using make_assemblers_cell_constant_data. Usually VariableInfo objects for all variables used in the class are initialized here. It assigns solution_index and fe_type.

Note
This function is called the first time the residual/matrix is computed. Note it cannot be called in initialize() since the FiniteElement objects have not yet been created.
This function is usually re-implemented in the derived equation classes based on the necessary solution variables.

Reimplemented in FuelCellShop::Equation::ThermalTransportEquation< dim >, FuelCellShop::Equation::NewFicksTransportEquation< dim >, FuelCellShop::Equation::LiquidSourceEquation< dim >, FuelCellShop::Equation::ProtonTransportEquation< dim >, FuelCellShop::Equation::ElectronTransportEquation< dim >, FuelCellShop::Equation::LiquidPressureEquation< dim >, FuelCellShop::Equation::LambdaTransportEquation< dim >, FuelCellShop::Equation::SaturationTransportEquation< dim >, FuelCellShop::Equation::SorptionSourceTerms< dim >, FuelCellShop::Equation::ReactionSourceTerms< dim >, FuelCellShop::Equation::ReactionSourceTermsKG< dim >, and FuelCellShop::Equation::CapillaryPressureEquation< dim >.

template<int dim>
virtual void FuelCellShop::Equation::EquationBase< dim >::make_boundary_types ( )
inlineprotectedvirtual
template<int dim>
virtual void FuelCellShop::Equation::EquationBase< dim >::make_component_boundaryID_value ( )
inlineprotectedvirtual

This function fills out component_boundaryID_value of a derived equation class.

template<int dim>
virtual void FuelCellShop::Equation::EquationBase< dim >::make_component_materialID_value ( )
inlineprotectedvirtual

This function fills out component_materialID_value of a derived equation class.

template<int dim>
virtual void FuelCellShop::Equation::EquationBase< dim >::make_internal_cell_couplings ( )
inlineprotectedvirtual
template<int dim>
virtual void FuelCellShop::Equation::EquationBase< dim >::make_internal_flux_couplings ( )
inlineprotectedvirtual

This function fills out internal_flux_couplings (DG FEM only) of a derived equation class.

template<int dim>
virtual void FuelCellShop::Equation::EquationBase< dim >::make_matrix_block_indices ( )
inlineprotectedvirtual

This function is only needed to provide the last argument to dealII_to_appframe.

It only need to be called once during initialization.

This function fills out matrix_block_indices of a derived equation class.

Reimplemented in FuelCellShop::Equation::CompressibleMultiComponentKGEquationsCoupled< dim >, FuelCellShop::Equation::IncompressibleSingleComponentNSEquations< dim >, and FuelCellShop::Equation::ReactionSourceTermsKG< dim >.

template<int dim>
virtual void FuelCellShop::Equation::EquationBase< dim >::make_multi_boundary_types ( )
inlineprotectedvirtual

This function fills out multi_boundary_types of a derived equation class.

template<int dim>
virtual void FuelCellShop::Equation::EquationBase< dim >::make_multi_output_types ( )
inlineprotectedvirtual

This function fills out multi_output_types of a derived equation class.

template<int dim>
virtual void FuelCellShop::Equation::EquationBase< dim >::make_output_types ( )
inlineprotectedvirtual
template<int dim>
virtual void FuelCellShop::Equation::EquationBase< dim >::make_residual_indices ( )
inlineprotectedvirtual

This function is only needed to provide the last argument to dealII_to_appframe.

It only need to be called once during initialization.

This function fills out residual_indices of a derived equation class.

Reimplemented in FuelCellShop::Equation::CompressibleMultiComponentKGEquationsCoupled< dim >, FuelCellShop::Equation::IncompressibleSingleComponentNSEquations< dim >, and FuelCellShop::Equation::ReactionSourceTermsKG< dim >.

template<int dim>
void FuelCellShop::Equation::EquationBase< dim >::print_caller_name ( const std::string &  caller_name) const
protected

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.

template<int dim>
virtual void FuelCellShop::Equation::EquationBase< dim >::print_equation_info ( ) const
inlinevirtual
template<int dim>
void FuelCellShop::Equation::EquationBase< dim >::select_cell_assemblers ( const typename FuelCell::ApplicationCore::DoFApplication< dim >::CellInfo &  cell_info,
FuelCellShop::Layer::BaseLayer< dim > *const  layer 
)
inlineprotected

This routine is used to select the make_assembly routines that need to be called inside assemble_cell_matrix to compute.

  • generic
  • constant_cell_data
  • variable_cell_data and resets the flags as appropriate
template<int dim>
virtual void FuelCellShop::Equation::EquationBase< dim >::set_parameters ( const std::vector< std::string > &  name_dvar,
const std::vector< double > &  value_dvar,
ParameterHandler &  param 
)
inlineprotectedvirtual

Set parameters using the parameter file, in order to run parametric/optimization studies.

Reimplemented in FuelCellShop::Equation::ElectronTransportEquation< dim >.

template<int dim>
void FuelCellShop::Equation::EquationBase< dim >::standard_to_block_wise ( FullMatrix< double > &  target) const
protected

This function changes the order of dealii::FullMatrix<double> target from standard to block-wise.

template<int dim>
void FuelCellShop::Equation::EquationBase< dim >::standard_to_block_wise ( Vector< double > &  target) const
protected

This function changes the order of dealii::Vector<double> target from standard to block-wise.

Member Data Documentation

template<int dim>
EquationFlags FuelCellShop::Equation::EquationBase< dim >::assemble_flags
protected

This vector contains a collection of internal flags to tell derived equation classes what needs to be re-computed.

EquationFlags can be easily extended.

Note
replaces counter
template<int dim>
DoFHandler<dim>::active_face_iterator FuelCellShop::Equation::EquationBase< dim >::bdry
protected

Currently active DoFHandler<dim> active boundary iterator.

template<int dim>
std::vector< BoundaryType > FuelCellShop::Equation::EquationBase< dim >::boundary_types
protected

The list of boundary types of a derived equation class.

template<int dim>
DoFHandler<dim>::active_cell_iterator FuelCellShop::Equation::EquationBase< dim >::cell
protected

Currently active DoFHandler<dim> active cell iterator.

template<int dim>
component_boundaryID_value_map FuelCellShop::Equation::EquationBase< dim >::component_boundaryID_value
protected

This object reflects the following structure (see FuelCell::InitialAndBoundaryData namespace docs):

  • first argument : name of the solution component,
  • second argument : boundary id,
  • third argument : value of the solution component.
template<int dim>
component_materialID_value_map FuelCellShop::Equation::EquationBase< dim >::component_materialID_value
protected

This object reflects the following structure (see FuelCell::InitialAndBoundaryData namespace docs):

  • first argument : name of the solution component,
  • second argument : material id,
  • third argument : value of the solution component.
template<int dim>
std::vector<bool> FuelCellShop::Equation::EquationBase< dim >::counter
protected

This vector contains the collection of internal "counters" used by the derived equation classes.

Deprecated:
use assemble_flags instead
template<int dim>
boost::shared_ptr< FuelCell::ApplicationCore::ApplicationData > FuelCellShop::Equation::EquationBase< dim >::data
protected

Data object for the application data to be passed to the equation classes.

template<int dim>
unsigned int FuelCellShop::Equation::EquationBase< dim >::dofs_per_cell
protected

Number of degrees of freedom per cell.

template<int dim>
std::string FuelCellShop::Equation::EquationBase< dim >::equation_name
protected
template<int dim>
couplings_map FuelCellShop::Equation::EquationBase< dim >::internal_cell_couplings
protected

This object contains the info on how the equations and solution variables of a derived equation class are coupled.

template<int dim>
couplings_map FuelCellShop::Equation::EquationBase< dim >::internal_flux_couplings
protected

This object contains the info on how the "X" and "Y" of a derived equation class are coupled (DG FEM only).

template<int dim>
std::vector<double> FuelCellShop::Equation::EquationBase< dim >::JxW_bdry
protected

Jacobian of mapping by Weight in the quadrature points of a boundary.

template<int dim>
std::vector<double> FuelCellShop::Equation::EquationBase< dim >::JxW_cell
protected

Jacobian of mapping by Weight in the quadrature points of a cell.

template<int dim>
std::vector<unsigned int> FuelCellShop::Equation::EquationBase< dim >::matrix_block_indices
protected

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).

template<int dim>
std::vector< std::vector< BoundaryType > > FuelCellShop::Equation::EquationBase< dim >::multi_boundary_types
protected

The list of multiple boundary types of a derived equation class.

template<int dim>
std::vector< std::vector< OutputType > > FuelCellShop::Equation::EquationBase< dim >::multi_output_types
protected

The list of multiple output types of a derived equation class.

template<int dim>
unsigned int FuelCellShop::Equation::EquationBase< dim >::n_q_points_bdry
protected

Number of quadrature points per boundary.

template<int dim>
unsigned int FuelCellShop::Equation::EquationBase< dim >::n_q_points_cell
protected
template<int dim>
std::string FuelCellShop::Equation::EquationBase< dim >::name_base_variable
protected

Const std::string member storing name of the base solution variable corresponding to the equation represented by this class.

Referenced by FuelCellShop::Equation::NewFicksTransportEquation< dim >::set_solute_and_solvent(), and FuelCellShop::Equation::FicksTransportEquation< dim >::set_solute_and_solvent().

template<int dim>
std::vector< Point<dim> > FuelCellShop::Equation::EquationBase< dim >::normal_vectors
protected

Normal vectors in the quadrature points of a boundary.

template<int dim>
std::vector< OutputType > FuelCellShop::Equation::EquationBase< dim >::output_types
protected

The list of output types of a derived equation class.

template<int dim>
std::vector<unsigned int> FuelCellShop::Equation::EquationBase< dim >::residual_indices
protected

The residual indices (a derived equation class) drawn from the global structure (a derived equation class + other active equation classes included into the computation).

template<int dim>
std::string FuelCellShop::Equation::EquationBase< dim >::residual_vector_name
protected

The name of the residual vector name in FEVectors.

template<int dim>
std::string FuelCellShop::Equation::EquationBase< dim >::solution_vector_name
protected

The name of the solution vector in FEVectors.

template<int dim>
FuelCell::SystemManagement* FuelCellShop::Equation::EquationBase< dim >::system_management
protected

Pointer to the external YourApplication<dim>::system_management object.

template<int dim>
std::vector< std::vector< Point<dim> > > FuelCellShop::Equation::EquationBase< dim >::tangential_vectors
protected

Tangential vectors in the quadrature points of a boundary.

template<int dim>
bool FuelCellShop::Equation::EquationBase< dim >::variable_boundary_data

true, if variable Dirichlet boundary conditions are prescribed on a part of the boundary.

false, otherwise.

template<int dim>
bool FuelCellShop::Equation::EquationBase< dim >::variable_initial_data

true, if variable initial data is prescribed on a part of the domain.

false, otherwise.


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