OpenFCST: The open-source Fuel Cell Simulation Toolbox
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
base_kinetics.h
Go to the documentation of this file.
1 //---------------------------------------------------------------------------
2 //
3 // FCST: Fuel Cell Simulation Toolbox
4 //
5 // Copyright (C) 2013 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: base_kinetics.h
11 // - Description: Base Kinetics. It implements the interface for other kinetics classes
12 // and some common methods.
13 // - Developers: M. Moore, M. Bhaiya and M. Secanell
14 // - $Id: base_kinetics.h 1459 2013-08-30 22:27:10Z secanell $
15 //
16 //---------------------------------------------------------------------------
17 
18 #ifndef _FUELCELLSHOP__BASE__KINETICS_H
19 #define _FUELCELLSHOP__BASE__KINETICS_H
20 
21 // Include deal.II classes
22 #include <base/parameter_handler.h>
23 #include <base/point.h>
24 #include <base/function.h>
25 #include <lac/vector.h>
26 #include <fe/fe_values.h>
27 
28 // Include STL
29 #include <cmath>
30 #include <iostream>
31 #include <map>
32 #include <algorithm>
33 
34 // Include OpenFCST routines:
35 #include "base_layer.h"
36 #include "catalyst_base.h"
38 #include "fcst_constants.h"
39 #include "system_management.h"
40 
41 using namespace dealii;
42 
43 namespace FuelCellShop
44 {
45  namespace Kinetics
46  {
113  {
114  public:
115 
117 
122  static void declare_Kinetics_parameters (ParameterHandler &param)
123  {
124 
125  for (typename FuelCellShop::Kinetics::BaseKinetics::_mapFactory::iterator iterator = FuelCellShop::Kinetics::BaseKinetics::get_mapFactory()->begin();
127  iterator++)
128  {
129  iterator->second->declare_parameters(param);
130  }
131  }
135  static void set_Kinetics_parameters (const std::vector<std::string>& name_dvar,
136  const std::vector<double>& value_dvar,
137  ParameterHandler &param)
138  {
139  for (typename FuelCellShop::Kinetics::BaseKinetics::_mapFactory::iterator iterator = FuelCellShop::Kinetics::BaseKinetics::get_mapFactory()->begin();
141  iterator++)
142  {
143  iterator->second->set_parameters(name_dvar,value_dvar, param);
144  }
145  }
164  static boost::shared_ptr<FuelCellShop::Kinetics::BaseKinetics > create_Kinetics (ParameterHandler &param,
165  std::string kinetics_name)
166  {
167 
168  boost::shared_ptr<FuelCellShop::Kinetics::BaseKinetics > pointer;
169 
170  typename FuelCellShop::Kinetics::BaseKinetics::_mapFactory::iterator iterator = FuelCellShop::Kinetics::BaseKinetics::get_mapFactory()->find(kinetics_name);
171 
173  {
174  if (iterator->second)
175  {
176  pointer = iterator->second->create_replica();
177  }
178  else
179  {
180  deallog<<"Pointer not initialized"<<std::endl;
181  abort();
182  }
183  }
184  else
185  {
186  deallog<<"Concrete name does not exist"<<std::endl;
187  abort();
188  }
189 
190  pointer->initialize(param);
191 
192  return pointer;
193  }
195 
197 
198 
203  void set_electrolyte_potential(const SolutionVariable& phi)
204  {
205  Assert( phi.get_variablename() == protonic_electrical_potential, ExcMessage("Wrong solution variable passed in BaseKinetics::set_electrolyte_potential.") );
206  phi_m = phi;
207  }
211  void set_solid_potential(const SolutionVariable& phi)
212  {
213  Assert( phi.get_variablename() == electronic_electrical_potential, ExcMessage("Wrong solution variable passed in BaseKinetics::set_solid_potential.") );
214  phi_s = phi;
215  }
219  void set_temperature(const SolutionVariable& temperature)
220  {
221  Assert( temperature.get_variablename() == temperature_of_REV, ExcMessage("Wrong solution variable passed in BaseKinetics::set_temperature") );
222  T = temperature;
223  }
227  void set_reactant_concentrations(const std::vector<SolutionVariable>& conc_vec)
228  {
229  Assert( conc_vec.size() > 0, ExcMessage("Atleast one reactant concentration should be passed inside BaseKinetics::set_reactant_concentrations. !!") );
230 
231  for ( unsigned int i=0; i<conc_vec.size(); ++i )
232  {
233  Assert( (conc_vec[i].get_variablename() != hydrogen_molar_fraction &&
234  conc_vec[i].get_variablename() != oxygen_molar_fraction &&
235  conc_vec[i].get_variablename() != water_molar_fraction), ExcMessage("Molar fractions can't be setup using BaseKinetics::set_reactant_concentrations.") );
236 
237  reactants_map[ conc_vec[i].get_variablename() ] = conc_vec[i];
238  }
239  }
240 
244  void set_derivative_flags(const std::vector<VariableNames>& flags)
245  {
246  derivative_flags = flags;
247  };
248 
253  void set_catalyst(FuelCellShop::Material::CatalystBase* cat_in)
254  {
255  catalyst = cat_in;
256  }
257 
263  void set_electrolyte(FuelCellShop::Material::PolymerElectrolyteBase* electrolyte_in)
264  {
265  electrolyte = electrolyte_in;
266  }
267 
273  virtual void set_reaction_kinetics(const std::string& name)
274  {
275  name_reaction_kinetics = name;
276  }
277 
279  void set_p_t(const double& P_Tot)
280  {
281  p_total = P_Tot;
282  }
284 
285 
286 
289  inline std::string get_reaction_name() const
290  {
291  return name_reaction_kinetics;
292  }
293 
296  {
297  return catalyst;
298 
299  }
300 
304  virtual void current_density (std::vector<double>&)
305  {
306  const std::type_info& info = typeid(*this);
307  deallog << "Pure function " << __FUNCTION__
308  << " called in Class "
309  << info.name() << std::endl;
310 
311  };
312 
319  virtual void derivative_current (std::map< VariableNames, std::vector<double> >&)
320  {
321  const std::type_info& info = typeid(*this);
322  deallog << "Pure function " << __FUNCTION__
323  << " called in Class "
324  << info.name() << std::endl;
325  }
326 
330  virtual void compute_coverages(const std::string& name_species,
331  std::vector<double>& coverage) const
332  {
333  const std::type_info& info = typeid(*this);
334  deallog << "Pure function " << __FUNCTION__
335  << " called in Class "
336  << info.name() << std::endl;
337  };
343  virtual void OH_coverage(std::vector<double>& ) const
344  {
345  const std::type_info& info = typeid(*this);
346  deallog << "Pure function " << __FUNCTION__
347  << " called in Class "
348  << info.name() << std::endl;
349  };
350 
356  virtual void O_coverage(std::vector<double>& ) const
357  {
358  const std::type_info& info = typeid(*this);
359  deallog << "Pure function " << __FUNCTION__
360  << " called in Class "
361  << info.name() << std::endl;
362  }
363 
365 
366  protected:
367 
369 
370 
373  typedef std::map< std::string, BaseKinetics* > _mapFactory;
375 
377 
378 
381  static _mapFactory * get_mapFactory()
382  {
383  static _mapFactory mapFactory;
384  return &mapFactory;
385  }
387 
389 
390 
395  virtual boost::shared_ptr<FuelCellShop::Kinetics::BaseKinetics > create_replica ()
396  {
397  const std::type_info& info = typeid(*this);
398  deallog << "Pure function " << __FUNCTION__
399  << " called in Class "
400  << info.name() << std::endl;
401  }
403 
405 
406 
410  {
411  //deallog << "Creating kinetics object of type ";
412 
413  catalyst = NULL;
414  electrolyte = NULL;
415 
416  // Universal constants
417  R = Constants::R();
418  F = Constants::F();
419  K = Constants::K();
420 
421  kin_param_initialized = false;
422  };
423 
427  virtual ~BaseKinetics()
428  {};
429 
433  virtual void declare_parameters (ParameterHandler&) const
434  {
435  const std::type_info& info = typeid(*this);
436  deallog << "Pure function " << __FUNCTION__
437  << " called in Class "
438  << info.name() << std::endl;
439 
440  };
441 
442  virtual void set_parameters(const std::vector<std::string>& name_dvar,
443  const std::vector<double>& value_dvar,
444  ParameterHandler& param) const
445  {
446  const std::type_info& info = typeid(*this);
447  deallog << "Pure function " << __FUNCTION__
448  << " called in Class "
449  << info.name() << std::endl;
450  }
451 
456  virtual void initialize (ParameterHandler&)
457  {
458  const std::type_info& info = typeid(*this);
459  deallog << "Pure function " << __FUNCTION__
460  << " called in Class "
461  << info.name() << std::endl;
462  };
464 
466 
467 
472  virtual void init_kin_param() = 0;
473 
478 
483 
488  std::vector<VariableNames> derivative_flags;
489 
493  unsigned int n_quad;
494 
506 
508 
509 
510  double R;
512  double F;
513  //Boltzmann constant
515  double K;
517 
519 
520 
521  double p_total;
522 
527  std::map< VariableNames, SolutionVariable > reactants_map;
528 
531 
533  SolutionVariable phi_s
534  ;
538  }; // BaseKinetics
539  } // Kinetics
540 } // FuelCellShop
541 
542 #endif