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 and M. Secanell
13 // - $Id: tafel_kinetics.h 1373 2013-08-21 00:34:53Z madhur $
14 //
15 //---------------------------------------------------------------------------
16 
17 #ifndef _FUELCELLSHOP__TAFEL_KINETICS_H
18 #define _FUELCELLSHOP__TAFEL_KINETICS_H
19 
20 // Include openFCST routines:
21 #include "base_kinetics.h"
22 #include "fcst_constants.h"
23 
24 // Include deal.II classes
25 #include <base/parameter_handler.h>
26 #include <base/point.h>
27 #include <base/function.h>
28 #include <lac/vector.h>
29 #include <fe/fe_values.h>
30 
31 //Include STL
32 #include <cmath>
33 #include <iostream>
34 
35 using namespace dealii;
36 
37 namespace FuelCellShop
38 {
39  namespace Kinetics
40  {
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&) const{};
95 
100  virtual void set_parameters(const std::vector<std::string>&,
101  const std::vector<double>&,
102  ParameterHandler&) const{};
103 
109  virtual void initialize(ParameterHandler&){};
110 
112 
114 
115 
126  static const std::string concrete_name;
128 
130 
131 
134  virtual void current_density (std::vector<double> &);
135 
142  virtual void derivative_current (std::map< VariableNames, std::vector<double> > &);
143 
145 
146  protected:
148 
149 
154  virtual boost::shared_ptr<FuelCellShop::Kinetics::BaseKinetics > create_replica ()
155  {
156  return boost::shared_ptr<FuelCellShop::Kinetics::BaseKinetics > (new FuelCellShop::Kinetics::TafelKinetics ( ));
157  }
161  static TafelKinetics const* PROTOTYPE;
163 
165 
166 
170  virtual void init_kin_param()
171  {
172  Assert( !kin_param_initialized, ExcInternalError() );
173  Assert( catalyst != NULL, ExcMessage("Catalyst object not initialized in the TafelKinetics object.") );
174  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.") );
175  Assert( reactants_map.size()>0, ExcMessage("Atleast one reactant should be set in the TafelKinetics object, using set_reactant_concentrations method.") );
176 
177  std::vector<VariableNames> names;
178  for ( std::map< VariableNames, SolutionVariable >::const_iterator iter=reactants_map.begin(); iter!=reactants_map.end(); ++iter )
179  names.push_back(iter->first);
180 
181  catalyst->reference_concentration(names, ref_conc);
182  catalyst->reaction_order(names, gamma);
183  catalyst->alpha_cathodic(alpha_c);
184 
185  kin_param_initialized = true;
186  }
187 
191  double alpha_c;
192 
196  std::map< VariableNames, double > ref_conc;
197 
201  std::map< VariableNames, double > gamma;
203 
204  };
205 
206  } //Kinetics
207 
208 } //FuelCellShop
209 
210 #endif //_FUELCELLSHOP__TAFEL_KINETICS_H