OpenFCST: The open-source Fuel Cell Simulation Toolbox
|
This class deals with Thermal Transport Equation. More...
#include <thermal_transport_equation.h>
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 ¶m) const |
Declare parameters. More... | |
virtual void | initialize (ParameterHandler ¶m) |
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_map & | get_internal_cell_couplings () const |
This function returns internal_cell_couplings of a derived equation class. More... | |
const couplings_map & | get_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 ¶m) |
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 ¶m_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/ )], 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, dim > | sigmaSeff_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 |
shape functions. More... | |
std::vector< std::vector < Tensor< 1, dim > > > | grad_phi_T_cell |
shape function gradients. More... | |
std::vector< std::vector < Tensor< 1, dim > > > | grad_phi_phiS_cell |
shape function gradients. More... | |
std::vector< std::vector < Tensor< 1, dim > > > | grad_phi_phiM_cell |
shape function gradients. More... | |
std::vector< std::vector < double > > | phi_lambda_cell |
shape functions. More... | |
std::vector< std::vector < Tensor< 1, dim > > > | grad_phi_lambda_cell |
shape function gradients. More... | |
std::vector< std::vector < double > > | phi_s_cell |
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 [W/ ], 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 [W/ ], 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 [W/ ], 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 shape function gradients. More... | |
std::vector< double > | electroosmotic_dhdT |
, at all quadrature points in the cell. More... | |
std::vector< double > | delectroosmotic_dhdT_dlambda |
, at all quadrature points in the cell. More... | |
std::vector< double > | delectroosmotic_dhdT_dT |
, at all quadrature points in the cell. More... | |
std::vector< double > | backdiff_dhdT |
, at all quadrature points in the cell. More... | |
std::vector< double > | dbackdiff_dhdT_dlambda |
, at all quadrature points in the cell. More... | |
std::vector< double > | dbackdiff_dhdT_dT |
, at all quadrature points in the cell. More... | |
std::vector< double > | thermoosmotic_dhdT |
, at all quadrature points in the cell. More... | |
std::vector< double > | dthermoosmotic_dhdT_dT |
, 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 |
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::SystemManagement * | system_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< BoundaryType > | boundary_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< OutputType > | output_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... | |
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.
(temperature_of_REV
)where, is in Kelvins.
This equation can be written as:
where
W/
(cm-K
)
], which can be a function of other variables, e.g. (Temperature).S/cm
].S/cm
].W/
(cm^3
)
].J/
(mol
)
] of the diffusing specie "i"
.mol/
(cm^2-s
)
] of the diffusing specie "i"
, for gaseous species it is given as: , where:cm^2/s
] of the diffusing gas "i"
.mol/cm^3
] of all the gaseous species involved in the cell."i"
."i"
can also be sorbed water, .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:
W/cm^2
]
is provided at a boundary. This value is provided using parameter file. If the value is positive, it implies that heat flux is leaving out of the boundary, otherwise negative implies entering into the boundary.These are specified in the parameter file under subsection "Boundary conditions", as a list of comma-separated map-like
values.
e.g.
where, boundary_id
"1" has constant heat flux value of 34.5 [W/
(coming out of the boundary) and cm^2
]boundary_id
"4" has constant heat flux value of -23.5 [W/
(going into the boundary).cm^2
]
heat
transfer
coefficient
, [W/
(cm^2-K
)] and ambient
temperature
, [K
] at the boundary.These are specified in the parameter file under subsection "Boundary conditions", as a list of comma-separated map-like
values.
e.g.
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., and at a particular boundary.
Neumann
boundary condition.electronic
and protonic
) can be turned ON or OFF, using flags defined under subsection "Boolean flags in the parameter file.Volts
], should be one of the solution variables, otherwise this class will throw an exception. Similarly, for protonic ohmic heating, [ Volts
], should be one of the solution variables.TRUE
, under subsection "Boolean flags : TRUE
, under subsection "Boolean flags : membrane_water_content
, as one of the solution variables in the application.protonic_electrical_potential
, is one of the solution variables and thermo-osmosis if it is set to TRUE
in parameter entry for Membrane Water Content Transport Equation.No
Heat
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.
FuelCellShop::Equation::ThermalTransportEquation< dim >::ThermalTransportEquation | ( | FuelCell::SystemManagement & | system_management, |
boost::shared_ptr< FuelCell::ApplicationCore::ApplicationData > | data = boost::shared_ptr< FuelCell::ApplicationCore::ApplicationData >() |
||
) |
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.
Reimplemented from FuelCellShop::Equation::EquationBase< dim >.
|
virtual |
Assemble local cell residual.
Reimplemented from FuelCellShop::Equation::EquationBase< dim >.
|
virtual |
Declare parameters.
Reimplemented from FuelCellShop::Equation::EquationBase< dim >.
|
protected |
This function computes Enthalpy transport assemblers - variable data (cell) corresponding to enthalpy transport via Fickian diffusion.
|
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.
|
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.
|
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.
|
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.
|
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.
|
virtual |
Initialize parameters.
This class will call make_internal_cell_couplings.
Reimplemented from FuelCellShop::Equation::EquationBase< dim >.
|
protected |
This function computes Enthalpy transport assemblers - variable data (cell) corresponding to enthalpy transport associated with lambda transport.
|
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 internal_cell_couplings
.
Reimplemented from FuelCellShop::Equation::EquationBase< dim >.
|
virtual |
This function prints out the equations info.
Reimplemented from FuelCellShop::Equation::EquationBase< dim >.
|
protected |
, at all quadrature points in the cell.
|
protected |
Counter set to TRUE when bdry_residual
is being assembled.
This ensures that only data required for bdry_residual is computed. (improves speed)
|
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)
|
protected |
Value
in the map represents [W/
], at previous iteration step at all quadrature points in the cell; corresponding to string (cm-K
)
Key
representing for e.g "oxygen_molar_fraction".
|
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/
. Positive heat flux represents heat leaving out of the boundary and Negative coming into the boundary. cm^2
]
|
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.
|
protected |
, at all quadrature points in the cell.
|
protected |
, at all quadrature points in the cell.
|
protected |
, at all quadrature points in the cell.
|
protected |
, at all quadrature points in the cell.
|
protected |
Derivative of effective thermal conductivity w.r.t "temperature_of_REV"
, at all quadrature points of the boundary.
|
protected |
Derivative of effective thermal conductivity w.r.t "temperature_of_REV"
, at all quadrature points of the cell.
|
protected |
|
protected |
Value
in the map represents [W/
], at previous iteration step at all quadrature points in the cell; corresponding to string (cm-K
)
Key
representing for e.g "oxygen_molar_fraction".
|
protected |
Derivative of effective protonic conductivity w.r.t "membrane_water_content"
, at all quadrature points in the cell.
|
protected |
Derivative of effective protonic conductivity w.r.t "temperature_of_REV"
, at all quadrature points in the cell.
|
protected |
Value
in the map represents [W/
], at previous iteration step at all quadrature points in the cell; corresponding to string (cm-K^2
)
Key
representing for e.g "oxygen_molar_fraction".
|
protected |
, at all quadrature points in the cell.
|
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().
|
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().
|
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().
|
protected |
, at all quadrature points in the cell.
|
protected |
This boolean data member indicates that the enthalpy transport via Fickian diffusion is ON
or OFF
.
|
protected |
This boolean data member indicates that the enthalpy transport associated with lambda (sorbed water) transport is ON
or OFF
.
|
protected |
This flag indicates that lambda transport due to thermo-osmotic diffusion is ON
or OFF
.
|
protected |
Gradient of "membrane_water_content"
at a previous Newton iteration, at all quadrature points in the cell.
|
protected |
shape function gradients.
grad_phi_lambda_cell
[
q
]
[
k
] denotes -th shape function gradient computed in -th quadrature point of the cell.
|
protected |
shape function gradients.
grad_phi_phiM_cell
[
q
]
[
k
] denotes -th shape function gradient computed in -th quadrature point of the cell.
|
protected |
shape function gradients.
grad_phi_phiS_cell
[
q
]
[
k
] denotes -th shape function gradient computed in -th quadrature point of the cell.
|
protected |
shape function gradients.
grad_phi_T_cell
[
q
]
[
k
] denotes -th shape function gradient computed in -th quadrature point of the cell.
|
protected |
Map of shape function gradients.
Value
in the grad_phi_xi_map
[
q
]
[
k
] denotes -th shape function gradient computed in -th quadrature point of the cell, corresponding to string
Key
representing for e.g "oxygen_molar_fraction".
|
protected |
Gradient of "protonic_electrical_potential"
at a previous Newton iteration, at all quadrature points in the cell.
|
protected |
Gradient of "electronic_electrical_potential"
at a previous Newton iteration, at all quadrature points in the cell.
|
protected |
Effective thermal conductivity, [W/
(cm-k
)], at all quadrature points of the cell.
|
protected |
VariableInfo structure corresponding to "membrane_water_content"
.
|
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 |
|
protected |
shape functions.
phi_lambda_cell
[
q
]
[
k
] denotes -th shape function computed in -th quadrature point of the cell.
|
protected |
VariableInfo structure corresponding to "protonic_electrical_potential"
.
|
protected |
|
protected |
VariableInfo structure corresponding to "electronic_electrical_potential"
.
|
protected |
shape functions.
phi_s_cell
[
q
]
[
k
] denotes -th shape function computed in -th quadrature point of the cell.
|
protected |
shape functions.
phi_T_bdry
[
q
]
[
k
] denotes -th shape function computed in -th quadrature point of the boundary.
|
protected |
shape functions.
phi_T_cell
[
q
]
[
k
] denotes -th shape function computed in -th quadrature point of the cell.
|
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().
|
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().
|
protected |
VariableInfo structure corresponding to "liquid_water_saturation"
.
|
protected |
Effective protonic conductivity, [S/cm
], at all quadrature points of the cell.
|
protected |
Effective electronic conductivity, [S/cm
], of the cell.
|
protected |
VariableInfo structure corresponding to base variable of this equation class, "temperature_of_REV"
.
|
protected |
, at all quadrature points in the cell.
|
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"
.