OpenFCST: The open-source Fuel Cell Simulation Toolbox
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
List of all members
FuelCellShop::MicroScale::DualPSD< dim > Class Template Reference

Dual Pore Size Distribution. More...

#include <PSD_dual.h>

Inheritance diagram for FuelCellShop::MicroScale::DualPSD< dim >:
Inheritance graph
[legend]
Collaboration diagram for FuelCellShop::MicroScale::DualPSD< dim >:
Collaboration graph
[legend]

Public Member Functions

Constructors, declarations and Initalization
 DualPSD ()
 Consturctor. More...
 
 DualPSD (std::string name)
 Constructor. More...
 
virtual ~DualPSD ()
 Destructor. More...
 
void declare_parameters (ParameterHandler &param) const
 Declare parameters for a parameter file. More...
 
void initialize (ParameterHandler &param)
 Member function used to read in data and initialize the necessary data to compute the coefficients. 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_capillary_pressure (const SolutionVariable &C_in)
 Member function used to set the capillary pressure [psi] at every quadrature point inside the cell. More...
 
void set_critical_radius ()
 Member function used to set the critical radius [nm] at every quadrature point inside the cell. More...
 
void set_saturation ()
 Member function used to set the saturation at every quadrature point inside the cell. More...
 
Accessors and info
const std::type_info & get_base_type () const
 This member function returns a type_info object with the name of the base layer type the inherited class belongs to, i.e. More...
 
virtual void get_saturation (std::vector< double > &S) const
 This function is used to compute the saturation by using PSD The saturation iof a porous medium can be obtained by integrating the change of cumulative pore volume fraction as a function of effective pore radius over the whole pore size domain. More...
 
virtual void get_derivative_saturation (std::vector< double > &S) const
 
virtual void get_global_saturated_permeability (double &saturated_permeability) const
 This function is used to compute the saturated_permeability by using PSD. More...
 
virtual void get_relative_liquid_permeability (std::vector< double > &liquid_permeability) const
 This function is used to compute the liquid_permeability by using PSD. More...
 
virtual void get_derivative_relative_liquid_permeability (std::vector< double > &liquid_permeability) const
 This function is used to compute the liquid_permeability by using PSD. More...
 
virtual void get_relative_gas_permeability (std::vector< double > &gas_permeability) const
 This function is used to compute the gas_permeability by using PSD. More...
 
virtual void get_liquid_gas_interfacial_surface (std::vector< double > &liquid_gas_interfacial_surface) const
 This function is used to compute the liquid_gas_interfacial_surface by using PSD. More...
 
virtual void get_derivative_liquid_gas_interfacial_surface (std::vector< double > &) const
 This function is used to compute the liquid_gas_interfacial_surface by using PSD. More...
 
virtual void get_wetted_wall_surface_area (std::vector< double > &wetted_wall_surface_area) const
 This function is used to compute the pore_wetted_wall by using PSD. More...
 
virtual void get_knudsen_radius (std::vector< double > &knudsen_radius) const
 This function is used to compute the knudsen_radius by using PSD. More...
 
virtual void get_diffusivity () const
 This function is used to compute the diffusivity by using PSD. More...
 
virtual void get_maximum_cross_sectional_areas (double &) const
 
virtual void get_PSD_plot (const std::vector< double >, std::vector< double > &, std::vector< double > &, std::vector< double > &) const
 This function is used to plot PSD configuration. More...
 
- Public Member Functions inherited from FuelCellShop::MicroScale::BasePSD< dim >
virtual ~BasePSD ()
 Destructor. More...
 
void set_porosity (double porosity)
 
double get_porosity () const
 
void set_derivative_flags (const std::vector< VariableNames > &flags)
 Set the names of FCST solution variables with respect to which you would like to compute the derivatives of material properties. 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 psd depend on the solution, the solution for a given cell should be passed to the class using this member function. More...
 
