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::ThermalTransportEquation< dim > Class Template Reference

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

#include <thermal_transport_equation.h>

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

Public Member Functions

Constructors, destructor, and initalization
 ThermalTransportEquation (FuelCell::SystemManagement &system_management, boost::shared_ptr< FuelCell::ApplicationCore::ApplicationData > data=boost::shared_ptr< FuelCell::ApplicationCore::ApplicationData >())
 Constructor. More...
 
virtual ~ThermalTransportEquation ()
 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
virtual void print_equation_info () const
 This function prints out the equations info. More...
 
bool get_electron_ohmic_heat_gdl () const
 Inline method to get whether the electronic ohmic heating in Gas diffusion layer is ON or OFF (electron_ohmic_heat_gdl). More...
 
bool get_electron_ohmic_heat_mpl () const
 Inline method to get whether the electronic ohmic heating in Microporous layer is ON or OFF (electron_ohmic_heat_mpl). More...
 
bool get_electron_ohmic_heat_cl () const
 Inline method to get whether the electronic ohmic heating in Catalyst layer is ON or OFF (electron_ohmic_heat_cl). More...
 
bool get_proton_ohmic_heat_cl () const
 Inline method to get whether the protonic ohmic heating in Catalyst layer is ON or OFF (proton_ohmic_heat_cl). More...
 
bool get_proton_ohmic_heat_ml () const
 Inline method to get whether the protonic ohmic heating in Membrane layer is ON or OFF (proton_ohmic_heat_ml). 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

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). More...
 
virtual void make_assemblers_cell_constant_data (const typename FuelCell::ApplicationCore::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) . More...
 
virtual void make_assemblers_bdry_constant_data (const typename FuelCell::ApplicationCore::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) . More...
 
virtual void make_assemblers_cell_variable_data (const typename FuelCell::ApplicationCore::DoFApplication< dim >::CellInfo &cell_info, FuelCellShop::Layer::BaseLayer< dim > *const layer)
 This function computes Local CG FEM based assemblers - variable data (cell) . More...
 
virtual void make_assemblers_bdry_variable_data (const typename FuelCell::ApplicationCore::DoFApplication< dim >::FaceInfo &bdry_info, FuelCellShop::Layer::BaseLayer< dim > *const layer)
 This function computes Local CG FEM based assemblers - variable data (boundary) . More...
 
template<typename porelayer >
void gas_enthalpy_transport_assemblers_cell_data (const typename FuelCell::ApplicationCore::DoFApplication< dim >::CellInfo &cell_info, porelayer *const layer)
 This function computes Enthalpy transport assemblers - variable data (cell) corresponding to enthalpy transport via Fickian diffusion. More...
 
template<typename polymer >
void lambda_enthalpy_transport_assemblers_cell_data (const typename FuelCell::ApplicationCore::DoFApplication< dim >::CellInfo &cell_info, const polymer *const layer)
 This function computes Enthalpy transport assemblers - variable data (cell) corresponding to enthalpy transport associated with lambda transport. More...
 
Other make_ functions
virtual void make_internal_cell_couplings ()
 This function fills out internal_cell_couplings. 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...
 
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 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...
 
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

bool cell_residual_counter
 Counter set to TRUE when cell_residual is being assembled. More...
 
bool bdry_residual_counter
 Counter set to TRUE when bdry_residual is being assembled. More...
 
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. More...
 
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. More...
 
Boolean flags
bool electron_ohmic_heat_gdl
 This boolean data member indicates that the electronic ohmic heating in Gas diffusion layer is ON or OFF. More...
 
bool electron_ohmic_heat_mpl
 This boolean data member indicates that the electronic ohmic heating in Microporous layer is ON or OFF. More...
 
bool electron_ohmic_heat_cl
 This boolean data member indicates that the electronic ohmic heating in Catalyst layer is ON or OFF. More...
 
bool proton_ohmic_heat_cl
 This boolean data member indicates that the protonic ohmic heating in Catalyst layer is ON or OFF. More...
 
bool proton_ohmic_heat_ml
 This boolean data member indicates that the protonic ohmic heating in Membrane layer is ON or OFF. More...
 
bool enthalpy_fickian_transport
 This boolean data member indicates that the enthalpy transport via Fickian diffusion is ON or OFF. More...
 
bool enthalpy_lambda_transport
 This boolean data member indicates that the enthalpy transport associated with lambda (sorbed water) transport is ON or OFF. More...
 
bool flag_thermoosmosis
 This flag indicates that lambda transport due to thermo-osmotic diffusion is ON or OFF. More...
 
