OpenFCST: The open-source Fuel Cell Simulation Toolbox
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Private Types | Private Member Functions | Private Attributes | List of all members
FuelCell::OperatingConditions Class Reference

Class used to store, read from file and define the operating conditions for a fuel cell. More...

#include <operating_conditions.h>

Collaboration diagram for FuelCell::OperatingConditions:
Collaboration graph
[legend]

Public Member Functions

 OperatingConditions ()
 Constructor. More...
 
 ~OperatingConditions ()
 Destructor. More...
 
void declare_parameters (ParameterHandler &param) const
 Declare all necessary parameters in order to compute the coefficients. More...
 
void initialize (ParameterHandler &param)
 Class used to read in data and initialize the necessary data to compute the coefficients. More...
 
void set_gas_mixtures (FuelCellShop::Material::GasMixture &c_mix, FuelCellShop::Material::GasMixture &a_mix)
 Initialize the temperature and pressure for anode gas mixture and store a copy of GasMixture class inside operating conditions. More...
 
std::vector< double > get_rho_anode () const
 Compute the density of each species in the anode based on GasMixture specified in set_gas_mixture_anode(). More...
 
std::vector< double > get_rho_cathode () const
 Compute the density of each species in the anode based on GasMixture specified in set_gas_mixture_cathode(). More...
 
double saturation_pressure () const
 Get the water vapour saturation pressure in atmospheres (atm) using cell temperature. More...
 
double get_x_wv (int) const
 Get the mole fraction of water at the selected electrode channel from pressure, temperature, and relative humidity. More...
 
double get_x_wv () const
 Get the mole fraction of water at cathode channel from pressure, temperature, and relative humidity. More...
 
double get_x_wv_anode () const
 Get the mole fraction of water at anode channel from pressure, temperature, and relative humidity. More...
 
double get_x_o2 () const
 Get the mole fraction of oxygen at the cathode channel from pressure, temperature, and relative humidity. More...
 
double get_x_h2 () const
 Get the mole fraction of hydrogen at the anode channel from pressure, temperature, and relative humidity. More...
 
double get_x_n2 (int) const
 Get the mole fraction of nitrogen at the anode channel from pressure, temperature, and relative humidity. More...
 
double get_rho_wv (int) const
 Get the density of water vapour (wv) at the selected electrode channel from pressure, temperature, mole fraction, and total concentration. More...
 
double get_rho_wv () const
 Get the density of water vapour (wv) at the cathode channel from pressure, temperature, mole fraction, and total concentration. More...
 
double get_rho_wv_anode () const
 Get the density of water vapour (wv) at the anode channel from pressure, temperature, mole fraction, and total concentration. More...
 
double get_rho_o2 () const
 Get the density of oxygen (O2) at the cathode channel from pressure, temperature, mole fraction, and total concentration. More...
 
double get_rho_h2 () const
 Get the density of hydrogen (O2) at the anode channel from pressure, temperature, mole fraction, and total concentration. More...
 
double get_rho_n2 (int) const
 Get the density of nitrogen (N2) depending on which electrode we are solving. More...
 
double get_rho_n2 () const
 Get the density of nitrogen (N2) at the cathode channel from pressure, temperature, mole fraction, and total concentration. More...
 
double get_rho_n2_anode () const
 Get the density of nitrogen (N2) at the anode channel from pressure, temperature, mole fraction, and total concentration. More...
 
void set_gas_map (std::map< unsigned int, std::string > tmp)
 This function is meant to be called by an application to give Operating Conditions the map relating species number (key) to gas species (value). More...
 
std::map< unsigned int,
std::string > 
get_gas_map () const
 Get the species map relating species number (key) to species material (value). More...
 
double get_density_from_map (std::string species) const
 Submit a string of an unsigned integer value and it will determine, what the respective species material is from the gasSpeciesMap and return the density at Operating Conditions. More...
 
