OpenFCST: The open-source Fuel Cell Simulation Toolbox
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
PSD_dual.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 #ifndef _FUELCELLSHOP__DUAL__PSD_H
17 #define _FUELCELLSHOP__DUAL__PSD_H
18 
19 // Include OpenFCST routines:
20 #include <microscale/PSD_HI.h>
21 #include <microscale/PSD_HO.h>
22 
23 
24 
25 using namespace dealii;
26 
27 namespace FuelCellShop
28 {
29 
30 
31  namespace MicroScale
32  {
96  template <int dim>
97  class DualPSD : public BasePSD<dim>
98  {
99  public:
100 
102 
103 
106  DualPSD();
107 
111  DualPSD (std::string name);
112 
116  virtual ~DualPSD();
117 
121  void declare_parameters (ParameterHandler &param) const;
126  void initialize (ParameterHandler &param );
127 
132  inline void set_temperature (const SolutionVariable& T_in)
133  {
134  Assert( T_in.get_variablename() == temperature_of_REV, ExcMessage("Wrong solution variable passed in PSD::set_temperature.") );
135  this->T_vector = T_in;
136  }
137 
142  inline void set_capillary_pressure (const SolutionVariable& C_in)
143  {
144  Assert( C_in.get_variablename() == capillary_pressure, ExcMessage("Wrong solution variable passed in PSD::capillary_pressure.") );
145  this->Capillary_pressure_vector = C_in;
146 
147  psd_hi.set_capillary_pressure(C_in);
148  psd_ho.set_capillary_pressure(C_in);
149  }
150 
157  inline void set_critical_radius()
158  {
159  psd_hi.set_critical_radius();
160  psd_ho.set_critical_radius();
161 
162  }
163 
170  inline void set_saturation()
171  {
172  psd_hi.set_saturation();
173  psd_ho.set_saturation();
174  }
175 
176 
178 
180 
181 
195  static const std::string concrete_name;
196 
198 
200 
201 
214  const std::type_info& get_base_type() const
215  {
216  return typeid(DualPSD<dim>);
217  }
218 
228  virtual void get_saturation(std::vector<double>& S) const ;
229  virtual void get_derivative_saturation(std::vector<double>& S) const ;
236  virtual void get_global_saturated_permeability(double& saturated_permeability) const ;
237 
244  virtual void get_relative_liquid_permeability(std::vector<double>& liquid_permeability) const ;
245  virtual void get_derivative_relative_liquid_permeability(std::vector<double>& liquid_permeability) const ;
252  virtual void get_relative_gas_permeability(std::vector<double>& gas_permeability) const ;
253 
260  virtual void get_liquid_gas_interfacial_surface(std::vector<double>& liquid_gas_interfacial_surface) const ;
261 
268  virtual void get_derivative_liquid_gas_interfacial_surface(std::vector<double>&) const ;
269 
276  virtual void get_wetted_wall_surface_area(std::vector<double>& wetted_wall_surface_area) const ;
277 
284  virtual void get_knudsen_radius(std::vector<double>& knudsen_radius) const ;
285 
289  virtual void get_diffusivity() const ;
290 
291  virtual void get_maximum_cross_sectional_areas(double&) const;
292 
296  virtual void get_PSD_plot(const std::vector<double> ,std::vector<double>&,std::vector<double>&,std::vector<double>&) const;
297 
298 
300 
301  protected:
303 
304 
309  virtual boost::shared_ptr<FuelCellShop::MicroScale::BasePSD<dim> > create_replica (const std::string &psd_section_name)
310  {
311  return boost::shared_ptr<FuelCellShop::MicroScale::BasePSD<dim> > (new FuelCellShop::MicroScale::DualPSD<dim> (psd_section_name));
312  }
314 
316 
317 
320  static DualPSD<dim> const* PROTOTYPE;
322 
324  // DATA //
326 
328 
329 
333 
337 
338  double pressure_c;
339 
342 
345 
347  };
348 
349  } // PSD
350 
351 } // FuelCellShop
352 
353 #endif
Definition: system_management.h:78
Hydrophobic Pore Size Distribution.
Definition: PSD_HO.h:91
VariableNames get_variablename() const
Function to get the VariableNames enumeration corresponding to this struct.
Definition: fcst_variables.h:170
SolutionVariable Capillary_pressure_vector
Capillary pressure at every quadrature point inside the cell.
Definition: PSD_dual.h:344
static const std::string concrete_name
Concrete name used for objects of this class.
Definition: PSD_dual.h:195
const std::type_info & get_base_type() const
This member function returns a type_info object with the name of the base layer type the inherited cl...
Definition: PSD_dual.h:214
This structure is used to encapsulate data from constant values and variable solution data that is us...
Definition: fcst_variables.h:86
Definition: system_management.h:76
SolutionVariable T_vector
Temperature at every quadrature point inside the cell.
Definition: PSD_dual.h:341
void set_critical_radius()
Member function used to set the critical radius [nm] at every quadrature point inside the cell...
Definition: PSD_dual.h:157
HIPSD< dim > psd_hi
Creating a pointer of PSD hydrophilic class.
Definition: PSD_dual.h:332
Hydrophilic Pore Size Distribution.
Definition: PSD_HI.h:89
HOPSD< dim > psd_ho
Creating a pointer of PSD hydrophobic class.
Definition: PSD_dual.h:336
void set_temperature(const SolutionVariable &T_in)
Member function used to set the temperature [Kelvin] at every quadrature point inside the cell...
Definition: PSD_dual.h:132
Dual Pore Size Distribution.
Definition: PSD_dual.h:97
Pore Size Distribution.
Definition: PSD_base.h:129
virtual boost::shared_ptr< FuelCellShop::MicroScale::BasePSD< dim > > create_replica(const std::string &psd_section_name)
This member function is used to create an object of type psd.
Definition: PSD_dual.h:309
static DualPSD< dim > const * PROTOTYPE
PROTOTYPE is the pointer is the dynamic pointer pointing to the DualPSD class itself.
Definition: PSD_dual.h:320
double pressure_c
Definition: PSD_dual.h:338
void set_saturation()
Member function used to set the saturation at every quadrature point inside the cell.
Definition: PSD_dual.h:170
void set_capillary_pressure(const SolutionVariable &C_in)
Member function used to set the capillary pressure [psi] at every quadrature point inside the cell...
Definition: PSD_dual.h:142