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

This class deals with Proton Transport Equation. More...

#include <proton_transport_equation.h>

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

Public Member Functions

Constructors, destructor, and initalization
 ProtonTransportEquation (FuelCell::SystemManagement &system_management)
 Constructor.
 
virtual ~ProtonTransportEquation ()
 Destructor.
 
virtual void declare_parameters (ParameterHandler &param) const
 Declare parameters.
 
virtual void initialize (ParameterHandler &param)
 Initialize parameters.
 
Local CG FEM based assemblers
virtual void assemble_cell_matrix (AppFrame::MatrixVector &cell_matrices, const typename AppFrame::DoFApplication< dim >::CellInfo &cell_info, FuelCellShop::Layer::BaseLayer< dim > *const layer)
 Assemble local cell matrix.
 
virtual void assemble_cell_residual (AppFrame::FEVector &cell_rhs, const typename AppFrame::DoFApplication< dim >::CellInfo &cell_info, FuelCellShop::Layer::BaseLayer< dim > *const layer)
 Assemble local cell residual.
 
virtual void assemble_bdry_matrix (AppFrame::MatrixVector &bdry_matrices, const typename AppFrame::DoFApplication< dim >::FaceInfo &bdry_info, FuelCellShop::Layer::BaseLayer< dim > *const layer)
 Assemble local boundary matrix.
 
virtual void assemble_bdry_residual (AppFrame::FEVector &bdry_rhs, const typename AppFrame::DoFApplication< dim >::FaceInfo &bdry_info, FuelCellShop::Layer::BaseLayer< dim > *const layer)
 Assemble local boundary residual.
 
Accessors and info
virtual void print_equation_info () const
 The function prints out the equation's info.
 
std::map< unsigned int, double > get_dirichlet_bdry_map () const
 Method to provide access to Dirichlet boundary conditions map filled using the parameter file.
 
void class_test ()
 This member function creates an object of its own type and runs test to diagnose if there are any problems with the routines that you have created.
 
- 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.
 
const couplings_mapget_internal_flux_couplings () const
 This function returns internal_flux_couplings (DG FEM only) of a derived equation class.
 
const
component_materialID_value_map
get_component_materialID_value () const
 This function returns component_materialID_value of a derived equation class.
 
const
component_boundaryID_value_map
get_component_boundaryID_value () const
 This function returns component_boundaryID_value of a derived equation class.
 
const std::vector< BoundaryType > & get_boundary_types () const
 This function returns boundary_types of a derived equation class.
 
const std::vector< std::vector
< BoundaryType > > & 
get_multi_boundary_types () const
 This function returns multi_boundary_types of a derived equation class.
 
const std::vector< OutputType > & get_output_types () const
 This function returns output_types of a derived equation class.
 
const std::vector< std::vector
< OutputType > > & 
get_multi_output_types () const
 This function returns multi_output_types of a derived equation class.
 
const std::string & get_equation_name () const
 This function returns equation_name of a derived equation class.
 
const std::vector< unsigned int > & get_matrix_block_indices () const
 This function returns matrix_block_indices of a derived equation class.
 
const std::vector< unsigned int > & get_residual_indices () const
 This function returns residual_indices of a derived equation class.
 

Protected Member Functions

Local CG FEM based assemblers - make_ functions
virtual void make_assemblers_generic_constant_data ()
 This function computes Local CG FEM based assemblers - constant data (generic).
 
virtual void make_assemblers_cell_constant_data (const typename AppFrame::DoFApplication< dim >::CellInfo &cell_info)
 This function computes Local CG FEM based assemblers - constant data (cell) and allocates the memory for shape functions, shape function gradients, and JxW_cell in Local CG FEM based assemblers - variable data (cell) .
 
virtual void make_assemblers_bdry_constant_data (const typename AppFrame::DoFApplication< dim >::FaceInfo &bdry_info)
 This function computes Local CG FEM based assemblers - constant data (boundary) and allocates the memory for shape functions, normal_vectors, and JxW_bdry in Local CG FEM based assemblers - variable data (boundary) .
 
virtual void make_assemblers_cell_variable_data (const typename AppFrame::DoFApplication< dim >::CellInfo &cell_info, FuelCellShop::Layer::BaseLayer< dim > *const layer)
 This function computes Local CG FEM based assemblers - variable data (cell) .
 