double voltage_cell_th ()
 NOTE: This function is redefined in base_kinetics class, considering variable temperature and gas pressures. More...
 
void print_operating_conditions () const
 Output operating conditions to screen. More...
 
void adjust_initial_solution (std::vector< component_materialID_value_map > &maps, const boost::shared_ptr< FuelCellShop::Geometry::GridBase< dim > > grid) const
 
void adjust_boundary_conditions (std::vector< component_boundaryID_value_map > &maps, const boost::shared_ptr< FuelCellShop::Geometry::GridBase< dim > > grid) const
 Routine used in order to adjust boundary conditions for O2 and cell voltage. More...
 
double get_c_c () const
 Get the total gas concentration in the cathode. More...
 
double get_T () const
 Return cell temperture as input in Operating Conditions subsection. More...
 
double get_V () const
 Return cell voltage as input in Operating Conditions subsection. More...
 
double get_dV_a () const
 Return the voltage drop in the anode. More...
 
double get_pc_Pa () const
 Return cathode pressure as input in Operating Conditions subsection. More...
 
double get_pc_atm () const
 Return cathode pressure as input in Operating Conditions subsection. More...
 
double get_c_a () const
 Get the total gas concentration in the anode. More...
 
double get_pa_Pa () const
 Return anode pressure as input in Operating Conditions subsection. More...
 
double get_pa_atm () const
 Return anode pressure as input in Operating Conditions subsection. More...
 
double get_RH_a () const
 Return anode relative humidity as input in Operating Conditions subsection. More...
 
double get_RH_c () const
 Return cathode relative humidity as input in Operating Conditions subsection. More...
 
double get_OCV () const
 Get the open circuit voltage for the cell. More...
 

Private Types

enum  Electrode { noElectrode, ANODE, CATHODE, MEMBRANE }
 Variable to switch between electrodes in molar fraction functions: More...
 

Private Member Functions

unsigned int get_species_index (std::string name) const
 Private member function used to extract the species index value from name. More...
 
double get_anode_gas_density (unsigned int species_index) const
 Private function used in adjust_initial_solution and adjust_boundary_conditions that estimates the density of each species in the anode compartment. More...
 
double get_cathode_gas_density (unsigned int species_index) const
 Private function used in adjust_initial_solution and adjust_boundary_conditions that estimates the density of each species in the cathode compartment of the fuel cell. More...
 

Private Attributes

bool adjust_BC
 Bool set to true if you want to modify boundary conditions. More...
 
double R
 
double channel_oxygen_mole_fraction
 Initial amount of oxygen in channel prior to humidification. More...
 
double channel_dry_hydrogen_mole_fraction
 Initial amount of hydrogen in channel prior to humidification. More...
 
double T_cell
 Operating temperature of the cell. More...
 
double V_cell
 Operating voltage of the cell. More...
 
double dV_a
 Voltage drop in the anode. More...
 
double OCV
 Open circuit voltage for the cell. More...
 
double E_th
 Theoretical voltage for the cell. More...
 
double p_a
 Pressure of the gas mixture in the anode B.C. More...
 
double c_a
 Concentration of the gas mixture in the anode B.C. More...
 
double RH_a
 Relative humidity of the gas mixture in the anode B.C. More...
 
double p_c
 Pressure of the gas mixture in the cathode B.C. More...
 
double c_c
 Concentration of the gas mixture in the cathode B.C. More...
 
double RH_c
 Relative humidity of the gas mixture in the anode B.C. More...
 
FuelCellShop::Material::WaterVapor water
 Get molar mass of water vapour. More...
 
const double M_water = water.get_molar_mass() * 1000
 
FuelCellShop::Material::Oxygen oxygen
 Get molar mass of oxygen. More...
 
const double M_oxygen = oxygen.get_molar_mass() * 1000
 
FuelCellShop::Material::Nitrogen nitrogen
 Get molar mass of nitrogen. More...
 
