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

This class implements a mass conservation equation for liquid transport in a fuel cell porous media. More...

#include <liquid_transport_equation.h>

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

Public Member Functions

Constructors, destructor, and initalization
 LiquidPressureEquation (FuelCell::SystemManagement &system_management, boost::shared_ptr< FuelCell::ApplicationCore::ApplicationData > data=boost::shared_ptr< FuelCell::ApplicationCore::ApplicationData >())
 Constructor. More...
 
virtual ~LiquidPressureEquation ()
 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_rhs, 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_rhs, 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
 The function prints out the equation's info. 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...
 
Other make_ functions
virtual void make_internal_cell_couplings ()
 This function fills out internal_cell_couplings. More...
 
virtual void make_boundary_types ()
 This function fills out boundary_types. More...
 
virtual void make_output_types ()
 This function fills out output_types. 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_multi_boundary_types ()
 This function fills out multi_boundary_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...
 
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...
 
bool breakthrough_phenomenon
 This data member is set to true when the following boundary condition is to be applied at the boundary $ N = k (p_{c(threshold)} - p_c) $. More...
 
std::vector< bool > reached_breakthrough_pressure
 Once the liquid pressure at the GDL and channel interface reaches the breakthrouhg pressure the following boundary condition is to set to be applied $ N = k (p_{c(threshold)} - p_c) $. More...
 
std::map< unsigned int, double > capillary_pressure_current_flux_map
 The map is used as the $ k_1 $ parameter in $ N = k_1 (p_{c(threshold)} - p_c) $. More...
 
double p_channel_GDL
 The boundary condition term which indicates the liquid pressure (Pa) $ p_{c(threshold1)} $ in $ N = k_1 (p_{c(threshold1)} - p_c) $. More...
 
double p_channel_GDL_second
 The boundary condition term which indicates the breakthrough liquid pressure (Pa). More...
 
double temperature
 Temperature [in K] of the layer. More...
 
Generic Constant Data
VariableInfo p_liquid_water
 VariableInfo structure corresponding to "liquid_pressure". More...
 
VariableInfo x_water
 VariableInfo structure corresponding to "water_molar_fraction". More...
 
VariableInfo t_rev
 VariableInfo structure corresponding to "temperature_of_REV". More...
 
unsigned int fickswater_blockindex_xwater
 Block index for "water_molar_fraction" corresponding to "Ficks Transport Equation - water". More...
 
unsigned int fickswater_blockindex_pliquid
 Block index for "liquid_pressure" corresponding to "Ficks Transport Equation - water". More...
 
unsigned int fickswater_blockindex_trev
 Block index for "temperature_of_REV" corresponding to "Ficks Transport Equation - water". More...
 
unsigned int thermal_blockindex_xwater
 Block index for "water_molar_fraction" corresponding to "Thermal Transport Equation". More...
 
unsigned int thermal_blockindex_pliquid
 Block index for "liquid_pressure" corresponding to "Thermal Transport Equation". More...
 
unsigned int thermal_blockindex_trev
 Block index for "temperature_of_REV" corresponding to "Thermal Transport Equation". More...
 
double M_water
 Molar weight of water in grams/mole. More...
 
double rho_l
 Density of liquid water, $ \rho_l ~\left[ \frac{g}{cm^3} \right] $. More...
 
double k_sat
 Saturated permeability inside the porous media, $ k_{sat} ~\left[ cm^2 \right] $. More...
 
double k_sat_bdry
 Saturated permeability at the boundary of the domain, $ k_{satbdry} ~\left[ cm^2 \right] $. More...
 
Local CG FEM based assemblers - variable data (cell)
double p_cell
 Total pressure $ \left[ Pa \right] $ in the cell. More...
 
std::vector< Tensor< 2, dim > > rhok_mu_cell
 $ \frac{\rho_l \hat{k}_l}{\mu_l} $, at all quadrature points in the cell. More...
 
std::vector< Tensor< 2, dim > > rhok_mu_cell_bdry
 $ \frac{\rho_l \hat{k}_l}{\mu_l} $, at all quadrature points in the boundary. More...
 
std::vector< Tensor< 2, dim > > rho_mu_dk_dp_cell
 Computing the derivitive with respect to $ p_c $ in the function $ \frac{\rho_l \hat{k}_l(p_c)}{\mu_l(T)} $, at all quadrature points in the cell. More...
 
std::vector< Tensor< 2, dim > > dT_rhok_mu_cell
 Computing the derivitive with respect to $ T $ in the function $ \frac{\rho_l \hat{k}_l(p_c)}{\mu_l(T)} $, at all quadrature points in the cell. More...
 
