17 #ifndef _FUELCELLSHOP__TAFEL_KINETICS_H
18 #define _FUELCELLSHOP__TAFEL_KINETICS_H
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>
35 using namespace dealii;
37 namespace FuelCellShop
95 virtual void declare_parameters(ParameterHandler& param)
const;
102 virtual void initialize(ParameterHandler& param);
135 virtual void derivative_current (std::map<
VariableNames, std::vector<double> > &);
147 virtual boost::shared_ptr<FuelCellShop::Kinetics::BaseKinetics >
create_replica ()
165 Assert( !kin_param_initialized, ExcInternalError() );
166 Assert( catalyst != NULL, ExcMessage(
"Catalyst object not initialized in the TafelKinetics object.") );
167 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.") );
168 Assert( reactants_map.size()>0, ExcMessage(
"Atleast one reactant should be set in the TafelKinetics object, using set_reactant_concentrations method.") );
170 std::vector<VariableNames> names;
171 for ( std::map< VariableNames, SolutionVariable >::const_iterator iter=reactants_map.begin(); iter!=reactants_map.end(); ++iter )
172 names.push_back(iter->first);
174 catalyst->reference_concentration(names, ref_conc);
175 catalyst->reaction_order(names, gamma);
176 catalyst->alpha_cathodic(alpha_c);
178 kin_param_initialized =
true;
194 std::map< VariableNames, double >
gamma;
203 #endif //_FUELCELLSHOP__TAFEL_KINETICS_H
std::map< VariableNames, double > gamma
Map of reaction orders.
Definition: tafel_kinetics.h:194
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:147
VariableNames
The enumeration containing the names of some of the available FCST solution variables and their deriv...
Definition: system_management.h:62
std::map< VariableNames, double > ref_conc
Map of reference concentrations.
Definition: tafel_kinetics.h:189
This class defines a simple Tafel kinetic model.
Definition: tafel_kinetics.h:67
virtual void init_kin_param()
Method used to initialize reference concentrations, reactions orders and cathodic transfer coefficien...
Definition: tafel_kinetics.h:163
Definition: system_management.h:79
Virtual class used to provide the interface for all kinetic/reaction children.
Definition: base_kinetics.h:107
static TafelKinetics const * PROTOTYPE
Create prototype for the layer.
Definition: tafel_kinetics.h:154
double alpha_c
Cathodic transfer coefficient.
Definition: tafel_kinetics.h:184
static const std::string concrete_name
Concrete name used for objects of this class.
Definition: tafel_kinetics.h:119