const double M_nitrogen = nitrogen.get_molar_mass() * 1000
 
FuelCellShop::Material::Hydrogen hydrogen
 Get molar mass of nitrogen. More...
 
const double M_hydrogen = hydrogen.get_molar_mass() * 1000
 
FuelCellShop::Material::GasMixtureanode_mix
 Container storing the gas mixture in the anode compartment of the cell. More...
 
FuelCellShop::Material::GasMixturecathode_mix
 Container storing the gas mixture in the anode compartment of the cell. More...
 
std::map< unsigned int,
std::string > 
gasSpeciesMap
 map relating species number (key) to species material (value) More...
 
bool anode = false
 Run an anode instead of a cathode: More...
 
Electrode compartment
 

Detailed Description

Class used to store, read from file and define the operating conditions for a fuel cell.

It stores the following:

This information, in conjunction with the water saturation equation is used to compute the mole fractions for each component in the mixture.

If the inlet is humidified air, please set the oxygen mole fraction to 0.21 (default).

In the input file, the following parameters can be specified (see declare_parameters ):

* subsection Fuel cell data
* (...)
* subsection Operating conditions
* set Adjust initial solution and boundary conditions = false
* set Temperature cell [K] = 353
* set Cathode pressure [Pa] = 101325
* set Cathode initial oxygen mole fraction (prior to humidification) = 0.21
* set Cathode relative humidity = 0.7
* set Anode pressure [Pa] = 101325
* set Anode relative humidity = 0.7
* set Voltage cell [V] = 0.6
* set Voltage drop in the anode [V] = 0.015
* set Open circuit voltage [V] = 1.23
* end
* end
* subsection Application
* subsection Electrode KG
* set Anode simulation = false
* end
* end
*

Usage Details:

In order to use this class, first an object of the class needs to be created. Usually, one such objects exists in every application. To create the object, include the .h file in the include application file and in the application data member section add the object. For example:

* #include "operating_conditions.h"
*
* // Then in the data member declaration (usually a private member)
*

Once the object is created, the section where the input data will be specified in the input file needs to be declared. To do so, in the declare_parameters section of your application call the following:

* //--------- IN DECLARE_PARAMETERS ------------------------------------------------------
* template <int dim>
* void
* NAME::AppCathode<dim>::declare_parameters(ParameterHandler& param)
* {
* (...)
* OC.declare_parameters(param);
* (...)
* }
*

Finally, once the input file has been read by our application, your class needs to be initialized. This is achieved using the function initialize()

* //--------- IN INITIALIZE ------------------------------------------------------
* template <int dim>
* void
* NAME::AppCathode<dim>::_initialize(ParameterHandler& param)
* {
* (...)
* OC.initialize(param);
* }
*

You are now ready to use your OperatingConditions object!.

Author
M. Secanell, 2009-2013

Member Enumeration Documentation

Variable to switch between electrodes in molar fraction functions:

Enumerator
noElectrode 
ANODE 
CATHODE 
MEMBRANE 

Constructor & Destructor Documentation

FuelCell::OperatingConditions::OperatingConditions ( )

Constructor.

FuelCell::OperatingConditions::~OperatingConditions ( )

Destructor.

Member Function Documentation

void FuelCell::OperatingConditions::adjust_boundary_conditions ( std::vector< component_boundaryID_value_map > &  maps,
const boost::shared_ptr< FuelCellShop::Geometry::GridBase< dim > >  grid 
) const

Routine used in order to adjust boundary conditions for O2 and cell voltage.

This routine should be called after the component_boundaryID_value_map has been initialized.

In order for the application to work well, the boundary IDs for the channel/GDL and land/GDL interfaces need to be properly identified in the grid section of the input file, i.e.