std::vector< std::vector
< double > > 
phi_p_cell
 $ \mathbf{p_c} $ shape functions. More...
 
std::vector< std::vector
< double > > 
phi_p_cell_bdry
 $ \mathbf{p_c} $ shape functions at the boundary. More...
 
std::vector< std::vector
< double > > 
phi_xwater_cell
 $ \mathbf{x_{H_2O}} $ shape functions. More...
 
std::vector< std::vector
< double > > 
phi_T_cell
 $ \mathbf{T} $ shape functions. More...
 
std::vector< std::vector
< Tensor< 1, dim > > > 
grad_phi_p_cell
 $ \mathbf{p_c} $ shape function gradients. More...
 
std::vector< std::vector
< Tensor< 1, dim > > > 
grad_phi_T_cell
 $ \mathbf{T} $ shape function gradients. 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::LiquidPressureEquation< dim >

This class implements a mass conservation equation for liquid transport in a fuel cell porous media.

The liquid water velocity is based on Darcy's law.

Governing equation

The governing equation is based on a mass balance combined with Darcy's law. The equation is given by:

$ \qquad -\mathbf{\nabla} \cdot \left[ \frac{\rho_l \hat{k}_l}{\mu_l} \mathbf{\nabla} p_c \right] = 0 \quad \in \quad \Omega $ where

The equation is solved with respect to:

The relative liquid permeability tensor is given by the layer classes. Currently the equation can be used with the following layer classes and its children:

Viscosity is a function of temperature. Therefore, it is suggested be solved this equation with the Thermal Transport Equation even though a constant temperature can also be provided.

The class currently implements the weak linearized form of the equation above. Therefore, cell_matrix returns the Jacobian of the equation and cell_residual contains the weak form of the equation using the solution at the previous Newton step.

Usage Details:

For users: An example of how the equation is setup in the input file is (assuming there is only one equation being solved in this case):

* subsection System management
* set Number of solution variables = 1
* subsection Solution variables
* set Solution variable 1 = liquid_pressure
* end
* subsection Equations
* set Equation 1 = Liquid Water Transport Equation
* end
* end
*

Boundary conditions

To be well-posed, these equations are equipped with the appropriate boundary conditions. All the boundary conditions can be described by boundary_id (p_c) and boundary_type.

Some boundary types require additional information, which can also be provided by the parameter file. We consider following types of boundary conditions:

Remarks
  • There is no provision to specify boundary indicators for No liquid water flux or Symmetric boundary conditions, as FEM formulation automatically implies a particular boundary is one of these cases, by default.

Usage Details:

For users:

* subsection Equations
* subsection Liquid Water Transport Equation
* set Liquid pressure at GDL Channel interface, [Pa] = 100.0
* set Breakthrough liquid pressure at GDL Channel interface, [Pa] = 0.0
* set Breakthrough phenomenon at GDL/channel interface = false
* subsection Boundary data
* set liquid_pressure = 4: -5e7
* end
* subsection Boundary conditions
* set Constant liquid pressure current flux boundary conditions = 2: 1e-3
* end
* end
*

Note that the breakthrough pressure is measured experimentally. If the breakthrough phenomenon has been turned off, a fixed pressure will be imposed on the GDL channel interface. The current flux of pressure on the boundary is the a parameter which controls how much liquid flux can go from the boundary.

Usage Details:

* // Creating Equation object (in Application Header file)
*
* // Declare parameters in application
* liquid_transport.declare_parameters(param);
*
* // Initialize in application
* liquid_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( liquid_transport.get_internal_cell_couplings() );
*
* // Making cell couplings using SystemManagement object created in the application
*
* // cell_matrix in application
* // Do a check against layer and it should match with the layers currently working for this equation class.
* // for eg: CCL is FuelCellShop::Layer::HomogeneousCL<dim> object.
* liquid_transport.assemble_cell_matrix(cell_matrices, cell_info, &CCL);
*
* // cell_residual in application
* liquid_transport.assemble_cell_residual(cell_vector, cell_info, &CCL);
*
Note
This class will work if thermal transport equation and fick's diffusion equation for water vapour is not solved for in the application. But it is always suggested to solve with those two equations together.This class does not assemble the phase change source/sink term for the water vapour transport inside the layers.
Author
J. Zhou, Marc Secanell 2015

Constructor & Destructor Documentation

Constructor.

Destructor.

Member Function Documentation