Boundary conditions
std::map< unsigned int, double > const_heat_flux_map
 std::map< unsigned int, double > container for details regarding Constant Heat Flux boundaries. More...
 
std::map< unsigned int,
std::vector< double > > 
conv_heat_flux_map
 Container for details regarding Convective Heat Flux boundaries. More...
 
Generic Constant Data
VariableInfo t_rev
 VariableInfo structure corresponding to base variable of this equation class, "temperature_of_REV". More...
 
VariableInfo phi_s
 VariableInfo structure corresponding to "electronic_electrical_potential". More...
 
VariableInfo phi_m
 VariableInfo structure corresponding to "protonic_electrical_potential". More...
 
VariableInfo lambda
 VariableInfo structure corresponding to "membrane_water_content". More...
 
VariableInfo s_liquid_water
 VariableInfo structure corresponding to "liquid_water_saturation". More...
 
VariableInfo p_liquid_water
 
std::map< std::string,
VariableInfo
xi_map
 Map of VariableInfo structures of various gaseous species, whose diffusion is being considered in the application. More...
 
Local CG FEM based assemblers - variable data (cell)
std::vector< Tensor< 2, dim > > keff_cell
 Effective thermal conductivity, [W/(cm-k )], at all quadrature points of the cell. More...
 
std::vector< Tensor< 2, dim > > dkeff_dT_cell
 Derivative of effective thermal conductivity w.r.t "temperature_of_REV", at all quadrature points of the cell. More...
 
Tensor< 2, dimsigmaSeff_cell
 Effective electronic conductivity, [S/cm], of the cell. More...
 
std::vector< double > sigmaMeff_cell
 Effective protonic conductivity, [S/cm], at all quadrature points of the cell. More...
 
std::vector< double > dsigmaMeff_dT_cell
 Derivative of effective protonic conductivity w.r.t "temperature_of_REV", at all quadrature points in the cell. More...
 
std::vector< double > dsigmaMeff_dlambda_cell
 Derivative of effective protonic conductivity w.r.t "membrane_water_content", at all quadrature points in the cell. More...
 
std::vector< std::vector
< double > > 
phi_T_cell
 $ \mathbf{T} $ shape functions. More...
 
std::vector< std::vector
< Tensor< 1, dim > > > 
grad_phi_T_cell
 $ \mathbf{T} $ shape function gradients. More...
 
std::vector< std::vector
< Tensor< 1, dim > > > 
grad_phi_phiS_cell
 $ \mathbf{\phi_s} $ shape function gradients. More...
 
std::vector< std::vector
< Tensor< 1, dim > > > 
grad_phi_phiM_cell
 $ \mathbf{\phi_m} $ shape function gradients. More...
 
std::vector< std::vector
< double > > 
phi_lambda_cell
 $ \mathbf{\lambda} $ shape functions. More...
 
std::vector< std::vector
< Tensor< 1, dim > > > 
grad_phi_lambda_cell
 $ \mathbf{\lambda} $ shape function gradients. More...
 
std::vector< std::vector
< double > > 
phi_s_cell
 $ \mathbf{s} $ shape functions. More...
 
std::vector< std::vector
< double > > 
phi_p_cell
 
std::vector< Tensor< 1, dim > > grad_phiM_cell_old
 Gradient of "protonic_electrical_potential" at a previous Newton iteration, at all quadrature points in the cell. More...
 
std::vector< Tensor< 1, dim > > grad_phiS_cell_old
 Gradient of "electronic_electrical_potential" at a previous Newton iteration, at all quadrature points in the cell. More...
 
std::vector< Tensor< 1, dim > > grad_lambda_cell_old
 Gradient of "membrane_water_content" at a previous Newton iteration, at all quadrature points in the cell. More...
 
Enthalpy transport assemblers - variable data (cell)
std::map< std::string,
std::vector< Tensor< 2, dim > > > 
conc_Deff_dHdT_map
 Value in the map represents $ \left( c_T \cdot \hat{D}_{i,eff} \cdot \frac{\partial \bar{H}_i}{\partial T} \right) ~$ [W/(cm-K)], at previous iteration step at all quadrature points in the cell; corresponding to string Key representing for e.g "oxygen_molar_fraction". More...
 
std::map< std::string,
std::vector< Tensor< 2, dim > > > 
dT_concDeffdHdT_map
 Value in the map represents $ \frac{\partial \left( c_T \cdot \hat{D}_{i,eff} \cdot \frac{\partial \bar{H}_i}{\partial T} \right) }{\partial T} ~$ [W/(cm-K^2)], at previous iteration step at all quadrature points in the cell; corresponding to string Key representing for e.g "oxygen_molar_fraction". More...
 