* subsection Grid generation
* subsection Internal mesh generator parameters
* subsection Boundary ID
* set c_Ch/GDL = 100
* set c_BPP/GDL = 200
* set a_Ch/GDL = 100
* set a_BPP/GDL = 200
* end
* end
* end
*
void FuelCell::OperatingConditions::adjust_initial_solution ( std::vector< component_materialID_value_map > &  maps,
const boost::shared_ptr< FuelCellShop::Geometry::GridBase< dim > >  grid 
) const
void FuelCell::OperatingConditions::declare_parameters ( ParameterHandler &  param) const

Declare all necessary parameters in order to compute the coefficients.

The parameters that can be specified in the input file are as follows:

* subsection Fuel cell data
* (...)
* subsection Operating conditions
* set Adjust initial solution and boundary conditions = false
* set Temperature cell = 353
* set Cathode pressure = 101325
* set Cathode initial oxygen mole fraction (prior to humidification) = 0.21
* set Cathode relative humidity = 0.7
* set Anode pressure = 101325
* set Anode relative humidity = 0.7
* set Voltage cell = 0.6
* set Voltage drop in the anode = 0.015
* set Open circuit voltage = 1.23
* end
* end
*
double FuelCell::OperatingConditions::get_anode_gas_density ( unsigned int  species_index) const
inlineprivate

Private function used in adjust_initial_solution and adjust_boundary_conditions that estimates the density of each species in the anode compartment.

The parameter species_index is the index obtained using get_species_index. Note that for the anode, the index is obtained as get_species_index - n_gases in cathode mix.

References Units::ATM_to_PA, Units::convert(), and Constants::R().

Here is the call graph for this function:

double FuelCell::OperatingConditions::get_c_a ( ) const
inline

Get the total gas concentration in the anode.

double FuelCell::OperatingConditions::get_c_c ( ) const
inline

Get the total gas concentration in the cathode.

NOTE: This is a constant value. Use only if the total gas concentration is assumed to be constant

double FuelCell::OperatingConditions::get_cathode_gas_density ( unsigned int  species_index) const
inlineprivate

Private function used in adjust_initial_solution and adjust_boundary_conditions that estimates the density of each species in the cathode compartment of the fuel cell.

The parameter species_index is the index obtained using get_species_index. Note that for the anode, the index is obtained as get_species_index in cathode mix.

References Units::ATM_to_PA, Units::convert(), and Constants::R().

Here is the call graph for this function:

double FuelCell::OperatingConditions::get_density_from_map ( std::string  species) const

Submit a string of an unsigned integer value and it will determine, what the respective species material is from the gasSpeciesMap and return the density at Operating Conditions.

Currently, it can only return density for oxygen, water, and nitrogen.

double FuelCell::OperatingConditions::get_dV_a ( ) const
inline

Return the voltage drop in the anode.

Note
Can be used as boundary condition for anode model or initial condition in a full MEA model
std::map<unsigned int, std::string> FuelCell::OperatingConditions::get_gas_map ( ) const

Get the species map relating species number (key) to species material (value).

double FuelCell::OperatingConditions::get_OCV ( ) const
inline

Get the open circuit voltage for the cell.

double FuelCell::OperatingConditions::get_pa_atm ( ) const
inline

Return anode pressure as input in Operating Conditions subsection.

References Units::ATM_to_PA, and Units::convert().

Here is the call graph for this function:

double FuelCell::OperatingConditions::get_pa_Pa ( ) const
inline

Return anode pressure as input in Operating Conditions subsection.

double FuelCell::OperatingConditions::get_pc_atm ( ) const
inline

Return cathode pressure as input in Operating Conditions subsection.

References Units::ATM_to_PA, and Units::convert().

Here is the call graph for this function:

double FuelCell::OperatingConditions::get_pc_Pa ( ) const
inline

Return cathode pressure as input in Operating Conditions subsection.

double FuelCell::OperatingConditions::get_RH_a ( ) const
inline

Return anode relative humidity as input in Operating Conditions subsection.

double FuelCell::OperatingConditions::get_RH_c ( ) const
inline

