16 #ifndef _FUELCELLSHOP__TAFEL_KINETICS_H
17 #define _FUELCELLSHOP__TAFEL_KINETICS_H
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>
34 using namespace dealii;
36 namespace FuelCellShop
94 virtual void declare_parameters(ParameterHandler& param)
const;
101 virtual void initialize(ParameterHandler& param);
134 virtual void derivative_current (std::map<
VariableNames, std::vector<double> > &);
146 virtual boost::shared_ptr<FuelCellShop::Kinetics::BaseKinetics >
create_replica ()
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.") );
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);
173 catalyst->reference_concentration(names, ref_conc);
174 catalyst->reaction_order(names, gamma);
175 catalyst->alpha_cathodic(alpha_c);
177 kin_param_initialized =
true;
186 double correction(1.0);
188 if (concentration < 0.0)
192 else if (derivative && (gamma == 1.0 || gamma == 3.0 || gamma == 5.0) )
194 else if (derivative && ( (gamma != 1.0) || (gamma != 3.0) || (gamma != 5.0) ) )
213 std::map< VariableNames, double >
gamma;
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