std::map< std::string,
std::vector< Tensor< 2, dim > > > 
ds_concDeffdHdT_map
 Value in the map represents $ \frac{\partial \left( c_T \cdot \hat{D}_{i,eff} \cdot \frac{\partial \bar{H}_i}{\partial T} \right) }{\partial s} ~$ [W/(cm-K)], at previous iteration step at all quadrature points in the cell; corresponding to string Key representing for e.g "oxygen_molar_fraction". More...
 
std::map< std::string,
std::vector< Tensor< 2, dim > > > 
dp_concDeffdHdT_map
 
std::map< std::string,
std::vector< std::vector
< Tensor< 1, dim > > > > 
grad_phi_xi_map
 Map of $ \mathbf{x_i} $ shape function gradients. More...
 
std::vector< double > electroosmotic_dhdT
 $ \left(\frac{n_d \sigma_{m,eff}}{F} \right) \frac{\partial \bar{H}_{\lambda}}{\partial T} $, at all quadrature points in the cell. More...
 
std::vector< double > delectroosmotic_dhdT_dlambda
 $ \frac{ \partial \left[ \left(\frac{n_d \sigma_{m,eff}}{F} \right) \frac{\partial \bar{H}_{\lambda}}{\partial T} \right]}{\partial \lambda} $, at all quadrature points in the cell. More...
 
std::vector< double > delectroosmotic_dhdT_dT
 $ \frac{ \partial \left[ \left(\frac{n_d \sigma_{m,eff}}{F} \right) \frac{\partial \bar{H}_{\lambda}}{\partial T} \right]}{\partial T} $, at all quadrature points in the cell. More...
 
std::vector< double > backdiff_dhdT
 $ \left(\frac{\rho_{dry}}{EW} \right) D_{\lambda ,eff} \frac{\partial \bar{H}_{\lambda}}{\partial T} $, at all quadrature points in the cell. More...
 
std::vector< double > dbackdiff_dhdT_dlambda
 $ \frac{ \partial \left[ \left(\frac{\rho_{dry}}{EW} \right) D_{\lambda ,eff} \frac{\partial \bar{H}_{\lambda}}{\partial T} \right]}{\partial \lambda} $, at all quadrature points in the cell. More...
 
std::vector< double > dbackdiff_dhdT_dT
 $ \frac{ \partial \left[ \left(\frac{\rho_{dry}}{EW} \right) D_{\lambda ,eff} \frac{\partial \bar{H}_{\lambda}}{\partial T} \right]}{\partial T} $, at all quadrature points in the cell. More...
 
std::vector< double > thermoosmotic_dhdT
 $ \left(\frac{D_{T,eff}}{M_{water}} \right) \frac{\partial \bar{H}_{\lambda}}{\partial T} $, at all quadrature points in the cell. More...
 
std::vector< double > dthermoosmotic_dhdT_dT
 $ \frac{ \partial \left[ \left(\frac{D_{T,eff}}{M_{water}} \right) \frac{\partial \bar{H}_{\lambda}}{\partial T} \right]}{\partial T} $, at all quadrature points in the cell. More...
 
Local CG FEM based assemblers - variable data (boundary)
std::vector< Tensor< 2, dim > > dkeff_dT_bdry
 Derivative of effective thermal conductivity w.r.t "temperature_of_REV", at all quadrature points of the boundary. More...
 
std::vector< std::vector
< double > > 
phi_T_bdry
 $ \mathbf{T} $ shape functions. More...
 
- 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::ThermalTransportEquation< dim >

This class deals with Thermal Transport Equation.