virtual void get_PSD_plot ()
 This function is used to create PSD configuration plot by outputing all the numbers. More...
 
const std::string & name_psd () const
 Return the name of the PSD. More...
 
virtual void print_psd_properties () const
 This function prints out the psd properties. More...
 

Static Public Attributes

concrete_name
static const std::string concrete_name
 Concrete name used for objects of this class. More...
 

Protected Member Functions

Instance Delivery (Replica creator)
virtual boost::shared_ptr
< FuelCellShop::MicroScale::BasePSD
< dim > > 
create_replica (const std::string &psd_section_name)
 This member function is used to create an object of type psd. More...
 
- Protected Member Functions inherited from FuelCellShop::MicroScale::BasePSD< dim >
 BasePSD ()
 Constructor. More...
 
 BasePSD (const std::string &name)
 Constructor. More...
 
void _initialize (ParameterHandler &param)
 Member function used to read in data and initialize the necessary data to compute the coefficients. More...
 

Protected Attributes

PSD properties
HIPSD< dimpsd_hi
 Creating a pointer of PSD hydrophilic class. More...
 
HOPSD< dimpsd_ho
 Creating a pointer of PSD hydrophobic class. More...
 
double pressure_c
 
SolutionVariable T_vector
 Temperature at every quadrature point inside the cell. More...
 
SolutionVariable Capillary_pressure_vector
 Capillary pressure at every quadrature point inside the cell. More...
 
- Protected Attributes inherited from FuelCellShop::MicroScale::BasePSD< dim >
const std::string name
 Name of the psd. More...
 
std::vector< VariableNamesderivative_flags
 Flags for derivatives: These flags are used to request derivatives of material properties. More...
 
std::map< VariableNames, double > constant_solutions
 Map storing values of solution variables constant in a particular application. More...
 
double gamma
 The gamma is the Water-air interface surface tension. More...
 
double contact_angle
 The contact_angle is the contact angle between air and water. More...
 
double lamda
 The lamda is measured from the mercury intrusion experiment. More...
 
double F_HI
 The F_HI is the fraction of the total volume corresponding to the hydrophilic pores. More...
 
double F_HO
 The F_HO is the fraction of the total volume corresponding to the hydrophobic pores. More...
 
std::vector< double > f_k
 The f_k is the contribution of the log-normal distribution k to the total PSD. More...
 
std::vector< double > r_k
 The r_k is the characteristic pore size of the distribution k. More...
 
std::vector< double > s_k
 The s_k is the spread of the distribution k. More...
 
double por
 The por is the porosity inherited from the layer class. More...
 

Static Protected Attributes

Instance Delivery (Prototype)
static DualPSD< dim > const * PROTOTYPE
 PROTOTYPE is the pointer is the dynamic pointer pointing to the DualPSD class itself. More...
 

Additional Inherited Members

- Static Public Member Functions inherited from FuelCellShop::MicroScale::BasePSD< dim >
static void declare_PSD_parameters (ParameterHandler &param)
 Function used to declare all the data necessary in the parameter files for all BasePSD children. More...
 
static boost::shared_ptr
< FuelCellShop::MicroScale::BasePSD
< dim > > 
create_PSD (const std::string &psd_section_name, ParameterHandler &param)
 Function used to select the appropriate CatalystLayer type as specified in the ParameterHandler under line. More...
 
- Protected Types inherited from FuelCellShop::MicroScale::BasePSD< dim >
typedef std::map< std::string,
BasePSD< dim > * > 
_mapFactory
 This object is used to store all objects of type psd. More...
 
- Static Protected Member Functions inherited from FuelCellShop::MicroScale::BasePSD< dim >
static _mapFactoryget_mapFactory ()
 Return the map library that stores all childrens of this class. More...
 

Detailed Description

template<int dim>
class FuelCellShop::MicroScale::DualPSD< dim >

Dual Pore Size Distribution.

