OpenFCST: The open-source Fuel Cell Simulation Toolbox
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
double_trap_kinetics.h
Go to the documentation of this file.
1 //---------------------------------------------------------------------------
2 //
3 // FCST: Fuel Cell Simulation Toolbox
4 //
5 // Copyright (C) 2010-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: double_trap_kinetics.h
11 // - Description: DoubleTrap Kinetics implements the double trap model proposed by
12 // Wang et al. and re-developed by Moore et al.
13 // - Developers: M. Moore, M. Bhaiya, M. Secanell, P. Wardlaw
14 // - $Id: double_trap_kinetics.h 2605 2014-08-15 03:36:44Z secanell $
15 //
16 //---------------------------------------------------------------------------
17 
18 #ifndef _FUELCELLSHOP__DOUBLETRAP_KINETICS_H
19 #define _FUELCELLSHOP__DOUBLETRAP_KINETICS_H
20 
21 // Include openFCST routines:
22 #include "base_kinetics.h"
23 #include "fcst_constants.h"
24 
25 // Include deal.II classes
26 #include <base/parameter_handler.h>
27 #include <base/point.h>
28 #include <base/function.h>
29 #include <lac/vector.h>
30 #include <fe/fe_values.h>
31 
32 //Include STL
33 #include <cmath>
34 #include <iostream>
35 
36 using namespace dealii;
37 
38 namespace FuelCellShop
39 {
40  namespace Kinetics
41  {
67  public BaseKinetics
68  {
69  public:
71 
72 
83  static const std::string concrete_name;
85 
87 
88 
93  DoubleTrapKinetics(const bool);
105  virtual void declare_parameters(ParameterHandler &param) const;
106 
111  virtual void initialize(ParameterHandler &param);
112 
114 
116 
117 
120  virtual void current_density (std::vector<double>& );
121 
128  virtual void derivative_current (std::map< VariableNames, std::vector<double> >& );
130 
132 
133 
134  inline void OH_coverage(std::vector<double>& temp) const
135  {
136  temp = theta_OH;
137  }
138 
140  inline void O_coverage(std::vector<double>& temp) const
141  {
142  temp = theta_O;
143  }
144 
146  inline void G_DA_rate(std::vector<double>& temp) const
147  {
148  double C;
149  if (reactants_map.at(oxygen_concentration)[0] < 0.0)
150  C = 1.0e-12;
151  else
152  C = pow(((reactants_map.at(oxygen_concentration)[0])/c_ref_oxygen),0.5);
153 
154  std::vector<double> t;
155 
156  for(unsigned int i =0; i < intrinsic_rates[0].size(); i++)
157  t.push_back(C*intrinsic_rates[0][i]);
158 
159  temp = t;//free_energies[0];
160  }
161 
163  inline void g_da_rate(std::vector<double>& temp) const
164  {
165  std::vector<double> t;
166 
167  for(unsigned int i =0; i < intrinsic_rates[1].size(); i++)
168  t.push_back(intrinsic_rates[1][i]*theta_O[i]);
169 
170  temp = t;//free_energies[1];
171  }
172 
174  inline void G_RA_rate(std::vector<double>& temp) const
175  {
176  double C;
177  if (reactants_map.at(oxygen_concentration)[0] < 0.0)
178  C = 1.0e-12;
179  else
180  C = pow(((reactants_map.at(oxygen_concentration)[0])/c_ref_oxygen),0.5);
181 
182  double Ch;
183  if(reactants_map.count(proton_concentration)){
184  Ch=(reactants_map.at(proton_concentration)[0])/c_ref_protons;
185  }
186  else
187  Ch = 1.0;
188 
189 
190  std::vector<double> t;
191  for(unsigned int i =0; i < intrinsic_rates[2].size(); i++)
192  t.push_back(C*Ch*intrinsic_rates[2][i]);
193 
194  temp = t;//free_energies[2];
195  }
196 
198  inline void g_ra_rate(std::vector<double>& temp) const
199  {
200  std::vector<double> t;
201  for(unsigned int i =0; i < intrinsic_rates[3].size(); i++)
202  t.push_back(theta_OH[i]*intrinsic_rates[3][i]);//free_energies[3];
203 
204  temp = t;
205  }
206 
208  inline void G_RT_rate(std::vector<double>& temp) const
209  {
210 
211  double Ch;
212  if(reactants_map.count(proton_concentration)){
213  Ch=(reactants_map.at(proton_concentration)[0])/c_ref_protons;
214  }
215  else
216  Ch = 1.0;
217 
218  std::vector<double> t;
219  for(unsigned int i =0; i < intrinsic_rates[4].size(); i++)
220  t.push_back(theta_O[i]*Ch*intrinsic_rates[4][i]);//free_energies[4];
221 
222  temp = t;
223  }
224 
226  inline void g_rt_rate(std::vector<double>& temp) const
227  {
228  std::vector<double> t;
229  for(unsigned int i =0; i < intrinsic_rates[5].size(); i++)
230  t.push_back(theta_OH[i]*intrinsic_rates[5][i]);
231 
232  temp = t;//free_energies[5];
233  }
234 
236  inline void G_RD_rate(std::vector<double>& temp) const
237  {
238  double Ch;
239  if(reactants_map.count(proton_concentration)){
240  Ch=(reactants_map.at(proton_concentration)[0])/c_ref_protons;
241  }
242  else
243  Ch = 1.0;
244 
245  std::vector<double> t;
246  for(unsigned int i =0; i < intrinsic_rates[6].size(); i++)
247  t.push_back(Ch*theta_OH[i]*intrinsic_rates[6][i]);//free_energies[6];
248 
249  temp = t;
250  }
251 
253  inline void g_rd_rate(std::vector<double>& temp) const
254  {
255  std::vector<double> t;
256  for(unsigned int i =0; i < intrinsic_rates[5].size(); i++)
257  t.push_back(intrinsic_rates[7][i]*(1.0 -theta_OH[i] -theta_O[i]));
258 
259  temp = t;//free_energies[7];
260  }
261 
268  void compute_energies();
270 
272 
273 
277  virtual void set_reaction_kinetics(const ReactionNames & name)
278  {
279  if (name == ORR)
280  {
281  name_reaction_kinetics = name;
282  }
283  else
284  {
285  const std::type_info& info = typeid(*this);
286  FcstUtilities::log << "Only ORR reaction is to be implemented in " << __FUNCTION__ << " called in Class " << info.name() << std::endl;
287  exit(1);
288  }
289  }
291 
292 
293  virtual bool has_coverage(const VariableNames& type){
294 
295  if((type == VariableNames::OH_coverage)or (type == VariableNames::O_coverage))
296  return true;
297 
298  return false;
299  }
300 
301  protected:
302 
304 
305 
310  virtual boost::shared_ptr<FuelCellShop::Kinetics::BaseKinetics > create_replica ()
311  {
312  return boost::shared_ptr<FuelCellShop::Kinetics::BaseKinetics > (new FuelCellShop::Kinetics::DoubleTrapKinetics () );
313  }
319 
320  private:
322 
323 
328  inline double negative_concentration_correction(unsigned int index) const
329  {
330  double correction(1.0);
331 
332  if (reactants_map.at(oxygen_concentration)[index] <= 0.0)
333  correction = 0.0;
334 
335  return correction;
336  }
337 
345  inline double oxygen_concentration_ratio(unsigned int index) const
346  {
347  double ratio;
348  double exponential = 0.5;
349 
350  if (reactants_map.at(oxygen_concentration)[index] < 0.0)
351  ratio = 1.0e-12;
352  else
353  ratio = pow((reactants_map.at(oxygen_concentration)[index]/c_ref_oxygen), exponential);
354 
355  return ratio;
356  }
357 
365  inline double proton_concentration_ratio(unsigned int index) const
366  {
367  double ratio = 1.0;
368 
369  if(reactants_map.count(proton_concentration))
370  ratio =(reactants_map.at(proton_concentration)[index])/c_ref_protons;
371 
372  return ratio;
373  }
377  virtual void init_kin_param()
378  {
379  Assert( !kin_param_initialized, ExcInternalError() );
380  Assert( catalyst != NULL, ExcMessage("Catalyst object not initialized in the DoubleTrapKinetics object.") );
381  Assert( catalyst->get_reaction_name() == ORR, ExcMessage("Catalyst object in the DoubleTrapKinetics not set to ORR reaction name.") );
382  Assert( phi_m.is_initialized() && phi_s.is_initialized() && T.is_initialized(), ExcMessage("Either phi_m/phi_s/T is not set in the DoubleTrapKinetics object.") );
383  Assert( reactants_map.find(oxygen_concentration) != reactants_map.end(), ExcMessage("Oxygen concentration is not set in the DoubleTrapKinetics object.") );
384 
385  std::vector<VariableNames> names(1, oxygen_concentration);
386  names.push_back(proton_concentration);
387  std::map< VariableNames, double > cref_map;
388  catalyst->reference_concentration(names, cref_map);
389 
390  c_ref_oxygen = cref_map[oxygen_concentration];
391  c_ref_protons = cref_map[proton_concentration];
392  kin_param_initialized = true;
393  }
394 
396  void d_OH_coverage_du (std::map< VariableNames, std::vector<double> >& ) const;
397 
399  void d_O_coverage_du (std::map< VariableNames, std::vector<double> >& ) const;
400 
404  void drate_du(std::vector<std::vector<double> >&, const VariableNames& ) const;
405 
417  std::vector<std::vector<double> > free_energies;
418 
429  double G_DA_0, G_RA_0, G_RT_0, G_RD_0;
430 
440  double G_O, G_OH;
441 
446  std::vector<std::vector<double> > intrinsic_rates;
447 
449  std::vector<double> theta_OH, theta_O;
450 
455  double prefactor;
456 
458  double c_ref_oxygen;
459 
462 
464  double alpha;
466  };
467  } //Kinetics
468 } //FuelCellShop
469 
470 #endif //_FUELCELLSHOP__DOUBLETRAP_KINETICS_H
Definition: system_management.h:85
double proton_concentration_ratio(unsigned int index) const
If needed, compute the proton concentration ratio, i.e., where is the variable exponential.
Definition: double_trap_kinetics.h:365
void G_RA_rate(std::vector< double > &temp) const
Returns the rate of the forward RA step.
Definition: double_trap_kinetics.h:174
void g_da_rate(std::vector< double > &temp) const
Returns the rate of the backward DA step.
Definition: double_trap_kinetics.h:163
virtual void set_reaction_kinetics(const ReactionNames &name)
Member function used to set the reaction name in the DoubleTrap kinetics object.
Definition: double_trap_kinetics.h:277
virtual bool has_coverage(const VariableNames &type)
Definition: double_trap_kinetics.h:293
void g_rt_rate(std::vector< double > &temp) const
Returns the rate of the backward RT step.
Definition: double_trap_kinetics.h:226
double negative_concentration_correction(unsigned int index) const
If any of the oxygen concentrations are negative, make the ORR reaction negative. ...
Definition: double_trap_kinetics.h:328
VariableNames
The enumeration containing the names of some of the available FCST solution variables and their deriv...
Definition: system_management.h:62
double alpha
Stores the transfer coefficient.
Definition: double_trap_kinetics.h:464
std::vector< std::vector< double > > intrinsic_rates
Values of the intrinsic rates for each reaction ( ).
Definition: double_trap_kinetics.h:446
std::vector< std::vector< double > > free_energies
Values of the free activation energies for each reaction ( ).
Definition: double_trap_kinetics.h:417
static DoubleTrapKinetics const * PROTOTYPE
Create prototype for the layer.
Definition: double_trap_kinetics.h:317
double c_ref_oxygen
Stores the reference concentration of oxygen.
Definition: double_trap_kinetics.h:458
Definition: system_management.h:83
Definition: system_management.h:82
std::vector< double > theta_OH
Stores the value of the coverage of the OH and O intermediates.
Definition: double_trap_kinetics.h:449
Definition: system_management.h:151
virtual void init_kin_param()
Method used to initialize reference concentration of oxygen for the reaction, and number of quadratur...
Definition: double_trap_kinetics.h:377
void O_coverage(std::vector< double > &temp) const
Returns the value of the coverage of the O intermediate.
Definition: double_trap_kinetics.h:140
void g_ra_rate(std::vector< double > &temp) const
Returns the rate of the backward RA step.
Definition: double_trap_kinetics.h:198
void G_DA_rate(std::vector< double > &temp) const
Returns the rate of the forward DA step.
Definition: double_trap_kinetics.h:146
double G_OH
Definition: double_trap_kinetics.h:440
FCSTLogStream log
Object used to output data to file and, if file attached recorded to a file as well.
ReactionNames
Definition: system_management.h:148
Definition: system_management.h:79
Virtual class used to provide the interface for all kinetic/reaction children.
Definition: base_kinetics.h:107
double G_RT_0
Definition: double_trap_kinetics.h:429
void G_RT_rate(std::vector< double > &temp) const
Returns the rate of the forward RT step.
Definition: double_trap_kinetics.h:208
void OH_coverage(std::vector< double > &temp) const
Returns the value of the coverage of the OH intermediate.
Definition: double_trap_kinetics.h:134
static const std::string concrete_name
Concrete name used for objects of this class.
Definition: double_trap_kinetics.h:83
double prefactor
prefactor is the reference prefactor that sets the scale for the current produced.
Definition: double_trap_kinetics.h:455
void G_RD_rate(std::vector< double > &temp) const
Returns the rate of the forward RD step.
Definition: double_trap_kinetics.h:236
void g_rd_rate(std::vector< double > &temp) const
Returns the rate of the backward RD step.
Definition: double_trap_kinetics.h:253
Definition: system_management.h:89
virtual boost::shared_ptr< FuelCellShop::Kinetics::BaseKinetics > create_replica()
This member function is used to create an object of type gas diffusion layer.
Definition: double_trap_kinetics.h:310
double c_ref_protons
Stores the reference concentration of protons.
Definition: double_trap_kinetics.h:461
double oxygen_concentration_ratio(unsigned int index) const
Compute the oxygen concentration ratio, i.e., where is the variable exponential.
Definition: double_trap_kinetics.h:345
This class contains the implementation of the double trap kinetic model as described in the following...
Definition: double_trap_kinetics.h:66