| OpenFCST: The open-source Fuel Cell Simulation Toolbox
    | 
This class deals with Electron Transport Equation. More...
#include <electron_transport_equation.h>


| Public Member Functions | |
| Constructors, destructor, and initalization | |
| ElectronTransportEquation (FuelCell::SystemManagement &system_management) | |
| Constructor. | |
| virtual | ~ElectronTransportEquation () | 
| Destructor. | |
| virtual void | declare_parameters (ParameterHandler ¶m) const | 
| Declare parameters. | |
| virtual void | initialize (ParameterHandler ¶m) | 
| 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 printing out the equations 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_map & | get_internal_cell_couplings () const | 
| This function returns internal_cell_couplingsof a derived equation class. | |
| const couplings_map & | get_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_valueof a derived equation class. | |
| const component_boundaryID_value_map & | get_component_boundaryID_value () const | 
| This function returns component_boundaryID_valueof a derived equation class. | |
| const std::vector< BoundaryType > & | get_boundary_types () const | 
| This function returns boundary_typesof a derived equation class. | |
| const std::vector< std::vector < BoundaryType > > & | get_multi_boundary_types () const | 
| This function returns multi_boundary_typesof a derived equation class. | |
| const std::vector< OutputType > & | get_output_types () const | 
| This function returns output_typesof a derived equation class. | |
| const std::vector< std::vector < OutputType > > & | get_multi_output_types () const | 
| This function returns multi_output_typesof a derived equation class. | |
| const std::string & | get_equation_name () const | 
| This function returns equation_nameof a derived equation class. | |
| const std::vector< unsigned int > & | get_matrix_block_indices () const | 
| This function returns matrix_block_indicesof a derived equation class. | |
| const std::vector< unsigned int > & | get_residual_indices () const | 
| This function returns residual_indicesof 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 shapefunctions,shapefunctiongradients, 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 shapefunctions, 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_valueof a derived equation class. | |
| virtual void | make_component_boundaryID_value () | 
| This function fills out component_boundaryID_valueof a derived equation class. | |
| virtual void | make_multi_boundary_types () | 
| This function fills out multi_boundary_typesof a derived equation class. | |
| virtual void | make_multi_output_types () | 
| This function fills out multi_output_typesof a derived equation class. | |
| virtual void | make_matrix_block_indices () | 
| This function fills out matrix_block_indicesof a derived equation class. | |
| virtual void | make_residual_indices () | 
| This function fills out residual_indicesof a derived equation class. | |
| void | standard_to_block_wise (FullMatrix< double > &target) const | 
| This function changes the order of dealii::FullMatrix<double> targetfrom standard to block-wise. | |
| void | standard_to_block_wise (Vector< double > &target) const | 
| This function changes the order of dealii::Vector<double> targetfrom 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> srcinto the block-wise ordered structure AppFrame::MatrixVectordst. | |
| 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> srcinto the block-wise ordered structure AppFrame::FEVectordst. | |
| bool | belongs_to_boundary (const unsigned int &tria_boundary_id, const unsigned int ¶m_boundary_id) const | 
| This function returns trueif 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 | |
| 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 | 
| Container of boundary_id(s) for Dirichlet boundary conditions Here,Key(unsignedint) represents theboundary_idandValue(double) represents the potential[V]. | |
| std::map< unsigned int, double > | electron_current_flux_map | 
| std::map<unsignedint,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_s | 
| VariableInfo structure corresponding to "electronic_electrical_potential". | |
| Local CG FEM based assemblers - variable data (cell) | |
| Tensor< 2, dim > | sigmaSeff_cell | 
| Effective electronic conductivity, [ S/cm], of the cell. | |
| std::vector< std::vector < Tensor< 1, dim > > > | grad_phi_phiS_cell | 
| \( \mathbf{\phi_s} \) shape function gradients. | |
| Local CG FEM based assemblers - variable data (boundary) | |
| std::vector< std::vector < double > > | phi_phiS_bdry | 
| \( \mathbf{\phi_s} \) 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::SystemManagement * | system_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< BoundaryType > | boundary_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< OutputType > | output_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. | |
This class deals with Electron Transport Equation.
This equation solves a steady-state Ohm's law for electron transport inside the layers.
It is solved with respect to:
where, \( \mathbf{\phi_{s}} \) is in Voltages.
This equation can be written as:
\( \qquad \mathbf{\nabla} \cdot \left( \hat{\sigma}_{s,eff} \cdot \mathbf{\nabla} \phi_s \right) = 0 \quad \in \quad \Omega \)
S/cm].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:
Dirichlet: At this boundary, we specify the \( \mathbf{\phi_{s}} \) values using parameter file. It is to be noted that these dirichlet values at the boundary should be supplied to Initial Solution methods used in a particular application. These are specified in the parameter file under subsection "Dirichlet Boundary Conditions", as a list of comma-separated map-like values.
e.g.
Galvanostatic (Constant Electron Current Flux): A constant electronic current flux, \( C ~ \) [A/cm^2] is provided at the boundary. This value is provided using parameter file. If the value is positive, it implies that electronic current flux is leaving out of the boundary, otherwise negative implies entering into the boundary.
\( \qquad \mathbf{n} \cdot \left( \hat{\sigma}_{s,eff} \cdot \mathbf{\nabla} \phi_s \right) = C \)
These are specified in the parameter file under subsection "Constant Electron Current Flux Boundary Conditions", as a list of comma-separated map-like values.
e.g.
 where, boundary_id "1" has constant electron current flux value of 2 [A/cm^2] (leaving out of the boundary) and boundary_id "4" has constant electron current flux value of -0.5 [A/cm^2] (entering into the boundary).