virtual void make_assemblers_bdry_variable_data (const typename AppFrame::DoFApplication< dim >::FaceInfo &bdry_info, FuelCellShop::Layer::BaseLayer< dim > *const layer)
 This function computes Local CG FEM based assemblers - variable data (boundary) .
 
Other make_ functions
virtual void make_internal_cell_couplings ()
 This function fills out internal_cell_couplings.
 
virtual void make_boundary_types ()
 This function fills out boundary_types.
 
virtual void make_output_types ()
 This function fills out output_types.
 
- Protected Member Functions inherited from FuelCellShop::Equation::EquationBase< dim >
 EquationBase (FuelCell::SystemManagement &system_management)
 Constructor.
 
virtual ~EquationBase ()
 Destructor.
 
virtual void make_internal_flux_couplings ()
 This function fills out internal_flux_couplings (DG FEM only) of a derived equation class.
 
virtual void make_component_materialID_value ()
 This function fills out component_materialID_value of a derived equation class.
 
virtual void make_component_boundaryID_value ()
 This function fills out component_boundaryID_value of a derived equation class.
 
virtual void make_multi_boundary_types ()
 This function fills out multi_boundary_types of a derived equation class.
 
virtual void make_multi_output_types ()
 This function fills out multi_output_types of a derived equation class.
 
virtual void make_matrix_block_indices ()
 This function fills out matrix_block_indices of a derived equation class.
 
virtual void make_residual_indices ()
 This function fills out residual_indices of a derived equation class.
 
void standard_to_block_wise (FullMatrix< double > &target) const
 This function changes the order of dealii::FullMatrix<double> target from standard to block-wise.
 
void standard_to_block_wise (Vector< double > &target) const
 This function changes the order of dealii::Vector<double> target from standard to block-wise.
 