This equation solves a steady-state, heat transfer equation (Fourier's law of conduction), alongwith enthalpy transport due to species diffusion in multi-component mixtures. Heat transfer via convection is generally very little (roughly 6 orders of magnitude smaller) as compared to conductive heat transport in a PEMFC, hence neglected here. Also currently, Fick's law of diffusion is being considered to account for enthalpy transport due to inter-diffusion amongst species.

where, $ \mathbf{T} $ is in Kelvins.

This equation can be written as:

$ \qquad -\mathbf{\nabla} \cdot \left( \hat{k}_{eff} \cdot \mathbf{\nabla} T \right) + \mathbf{\nabla} \cdot \left( \bar{H}_i \vec{J}_i \right) = \mathbf{S_{thermal}} \quad \in \quad \Omega $

where

$ \qquad \mathbf{S_{thermal}} = \begin{cases} \hat{\sigma}_{s,eff} \cdot \left( \mathbf{\nabla} \phi_s \otimes \mathbf{\nabla} \phi_s \right) \quad - \quad \text{Electronic ohmic heating in all layers except Membrane layer} \\ ~ \\ \hat{\sigma}_{m,eff} \cdot \left( \mathbf{\nabla} \phi_m \otimes \mathbf{\nabla} \phi_m \right) \quad - \quad \text{Protonic ohmic heating in Catalyst layers and Membrane layer} \end{cases} $

To be well-posed, these equations are equipped with the appropriate boundary conditions. Dirichlet boundary conditions can be mentioned in the subsection "Boundary data". However, some other boundary types require additional information, which can also be provided by the parameter file. We consider several types of boundary conditions:

$ \qquad -\mathbf{n} \cdot \left( \hat{k}_{eff} \cdot \mathbf{\nabla} T \right) = C $

These are specified in the parameter file under subsection "Boundary conditions", as a list of comma-separated map-like values.

e.g.

set Constant Heat Flux Boundary Conditions = 1: 34.5 , 4: -23.5

where, boundary_id "1" has constant heat flux value of 34.5 [W/cm^2] (coming out of the boundary) and boundary_id "4" has constant heat flux value of -23.5 [W/cm^2] (going into the boundary).

$ \qquad -\mathbf{n} \cdot \left( \hat{k}_{eff} \cdot \mathbf{\nabla} T \right) = h(T - T_{\infty}) $

These are specified in the parameter file under subsection "Boundary conditions", as a list of comma-separated map-like values.

e.g.

set Convective Heat Flux Boundary Conditions = 1: 1.5;340, 3: 0.4;230

where, boundary_id "1" has convective heat transfer coefficient of 1.5 [W/(cm^2-K )] and ambient temperature of 340 [K]; and boundary_id "3" has convective heat transfer coefficient of 0.4 [K] and ambient temperature of 230 [K]. Please note that ";" is necessary to provide both of the required parameters, viz., $ h $ and $ T_{\infty} $ at a particular boundary.

$ \qquad -\mathbf{n} \cdot \left( \hat{k}_{eff} \cdot \mathbf{\nabla} T \right) = 0 $

Remarks
  • Ohmic heating source terms (electronic and protonic) can be turned ON or OFF, using flags defined under subsection "Boolean flags in the parameter file.
  • In order to have electronic ohmic heating source term, $ \phi_s $ [Volts ], should be one of the solution variables, otherwise this class will throw an exception. Similarly, for protonic ohmic heating, $ \phi_m $ [ Volts ], should be one of the solution variables.
  • Enthalpy tranport due to Fickian diffusion of species is defaulted to OFF. It can be turned ON by setting following paramter entry to TRUE, under subsection "Boolean flags :
    set Enthalpy transport due to fickian diffusion of gases = true
  • Please note that while considering enthalpy transport due to fickian diffusion, it is very necessary to set the solvent gas for a particular layer as the last entry in the input vector to FuelCellShop::Layer::PorousLayer::set_gases method, in the initialization section of the application.
  • Enthalpy transport due to transport of sorbed water is defaulted to OFF. It can be turned ON by setting following parameter entry to TRUE, under subsection "Boolean flags :
    set Enthalpy transport associated with lambda transport = true
  • Please note that while considering enthalpy transport due to lambda transport, it is very necessary to include membrane_water_content, $ \lambda $ as one of the solution variables in the application.
  • Lambda transport accounts for diffusion, electro-osmosis if protonic_electrical_potential, $ \phi_m $ is one of the solution variables and thermo-osmosis if it is set to TRUE in parameter entry for Membrane Water Content Transport Equation.
  • There is no provision to specify boundary indicators for No Heat Flux or Symmetric boundary conditions, as FEM formulation automatically implies a particular boundary is one of these cases, by default.
  • This class currently works with the following layer classes:
    • FuelCellShop::Layer::GasDiffusionLayer<dim>
    • FuelCellShop::Layer::MicroPorousLayer<dim>
    • 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.

Usage Details:

* // Creating Equation object (in Application Header file)
*
* // Declare parameters in application
* thermal_transport.declare_parameters(param);
*
* // Initialize in application
* thermal_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( thermal_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: CGDL is FuelCellShop::Layer::DesignFibrousGDL<dim> object.
* thermal_transport.assemble_cell_matrix(cell_matrices, cell_info, &CGDL);
*
* // cell_residual in application
* thermal_transport.assemble_cell_residual(cell_vector, cell_info, &CGDL);
*
* // For bdry_matrices and bdry_vector; check for the layer present at that boundary. For eg: CGDL
* thermal_transport.assemble_bdry_matrix(bdry_matrices, bdry_info, &CGDL);
* thermal_transport.assemble_bdry_residual(bdry_vector, bdry_info, &CGDL);
*
Note
This class doesn't assemble for reaction heat 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 heat source terms due to reaction are being considered, it's very important to use ReactionSourceTerms::adjust_internal_cell_couplings method, before using FuelCell::SystemManagement::make_cell_couplings at the application level.
Author
Madhur Bhaiya, 2013

Constructor & Destructor Documentation

Constructor.

Destructor.

Member Function Documentation

template<int dim>
virtual void FuelCellShop::Equation::ThermalTransportEquation< 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::ThermalTransportEquation< 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::ThermalTransportEquation< 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::ThermalTransportEquation< 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::ThermalTransportEquation< dim >::declare_parameters ( ParameterHandler &  param) const
virtual

Declare parameters.

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

template<int dim>
template<typename porelayer >
void FuelCellShop::Equation::ThermalTransportEquation< dim >::gas_enthalpy_transport_assemblers_cell_data ( const typename FuelCell::ApplicationCore::DoFApplication< dim >::CellInfo &  cell_info,
porelayer *const  layer 
)
protected

This function computes Enthalpy transport assemblers - variable data (cell) corresponding to enthalpy transport via Fickian diffusion.

template<int dim>
bool FuelCellShop::Equation::ThermalTransportEquation< dim >::get_electron_ohmic_heat_cl ( ) const
inline

Inline method to get whether the electronic ohmic heating in Catalyst layer is ON or OFF (electron_ohmic_heat_cl).

The typical use of this method is for the post-processing routines in the application.

References FuelCellShop::Equation::ThermalTransportEquation< dim >::electron_ohmic_heat_cl.

template<int dim>
bool FuelCellShop::Equation::ThermalTransportEquation< dim >::get_electron_ohmic_heat_gdl ( ) const
inline

Inline method to get whether the electronic ohmic heating in Gas diffusion layer is ON or OFF (electron_ohmic_heat_gdl).

The typical use of this method is for the post-processing routines in the application.

References FuelCellShop::Equation::ThermalTransportEquation< dim >::electron_ohmic_heat_gdl.

template<int dim>
bool FuelCellShop::Equation::ThermalTransportEquation< dim >::get_electron_ohmic_heat_mpl ( ) const
inline

Inline method to get whether the electronic ohmic heating in Microporous layer is ON or OFF (electron_ohmic_heat_mpl).

The typical use of this method is for the post-processing routines in the application.

References FuelCellShop::Equation::ThermalTransportEquation< dim >::electron_ohmic_heat_mpl.

template<int dim>
bool FuelCellShop::Equation::ThermalTransportEquation< dim >::get_proton_ohmic_heat_cl ( ) const
inline

Inline method to get whether the protonic ohmic heating in Catalyst layer is ON or OFF (proton_ohmic_heat_cl).

The typical use of this method is for the post-processing routines in the application.

References FuelCellShop::Equation::ThermalTransportEquation< dim >::proton_ohmic_heat_cl.

template<int dim>
bool FuelCellShop::Equation::ThermalTransportEquation< dim >::get_proton_ohmic_heat_ml ( ) const
inline

Inline method to get whether the protonic ohmic heating in Membrane layer is ON or OFF (proton_ohmic_heat_ml).

The typical use of this method is for the post-processing routines in the application.

References FuelCellShop::Equation::ThermalTransportEquation< dim >::proton_ohmic_heat_ml.

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

Initialize parameters.

This class will call make_internal_cell_couplings.

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

template<int dim>
template<typename polymer >
void FuelCellShop::Equation::ThermalTransportEquation< dim >::lambda_enthalpy_transport_assemblers_cell_data ( const typename FuelCell::ApplicationCore::DoFApplication< dim >::CellInfo &  cell_info,
const polymer *const  layer 
)
protected

This function computes Enthalpy transport assemblers - variable data (cell) corresponding to enthalpy transport associated with lambda transport.

template<int dim>
virtual void FuelCellShop::Equation::ThermalTransportEquation< dim >::make_assemblers_bdry_constant_data ( const typename FuelCell::ApplicationCore::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::ThermalTransportEquation< dim >::make_assemblers_bdry_variable_data ( const typename FuelCell::ApplicationCore::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::ThermalTransportEquation< dim >::make_assemblers_cell_constant_data ( const typename FuelCell::ApplicationCore::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::ThermalTransportEquation< dim >::make_assemblers_cell_variable_data ( const typename FuelCell::ApplicationCore::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::ThermalTransportEquation< 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::ThermalTransportEquation< 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::ThermalTransportEquation< 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::vector<double> FuelCellShop::Equation::ThermalTransportEquation< dim >::backdiff_dhdT
protected

$ \left(\frac{\rho_{dry}}{EW} \right) D_{\lambda ,eff} \frac{\partial \bar{H}_{\lambda}}{\partial T} $, at all quadrature points in the cell.

template<int dim>
bool FuelCellShop::Equation::ThermalTransportEquation< 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::ThermalTransportEquation< dim >::cell_residual_counter
protected

Counter set to TRUE when cell_residual is being assembled.

This ensures that only effective transport properties are calculated, not their derivatives. (improves speed)

template<int dim>
std::map< std::string, std::vector< Tensor<2,dim> > > FuelCellShop::Equation::ThermalTransportEquation< dim >::conc_Deff_dHdT_map
protected

Value in the map represents $ \left( c_T \cdot \hat{D}_{i,eff} \cdot \frac{\partial \bar{H}_i}{\partial T} \right) ~$ [W/(cm-K)], at previous iteration step at all quadrature points in the cell; corresponding to string Key representing for e.g "oxygen_molar_fraction".

template<int dim>
std::map<unsigned int, double> FuelCellShop::Equation::ThermalTransportEquation< dim >::const_heat_flux_map
protected

std::map< unsigned int, double > container for details regarding Constant Heat Flux boundaries.

Here, Key (unsigned int) represents the boundary_id and Value (double) represents the constant heat flux values [W/cm^2]. Positive heat flux represents heat leaving out of the boundary and Negative coming into the boundary.

template<int dim>
std::map<unsigned int, std::vector<double> > FuelCellShop::Equation::ThermalTransportEquation< dim >::conv_heat_flux_map
protected

Container for details regarding Convective Heat Flux boundaries.

Here, Key represents the boundary_id, and Value vector stores the heat transfer coefficient [W/(cm^2-K )] and ambient temperature [K] respectively.

template<int dim>
std::vector<double> FuelCellShop::Equation::ThermalTransportEquation< dim >::dbackdiff_dhdT_dlambda
protected

$ \frac{ \partial \left[ \left(\frac{\rho_{dry}}{EW} \right) D_{\lambda ,eff} \frac{\partial \bar{H}_{\lambda}}{\partial T} \right]}{\partial \lambda} $, at all quadrature points in the cell.

template<int dim>
std::vector<double> FuelCellShop::Equation::ThermalTransportEquation< dim >::dbackdiff_dhdT_dT
protected

$ \frac{ \partial \left[ \left(\frac{\rho_{dry}}{EW} \right) D_{\lambda ,eff} \frac{\partial \bar{H}_{\lambda}}{\partial T} \right]}{\partial T} $, at all quadrature points in the cell.

template<int dim>
std::vector<double> FuelCellShop::Equation::ThermalTransportEquation< dim >::delectroosmotic_dhdT_dlambda
protected

$ \frac{ \partial \left[ \left(\frac{n_d \sigma_{m,eff}}{F} \right) \frac{\partial \bar{H}_{\lambda}}{\partial T} \right]}{\partial \lambda} $, at all quadrature points in the cell.

template<int dim>
std::vector<double> FuelCellShop::Equation::ThermalTransportEquation< dim >::delectroosmotic_dhdT_dT
protected

$ \frac{ \partial \left[ \left(\frac{n_d \sigma_{m,eff}}{F} \right) \frac{\partial \bar{H}_{\lambda}}{\partial T} \right]}{\partial T} $, at all quadrature points in the cell.

template<int dim>
std::vector< Tensor<2,dim> > FuelCellShop::Equation::ThermalTransportEquation< dim >::dkeff_dT_bdry
protected

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

template<int dim>
std::vector< Tensor<2,dim> > FuelCellShop::Equation::ThermalTransportEquation< dim >::dkeff_dT_cell
protected

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

template<int dim>
std::map< std::string, std::vector< Tensor<2,dim> > > FuelCellShop::Equation::ThermalTransportEquation< dim >::dp_concDeffdHdT_map
protected
template<int dim>
std::map< std::string, std::vector< Tensor<2,dim> > > FuelCellShop::Equation::ThermalTransportEquation< dim >::ds_concDeffdHdT_map
protected

Value in the map represents $ \frac{\partial \left( c_T \cdot \hat{D}_{i,eff} \cdot \frac{\partial \bar{H}_i}{\partial T} \right) }{\partial s} ~$ [W/(cm-K)], at previous iteration step at all quadrature points in the cell; corresponding to string Key representing for e.g "oxygen_molar_fraction".

template<int dim>
std::vector<double> FuelCellShop::Equation::ThermalTransportEquation< 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::ThermalTransportEquation< 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::map< std::string, std::vector< Tensor<2,dim> > > FuelCellShop::Equation::ThermalTransportEquation< dim >::dT_concDeffdHdT_map
protected

Value in the map represents $ \frac{\partial \left( c_T \cdot \hat{D}_{i,eff} \cdot \frac{\partial \bar{H}_i}{\partial T} \right) }{\partial T} ~$ [W/(cm-K^2)], at previous iteration step at all quadrature points in the cell; corresponding to string Key representing for e.g "oxygen_molar_fraction".

template<int dim>
std::vector<double> FuelCellShop::Equation::ThermalTransportEquation< dim >::dthermoosmotic_dhdT_dT
protected

$ \frac{ \partial \left[ \left(\frac{D_{T,eff}}{M_{water}} \right) \frac{\partial \bar{H}_{\lambda}}{\partial T} \right]}{\partial T} $, at all quadrature points in the cell.

template<int dim>
bool FuelCellShop::Equation::ThermalTransportEquation< dim >::electron_ohmic_heat_cl
protected

This boolean data member indicates that the electronic ohmic heating in Catalyst layer is ON or OFF.

Referenced by FuelCellShop::Equation::ThermalTransportEquation< dim >::get_electron_ohmic_heat_cl().

template<int dim>
bool FuelCellShop::Equation::ThermalTransportEquation< dim >::electron_ohmic_heat_gdl
protected

This boolean data member indicates that the electronic ohmic heating in Gas diffusion layer is ON or OFF.

Referenced by FuelCellShop::Equation::ThermalTransportEquation< dim >::get_electron_ohmic_heat_gdl().

template<int dim>
bool FuelCellShop::Equation::ThermalTransportEquation< dim >::electron_ohmic_heat_mpl
protected

This boolean data member indicates that the electronic ohmic heating in Microporous layer is ON or OFF.

Referenced by FuelCellShop::Equation::ThermalTransportEquation< dim >::get_electron_ohmic_heat_mpl().

template<int dim>
std::vector<double> FuelCellShop::Equation::ThermalTransportEquation< dim >::electroosmotic_dhdT
protected

$ \left(\frac{n_d \sigma_{m,eff}}{F} \right) \frac{\partial \bar{H}_{\lambda}}{\partial T} $, at all quadrature points in the cell.

template<int dim>
bool FuelCellShop::Equation::ThermalTransportEquation< dim >::enthalpy_fickian_transport
protected

This boolean data member indicates that the enthalpy transport via Fickian diffusion is ON or OFF.

template<int dim>
bool FuelCellShop::Equation::ThermalTransportEquation< dim >::enthalpy_lambda_transport
protected

This boolean data member indicates that the enthalpy transport associated with lambda (sorbed water) transport is ON or OFF.

template<int dim>
bool FuelCellShop::Equation::ThermalTransportEquation< dim >::flag_thermoosmosis
protected

This flag indicates that lambda transport due to thermo-osmotic diffusion is ON or OFF.

template<int dim>
std::vector< Tensor<1,dim> > FuelCellShop::Equation::ThermalTransportEquation< dim >::grad_lambda_cell_old
protected

Gradient of "membrane_water_content" at a previous Newton iteration, at all quadrature points in the cell.

template<int dim>
std::vector< std::vector< Tensor<1,dim> > > FuelCellShop::Equation::ThermalTransportEquation< 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::ThermalTransportEquation< 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::ThermalTransportEquation< dim >::grad_phi_phiS_cell
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.

template<int dim>
std::vector< std::vector< Tensor<1,dim> > > FuelCellShop::Equation::ThermalTransportEquation< 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>
std::map< std::string, std::vector< std::vector< Tensor<1,dim> > > > FuelCellShop::Equation::ThermalTransportEquation< dim >::grad_phi_xi_map
protected

Map of $ \mathbf{x_i} $ shape function gradients.

Value in the grad_phi_xi_map [ q ] [ k ] denotes $ k $-th $ \mathbf{x_i} $ shape function gradient computed in $ q $-th quadrature point of the cell, corresponding to string Key representing for e.g "oxygen_molar_fraction".

template<int dim>
std::vector< Tensor<1,dim> > FuelCellShop::Equation::ThermalTransportEquation< dim >::grad_phiM_cell_old
protected

Gradient of "protonic_electrical_potential" at a previous Newton iteration, at all quadrature points in the cell.

template<int dim>
std::vector< Tensor<1,dim> > FuelCellShop::Equation::ThermalTransportEquation< dim >::grad_phiS_cell_old
protected

Gradient of "electronic_electrical_potential" at a previous Newton iteration, at all quadrature points in the cell.

template<int dim>
std::vector< Tensor<2,dim> > FuelCellShop::Equation::ThermalTransportEquation< dim >::keff_cell
protected

Effective thermal conductivity, [W/(cm-k )], at all quadrature points of the cell.

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

VariableInfo structure corresponding to "membrane_water_content".

template<int dim>
unsigned int FuelCellShop::Equation::ThermalTransportEquation< 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::ThermalTransportEquation< 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>
VariableInfo FuelCellShop::Equation::ThermalTransportEquation< dim >::p_liquid_water
protected
template<int dim>
std::vector< std::vector<double> > FuelCellShop::Equation::ThermalTransportEquation< dim >::phi_lambda_cell
protected

$ \mathbf{\lambda} $ shape functions.

phi_lambda_cell [ q ] [ k ] denotes $ k $-th $ \mathbf{\lambda} $ shape function computed in $ q $-th quadrature point of the cell.

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

VariableInfo structure corresponding to "protonic_electrical_potential".

template<int dim>
std::vector< std::vector<double> > FuelCellShop::Equation::ThermalTransportEquation< dim >::phi_p_cell
protected
template<int dim>
VariableInfo FuelCellShop::Equation::ThermalTransportEquation< dim >::phi_s
protected

VariableInfo structure corresponding to "electronic_electrical_potential".

template<int dim>
std::vector< std::vector<double> > FuelCellShop::Equation::ThermalTransportEquation< dim >::phi_s_cell
protected

$ \mathbf{s} $ shape functions.

phi_s_cell [ q ] [ k ] denotes $ k $-th $ \mathbf{s} $ shape function computed in $ q $-th quadrature point of the cell.

template<int dim>
std::vector< std::vector<double> > FuelCellShop::Equation::ThermalTransportEquation< 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::vector< std::vector<double> > FuelCellShop::Equation::ThermalTransportEquation< dim >::phi_T_cell
protected

$ \mathbf{T} $ shape functions.

phi_T_cell [ q ] [ k ] denotes $ k $-th $ \mathbf{T} $ shape function computed in $ q $-th quadrature point of the cell.

template<int dim>
bool FuelCellShop::Equation::ThermalTransportEquation< dim >::proton_ohmic_heat_cl
protected

This boolean data member indicates that the protonic ohmic heating in Catalyst layer is ON or OFF.

Referenced by FuelCellShop::Equation::ThermalTransportEquation< dim >::get_proton_ohmic_heat_cl().

template<int dim>
bool FuelCellShop::Equation::ThermalTransportEquation< dim >::proton_ohmic_heat_ml
protected

This boolean data member indicates that the protonic ohmic heating in Membrane layer is ON or OFF.

Referenced by FuelCellShop::Equation::ThermalTransportEquation< dim >::get_proton_ohmic_heat_ml().

template<int dim>
VariableInfo FuelCellShop::Equation::ThermalTransportEquation< dim >::s_liquid_water
protected

VariableInfo structure corresponding to "liquid_water_saturation".

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

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

template<int dim>
Tensor<2,dim> FuelCellShop::Equation::ThermalTransportEquation< dim >::sigmaSeff_cell
protected

Effective electronic conductivity, [S/cm], of the cell.

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

VariableInfo structure corresponding to base variable of this equation class, "temperature_of_REV".

template<int dim>
std::vector<double> FuelCellShop::Equation::ThermalTransportEquation< dim >::thermoosmotic_dhdT
protected

$ \left(\frac{D_{T,eff}}{M_{water}} \right) \frac{\partial \bar{H}_{\lambda}}{\partial T} $, at all quadrature points in the cell.

template<int dim>
std::map< std::string, VariableInfo > FuelCellShop::Equation::ThermalTransportEquation< dim >::xi_map
protected

Map of VariableInfo structures of various gaseous species, whose diffusion is being considered in the application.

For instance, in case of oxygen gas: String Key is "oxygen_molar_fraction", while Value is VariableInfo structure corresponding to "oxygen_molar_fraction".


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