Based on the results of the mecury intrusion experiment, this class calculates the coefficients such as relative liquid permeability, permeability and knudsen radius...

This class is the child class of the all the BasePSD class which reimplements all the functions that are decleared in the BasePSD class.

Input parameters

The input parameters for this class are:Mode probability global, Mode characteristic radius global, Mode width global, Volume fraction Hydrophobic, Volume fraction Hydrophilic... All data below is from the Hydrophilic unit test

* subsection PSD parameters
* subsection BasePSD
* set psd type = DualPSD
* set Mode probability global = 0.72, 0.28
* set Mode characteristic radius global = 34.0, 14.2
* set Mode width global = 0.35, 1.0
* set Gamma = 0.24
* set contact_angle = 1.396
* set lamda = 1.0
* set probability P_b = 1
* subsection HIPSD
* set Hydrophilic Mode probability global = 0.72, 0.28
* set Hydrophilic Mode characteristic radius global = 34.0, 14.2
* set Hydrophilic Mode width global = 0.35, 1.0
* end
* subsection HOPSD
* set Hydrophobic Mode probability global = 0.72, 0.28
* set Hydrophobic Mode characteristic radius global = 34.0, 14.2
* set Hydrophobic Mode width global = 0.35, 1.0
* end
* end
* end
*

Usage details

If you want to use the object of the DualPSD, use the following code, but you should not create an object of the DualPSD out of the PSD scope. If you want to use the PSD in the layer level, use the code in BasePSD.

* //Create an object of TemplateClass
*

References

[1] Pedro Abdiel Mateo Villanueva, A MIXED WETTABILITY PORE SIZE DISTRIBUTION MODEL FOR THE ANALYSIS OF WATER TRANSPORT IN PEMFC MATERIALS, M. Sc. thesis, University of Alberta, 2013

Author
J. Zhou

Marc Secanell

Date
2014

Constructor & Destructor Documentation

template<int dim>
FuelCellShop::MicroScale::DualPSD< dim >::DualPSD ( )

Consturctor.

template<int dim>
FuelCellShop::MicroScale::DualPSD< dim >::DualPSD ( std::string  name)

Constructor.

template<int dim>
virtual FuelCellShop::MicroScale::DualPSD< dim >::~DualPSD ( )
virtual

Destructor.

Member Function Documentation

template<int dim>
virtual boost::shared_ptr<FuelCellShop::MicroScale::BasePSD<dim> > FuelCellShop::MicroScale::DualPSD< dim >::create_replica ( const std::string &  psd_section_name)
inlineprotectedvirtual

This member function is used to create an object of type psd.

Warning
This class MUST be redeclared in every child.

Reimplemented from FuelCellShop::MicroScale::BasePSD< dim >.

template<int dim>
void FuelCellShop::MicroScale::DualPSD< dim >::declare_parameters ( ParameterHandler &  param) const
virtual

Declare parameters for a parameter file.

Reimplemented from FuelCellShop::MicroScale::BasePSD< dim >.

template<int dim>
const std::type_info& FuelCellShop::MicroScale::DualPSD< dim >::get_base_type ( ) const
inlinevirtual

This member function returns a type_info object with the name of the base layer type the inherited class belongs to, i.e.

HIPSD HOPSD DualPSD

Note that this is necessary if we want to find out not the name of the actual class which can be obtain using

const std::type_info& name = typeid(*this)

but the name of the parent class.

Note
Do not re-implement this class in children classes

Reimplemented from FuelCellShop::MicroScale::BasePSD< dim >.

template<int dim>
virtual void FuelCellShop::MicroScale::DualPSD< dim >::get_derivative_liquid_gas_interfacial_surface ( std::vector< double > &  ) const
virtual

This function is used to compute the liquid_gas_interfacial_surface by using PSD.

$ a_{LV} = \left[ \frac{a(r)_{HI}}{V_{T}} + \frac{a(r)_{HO}}{V_{T}} \right]\frac{1}{\lambda_{2}} $