No current flux (Insulated) / Symmetric: A particular case of Neumann boundary condition.
\( \qquad \mathbf{n} \cdot \left( \hat{\sigma}_{s,eff} \cdot \mathbf{\nabla} \phi_s \right) = 0 \)
No current flux or Symmetric boundary conditions, as FEM formulation automatically implies a particular boundary is one of these cases, by default.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 electron transport by Ohm's law,
ReactionSourceTerms class. Please read the documentation of ReactionSourceTerms class, for additional methods to be implemented in the application.adjust_internal_cell_couplings member function of ReactionSourceTerms class, before using make_cell_couplings of SystemManagement at the application level.| FuelCellShop::Equation::ElectronTransportEquation< dim >::ElectronTransportEquation | ( | FuelCell::SystemManagement & | system_management | ) | 
Constructor.
| 
 | virtual | 
Destructor.
| 
 | virtual | 
Assemble local boundary matrix.
Reimplemented from FuelCellShop::Equation::EquationBase< dim >.
| 
 | virtual | 
Assemble local boundary residual.
Reimplemented from FuelCellShop::Equation::EquationBase< dim >.
| 
 | virtual | 
Assemble local cell matrix.
Computes the contribution of Ohm's law to the cell_matrix:
\[ -\frac{\partial R}{\partial u} \delta u \]
Reimplemented from FuelCellShop::Equation::EquationBase< dim >.
| 
 | virtual | 
Assemble local cell residual.
Reimplemented from FuelCellShop::Equation::EquationBase< dim >.
| void FuelCellShop::Equation::ElectronTransportEquation< 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.
| 
 | virtual | 
Declare parameters.
Reimplemented from FuelCellShop::Equation::EquationBase< dim >.
| 
 | 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::ElectronTransportEquation< dim >::dirichlet_bdry_map.
| 
 | virtual | 
Initialize parameters.
Reimplemented from FuelCellShop::Equation::EquationBase< dim >.
| 
 | 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 >.
| 
 | protectedvirtual | 
This function computes Local CG FEM based assemblers - variable data (boundary) .
Reimplemented from FuelCellShop::Equation::EquationBase< dim >.
| 
 | 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 >.
| 
 | protectedvirtual | 
This function computes Local CG FEM based assemblers - variable data (cell) .
Reimplemented from FuelCellShop::Equation::EquationBase< dim >.
| 
 | protectedvirtual | 
This function computes Local CG FEM based assemblers - constant data (generic).
Reimplemented from FuelCellShop::Equation::EquationBase< dim >.
| 
 | protectedvirtual | 
This function fills out boundary_types. 
Reimplemented from FuelCellShop::Equation::EquationBase< dim >.
| 
 | protectedvirtual | 
This function fills out internal_cell_couplings. 
Reimplemented from FuelCellShop::Equation::EquationBase< dim >.
| 
 | inlineprotectedvirtual | 
This function fills out output_types. 
Reimplemented from FuelCellShop::Equation::EquationBase< dim >.
| 
 | virtual | 
The function printing out the equations info.
Reimplemented from FuelCellShop::Equation::EquationBase< dim >.
| 
 | protected | 
Container of boundary_id (s ) for Dirichlet boundary conditions Here, Key (unsigned int) represents the boundary_id and Value (double) represents the potential [V]. 
Referenced by FuelCellShop::Equation::ElectronTransportEquation< dim >::get_dirichlet_bdry_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 electronic current flux values [A/. cm^2]
| 
 | protected | 
\( \mathbf{\phi_s} \) shape function gradients.
grad_phi_phiS_cell [ q ] [ k ] denotes \( k \)-th \( \mathbf{\phi_s} \) shape function gradient computed in \( q \)-th quadrature point of the cell. 
| 
 | 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.
| 
 | 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.
| 
 | protected | 
Const std::string member storing name of the base solution variable corresponding to the equation represented by this class.
| 
 | protected | 
\( \mathbf{\phi_s} \) shape functions.
phi_phiS_bdry [ q ] [ k ] denotes \( k \)-th \( \mathbf{\phi_s} \) shape function computed in \( q \)-th quadrature point of the boundary. 
| 
 | protected | 
VariableInfo structure corresponding to "electronic_electrical_potential". 
| 
 | protected | 
Effective electronic conductivity, [S/cm], of the cell. 
 1.8.2
 1.8.2