template<int dim>
virtual void FuelCellShop::Equation::LiquidPressureEquation< 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.

Currently, NOT IMPLEMENTED.

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

template<int dim>
virtual void FuelCellShop::Equation::LiquidPressureEquation< dim >::assemble_bdry_residual ( FuelCell::ApplicationCore::FEVector bdry_rhs,
const typename FuelCell::ApplicationCore::DoFApplication< dim >::FaceInfo &  bdry_info,
FuelCellShop::Layer::BaseLayer< dim > *const  layer 
)
virtual

Assemble local boundary residual.

Currently, NOT IMPLEMENTED.

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

template<int dim>
virtual void FuelCellShop::Equation::LiquidPressureEquation< 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::LiquidPressureEquation< dim >::assemble_cell_residual ( FuelCell::ApplicationCore::FEVector cell_rhs,
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::LiquidPressureEquation< dim >::declare_parameters ( ParameterHandler &  param) const
virtual

Declare parameters.

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

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

Initialize parameters.

This class will call make_internal_cell_couplings and make_boundary_types.

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

template<int dim>
virtual void FuelCellShop::Equation::LiquidPressureEquation< 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) .

Currently, NOT IMPLEMENTED.

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

template<int dim>
virtual void FuelCellShop::Equation::LiquidPressureEquation< 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) .

Currently, NOT IMPLEMENTED.

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

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

This function fills out boundary_types.

Currently, NOT IMPLEMENTED.

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

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

This function fills out output_types.

Currently, NOT IMPLEMENTED.

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

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

The function prints out the equation's info.

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

Member Data Documentation

template<int dim>
bool FuelCellShop::Equation::LiquidPressureEquation< dim >::breakthrough_phenomenon
protected

This data member is set to true when the following boundary condition is to be applied at the boundary $ N = k (p_{c(threshold)} - p_c) $.

* subsection subsection Equations
* subsection Liquid Water liquid Transport Equation
* set Breakthrough phenomenon at GDL/channel interface = true
* end
* end
*
template<int dim>
std::map<unsigned int, double> FuelCellShop::Equation::LiquidPressureEquation< dim >::capillary_pressure_current_flux_map
protected

The map is used as the $ k_1 $ parameter in $ N = k_1 (p_{c(threshold)} - p_c) $.

* subsection subsection Equations
* subsection Liquid Water liquid Transport Equation
* set Constant liquid Pressure Current Flux Boundary Conditions = 2: 1e-3
* end
* end
*
template<int dim>
bool FuelCellShop::Equation::LiquidPressureEquation< 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::vector< Tensor<2,dim> > FuelCellShop::Equation::LiquidPressureEquation< dim >::dT_rhok_mu_cell
protected

Computing the derivitive with respect to $ T $ in the function $ \frac{\rho_l \hat{k}_l(p_c)}{\mu_l(T)} $, at all quadrature points in the cell.

template<int dim>
unsigned int FuelCellShop::Equation::LiquidPressureEquation< dim >::fickswater_blockindex_pliquid
protected

Block index for "liquid_pressure" corresponding to "Ficks Transport Equation - water".

template<int dim>
unsigned int FuelCellShop::Equation::LiquidPressureEquation< dim >::fickswater_blockindex_trev
protected

Block index for "temperature_of_REV" corresponding to "Ficks Transport Equation - water".

template<int dim>
unsigned int FuelCellShop::Equation::LiquidPressureEquation< dim >::fickswater_blockindex_xwater
protected

Block index for "water_molar_fraction" corresponding to "Ficks Transport Equation - water".

template<int dim>
std::vector< std::vector< Tensor<1,dim> > > FuelCellShop::Equation::LiquidPressureEquation< dim >::grad_phi_p_cell
protected

$ \mathbf{p_c} $ shape function gradients.

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

template<int dim>
std::vector< std::vector< Tensor<1,dim> > > FuelCellShop::Equation::LiquidPressureEquation< 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>
double FuelCellShop::Equation::LiquidPressureEquation< dim >::k_sat
protected

Saturated permeability inside the porous media, $ k_{sat} ~\left[ cm^2 \right] $.

template<int dim>
double FuelCellShop::Equation::LiquidPressureEquation< dim >::k_sat_bdry
protected

Saturated permeability at the boundary of the domain, $ k_{satbdry} ~\left[ cm^2 \right] $.

template<int dim>
unsigned int FuelCellShop::Equation::LiquidPressureEquation< 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::LiquidPressureEquation< 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>
double FuelCellShop::Equation::LiquidPressureEquation< dim >::M_water
protected