Return cathode relative humidity as input in Operating Conditions subsection.

std::vector<double> FuelCell::OperatingConditions::get_rho_anode ( ) const

Compute the density of each species in the anode based on GasMixture specified in set_gas_mixture_anode().

std::vector<double> FuelCell::OperatingConditions::get_rho_cathode ( ) const

Compute the density of each species in the anode based on GasMixture specified in set_gas_mixture_cathode().

double FuelCell::OperatingConditions::get_rho_h2 ( ) const

Get the density of hydrogen (O2) at the anode channel from pressure, temperature, mole fraction, and total concentration.

double FuelCell::OperatingConditions::get_rho_n2 ( int  ) const

Get the density of nitrogen (N2) depending on which electrode we are solving.

double FuelCell::OperatingConditions::get_rho_n2 ( ) const

Get the density of nitrogen (N2) at the cathode channel from pressure, temperature, mole fraction, and total concentration.

double FuelCell::OperatingConditions::get_rho_n2_anode ( ) const

Get the density of nitrogen (N2) at the anode channel from pressure, temperature, mole fraction, and total concentration.

double FuelCell::OperatingConditions::get_rho_o2 ( ) const

Get the density of oxygen (O2) at the cathode channel from pressure, temperature, mole fraction, and total concentration.

double FuelCell::OperatingConditions::get_rho_wv ( int  ) const

Get the density of water vapour (wv) at the selected electrode channel from pressure, temperature, mole fraction, and total concentration.

double FuelCell::OperatingConditions::get_rho_wv ( ) const

Get the density of water vapour (wv) at the cathode channel from pressure, temperature, mole fraction, and total concentration.

double FuelCell::OperatingConditions::get_rho_wv_anode ( ) const

Get the density of water vapour (wv) at the anode channel from pressure, temperature, mole fraction, and total concentration.

unsigned int FuelCell::OperatingConditions::get_species_index ( std::string  name) const
inlineprivate

Private member function used to extract the species index value from name.

Note that density_species_1 will have index 0.

double FuelCell::OperatingConditions::get_T ( ) const
inline

Return cell temperture as input in Operating Conditions subsection.

double FuelCell::OperatingConditions::get_V ( ) const
inline

Return cell voltage as input in Operating Conditions subsection.

double FuelCell::OperatingConditions::get_x_h2 ( ) const

Get the mole fraction of hydrogen at the anode channel from pressure, temperature, and relative humidity.

double FuelCell::OperatingConditions::get_x_n2 ( int  ) const

Get the mole fraction of nitrogen at the anode channel from pressure, temperature, and relative humidity.

double FuelCell::OperatingConditions::get_x_o2 ( ) const

Get the mole fraction of oxygen at the cathode channel from pressure, temperature, and relative humidity.

double FuelCell::OperatingConditions::get_x_wv ( int  ) const

Get the mole fraction of water at the selected electrode channel from pressure, temperature, and relative humidity.

double FuelCell::OperatingConditions::get_x_wv ( ) const

Get the mole fraction of water at cathode channel from pressure, temperature, and relative humidity.

double FuelCell::OperatingConditions::get_x_wv_anode ( ) const

Get the mole fraction of water at anode channel from pressure, temperature, and relative humidity.

void FuelCell::OperatingConditions::initialize ( ParameterHandler &  param)

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

void FuelCell::OperatingConditions::print_operating_conditions ( ) const

Output operating conditions to screen.

double FuelCell::OperatingConditions::saturation_pressure ( ) const

Get the water vapour saturation pressure in atmospheres (atm) using cell temperature.

void FuelCell::OperatingConditions::set_gas_map ( std::map< unsigned int, std::string >  tmp)

This function is meant to be called by an application to give Operating Conditions the map relating species number (key) to gas species (value).

The map is set to the private map gasSpeciesMap. At this time this is only needed by appCathodeKG as it uses density for fluid flow

