OpenFCST: The open-source Fuel Cell Simulation Toolbox
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
PSD_base.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: PSD_base.h
11 // - Description: Base class for pore size distribution model.
12 // - Developers: 2009-13 by Marc Secanell, University of Alberta
13 // 2013-14 by Jie Zhou, University of Alberta
14 // - $ $
15 //
16 //---------------------------------------------------------------------------
17 #ifndef _FUELCELLSHOP__BASE__PSD_H
18 #define _FUELCELLSHOP__BASE__PSD_H
19 
20 // Include deal.II classes
21 #include <base/parameter_handler.h>
22 #include <base/point.h>
23 #include <base/function.h>
24 #include <lac/vector.h>
25 #include <fe/fe_values.h>
26 
27 //Include STL
28 #include <cmath>
29 #include <iostream>
30 #include <math.h>
31 
32 // Include OpenFCST routines:
35 #include "utils/fcst_utilities.h"
36 #include "utils/fcst_constants.h"
37 #include "porous_layer.h"
38 
39 using namespace dealii;
40 
41 namespace FuelCellShop
42 {
43 
44 
45  namespace MicroScale
46  {
129  template <int dim>
130  class BasePSD : public Subscriptor
131  {
132  public:
133 
135 
136 
139  virtual ~BasePSD();
141 
143 
144 
163  static void declare_PSD_parameters (ParameterHandler &param)
164  {
167  iterator++)
168  {
169  iterator->second->declare_parameters(param);
170  }
171  }
172 
186  static boost::shared_ptr<FuelCellShop::MicroScale::BasePSD<dim> > create_PSD (const std::string& psd_section_name,
187  ParameterHandler &param)
188  {
189  boost::shared_ptr<FuelCellShop::MicroScale::BasePSD<dim> > pointer;
190 
191  std::string concrete_name;
192 
193  param.enter_subsection("PSD parameters");
194  {
195  param.enter_subsection("BasePSD");
196  {
197  concrete_name = param.get("psd type");
198  }
199  param.leave_subsection();
200  }
201  param.leave_subsection();
202 
204 
206  {
207  if (iterator->second)
208  {
209  pointer = iterator->second->create_replica(psd_section_name);
210  }
211  else
212  {
213  deallog<<"Pointer not initialized"<<std::endl;
214  abort();
215  }
216  }
217  else
218  {
219  deallog<<"Concrete name in FuelCellShop::MicroScale::BasePSD::create_psd does not exist"<<std::endl;
220  abort();
221  }
222 
223  pointer->initialize(param);
224 
225  return pointer;
226  }
228 
230 
231 
234  inline void set_porosity(double porosity )
235  {
236  this->por = porosity;
237  }
238 
242  inline double get_porosity() const
243  {
244  return this->por;
245  }
246 
252  inline void set_derivative_flags(const std::vector<VariableNames>& flags)
253  {
254  this->derivative_flags = flags;
255  }
256 
268  virtual void set_constant_solution(const double& value, const VariableNames& name)
269  {
270  constant_solutions[name] = value;
271  }
272 
282  virtual void set_solution(const std::vector< SolutionVariable >&)
283  {
284  const std::type_info& info = typeid(*this);
285  deallog << "Pure function " << __FUNCTION__
286  << " called in Class "
287  << info.name() << std::endl;
288  }
290 
292 
293 
296  inline const std::string& name_psd() const
297  {
298  return name;
299  }
300 
305  virtual void print_psd_properties() const
306  {
307  const std::type_info& info = typeid(*this);
308  deallog << "Pure function " << __FUNCTION__
309  << " called in Class "
310  << info.name() << std::endl;
311  }
312 
325  virtual void get_saturation(std::vector<double>& ) const = 0;
326 
332  virtual void get_global_saturated_permeability(double& ) const = 0;
333 
339  virtual void get_relative_liquid_permeability(std::vector<double>& ) const = 0;
340 
346  virtual void get_relative_gas_permeability(std::vector<double>& ) const = 0;
347 
353  virtual void get_liquid_gas_interfacial_surface(std::vector<double>& ) const = 0;
354 
360  virtual void get_wetted_wall_surface_area(std::vector<double>& ) const = 0;
361 
367  virtual void get_knudsen_radius(std::vector<double>& ) const = 0;
368 
370 
371  protected:
373 
374 
384  {};
385 
389  BasePSD(const std::string& name);
390 
398  void declare_parameters (ParameterHandler &param) const;
399 
404  void initialize (ParameterHandler &param) ;
405 
407 
409 
412  typedef std::map< std::string, BasePSD<dim>* > _mapFactory;
414 
416 
417 
465  {
466  static _mapFactory mapFactory;
467  return &mapFactory;
468  }
470 
471 
476  virtual boost::shared_ptr<FuelCellShop::MicroScale::BasePSD<dim> > create_replica (const std::string &name)
477  {
478  const std::type_info& info = typeid(*this);
479  deallog << "Pure function " << __FUNCTION__
480  << " called in Class "
481  << info.name() << std::endl;
482  }
484 
486 
487 
490  const std::string name;
491 
493  std::vector<VariableNames> derivative_flags;
494 
498  std::map< VariableNames, double > constant_solutions;
499 
503  double gamma;
504 
509 
513  double lamda;
514 
518  double P_b;
519 
523  double F_HI;
524 
528  double F_HO;
529 
533  std::vector<double> f_k;
534 
538  std::vector<double> r_k;
539 
543  std::vector<double> s_k;
544 
548  double por;
550  };
551 
552  } // PSD
553 
554 } // FuelCellShop
555 
556 #endif
virtual boost::shared_ptr< FuelCellShop::MicroScale::BasePSD< dim > > create_replica(const std::string &name)
This member function is used to create an object of type gas diffusion layer.
Definition: PSD_base.h:476
double get_porosity() const
Definition: PSD_base.h:242
VariableNames
The enumeration containing the names of some of the available FCST solution variables and their deriv...
Definition: system_management.h:62
std::vector< double > f_k
The f_k is the contribution of the log-normal distribution k to the total PSD.
Definition: PSD_base.h:533
const std::string & name_psd() const
Return the name of the PSD.
Definition: PSD_base.h:296
double P_b
The effects of pore interconnectivity is represented by probability P_b.
Definition: PSD_base.h:518
void set_porosity(double porosity)
Definition: PSD_base.h:234
static void declare_PSD_parameters(ParameterHandler &param)
Function used to declare all the data necessary in the parameter files for all BasePSD children...
Definition: PSD_base.h:163
double lamda
The lamda is measured from the mercury intrusion experiment.
Definition: PSD_base.h:513
virtual void set_constant_solution(const double &value, const VariableNames &name)
Set those solution variables which are constant in the particular application.
Definition: PSD_base.h:268
virtual void print_psd_properties() const
This function prints out the psd properties.
Definition: PSD_base.h:305
std::map< VariableNames, double > constant_solutions
Map storing values of solution variables constant in a particular application.
Definition: PSD_base.h:498
std::vector< VariableNames > derivative_flags
Flags for derivatives: These flags are used to request derivatives of material properties.
Definition: PSD_base.h:493
const std::string name
Name of the psd.
Definition: PSD_base.h:490
double por
The por is the porosity inherited from the layer class.
Definition: PSD_base.h:548
std::vector< double > r_k
The r_k is the characteristic pore size of the distribution k.
Definition: PSD_base.h:538
static _mapFactory * get_mapFactory()
Return the map library that stores all childrens of this class.
Definition: PSD_base.h:464
double contact_angle
The contact_angle is the contact angle between air and water.
Definition: PSD_base.h:508
virtual void set_solution(const std::vector< SolutionVariable > &)
If the effective properties in the psd depend on the solution, the solution for a given cell should b...
Definition: PSD_base.h:282
std::map< std::string, BasePSD< dim > * > _mapFactory
This object is used to store all objects of type psd.
Definition: PSD_base.h:412
double gamma
The gamma is the Water-air interface surface tension.
Definition: PSD_base.h:503
Pore Size Distribution.
Definition: PSD_base.h:130
static boost::shared_ptr< FuelCellShop::MicroScale::BasePSD< dim > > create_PSD(const std::string &psd_section_name, ParameterHandler &param)
Function used to select the appropriate CatalystLayer type as specified in the ParameterHandler under...
Definition: PSD_base.h:186
double F_HO
The F_HO is the fraction of the total volume corresponding to the hydrophobic pores.
Definition: PSD_base.h:528
double F_HI
The F_HI is the fraction of the total volume corresponding to the hydrophilic pores.
Definition: PSD_base.h:523
BasePSD()
Constructor.
Definition: PSD_base.h:383
void set_derivative_flags(const std::vector< VariableNames > &flags)
Set the names of FCST solution variables with respect to which you would like to compute the derivati...
Definition: PSD_base.h:252
std::vector< double > s_k
The s_k is the spread of the distribution k.
Definition: PSD_base.h:543