Implements FuelCellShop::MicroScale::BasePSD< dim >.

template<int dim>
virtual void FuelCellShop::MicroScale::DualPSD< dim >::get_derivative_relative_liquid_permeability ( std::vector< double > &  ) const
virtual

This function is used to compute the liquid_permeability by using PSD.

Warning
This function is a virtual function and it must be reimplemented in the child class

Implements FuelCellShop::MicroScale::BasePSD< dim >.

template<int dim>
virtual void FuelCellShop::MicroScale::DualPSD< dim >::get_derivative_saturation ( std::vector< double > &  S) const
virtual
template<int dim>
virtual void FuelCellShop::MicroScale::DualPSD< dim >::get_diffusivity ( ) const
virtual

This function is used to compute the diffusivity by using PSD.

template<int dim>
virtual void FuelCellShop::MicroScale::DualPSD< dim >::get_global_saturated_permeability ( double &  saturated_permeability) const
virtual

This function is used to compute the saturated_permeability by using PSD.

$ k_{sat} = \frac{1}{8}\left[\frac{\varepsilon_{o}}{\lambda}\right]^{2} \sum_{k} \,\exp{(-2\,s_{k}^{2})}\,r_{k}^{2}f_{k} $

Implements FuelCellShop::MicroScale::BasePSD< dim >.

template<int dim>
virtual void FuelCellShop::MicroScale::DualPSD< dim >::get_knudsen_radius ( std::vector< double > &  knudsen_radius) const
virtual

This function is used to compute the knudsen_radius by using PSD.

$ r_{Kn} = \frac{C_{1} + C_{2}}{C_{3} + C_{4}} $

Implements FuelCellShop::MicroScale::BasePSD< dim >.

template<int dim>
virtual void FuelCellShop::MicroScale::DualPSD< dim >::get_liquid_gas_interfacial_surface ( std::vector< double > &  liquid_gas_interfacial_surface) const
virtual

This function is used to compute the liquid_gas_interfacial_surface by using PSD.

$ a_{LV} = \left[ \frac{a(r)_{HI}}{V_{T}} + \frac{a(r)_{HO}}{V_{T}} \right]\frac{1}{\lambda_{2}} $

Implements FuelCellShop::MicroScale::BasePSD< dim >.

template<int dim>
virtual void FuelCellShop::MicroScale::DualPSD< dim >::get_maximum_cross_sectional_areas ( double &  ) const
virtual
template<int dim>
virtual void FuelCellShop::MicroScale::DualPSD< dim >::get_PSD_plot ( const std::vector< double >  ,
std::vector< double > &  ,
std::vector< double > &  ,
std::vector< double > &   
) const
virtual

This function is used to plot PSD configuration.

template<int dim>
virtual void FuelCellShop::MicroScale::DualPSD< dim >::get_relative_gas_permeability ( std::vector< double > &  gas_permeability) const
virtual

This function is used to compute the gas_permeability by using PSD.

$ k_{r,G} = \frac{k_{G,HI} + k_{G,HO}}{k_{sat}} $

Implements FuelCellShop::MicroScale::BasePSD< dim >.

template<int dim>
virtual void FuelCellShop::MicroScale::DualPSD< dim >::get_relative_liquid_permeability ( std::vector< double > &  liquid_permeability) const
virtual

This function is used to compute the liquid_permeability by using PSD.

$ k_{r,L} = \frac{k_{L,HI} + k_{L,HO}}{k_{sat}} $

Implements FuelCellShop::MicroScale::BasePSD< dim >.

template<int dim>
virtual void FuelCellShop::MicroScale::DualPSD< dim >::get_saturation ( std::vector< double > &  S) const
virtual

This function is used to compute the saturation by using PSD The saturation iof a porous medium can be obtained by integrating the change of cumulative pore volume fraction as a function of effective pore radius over the whole pore size domain.