void dealII_to_appframe (AppFrame::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 AppFrame::MatrixVector dst.
 
void dealII_to_appframe (AppFrame::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 AppFrame::FEVector dst.
 
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.
 
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.
 

Protected Attributes

bool cell_residual_counter
 Counter set to TRUE when cell_residual is being assembled.
 
bool bdry_residual_counter
 Counter set to TRUE when bdry_residual is being assembled.
 
unsigned int last_iter_cell
 Variable used to store the index in cell_info->global_data of the previous Newton solution The solution at the previous iteration is used to compute cell_matrix and cell_residual.
 
unsigned int last_iter_bdry
 Variable used to store the index in bdry_info->global_data of the previous Newton solution The solution at the previous iteration is used to compute bdry_matrix and bdry_residual.
 
Boundary conditions
std::map< unsigned int, double > dirichlet_bdry_map
 std::map< unsigned int, double > container for Dirichlet boundary conditions.
 
std::map< unsigned int, double > proton_current_flux_map
 std::map< unsigned int, double > container for details regarding Galvanostatic boundary conditions.
 
Generic Constant Data
const std::string name_base_variable
 Const std::string member storing name of the base solution variable corresponding to the equation represented by this class.
 
VariableInfo phi_m
 VariableInfo structure corresponding to "protonic_electrical_potential".
 
VariableInfo lambda
 VariableInfo structure corresponding to "membrane_water_content".
 
VariableInfo t_rev
 VariableInfo structure corresponding to "temperature_of_REV".
 
Local CG FEM based assemblers - variable data (cell)
std::vector< double > sigmaMeff_cell
 Effective proton conductivity, [S/cm], at all quadrature points of the cell.
 
std::vector< double > dsigmaMeff_dT_cell
 Derivative of effective protonic conductivity w.r.t "temperature_of_REV", at all quadrature points in the cell.
 
std::vector< double > dsigmaMeff_dlambda_cell
 Derivative of effective protonic conductivity w.r.t.
 
std::vector< std::vector
< double > > 
phi_phiM_cell
 \( \mathbf{\phi_{M}} \) shape functions.
 
std::vector< std::vector
< Tensor< 1, dim > > > 
grad_phi_phiM_cell
 \( \mathbf{\phi_m} \) shape function gradients.
 
std::vector< std::vector
< Tensor< 1, dim > > > 
grad_phi_T_cell
 \( \mathbf{T} \) shape function gradients.
 
std::vector< std::vector
< Tensor< 1, dim > > > 
grad_phi_lambda_cell
 \( \mathbf{\lambda} \) shape function gradients.
 
Local CG FEM based assemblers - variable data (boundary)
std::vector< double > dsigmaMeff_dT_bdry
 Derivative of effective protonic conductivity w.r.t "temperature_of_REV", at all quadrature points in the boundary.
 
std::vector< double > dsigmaMeff_dlambda_bdry
 Derivative of effective protonic conductivity w.r.t.
 
std::vector< std::vector
< double > > 
phi_phiM_bdry
 \( \mathbf{\phi_m} \) shape functions.
 
std::vector< std::vector
< double > > 
phi_lambda_bdry
 \( \mathbf{\lambda} \) shape functions.
 
std::vector< std::vector
< double > > 
phi_T_bdry
 \( \mathbf{T} \) shape functions.
 
- Protected Attributes inherited from FuelCellShop::Equation::EquationBase< dim >
unsigned int dofs_per_cell
 Number of degrees of freedom per cell.
 
unsigned int n_q_points_cell
 Number of quadrature points per cell.
 
unsigned int n_q_points_bdry
 Number of quadrature points per boundary.
 
DoFHandler< dim >
::active_cell_iterator 
cell
 Currently active DoFHandler<dim> active cell iterator.
 
DoFHandler< dim >
::active_face_iterator 
bdry
 Currently active DoFHandler<dim> active boundary iterator.
 
std::vector< double > JxW_cell
 Jacobian of mapping by Weight in the quadrature points of a cell.
 
std::vector< double > JxW_bdry
 Jacobian of mapping by Weight in the quadrature points of a boundary.
 
std::vector< Point< dim > > normal_vectors
 Normal vectors in the quadrature points of a boundary.
 
std::vector< std::vector
< Point< dim > > > 
tangential_vectors
 Tangential vectors in the quadrature points of a boundary.
 
FuelCell::SystemManagementsystem_management
 Pointer to the external YourApplication<dim>::system_management object.
 
couplings_map internal_cell_couplings
 This object contains the info on how the equations and solution variables of a derived equation class are coupled.
 
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).
 
component_materialID_value_map component_materialID_value
 This object reflects the following structure (see FuelCell::InitialAndBoundaryData namespace docs):
 
component_boundaryID_value_map component_boundaryID_value
 This object reflects the following structure (see FuelCell::InitialAndBoundaryData namespace docs):
 
std::vector< BoundaryTypeboundary_types
 The list of boundary types of a derived equation class.
 
std::vector< std::vector
< BoundaryType > > 
multi_boundary_types
 The list of multiple boundary types of a derived equation class.
 
std::vector< OutputTypeoutput_types
 The list of output types of a derived equation class.
 
std::vector< std::vector
< OutputType > > 
multi_output_types
 The list of multiple output types of a derived equation class.
 
std::string equation_name
 The name of a derived equation class.
 
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).
 
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).
 
std::vector< bool > counter
 This vector contains the collection of internal "counters" used by the derived equation classes.
 

Detailed Description

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

This class deals with Proton Transport Equation.

This equation solves a steady-state Ohm's law for proton transport inside the polymer electrolyte. Note that Ohm's law can be used in the polymer electrolyte because electro-neutrality is assumed due to the constant number of negative ionic charges that are attached to the electrolyte. Due to electro-neutrality, convective and diffusive current terms in the Nernst-Planck Equation become negligible, hence it simplifies to Ohm's Law.

It is solved with respect to:

where, \( \mathbf{\phi_{m}} \) is in Voltages.

This equation can be written as:

\( \qquad \mathbf{\nabla} \cdot \left( \hat{\sigma}_{m,eff} \cdot \mathbf{\nabla} \phi_m \right) = 0 \quad \in \quad \Omega \)

To be well-posed, these equations are equipped with the appropriate boundary conditions. All the boundary conditions can be described by boundary_id (s ) and boundary_type. Besides, this some boundary types require additional information, which can also be provided by the parameter file. We consider several types of boundary conditions:

These are specified in the parameter file under subsection "Constant Proton Current Flux Boundary Conditions", as a list of comma-separated map-like values.

e.g.

set Constant Proton Current Flux Boundary Conditions = 1: 2 , 4: -0.5

where, boundary_id "1" has constant proton current flux value of 2 [A/cm^2] (leaving out of the boundary) and boundary_id "4" has constant proton current flux value of -0.5 [A/cm^2] (entering into the boundary).

Remarks
  • There is no provision to specify boundary indicators for No current flux or Symmetric boundary conditions, as FEM formulation automatically implies a particular boundary is one of these cases, by default.
  • This class works with the following layer classes only:
    • FuelCellShop::Layer::CatalystLayer<dim>
    • FuelCellShop::Layer::MembraneLayer<dim>

We solve the whole problem by linearizing the governing equation at each Newton iteration with subsequent CG FEM discretization in space. The class contains the necessary class members to add the necessary contributions to cell_matrix and cell_residual to the governing equations used to analyze proton transport by Ohm's law,

Usage Details:

// Creating Equation object (in Application Header file)
// Declare parameters in application
proton_transport.declare_parameters(param);
// Initialize in application
proton_transport.initialize(param);
// Create a temporary vector in the application for storing couplings_map from all the equation used in the application.
std::vector<couplings_map> tmp;
... // other equations
tmp.push_back( proton_transport.get_internal_cell_couplings() );
// Look at ReactionSourceTerms class here, if source terms due to electrochemical reaction are to be considered.
// Making cell couplings using SystemManagement object created in the application
// cell_matrix in application
// Do a check against layer and it should match with the layers currently working for this equation class.
// for eg: CCL is FuelCellShop::Layer::HomogeneousCL<dim> object.
proton_transport.assemble_cell_matrix(cell_matrices, cell_info, &CCL);
// cell_residual in application
proton_transport.assemble_cell_residual(cell_vector, cell_info, &CCL);
Note
This class doesn't assemble for reaction source terms; that is taken care off by ReactionSourceTerms class. Please read the documentation of ReactionSourceTerms class, for additional methods to be implemented in the application.
Warning
If current source terms due to reaction are being considered, it's very important to use adjust_internal_cell_couplings member function of ReactionSourceTerms class, before using make_cell_couplings of SystemManagement at the application level.
Author
Madhur Bhaiya, 2013

Constructor & Destructor Documentation

Constructor.

Destructor.

Member Function Documentation

