OpenFCST: The open-source Fuel Cell Simulation Toolbox
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
tafel_kinetics.h
Go to the documentation of this file.
1 //---------------------------------------------------------------------------
2 //
3 // FCST: Fuel Cell Simulation Toolbox
4 //
5 // Copyright (C) 2006-13 by Energy Systems Design Laboratory, University of Alberta
6 //
7 // This software is distributed under the MIT License.
8 // For more information, see the README file in /doc/LICENSE
9 //
10 // - Class: tafel_kinetics.h
11 // - Description: Header file for Tafel Kinetics model class.
12 // - Developers: M. Moore, M. Bhaiya, M. Secanell and V. Zingan
13 //
14 //---------------------------------------------------------------------------
15 
16 #ifndef _FUELCELLSHOP__TAFEL_KINETICS_H
17 #define _FUELCELLSHOP__TAFEL_KINETICS_H
18 
19 // Include openFCST routines:
21 #include <utils/fcst_constants.h>
22 
23 // Include deal.II classes
24 #include <deal.II/base/parameter_handler.h>
25 #include <deal.II/base/point.h>
26 #include <deal.II/base/function.h>
27 #include <deal.II/lac/vector.h>
28 #include <deal.II/fe/fe_values.h>
29 
30 //Include STL
31 #include <cmath>
32 #include <iostream>
33 
34 using namespace dealii;
35 
36 namespace FuelCellShop
37 {
38  namespace Kinetics
39  {
67  :
68  public BaseKinetics
69  {
70  public:
71 
73 
74 
77  TafelKinetics();
78 
83  TafelKinetics(const bool);
84 
88  ~TafelKinetics();
89 
94  virtual void declare_parameters(ParameterHandler& param) const;
95 
101  virtual void initialize(ParameterHandler& param);
102 
104 
106 
107 
118  static const std::string concrete_name;
120 
122 
123 
126  virtual void current_density (std::vector<double> &);
127 
134  virtual void derivative_current (std::map< VariableNames, std::vector<double> > &);
135 
137 
138  protected:
140 
141 
146  virtual boost::shared_ptr<FuelCellShop::Kinetics::BaseKinetics > create_replica ()
147  {
148  return boost::shared_ptr<FuelCellShop::Kinetics::BaseKinetics > (new FuelCellShop::Kinetics::TafelKinetics ( ));
149  }
153  static TafelKinetics const* PROTOTYPE;
155 
157 
158 
162  virtual void init_kin_param()
163  {
164  Assert( !kin_param_initialized, ExcInternalError() );
165  Assert( catalyst != NULL, ExcMessage("Catalyst object not initialized in the TafelKinetics object.") );
166  Assert( phi_m.is_initialized() && phi_s.is_initialized() && T.is_initialized(), ExcMessage("Either phi_m/phi_s/T is not set in the TafelKinetics object.") );
167  Assert( reactants_map.size()>0, ExcMessage("At least one reactant should be set in the TafelKinetics object, using set_reactant_concentrations method.") );
168 
169  std::vector<VariableNames> names;
170  for ( std::map< VariableNames, SolutionVariable >::const_iterator iter=reactants_map.begin(); iter!=reactants_map.end(); ++iter )
171  names.push_back(iter->first);
172 
173  catalyst->reference_concentration(names, ref_conc);
174  catalyst->reaction_order(names, gamma);
175  catalyst->alpha_cathodic(alpha_c);
176 
177  kin_param_initialized = true;
178  }
184  inline double negative_concentration_correction(double concentration, double gamma, bool derivative = true) const
185  {
186  double correction(1.0);
187 
188  if (concentration < 0.0)
189  {
190  if ( !derivative )
191  correction = -1.0;
192  else if (derivative && (gamma == 1.0 || gamma == 3.0 || gamma == 5.0) )
193  correction = 1.0; // <-- Note: 1.0 works well with gamma = 1.0 since (gamma-1) is zero, i.e. sign of derivative independent of concentration.
194  else if (derivative && ( (gamma != 1.0) || (gamma != 3.0) || (gamma != 5.0) ) )
195  correction = -1.0;
196  }
197 
198  return correction;
199  }
203  double alpha_c;
204 
208  std::map< VariableNames, double > ref_conc;
209 
213  std::map< VariableNames, double > gamma;
215 
216  };
217 
218  } //Kinetics
219 
220 } //FuelCellShop
221 
222 #endif //_FUELCELLSHOP__TAFEL_KINETICS_H
std::map< VariableNames, double > gamma
Map of reaction orders.
Definition: tafel_kinetics.h:213
virtual boost::shared_ptr< FuelCellShop::Kinetics::BaseKinetics > create_replica()
This member function is used to create an object of type gas diffusion layer.
Definition: tafel_kinetics.h:146
VariableNames
The enumeration containing the names of some of the available FCST solution variables and their deriv...
Definition: system_management.h:63
std::map< VariableNames, double > ref_conc
Map of reference concentrations.
Definition: tafel_kinetics.h:208
This class defines a simple Tafel kinetic model.
Definition: tafel_kinetics.h:66
virtual void init_kin_param()
Method used to initialize reference concentrations, reactions orders and cathodic transfer coefficien...
Definition: tafel_kinetics.h:162
Definition: system_management.h:80
Virtual class used to provide the interface for all kinetic/reaction children.
Definition: base_kinetics.h:102
static TafelKinetics const * PROTOTYPE
Create prototype for the layer.
Definition: tafel_kinetics.h:153
double alpha_c
Cathodic transfer coefficient.
Definition: tafel_kinetics.h:203
double negative_concentration_correction(double concentration, double gamma, bool derivative=true) const
If any of the oxygen concentrations are negative, make the ORR reaction negative. ...
Definition: tafel_kinetics.h:184
static const std::string concrete_name
Concrete name used for objects of this class.
Definition: tafel_kinetics.h:118