$ S_{total} = S_{HI} + S_{HO} $

Implements FuelCellShop::MicroScale::BasePSD< dim >.

template<int dim>
virtual void FuelCellShop::MicroScale::DualPSD< dim >::get_wetted_wall_surface_area ( std::vector< double > &  wetted_wall_surface_area) const
virtual

This function is used to compute the pore_wetted_wall by using PSD.

$ a_{wall} = a_{wall,HI} + a_{wall,HO} $

Implements FuelCellShop::MicroScale::BasePSD< dim >.

template<int dim>
void FuelCellShop::MicroScale::DualPSD< dim >::initialize ( ParameterHandler &  param)
virtual

Member function used to read in data and initialize the necessary data to compute the coefficients.

Reimplemented from FuelCellShop::MicroScale::BasePSD< dim >.

template<int dim>
void FuelCellShop::MicroScale::DualPSD< dim >::set_capillary_pressure ( const SolutionVariable C_in)
inlinevirtual

Member function used to set the capillary pressure [psi] at every quadrature point inside the cell.

This function should particulary be used in the case of non-isothermal application.

Reimplemented from FuelCellShop::MicroScale::BasePSD< dim >.

References capillary_pressure, and FuelCellShop::SolutionVariable::get_variablename().

Here is the call graph for this function:

template<int dim>
void FuelCellShop::MicroScale::DualPSD< dim >::set_critical_radius ( )
inlinevirtual

Member function used to set the critical radius [nm] at every quadrature point inside the cell.

This function should particulary be used in the case of non-isothermal application. It needs to be implemented at each iteration to reduce the computational time.

Reimplemented from FuelCellShop::MicroScale::BasePSD< dim >.

template<int dim>
void FuelCellShop::MicroScale::DualPSD< dim >::set_saturation ( )
inlinevirtual

Member function used to set the saturation at every quadrature point inside the cell.

This function should particulary be used in the case of non-isothermal application. It needs to be implemented at each iteration to reduce the computational time.

Reimplemented from FuelCellShop::MicroScale::BasePSD< dim >.

template<int dim>
void FuelCellShop::MicroScale::DualPSD< dim >::set_temperature ( const SolutionVariable T_in)
inlinevirtual

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.

Reimplemented from FuelCellShop::MicroScale::BasePSD< dim >.

References FuelCellShop::SolutionVariable::get_variablename(), and temperature_of_REV.

Here is the call graph for this function:

Member Data Documentation

template<int dim>
SolutionVariable FuelCellShop::MicroScale::DualPSD< dim >::Capillary_pressure_vector
protected

Capillary pressure at every quadrature point inside the cell.

template<int dim>
const std::string FuelCellShop::MicroScale::DualPSD< dim >::concrete_name
static

Concrete name used for objects of this class.

This name is used when setting up the subsection where the data is stored in the input file.

The data will be store under

* subsection name_specified_in_constructor
* set psd type = DualPSD # <-here I select the type of object of type psd
* subsection DualPSD # <- this is the concrete_name for this class
* set all info relevant to this object
* end
* end
*
template<int dim>
double FuelCellShop::MicroScale::DualPSD< dim >::pressure_c
protected
template<int dim>
DualPSD<dim> const* FuelCellShop::MicroScale::DualPSD< dim >::PROTOTYPE
staticprotected

PROTOTYPE is the pointer is the dynamic pointer pointing to the DualPSD class itself.

template<int dim>
HIPSD<dim> FuelCellShop::MicroScale::DualPSD< dim >::psd_hi
protected

Creating a pointer of PSD hydrophilic class.

template<int dim>
HOPSD<dim> FuelCellShop::MicroScale::DualPSD< dim >::psd_ho
protected

Creating a pointer of PSD hydrophobic class.

template<int dim>
SolutionVariable FuelCellShop::MicroScale::DualPSD< dim >::T_vector
protected

Temperature at every quadrature point inside the cell.


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