template<int dim>
virtual void FuelCellShop::Equation::ProtonTransportEquation< dim >::assemble_bdry_matrix ( AppFrame::MatrixVector bdry_matrices,
const typename AppFrame::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::ProtonTransportEquation< dim >::assemble_bdry_residual ( AppFrame::FEVector bdry_rhs,
const typename AppFrame::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::ProtonTransportEquation< dim >::assemble_cell_matrix ( AppFrame::MatrixVector cell_matrices,
const typename AppFrame::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::ProtonTransportEquation< dim >::assemble_cell_residual ( AppFrame::FEVector cell_rhs,
const typename AppFrame::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>
void FuelCellShop::Equation::ProtonTransportEquation< dim >::class_test ( )

This member function creates an object of its own type and runs test to diagnose if there are any problems with the routines that you have created.

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

Declare parameters.

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

template<int dim>
std::map<unsigned int, double> FuelCellShop::Equation::ProtonTransportEquation< dim >::get_dirichlet_bdry_map ( ) const
inline

Method to provide access to Dirichlet boundary conditions map filled using the parameter file.

This function is useful in accessing the dircihlet values to be set in initial conditions methods of the application.

References FuelCellShop::Equation::ProtonTransportEquation< dim >::dirichlet_bdry_map.

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

Initialize parameters.

This class will call make_internal_cell_couplings and make_boundary_types.

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

template<int dim>
virtual void FuelCellShop::Equation::ProtonTransportEquation< dim >::make_assemblers_bdry_constant_data ( const typename AppFrame::DoFApplication< dim >::FaceInfo &  bdry_info)
protectedvirtual

This function computes Local CG FEM based assemblers - constant data (boundary) and allocates the memory for shape functions, normal_vectors, and JxW_bdry in Local CG FEM based assemblers - variable data (boundary) .

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

template<int dim>
virtual void FuelCellShop::Equation::ProtonTransportEquation< dim >::make_assemblers_bdry_variable_data ( const typename AppFrame::DoFApplication< dim >::FaceInfo &  bdry_info,
FuelCellShop::Layer::BaseLayer< dim > *const  layer 
)
protectedvirtual

This function computes Local CG FEM based assemblers - variable data (boundary) .

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

template<int dim>
virtual void FuelCellShop::Equation::ProtonTransportEquation< dim >::make_assemblers_cell_constant_data ( const typename AppFrame::DoFApplication< dim >::CellInfo &  cell_info)
protectedvirtual

This function computes Local CG FEM based assemblers - constant data (cell) and allocates the memory for shape functions, shape function gradients, and JxW_cell in Local CG FEM based assemblers - variable data (cell) .

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

template<int dim>
virtual void FuelCellShop::Equation::ProtonTransportEquation< dim >::make_assemblers_cell_variable_data ( const typename AppFrame::DoFApplication< dim >::CellInfo &  cell_info,
FuelCellShop::Layer::BaseLayer< dim > *const  layer 
)
protectedvirtual

This function computes Local CG FEM based assemblers - variable data (cell) .

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

template<int dim>
virtual void FuelCellShop::Equation::ProtonTransportEquation< dim >::make_assemblers_generic_constant_data ( )
protectedvirtual

This function computes Local CG FEM based assemblers - constant data (generic).

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

template<int dim>
virtual void FuelCellShop::Equation::ProtonTransportEquation< dim >::make_boundary_types ( )
protectedvirtual

This function fills out boundary_types.

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

template<int dim>
virtual void FuelCellShop::Equation::ProtonTransportEquation< 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::ProtonTransportEquation< dim >::make_output_types ( )
inlineprotectedvirtual

This function fills out output_types.

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

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

The function prints out the equation's info.

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

Member Data Documentation

template<int dim>
bool FuelCellShop::Equation::ProtonTransportEquation< dim >::bdry_residual_counter
protected

Counter set to TRUE when bdry_residual is being assembled.

This ensures that only data required for bdry_residual is computed. (improves speed)

template<int dim>
bool FuelCellShop::Equation::ProtonTransportEquation< dim >::cell_residual_counter
protected

Counter set to TRUE when cell_residual is being assembled.

This ensures that only effective protonic conductivity is being calculated, not its derivatives. (improves speed)

template<int dim>
std::map<unsigned int, double> FuelCellShop::Equation::ProtonTransportEquation< dim >::dirichlet_bdry_map
protected

std::map< unsigned int, double > container for Dirichlet boundary conditions.

Here, Key (unsigned int ) represents the boundary_id and Value (double) represents the potential [V] .

Referenced by FuelCellShop::Equation::ProtonTransportEquation< dim >::get_dirichlet_bdry_map().

template<int dim>
std::vector<double> FuelCellShop::Equation::ProtonTransportEquation< dim >::dsigmaMeff_dlambda_bdry
protected

Derivative of effective protonic conductivity w.r.t.

"membrane_water_content", at all quadrature points in the boundary.

template<int dim>
std::vector<double> FuelCellShop::Equation::ProtonTransportEquation< dim >::dsigmaMeff_dlambda_cell
protected

Derivative of effective protonic conductivity w.r.t.

"membrane_water_content", at all quadrature points in the cell.

template<int dim>
std::vector<double> FuelCellShop::Equation::ProtonTransportEquation< dim >::dsigmaMeff_dT_bdry
protected

Derivative of effective protonic conductivity w.r.t "temperature_of_REV", at all quadrature points in the boundary.

template<int dim>
std::vector<double> FuelCellShop::Equation::ProtonTransportEquation< dim >::dsigmaMeff_dT_cell
protected

Derivative of effective protonic conductivity w.r.t "temperature_of_REV", at all quadrature points in the cell.

template<int dim>
std::vector< std::vector< Tensor<1,dim> > > FuelCellShop::Equation::ProtonTransportEquation< dim >::grad_phi_lambda_cell
protected

\( \mathbf{\lambda} \) shape function gradients.

grad_phi_lambda_cell [ q ] [ k ] denotes \( k \)-th \( \mathbf{\lambda} \) shape function gradient computed in \( q \)-th quadrature point of the cell.

template<int dim>
std::vector< std::vector< Tensor<1,dim> > > FuelCellShop::Equation::ProtonTransportEquation< dim >::grad_phi_phiM_cell
protected

\( \mathbf{\phi_m} \) shape function gradients.

grad_phi_phiM_cell [ q ] [ k ] denotes \( k \)-th \( \mathbf{\phi_m} \) shape function gradient computed in \( q \)-th quadrature point of the cell.

template<int dim>
std::vector< std::vector< Tensor<1,dim> > > FuelCellShop::Equation::ProtonTransportEquation< dim >::grad_phi_T_cell
protected

\( \mathbf{T} \) shape function gradients.

grad_phi_T_cell [ q ] [ k ] denotes \( k \)-th \( \mathbf{T} \) shape function gradient computed in \( q \)-th quadrature point of the cell.

template<int dim>
VariableInfo FuelCellShop::Equation::ProtonTransportEquation< dim >::lambda
protected

VariableInfo structure corresponding to "membrane_water_content".

template<int dim>
unsigned int FuelCellShop::Equation::ProtonTransportEquation< dim >::last_iter_bdry
protected

Variable used to store the index in bdry_info->global_data of the previous Newton solution The solution at the previous iteration is used to compute bdry_matrix and bdry_residual.

template<int dim>
unsigned int FuelCellShop::Equation::ProtonTransportEquation< dim >::last_iter_cell
protected

Variable used to store the index in cell_info->global_data of the previous Newton solution The solution at the previous iteration is used to compute cell_matrix and cell_residual.

template<int dim>
const std::string FuelCellShop::Equation::ProtonTransportEquation< dim >::name_base_variable
protected

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

template<int dim>
std::vector< std::vector<double> > FuelCellShop::Equation::ProtonTransportEquation< dim >::phi_lambda_bdry
protected

\( \mathbf{\lambda} \) shape functions.

phi_lambda_bdry [ q ] [ k ] denotes \( k \)-th \( \mathbf{\lambda} \) shape function computed in \( q \)-th quadrature point of the boundary.

template<int dim>
VariableInfo FuelCellShop::Equation::ProtonTransportEquation< dim >::phi_m
protected

VariableInfo structure corresponding to "protonic_electrical_potential".

template<int dim>
std::vector< std::vector<double> > FuelCellShop::Equation::ProtonTransportEquation< dim >::phi_phiM_bdry
protected

\( \mathbf{\phi_m} \) shape functions.

phi_phiM_bdry [ q ] [ k ] denotes \( k \)-th \( \mathbf{\phi_m} \) shape function computed in \( q \)-th quadrature point of the boundary.

template<int dim>
std::vector< std::vector<double> > FuelCellShop::Equation::ProtonTransportEquation< dim >::phi_phiM_cell
protected

\( \mathbf{\phi_{M}} \) shape functions.

phi_phiM_cell [ q ] [ k ] denotes \( k \)-th \( \mathbf{\phi_{M}} \) shape function computed in \( q \)-th quadrature point of the cell.

template<int dim>
std::vector< std::vector<double> > FuelCellShop::Equation::ProtonTransportEquation< dim >::phi_T_bdry
protected

\( \mathbf{T} \) shape functions.

phi_T_bdry [ q ] [ k ] denotes \( k \)-th \( \mathbf{T} \) shape function computed in \( q \)-th quadrature point of the boundary.

template<int dim>
std::map<unsigned int, double> FuelCellShop::Equation::ProtonTransportEquation< dim >::proton_current_flux_map
protected

std::map< unsigned int, double > container for details regarding Galvanostatic boundary conditions.

Here, Key (unsigned int) represents the boundary_id and Value (double) represents the constant protonic current flux values [A/cm^2].

template<int dim>
std::vector<double> FuelCellShop::Equation::ProtonTransportEquation< dim >::sigmaMeff_cell
protected

Effective proton conductivity, [S/cm], at all quadrature points of the cell.

template<int dim>
VariableInfo FuelCellShop::Equation::ProtonTransportEquation< dim >::t_rev
protected

VariableInfo structure corresponding to "temperature_of_REV".


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