void FuelCell::OperatingConditions::set_gas_mixtures ( FuelCellShop::Material::GasMixture c_mix,
FuelCellShop::Material::GasMixture a_mix 
)
inline

Initialize the temperature and pressure for anode gas mixture and store a copy of GasMixture class inside operating conditions.

Note: The gas mixture species should already be specified.

References FuelCellShop::Material::GasMixture::set_temperature(), and FuelCellShop::Material::GasMixture::set_total_pressure().

Here is the call graph for this function:

double FuelCell::OperatingConditions::voltage_cell_th ( )

NOTE: This function is redefined in base_kinetics class, considering variable temperature and gas pressures.

The function here should only be used for defining Initial condition values. Get the theoretical cell voltage using the cell temperature, and the reactant gas pressures

Member Data Documentation

bool FuelCell::OperatingConditions::adjust_BC
private

Bool set to true if you want to modify boundary conditions.

bool FuelCell::OperatingConditions::anode = false
private

Run an anode instead of a cathode:

FuelCellShop::Material::GasMixture* FuelCell::OperatingConditions::anode_mix
private

Container storing the gas mixture in the anode compartment of the cell.

double FuelCell::OperatingConditions::c_a
private

Concentration of the gas mixture in the anode B.C.

double FuelCell::OperatingConditions::c_c
private

Concentration of the gas mixture in the cathode B.C.

FuelCellShop::Material::GasMixture* FuelCell::OperatingConditions::cathode_mix
private

Container storing the gas mixture in the anode compartment of the cell.

double FuelCell::OperatingConditions::channel_dry_hydrogen_mole_fraction
private

Initial amount of hydrogen in channel prior to humidification.

double FuelCell::OperatingConditions::channel_oxygen_mole_fraction
private

Initial amount of oxygen in channel prior to humidification.

Electrode FuelCell::OperatingConditions::compartment
private
double FuelCell::OperatingConditions::dV_a
private

Voltage drop in the anode.

double FuelCell::OperatingConditions::E_th
private

Theoretical voltage for the cell.

std::map<unsigned int, std::string> FuelCell::OperatingConditions::gasSpeciesMap
private

map relating species number (key) to species material (value)

FuelCellShop::Material::Hydrogen FuelCell::OperatingConditions::hydrogen
private

Get molar mass of nitrogen.

const double FuelCell::OperatingConditions::M_hydrogen = hydrogen.get_molar_mass() * 1000
private
const double FuelCell::OperatingConditions::M_nitrogen = nitrogen.get_molar_mass() * 1000
private
const double FuelCell::OperatingConditions::M_oxygen = oxygen.get_molar_mass() * 1000
private
const double FuelCell::OperatingConditions::M_water = water.get_molar_mass() * 1000
private
FuelCellShop::Material::Nitrogen FuelCell::OperatingConditions::nitrogen
private

Get molar mass of nitrogen.

double FuelCell::OperatingConditions::OCV
private

Open circuit voltage for the cell.

FuelCellShop::Material::Oxygen FuelCell::OperatingConditions::oxygen
private

Get molar mass of oxygen.

double FuelCell::OperatingConditions::p_a
private

Pressure of the gas mixture in the anode B.C.

double FuelCell::OperatingConditions::p_c
private

Pressure of the gas mixture in the cathode B.C.

double FuelCell::OperatingConditions::R
private
double FuelCell::OperatingConditions::RH_a
private

Relative humidity of the gas mixture in the anode B.C.

double FuelCell::OperatingConditions::RH_c
private

Relative humidity of the gas mixture in the anode B.C.

double FuelCell::OperatingConditions::T_cell
private

Operating temperature of the cell.

double FuelCell::OperatingConditions::V_cell
private

Operating voltage of the cell.

FuelCellShop::Material::WaterVapor FuelCell::OperatingConditions::water
private

Get molar mass of water vapour.


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