OpenFCST: The open-source Fuel Cell Simulation Toolbox
|
Virtual class used to implement properties that are characteristic of a porous layer. More...
#include <porous_layer.h>
Public Member Functions | |
Initalization | |
void | set_gases_and_compute (std::vector< FuelCellShop::Material::PureGas * > &gases_in, const double &pressure_in, const double &temperature_in) |
Member function used to store all the gases that are in the pore space in the gas diffusion layer as well as their temperature [Kelvin ] and total pressure [atm ]. More... | |
void | compute_gas_diffusion (FuelCellShop::Material::PureGas *solute_gas, FuelCellShop::Material::PureGas *solvent_gas) |
Member function used to compute bulk diffusion coefficients (and derivatives w.r.t temperature for non-isothermal case and store inside the layer). More... | |
void | set_gases (std::vector< FuelCellShop::Material::PureGas * > &gases_in, const double &pressure_in) |
Member function used to store all the gases that are in the pore space in the porous layer. More... | |
void | set_gas_mixture (FuelCellShop::Material::GasMixture &rgas_mixture) |
Set gas_mixture . More... | |
void | set_porosity_permeability_tortuosity_booleans (const bool &rporosity_is_constant, const bool &rpermeability_is_constant, const bool &rtortuosity_is_constant) |
Set. More... | |
void | set_pressure (const SolutionVariable &p_in) |
Member function used to set the temperature [Kelvin ] at every quadrature point inside the cell. More... | |
void | set_temperature (const SolutionVariable &T_in) |
Member function used to set the temperature [Kelvin ] at every quadrature point inside the cell. More... | |
void | set_saturation (const SolutionVariable &s_in) |
Member function used to set the liquid water saturation at every quadrature point inside the cell. More... | |
void | set_capillary_pressure (const SolutionVariable &p_in) |
Member function used to set the liquid water saturation at every quadrature point inside the cell. More... | |
Accessors and info | |
FuelCellShop::Material::PureGas * | get_gas_pointer (int index) const |
Return the FuelCellShop::Material::PureGas pointer that is stored inside the class in the ith position. More... | |
std::vector < FuelCellShop::Material::PureGas * > | get_gases () const |
Returns the vector of FuelCellShop::Material::PureGas pointers stored in the porous layer. More... | |
const FuelCellShop::Material::GasMixture *const | get_gas_mixture () const |
This function returns gas_mixture . More... | |
void | get_gas_index (FuelCellShop::Material::PureGas *gas_type, int &index) const |
Return the gas index in the GDL class. More... | |
void | get_T_and_p (double &T, double &p) const |
Return the constant temperature [Kelvin ] and constant pressure [atm ] inside the layer. More... | |
void | get_p (double &p) const |
Return the constant pressure [atm ] inside the layer. More... | |
const bool & | get_porosity_is_constant () const |
This function returns porosity_is_constant . More... | |
const bool & | get_permeability_is_constant () const |
This function returns permeability_is_constant . More... | |
const bool & | get_tortuosity_is_constant () const |
This function returns tortuosity_is_constant . More... | |
double | get_porosity () const |
This function computes constant porosity in quadrature points of a mesh entity. More... | |
void | get_porosity (std::vector< double > &dst) const |
This function computes constant porosity in quadrature points of a mesh entity. More... | |
void | get_porosity (std::vector< double > &dst, const std::vector< Point< dim > > &points) const |
This function computes variable porosity in quadrature points of a mesh entity. More... | |
void | get_permeability (std::vector< SymmetricTensor< 2, dim > > &dst) const |
This function computes constant permeability in quadrature points of a mesh entity. More... | |
void | get_permeability (std::vector< SymmetricTensor< 2, dim > > &dst, const std::vector< Point< dim > > &points) const |
This function computes variable permeability in quadrature points of a mesh entity. More... | |
void | get_SQRT_permeability (std::vector< SymmetricTensor< 2, dim > > &dst) const |
This function computes square root of constant permeability in quadrature points of a mesh entity. More... | |
void | get_SQRT_permeability (std::vector< SymmetricTensor< 2, dim > > &dst, const std::vector< Point< dim > > &points) const |
This function computes square root of variable permeability in quadrature points of a mesh entity. More... | |
void | get_permeability_INV (std::vector< SymmetricTensor< 2, dim > > &dst) const |
This function computes inverse of constant permeability in quadrature points of a mesh entity. More... | |
void | get_permeability_INV (std::vector< SymmetricTensor< 2, dim > > &dst, const std::vector< Point< dim > > &points) const |
This function computes inverse of variable permeability in quadrature points of a mesh entity. More... | |
void | get_SQRT_permeability_INV (std::vector< SymmetricTensor< 2, dim > > &dst) const |
This function computes inverse of square root of constant permeability in quadrature points of a mesh entity. More... | |
void | get_SQRT_permeability_INV (std::vector< SymmetricTensor< 2, dim > > &dst, const std::vector< Point< dim > > &points) const |
This function computes inverse of square root of variable permeability in quadrature points of a mesh entity. More... | |
void | get_Forchheimer_permeability (std::vector< SymmetricTensor< 2, dim > > &dst) const |
This function computes constant Forchheimer permeability in quadrature points of a mesh entity. More... | |
void | get_Forchheimer_permeability (std::vector< SymmetricTensor< 2, dim > > &dst, const std::vector< Point< dim > > &points) const |
This function computes variable Forchheimer permeability in quadrature points of a mesh entity. More... | |
void | get_tortuosity (std::vector< SymmetricTensor< 2, dim > > &dst) const |
This function computes constant tortuosity in quadrature points of a mesh entity. More... | |
void | get_tortuosity (std::vector< SymmetricTensor< 2, dim > > &dst, const std::vector< Point< dim > > &points) const |
This function computes variable tortuosity in quadrature points of a mesh entity. More... | |
virtual void | effective_gas_diffusivity (std::vector< Tensor< 2, dim > > &) const |
Return the effective diffusivity [m^2/s ] in the GDL. More... | |
virtual void | effective_gas_diffusivity (Table< 2, double > &D_eff) const |
Return the effective diffusivty in the GDL for all the gases assigned to the layer using set_gases_and_compute. More... | |
virtual void | effective_gas_diffusivity (Table< 2, Tensor< 2, dim > > &D_eff) const |
Return a tensor with the effective diffusivty in the GDL for all the gases assigned to the layer using set_gases_and_compute. More... | |
virtual void | derivative_effective_gas_diffusivity (std::map< VariableNames, std::vector< Tensor< 2, dim > > > &) const |
Return the derivative of effective diffusivity w.r.t solution variables/design parameters in the GDL. More... | |
virtual void | gas_diffusion_coefficient (std::vector< double > &D_b) const |
Member function used to compute diffusion for a solute_gas, solvent_gas combination at a given temperature and pressure. More... | |
virtual void | gas_diffusion_coefficient (std::vector< double > &D_b, std::vector< double > &dD_b_dT) const |
Member function used to compute diffusion for a solute_gas, solvent_gas combination at a given temperature and pressure. More... | |
void | molecular_gas_diffusion_coefficient (std::vector< double > &D_m) const |
Member function used to compute molecular diffusion for a solute_gas, solvent_gas combination at a given temperature and pressure. More... | |
void | molecular_gas_diffusion_coefficient (std::vector< double > &D_m, std::vector< double > &dD_m_dT) const |
Member function used to compute molecular diffusion for a solute_gas, solvent_gas combination at a given temperature and pressure. More... | |
void | Knudsen_diffusion (std::vector< double > &D) const |
Member function used to get the Knudsen diffusivity in the layer after calling compute_gas_diffusion. More... | |
void | Knudsen_diffusion (std::vector< double > &D, std::vector< double > &dD_dT) const |
Member function used to compute the Knudsen diffusivity in the layer.after calling compute_gas_diffusion. More... | |
virtual void | pcapillary (std::vector< double > &) const |
Compute , at all quadrature points in the cell. More... | |
virtual void | saturation_from_capillary_equation (std::vector< double > &) const |
virtual void | derivative_saturation_from_capillary_equation_PSD (std::vector< double > &) const |
virtual void | saturated_liquid_permeablity_PSD (double &) const |
virtual void | interfacial_surface_area_PSD (std::vector< double > &) const |
virtual void | relative_liquid_permeability_PSD (std::vector< Tensor< 2, dim > > &) const |
virtual void | derivative_relative_liquid_permeablity_PSD (std::vector< double > &) const |
virtual void | derivative_interfacial_surface_area_PSD (std::vector< double > &) const |
Computing member functions | |
void | compute_Knudsen_diffusion (const FuelCellShop::Material::PureGas *solute_gas, const SolutionVariable &T_in, std::vector< double > &D_k) const |
Member function used to compute the Knudsen diffusivity in the layer. More... | |
void | compute_Knudsen_diffusion (const FuelCellShop::Material::PureGas *solute_gas, const SolutionVariable &T_in, std::vector< double > &D_k, std::vector< double > &dD_k_dT) const |
Member function used to compute the Knudsen diffusivity in the layer. More... | |
Debugging and info | |
virtual void | print_layer_properties () const |
This member function is a virtual class that can be used to output to screen information from the layer. More... | |
Public Member Functions inherited from FuelCellShop::Layer::BaseLayer< dim > | |
virtual void | set_derivative_flags (const std::vector< VariableNames > &flags) |
Set the variables for which you would like to compute the derivatives. More... | |
void | set_position (const std::vector< Point< dim > > &p) |
Member function used by some applications such as dummyGDL in order to know which value to return. More... | |
virtual void | set_local_material_id (const unsigned int &id) |
Function for setting local material id, for unit testing purposes. More... | |
void | unset_local_material_id () |
Function for unsetting local material id, so that it isn't incorrectly used later Once the key is "unset" to some invalid value, an error will be thrown if the key is requested again without being set. More... | |
virtual void | set_constant_solution (const double &value, const VariableNames &name) |
Set those solution variables which are constant in the particular application. More... | |
virtual void | set_solution (const std::vector< SolutionVariable > &) |
If the effective properties in the layer depend on the solution, the solution for a given cell should be passed to the class using this member function. More... | |
bool | belongs_to_material (const unsigned int material_id) |
Check if a given cell belongs to the catalyst layer and assign. More... | |
const std::string & | name_layer () const |
Return the name of the layer. More... | |
virtual const std::type_info & | get_base_type () const |
This member function return the name of the type of layer, i.e. More... | |
virtual bool | test_layer () |
This virtual class should be used for any derived class to be able to test the functionality of the class. More... | |
std::vector< unsigned int > | get_material_ids () |
Return the local material id of the layer. More... | |
unsigned int | local_material_id () const |
Return the local material id of the layer, performs a check. More... | |
Protected Member Functions | |
Constructors, destructor, and initialization | |
PorousLayer (const std::string &name) | |
Constructor. More... | |
PorousLayer () | |
Constructor. More... | |
PorousLayer (const std::string &name, FuelCellShop::Material::GasMixture &gas_mixture) | |
Constructor. More... | |
virtual | ~PorousLayer () |
Destructor. More... | |
virtual void | declare_parameters (const std::string &name, ParameterHandler ¶m) const |
Declare parameters for a parameter file. More... | |
virtual void | declare_parameters (ParameterHandler ¶m) const |
Declare parameters for a parameter file. More... | |
virtual void | initialize (ParameterHandler ¶m) |
Member function used to read in data and initialize the necessary data to compute the coefficients. More... | |
Minor functions | |
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... | |
Internal accessors and info | |
virtual void | gas_diffusion_coefficients (Table< 2, double > &) const |
Return the molecular diffusivty all the gases assigned to the layer using set_gases_and_compute. More... | |
virtual void | derivative_gas_diffusion_coefficients (std::vector< Table< 2, double > > &) const |
Return the derivative of the molecular diffusion coefficient with respect to the derivative flags for all the gases assigned to the layer using set_gases_and_compute. More... | |
Protected Member Functions inherited from FuelCellShop::Layer::BaseLayer< dim > | |
BaseLayer () | |
Constructor. More... | |
BaseLayer (const std::string &name) | |
Constructor. More... | |
virtual | ~BaseLayer () |
Destructor. More... | |
virtual void | set_parameters (const std::string &object_name, const std::vector< std::string > &name_dvar, const std::vector< double > &value_dvar, ParameterHandler ¶m) |
Member function used to change the values in the parameter file for a given list of parameters. More... | |
virtual void | set_parameters (const std::vector< std::string > &name_dvar, const std::vector< double > &value_dvar, ParameterHandler ¶m) |
Set parameters in parameter file. More... | |
Protected Attributes | |
Layer properties | |
FuelCellShop::Material::GasMixture * | gas_mixture |
Gas mixture. More... | |
std::vector < FuelCellShop::Material::PureGas * > | gases |
Gases inside a porous layer. More... | |
bool | porosity_is_constant |
Variable defining if the porosity is constant. More... | |
bool | permeability_is_constant |
Variable defining if the permeability is constant. More... | |
bool | tortuosity_is_constant |
Variable defining if the tortuosity is constant. More... | |
double | porosity |
User defined constant porosity. More... | |
bool | use_Bosanquet |
Boolean flag that specifies if Knudsen effects should be accounted for. More... | |
double | Knudsen_radius |
Parameter used to define Knudsen pore radius. More... | |
SymmetricTensor< 2, dim > | permeability |
User defined constant permeability, m^2. More... | |
SymmetricTensor< 2, dim > | SQRT_permeability |
Square root of user defined constant permeability, m. More... | |
SymmetricTensor< 2, dim > | permeability_INV |
Inverse of user defined constant permeability, 1/m^2. More... | |
SymmetricTensor< 2, dim > | SQRT_permeability_INV |
Inverse of square root of user defined constant permeability, 1/m. More... | |
SymmetricTensor< 2, dim > | Forchheimer_permeability |
User defined constant Forchheimer permeability, 1/m. More... | |
SymmetricTensor< 2, dim > | tortuosity |
User defined constant tortuosity. More... | |
std::string | diffusion_species_name |
If GDL properties are stored inside the class (e.g DummyGDL) then, return the property stored under coefficient_name name. More... | |
double | temperature |
Temperature [K ] used to compute gas diffusivity. More... | |
double | pressure |
Total pressure [atm ] used to compute gas diffusivity. More... | |
SolutionVariable | p_vector |
Pressure at every quadrature point inside the cell in [Pa]. More... | |
SolutionVariable | T_vector |
Temperature at every quadrature point inside the cell in [K]. More... | |
SolutionVariable | s_vector |
Liquid water saturation at every quadrature point inside the cell [-]. More... | |
SolutionVariable | capillary_pressure_vector |
Liquid water capillary pressure at every quadrature point inside the cell in [Pa]. More... | |
Table< 2, double > | D_ECtheory |
Tensor of diffusion coefficients This are computed with setting up the gas so that they do not need to be recomputed all the time. More... | |
std::vector< Table< 2, double > > | dD_ECtheory_dx |
Vector of tensors for the derivative of the diffusion coefficients – This are computed with setting up the gas so that they do not need to be recomputed all the time. More... | |
std::vector< double > | D_molecular |
Vector of molecular diffusion coefficients at every quadrature point inside the cell in m^2/s. More... | |
std::vector< double > | dD_molecular_dT |
Vector of derivatives for molecular diffusion coefficients w.r.t temperature, at every quadrature in m^2/s. More... | |
std::vector< double > | D_k |
Vector of Knudsen diffusion coefficients at every quadrature point inside the cell in m^2/s. More... | |
std::vector< double > | dD_k_dT |
Vector of derivatives for Knudsen diffusion coefficients w.r.t temperature, at every quadrature in m^2/s. More... | |
std::vector< double > | D_bulk |
Vector of bulk diffusion coefficients at every quadrature point inside the cell. More... | |
std::vector< double > | dD_bulk_dT |
Vector of derivative of bulk diffusion coefficients w.r.t temperature, at every quadrature point inside the cell. More... | |
bool | PSD_is_used |
Boolean flag to specify if a PSD is to be used to estimate saturation, permeability, etc. More... | |
std::string | PSD_type |
PSD class type from input file. More... | |
boost::shared_ptr < FuelCellShop::MicroScale::BasePSD < dim > > | PSD |
Pointer to the PSD object. More... | |
FuelCellShop::MicroScale::BasePSD < dim > * | psd_pointer |
Pointer to the PSD object. More... | |
FuelCellShop::Material::PureGas * | solute_gas |
Pointer used to store the solute gas for computing binary diffusion coefficients. More... | |
FuelCellShop::Material::PureGas * | solvent_gas |
Pointer used to store the solute gas for computing binary diffusion coefficients. More... | |
Protected Attributes inherited from FuelCellShop::Layer::BaseLayer< dim > | |
const std::string | name |
Name of the layer. More... | |
std::vector< unsigned int > | material_ids |
List of material IDs that belong to the layer. More... | |
std::vector< Point< dim > > | point |
Coordinates of the point where we would like to compute the effective properties. More... | |
std::vector< VariableNames > | derivative_flags |
Flags for derivatives: These flags are used to request derivatives. More... | |
std::map< VariableNames, double > | constant_solutions |
Map storing values of solution variables constant in a particular application. More... | |
Virtual class used to implement properties that are characteristic of a porous layer.
This layer is used to setup information needed to compute molecular diffusivities and Knudsen diffusivities (pending implementation) for any gases. Therefore, it provides member functions to specify the gases that exist within the layer as well as their temperature and pressure.
It also inherits all properties of BaseLayer.
Note you cannot make an object of this class since it still lacks many member functions that are necessary for a full layer implementation.
With respect to gases: you can either use gases
or gas_mixture
. Whatever you prefer. However, gas_mixture
offers more methods in its structure. You cannot use both.
In order to use this class, create a child object of this class. Then, initialize the object. At every cell use set pressure and temperature using set_pressure and set_temperature to specify the operating conditions at that location. Then, use compute_gas_diffusion to obtain values for diffusivity.
A similar process would work for other functions.
|
inlineprotected |
Constructor.
References temperature_of_REV, and total_pressure.
|
inlineprotected |
Constructor.
|
inlineprotected |
Constructor.
References temperature_of_REV, and total_pressure.
|
inlineprotectedvirtual |
Destructor.
void FuelCellShop::Layer::PorousLayer< dim >::compute_gas_diffusion | ( | FuelCellShop::Material::PureGas * | solute_gas, |
FuelCellShop::Material::PureGas * | solvent_gas | ||
) |
Member function used to compute bulk diffusion coefficients (and derivatives w.r.t temperature for non-isothermal case and store inside the layer).
This function takes solute and solvent gas as an input argument (Fick's diffusion model). Before calling this function, pressure [atm
] and temperature [Kelvin
] must be set using set_pressure and set_temperature method.
void FuelCellShop::Layer::PorousLayer< dim >::compute_Knudsen_diffusion | ( | const FuelCellShop::Material::PureGas * | solute_gas, |
const SolutionVariable & | T_in, | ||
std::vector< double > & | D_k | ||
) | const |
Member function used to compute the Knudsen diffusivity in the layer.
void FuelCellShop::Layer::PorousLayer< dim >::compute_Knudsen_diffusion | ( | const FuelCellShop::Material::PureGas * | solute_gas, |
const SolutionVariable & | T_in, | ||
std::vector< double > & | D_k, | ||
std::vector< double > & | dD_k_dT | ||
) | const |
Member function used to compute the Knudsen diffusivity in the layer.
The equation used to compute Knudsen diffusivity is
where R is 8.3144 J/mol K in SI units, molecular weight M_{i} is expressed in units of kg/mol and temperature T has units of K.
|
protectedvirtual |
Declare parameters for a parameter file.
Reimplemented from FuelCellShop::Layer::BaseLayer< dim >.
Reimplemented in FuelCellShop::Layer::CatalystLayer< dim >, FuelCellShop::Layer::GasDiffusionLayer< dim >, FuelCellShop::Layer::MicroPorousLayer< dim >, FuelCellShop::Layer::ConventionalCL< dim >, FuelCellShop::Layer::DesignMPL< dim >, FuelCellShop::Layer::MultiScaleCL< dim >, FuelCellShop::Layer::MultiScaleCL< deal_II_dimension >, FuelCellShop::Layer::DesignFibrousGDL< dim >, FuelCellShop::Layer::DummyCL< dim >, FuelCellShop::Layer::DummyGDL< dim >, FuelCellShop::Layer::SGL24BA< dim >, and FuelCellShop::Layer::HomogeneousCL< dim >.
|
inlineprotectedvirtual |
Declare parameters for a parameter file.
Reimplemented from FuelCellShop::Layer::BaseLayer< dim >.
Reimplemented in FuelCellShop::Layer::GasDiffusionLayer< dim >, FuelCellShop::Layer::MultiScaleCL< dim >, FuelCellShop::Layer::MultiScaleCL< deal_II_dimension >, FuelCellShop::Layer::DesignMPL< dim >, FuelCellShop::Layer::SGL24BA< dim >, FuelCellShop::Layer::SGL24BC< dim >, FuelCellShop::Layer::DummyGDL< dim >, FuelCellShop::Layer::DummyCL< dim >, and FuelCellShop::Layer::DummyMPL< dim >.
|
inlinevirtual |
Return the derivative of effective diffusivity w.r.t solution variables/design parameters in the GDL.
It transforms bulk diffusion derivative properties computed using compute_gas_diffusion method and transforms it into an effective property, taking into account the porosity, saturation and GDL structure (Anisotropic case), at all quadrature points of the cell.
Reimplemented in FuelCellShop::Layer::CatalystLayer< dim >, FuelCellShop::Layer::GasDiffusionLayer< dim >, FuelCellShop::Layer::MicroPorousLayer< dim >, FuelCellShop::Layer::DesignMPL< dim >, FuelCellShop::Layer::ConventionalCL< dim >, FuelCellShop::Layer::DummyGDL< dim >, FuelCellShop::Layer::DesignFibrousGDL< dim >, FuelCellShop::Layer::DummyCL< dim >, and FuelCellShop::Layer::DummyMPL< dim >.
|
protectedvirtual |
Return the derivative of the molecular diffusion coefficient with respect to the derivative flags for all the gases assigned to the layer using set_gases_and_compute.
This routine uses FuelCellShop::Mixture::ChapmanEnskog to compute the binary diffusivity for each gas with respect to each other. If the layer contains three gases the vector returns D12, D13, D23. For 4 gases, it returns D12, D13, D14, D23, D24, D34 . Derivatives are provided in the vector w.r.t following parameters in the same order as below:
|
inlinevirtual |
Reimplemented in FuelCellShop::Layer::CatalystLayer< dim >, FuelCellShop::Layer::GasDiffusionLayer< dim >, FuelCellShop::Layer::MicroPorousLayer< dim >, FuelCellShop::Layer::ConventionalCL< dim >, FuelCellShop::Layer::DesignMPL< dim >, and FuelCellShop::Layer::DesignFibrousGDL< dim >.
References FcstUtilities::log.
|
inlinevirtual |
Reimplemented in FuelCellShop::Layer::CatalystLayer< dim >, FuelCellShop::Layer::GasDiffusionLayer< dim >, FuelCellShop::Layer::MicroPorousLayer< dim >, FuelCellShop::Layer::ConventionalCL< dim >, FuelCellShop::Layer::DesignMPL< dim >, and FuelCellShop::Layer::DesignFibrousGDL< dim >.
References FcstUtilities::log.
|
inlinevirtual |
Reimplemented in FuelCellShop::Layer::CatalystLayer< dim >, FuelCellShop::Layer::GasDiffusionLayer< dim >, FuelCellShop::Layer::MicroPorousLayer< dim >, FuelCellShop::Layer::ConventionalCL< dim >, FuelCellShop::Layer::DesignMPL< dim >, and FuelCellShop::Layer::DesignFibrousGDL< dim >.
References FcstUtilities::log.
|
inlinevirtual |
Return the effective diffusivity [m^2/s
] in the GDL.
It transforms bulk diffusivity, computed using compute_gas_diffusion method and transforms it into an effective property, taking into account the porosity, saturation and GDL structure (Anisotropic case), at all quadrature points of the cell.
Reimplemented in FuelCellShop::Layer::CatalystLayer< dim >, FuelCellShop::Layer::GasDiffusionLayer< dim >, FuelCellShop::Layer::MicroPorousLayer< dim >, FuelCellShop::Layer::DesignMPL< dim >, FuelCellShop::Layer::ConventionalCL< dim >, FuelCellShop::Layer::DummyCL< dim >, FuelCellShop::Layer::DummyGDL< dim >, FuelCellShop::Layer::DesignFibrousGDL< dim >, and FuelCellShop::Layer::DummyMPL< dim >.
|
inlinevirtual |
Return the effective diffusivty in the GDL for all the gases assigned to the layer using set_gases_and_compute.
This routine uses FuelCellShop::Mixture::ChapmanEnskog to compute the binary diffusivity for each gas with respect to each other. If the layer contains three gases the vector returns D12, D13, D23. For 4 gases, it returns D12, D13, D14, D23, D24, D34.
Reimplemented in FuelCellShop::Layer::GasDiffusionLayer< dim >, and FuelCellShop::Layer::DummyGDL< dim >.
References FcstUtilities::log.
|
inlinevirtual |
Return a tensor with the effective diffusivty in the GDL for all the gases assigned to the layer using set_gases_and_compute.
This routine uses FuelCellShop::Mixture::ChapmanEnskog to compute the binary diffusivity for each gas with respect to each other. If the layer contains three gases the vector returns D12, D13, D23. For 4 gases, it returns D12, D13, D14, D23, D24, D34 (Anisotropic case).
Reimplemented in FuelCellShop::Layer::CatalystLayer< dim >, FuelCellShop::Layer::GasDiffusionLayer< dim >, FuelCellShop::Layer::MicroPorousLayer< dim >, FuelCellShop::Layer::DesignMPL< dim >, FuelCellShop::Layer::ConventionalCL< dim >, FuelCellShop::Layer::DummyGDL< dim >, FuelCellShop::Layer::SGL24BA< dim >, FuelCellShop::Layer::DesignFibrousGDL< dim >, FuelCellShop::Layer::SGL24BC< dim >, FuelCellShop::Layer::DummyCL< dim >, and FuelCellShop::Layer::DummyMPL< dim >.
References FcstUtilities::log.
|
inlinevirtual |
Member function used to compute diffusion for a solute_gas, solvent_gas combination at a given temperature and pressure.
The diffusion coefficient is returned in SI units, i.e., m^2/s.
The effective diffusivity returned is either the molecular diffusivity or, if the flag "set Use Bosanquet approx." use set to true, using the Bousinesq approximation is used in order to compute the bulk combined diffusivity such that
where is the molecular diffusion coefficient obtained using effective_molecular_gas_diffusion_coefficient and is the Knudsen diffusivity computed using Knudsen_diffusion
By defaul the diffusion coefficient is computed using the pore radius in:
If you would like to neglect the diffusion coefficient, simply make this value very large (the default basically results on negligible Knudsen effects – the value is: 161.0949 (SI) compared to 1e-5 (SI) for molecular diff for gases – Note that you add the inverse of these values).
References
|
inlinevirtual |
Member function used to compute diffusion for a solute_gas, solvent_gas combination at a given temperature and pressure.
The diffusion coefficient is returned in SI units, i.e., m^2/s.
The effective diffusivity returned is either the molecular diffusivity or, if the flag "set Use Bosanquet approx." use set to true, using the Bousinesq approximation is used in order to compute the bulk combined diffusivity such that
where is the molecular diffusion coefficient obtained using effective_molecular_gas_diffusion_coefficient and is the Knudsen diffusivity computed using Knudsen_diffusion
References
|
protectedvirtual |
Return the molecular diffusivty all the gases assigned to the layer using set_gases_and_compute.
This routine uses FuelCellShop::Mixture::ChapmanEnskog to compute the binary diffusivity for each gas with respect to each other. If the layer contains three gases the vector returns D12, D13, D23. For 4 gases, it returns D12, D13, D14, D23, D24, D34.
void FuelCellShop::Layer::PorousLayer< dim >::get_Forchheimer_permeability | ( | std::vector< SymmetricTensor< 2, dim > > & | dst | ) | const |
This function computes constant Forchheimer permeability in quadrature points of a mesh entity.
void FuelCellShop::Layer::PorousLayer< dim >::get_Forchheimer_permeability | ( | std::vector< SymmetricTensor< 2, dim > > & | dst, |
const std::vector< Point< dim > > & | points | ||
) | const |
This function computes variable Forchheimer permeability in quadrature points of a mesh entity.
void FuelCellShop::Layer::PorousLayer< dim >::get_gas_index | ( | FuelCellShop::Material::PureGas * | gas_type, |
int & | index | ||
) | const |
Return the gas index in the GDL class.
If the gas does not exist in the GDL it will return -1 as the index.
|
inline |
This function returns gas_mixture
.
|
inline |
Return the FuelCellShop::Material::PureGas pointer that is stored inside the class in the ith position.
The diffusion coefficient D[i][j] would be the diffusion coefficient of gas type i with j.
|
inline |
Returns the vector of FuelCellShop::Material::PureGas pointers stored in the porous layer.
This method is useful to access all gases which are being transported inside a particular porous layer object.
|
inline |
Return the constant pressure [atm
] inside the layer.
void FuelCellShop::Layer::PorousLayer< dim >::get_permeability | ( | std::vector< SymmetricTensor< 2, dim > > & | dst | ) | const |
This function computes constant permeability in quadrature points of a mesh entity.
void FuelCellShop::Layer::PorousLayer< dim >::get_permeability | ( | std::vector< SymmetricTensor< 2, dim > > & | dst, |
const std::vector< Point< dim > > & | points | ||
) | const |
This function computes variable permeability in quadrature points of a mesh entity.
void FuelCellShop::Layer::PorousLayer< dim >::get_permeability_INV | ( | std::vector< SymmetricTensor< 2, dim > > & | dst | ) | const |
This function computes inverse of constant permeability in quadrature points of a mesh entity.
void FuelCellShop::Layer::PorousLayer< dim >::get_permeability_INV | ( | std::vector< SymmetricTensor< 2, dim > > & | dst, |
const std::vector< Point< dim > > & | points | ||
) | const |
This function computes inverse of variable permeability in quadrature points of a mesh entity.
|
inline |
This function returns permeability_is_constant
.
|
inline |
This function computes constant porosity in quadrature points of a mesh entity.
|
inline |
This function computes constant porosity in quadrature points of a mesh entity.
void FuelCellShop::Layer::PorousLayer< dim >::get_porosity | ( | std::vector< double > & | dst, |
const std::vector< Point< dim > > & | points | ||
) | const |
This function computes variable porosity in quadrature points of a mesh entity.
|
inline |
This function returns porosity_is_constant
.
void FuelCellShop::Layer::PorousLayer< dim >::get_SQRT_permeability | ( | std::vector< SymmetricTensor< 2, dim > > & | dst | ) | const |
This function computes square root of constant permeability in quadrature points of a mesh entity.
void FuelCellShop::Layer::PorousLayer< dim >::get_SQRT_permeability | ( | std::vector< SymmetricTensor< 2, dim > > & | dst, |
const std::vector< Point< dim > > & | points | ||
) | const |
This function computes square root of variable permeability in quadrature points of a mesh entity.
void FuelCellShop::Layer::PorousLayer< dim >::get_SQRT_permeability_INV | ( | std::vector< SymmetricTensor< 2, dim > > & | dst | ) | const |
This function computes inverse of square root of constant permeability in quadrature points of a mesh entity.
void FuelCellShop::Layer::PorousLayer< dim >::get_SQRT_permeability_INV | ( | std::vector< SymmetricTensor< 2, dim > > & | dst, |
const std::vector< Point< dim > > & | points | ||
) | const |
This function computes inverse of square root of variable permeability in quadrature points of a mesh entity.
|
inline |
Return the constant temperature [Kelvin
] and constant pressure [atm
] inside the layer.
void FuelCellShop::Layer::PorousLayer< dim >::get_tortuosity | ( | std::vector< SymmetricTensor< 2, dim > > & | dst | ) | const |
This function computes constant tortuosity in quadrature points of a mesh entity.
void FuelCellShop::Layer::PorousLayer< dim >::get_tortuosity | ( | std::vector< SymmetricTensor< 2, dim > > & | dst, |
const std::vector< Point< dim > > & | points | ||
) | const |
This function computes variable tortuosity in quadrature points of a mesh entity.
|
inline |
This function returns tortuosity_is_constant
.
|
protectedvirtual |
Member function used to read in data and initialize the necessary data to compute the coefficients.
Reimplemented from FuelCellShop::Layer::BaseLayer< dim >.
Reimplemented in FuelCellShop::Layer::CatalystLayer< dim >, FuelCellShop::Layer::GasDiffusionLayer< dim >, FuelCellShop::Layer::ConventionalCL< dim >, FuelCellShop::Layer::MultiScaleCL< dim >, FuelCellShop::Layer::MultiScaleCL< deal_II_dimension >, FuelCellShop::Layer::MicroPorousLayer< dim >, FuelCellShop::Layer::DesignMPL< dim >, FuelCellShop::Layer::HomogeneousCL< dim >, FuelCellShop::Layer::SGL24BA< dim >, FuelCellShop::Layer::SGL24BC< dim >, FuelCellShop::Layer::DummyGDL< dim >, FuelCellShop::Layer::DummyCL< dim >, FuelCellShop::Layer::DummyMPL< dim >, and FuelCellShop::Layer::DesignFibrousGDL< dim >.
|
inlinevirtual |
Reimplemented in FuelCellShop::Layer::CatalystLayer< dim >, FuelCellShop::Layer::GasDiffusionLayer< dim >, FuelCellShop::Layer::MicroPorousLayer< dim >, FuelCellShop::Layer::ConventionalCL< dim >, FuelCellShop::Layer::DesignMPL< dim >, and FuelCellShop::Layer::DesignFibrousGDL< dim >.
References FcstUtilities::log.
|
inline |
Member function used to get the Knudsen diffusivity in the layer after calling compute_gas_diffusion.
|
inline |
Member function used to compute the Knudsen diffusivity in the layer.after calling compute_gas_diffusion.
|
inline |
Member function used to compute molecular diffusion for a solute_gas, solvent_gas combination at a given temperature and pressure.
The molecular diffusivity is returned in units of m^2/s
|
inline |
Member function used to compute molecular diffusion for a solute_gas, solvent_gas combination at a given temperature and pressure.
The molecular diffusivity is returned in units of m^2/s
|
inlinevirtual |
Compute , at all quadrature points in the cell.
Reimplemented in FuelCellShop::Layer::CatalystLayer< dim >, FuelCellShop::Layer::GasDiffusionLayer< dim >, FuelCellShop::Layer::MicroPorousLayer< dim >, FuelCellShop::Layer::ConventionalCL< dim >, FuelCellShop::Layer::DesignMPL< dim >, and FuelCellShop::Layer::DesignFibrousGDL< dim >.
References FcstUtilities::log.
|
protected |
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.
|
virtual |
This member function is a virtual class that can be used to output to screen information from the layer.
This function should be re-implemented in each layer with data that is relevant in each case.
Reimplemented from FuelCellShop::Layer::BaseLayer< dim >.
Reimplemented in FuelCellShop::Layer::MultiScaleCL< dim >, FuelCellShop::Layer::MultiScaleCL< deal_II_dimension >, and FuelCellShop::Layer::ConventionalCL< dim >.
|
inlinevirtual |
Reimplemented in FuelCellShop::Layer::CatalystLayer< dim >, FuelCellShop::Layer::GasDiffusionLayer< dim >, FuelCellShop::Layer::MicroPorousLayer< dim >, FuelCellShop::Layer::ConventionalCL< dim >, FuelCellShop::Layer::DesignMPL< dim >, and FuelCellShop::Layer::DesignFibrousGDL< dim >.
References FcstUtilities::log.
|
inlinevirtual |
Reimplemented in FuelCellShop::Layer::CatalystLayer< dim >, FuelCellShop::Layer::GasDiffusionLayer< dim >, FuelCellShop::Layer::MicroPorousLayer< dim >, FuelCellShop::Layer::ConventionalCL< dim >, FuelCellShop::Layer::DesignMPL< dim >, and FuelCellShop::Layer::DesignFibrousGDL< dim >.
References FcstUtilities::log.
|
inlinevirtual |
Reimplemented in FuelCellShop::Layer::CatalystLayer< dim >, FuelCellShop::Layer::GasDiffusionLayer< dim >, FuelCellShop::Layer::MicroPorousLayer< dim >, FuelCellShop::Layer::ConventionalCL< dim >, FuelCellShop::Layer::DesignMPL< dim >, and FuelCellShop::Layer::DesignFibrousGDL< dim >.
References FcstUtilities::log.
|
inline |
Member function used to set the liquid water saturation at every quadrature point inside the cell.
This function should particulary be used in the case of multi-phase application.
References capillary_pressure, and FuelCellShop::SolutionVariable::get_variablename().
|
inline |
Set gas_mixture
.
|
inline |
Member function used to store all the gases that are in the pore space in the porous layer.
Besides this, it takes total pressure [atm
] as input. This method can be particulary used for non-isothermal applications when temperature changes in every cell, hence this method doesn't compute the diffusion coefficients but store the relevant constant data.
solvent
gas
should be the last object inside the input vector. On not ensuring this may not give any errors, as this function may also be used for multi-component gas transport but for fickian diffusion, it may give wrong results. void FuelCellShop::Layer::PorousLayer< dim >::set_gases_and_compute | ( | std::vector< FuelCellShop::Material::PureGas * > & | gases_in, |
const double & | pressure_in, | ||
const double & | temperature_in | ||
) |
Member function used to store all the gases that are in the pore space in the gas diffusion layer as well as their temperature [Kelvin
] and total pressure [atm
].
Then, when computing diffusion, the class will return each one of the diffusion coefficients for the gases. For example, if we pass a vector with for gases, say oxygen and nitrogen, effective_gas_diffusivity will return
|
inline |
Set.
porosity_is_constant
,permeability_is_constant
,tortuosity_is_constant
.
|
inline |
Member function used to set the temperature [Kelvin
] at every quadrature point inside the cell.
This function should particulary be used in the case of non-isothermal application.
References FuelCellShop::SolutionVariable::get_variablename(), and total_pressure.
|
inline |
Member function used to set the liquid water saturation at every quadrature point inside the cell.
This function should particulary be used in the case of multi-phase application.
References FuelCellShop::SolutionVariable::get_variablename(), and liquid_water_saturation.
|
inline |
Member function used to set the temperature [Kelvin
] at every quadrature point inside the cell.
This function should particulary be used in the case of non-isothermal application.
References FuelCellShop::SolutionVariable::get_variablename(), and temperature_of_REV.
|
protected |
Liquid water capillary pressure at every quadrature point inside the cell in [Pa].
|
protected |
Vector of bulk diffusion coefficients at every quadrature point inside the cell.
This might store the molecular diffusivity or the diffusivity computed using Bosanquet approx.
|
protected |
Tensor of diffusion coefficients This are computed with setting up the gas so that they do not need to be recomputed all the time.
|
protected |
Vector of Knudsen diffusion coefficients at every quadrature point inside the cell in m^2/s.
|
protected |
Vector of molecular diffusion coefficients at every quadrature point inside the cell in m^2/s.
|
protected |
Vector of derivative of bulk diffusion coefficients w.r.t temperature, at every quadrature point inside the cell.
These are computed using compute_gas_diffusion method.
This might store the molecular diffusivity or the diffusivity computed using Bosanquet approx.
|
protected |
Vector of tensors for the derivative of the diffusion coefficients – This are computed with setting up the gas so that they do not need to be recomputed all the time.
|
protected |
Vector of derivatives for Knudsen diffusion coefficients w.r.t temperature, at every quadrature in m^2/s.
|
protected |
Vector of derivatives for molecular diffusion coefficients w.r.t temperature, at every quadrature in m^2/s.
|
protected |
If GDL properties are stored inside the class (e.g DummyGDL) then, return the property stored under coefficient_name name.
|
protected |
User defined constant Forchheimer permeability, 1/m.
|
protected |
Gas mixture.
|
protected |
Gases inside a porous layer.
gas_mixture
instead.
|
protected |
Parameter used to define Knudsen pore radius.
|
protected |
Pressure at every quadrature point inside the cell in [Pa].
|
protected |
User defined constant permeability, m^2.
|
protected |
Inverse of user defined constant permeability, 1/m^2.
|
protected |
Variable defining if the permeability is constant.
|
protected |
User defined constant porosity.
|
protected |
Variable defining if the porosity is constant.
|
protected |
Total pressure [atm
] used to compute gas diffusivity.
. Use total_pressure
stored in gas_mixture
instead. Also you should be using p_vector.
|
protected |
Pointer to the PSD object.
|
protected |
Boolean flag to specify if a PSD is to be used to estimate saturation, permeability, etc.
|
protected |
Pointer to the PSD object.
|
protected |
PSD class type from input file.
|
protected |
Liquid water saturation at every quadrature point inside the cell [-].
|
protected |
Pointer used to store the solute gas for computing binary diffusion coefficients.
This value is used in child classes and is setup using compute_gas_diffusion
|
protected |
Pointer used to store the solute gas for computing binary diffusion coefficients.
This value is used in child classes and is setup using compute_gas_diffusion
|
protected |
Square root of user defined constant permeability, m.
|
protected |
Inverse of square root of user defined constant permeability, 1/m.
|
protected |
Temperature at every quadrature point inside the cell in [K].
|
protected |
Temperature [K
] used to compute gas diffusivity.
. Use temperature
stored in gas_mixture
instead. Also, you should be using T_vector
|
protected |
User defined constant tortuosity.
|
protected |
Variable defining if the tortuosity is constant.
|
protected |
Boolean flag that specifies if Knudsen effects should be accounted for.
This value is specified in Generic>>Use Bosanquet approx. Its default value is false.