Molar weight of water in grams/mole.

template<int dim>
double FuelCellShop::Equation::LiquidPressureEquation< dim >::p_cell
protected

Total pressure $ \left[ Pa \right] $ in the cell.

template<int dim>
double FuelCellShop::Equation::LiquidPressureEquation< dim >::p_channel_GDL
protected

The boundary condition term which indicates the liquid pressure (Pa) $ p_{c(threshold1)} $ in $ N = k_1 (p_{c(threshold1)} - p_c) $.

* subsection subsection Equations
* subsection Liquid Water liquid Transport Equation
* set liquid pressure at GDL Channel interface, [Pa] = 100.0
* end
* end
*
template<int dim>
double FuelCellShop::Equation::LiquidPressureEquation< dim >::p_channel_GDL_second
protected

The boundary condition term which indicates the breakthrough liquid pressure (Pa).

This is initialized

* subsection subsection Equations
* subsection Liquid Water liquid Transport Equation
* set Breakthrough liquid pressure at GDL Channel interface, [Pa] = 0.0
* end
* end
*
template<int dim>
VariableInfo FuelCellShop::Equation::LiquidPressureEquation< dim >::p_liquid_water
protected

VariableInfo structure corresponding to "liquid_pressure".

template<int dim>
std::vector< std::vector<double> > FuelCellShop::Equation::LiquidPressureEquation< dim >::phi_p_cell
protected

$ \mathbf{p_c} $ shape functions.

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

template<int dim>
std::vector< std::vector<double> > FuelCellShop::Equation::LiquidPressureEquation< dim >::phi_p_cell_bdry
protected

$ \mathbf{p_c} $ shape functions at the boundary.

phi_p_cell_bdry [ q ] [ k ] denotes $ k $-th $ \mathbf{p_c} $ shape function computed in $ q $-th quadrature point at the boundary.

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

$ \mathbf{x_{H_2O}} $ shape functions.

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

template<int dim>
std::vector<bool> FuelCellShop::Equation::LiquidPressureEquation< dim >::reached_breakthrough_pressure
protected

Once the liquid pressure at the GDL and channel interface reaches the breakthrouhg pressure the following boundary condition is to set to be applied $ N = k (p_{c(threshold)} - p_c) $.

And this vector of bool is set to be true to lock the boundary condition.

template<int dim>
double FuelCellShop::Equation::LiquidPressureEquation< dim >::rho_l
protected

Density of liquid water, $ \rho_l ~\left[ \frac{g}{cm^3} \right] $.

template<int dim>
std::vector< Tensor<2,dim> > FuelCellShop::Equation::LiquidPressureEquation< dim >::rho_mu_dk_dp_cell
protected

Computing the derivitive with respect to $ p_c $ in the function $ \frac{\rho_l \hat{k}_l(p_c)}{\mu_l(T)} $, at all quadrature points in the cell.

template<int dim>
std::vector< Tensor<2,dim> > FuelCellShop::Equation::LiquidPressureEquation< dim >::rhok_mu_cell
protected

$ \frac{\rho_l \hat{k}_l}{\mu_l} $, at all quadrature points in the cell.

template<int dim>
std::vector< Tensor<2,dim> > FuelCellShop::Equation::LiquidPressureEquation< dim >::rhok_mu_cell_bdry
protected

$ \frac{\rho_l \hat{k}_l}{\mu_l} $, at all quadrature points in the boundary.

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

VariableInfo structure corresponding to "temperature_of_REV".

template<int dim>
double FuelCellShop::Equation::LiquidPressureEquation< dim >::temperature
protected

Temperature [in K] of the layer.

template<int dim>
unsigned int FuelCellShop::Equation::LiquidPressureEquation< dim >::thermal_blockindex_pliquid
protected

Block index for "liquid_pressure" corresponding to "Thermal Transport Equation".

template<int dim>
unsigned int FuelCellShop::Equation::LiquidPressureEquation< dim >::thermal_blockindex_trev
protected

Block index for "temperature_of_REV" corresponding to "Thermal Transport Equation".

template<int dim>
unsigned int FuelCellShop::Equation::LiquidPressureEquation< dim >::thermal_blockindex_xwater
protected

Block index for "water_molar_fraction" corresponding to "Thermal Transport Equation".

template<int dim>
VariableInfo FuelCellShop::Equation::LiquidPressureEquation< dim >::x_water
protected

VariableInfo structure corresponding to "water_molar_fraction".


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