OpenFCST: The open-source Fuel Cell Simulation Toolbox
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
PSD_HI.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__HI__PSD_H
18 #define _FUELCELLSHOP__HI__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 
31 // Include OpenFCST routines:
34 #include "utils/fcst_utilities.h"
35 #include "utils/fcst_constants.h"
36 #include "PSD_base.h"
37 
38 using namespace dealii;
39 
40 namespace FuelCellShop
41 {
42 
43 
44  namespace MicroScale
45  {
104  template <int dim>
105  class HIPSD
106  :
107  public BasePSD<dim>
108  {
109  public:
111 
112 
119  HIPSD();
120 
124  HIPSD (std::string name);
125 
129  ~HIPSD() {}
130 
134  void declare_parameters (ParameterHandler &param) const;
135 
140  void initialize ( ParameterHandler &param);
141 
146  inline void set_temperature (const SolutionVariable& T_in)
147  {
148  Assert( T_in.get_variablename() == temperature_of_REV, ExcMessage("Wrong solution variable passed in PSD::set_temperature.") );
149  this->T_vector = T_in;
150  }
151 
156  inline void set_capillary_pressure (const SolutionVariable& C_in)
157  {
158  Assert( C_in.get_variablename() == capillary_pressure, ExcMessage("Wrong solution variable passed in PSD::capillary_pressure.") );
159  this->Capillary_pressure_vector = C_in;
160 
161  critical_radius_is_initialized = false;
162  saturation_is_initialized = false;
163  critical_radius_computed.clear();
164  saturation_computed.clear();
165 
166  }
167 
174  inline void set_critical_radius()
175  {
176  get_critical_radius(critical_radius_computed);
177 
178  critical_radius_is_initialized = true;
179  }
180 
187  inline void set_saturation()
188  {
189  get_saturation(saturation_computed);
190 
191  saturation_is_initialized = true;
192  }
206  static const std::string concrete_name;
207 
209 
211 
212 
221  virtual void get_saturation(std::vector<double>& S) const ;
222 
229  virtual void get_global_saturated_permeability(double& saturated_permeability) const ;
230 
237  virtual void get_pore_HI_liquid_saturated_permeability(std::vector<double>& liquid_permeability) const ;
238 
243  virtual void get_relative_liquid_permeability(std::vector<double>& liquid_permeability) const ;
244 
251  virtual void get_pore_HI_gas_saturated_permeability(std::vector<double>& saturated_HI_permeability) const ;
252 
257  virtual void get_relative_gas_permeability(std::vector<double>& gas_permeability) const ;
258 
265  virtual void get_liquid_gas_interfacial_surface(std::vector<double>& HI_liquid_gas_interfacial_surface) const ;
266 
272  virtual void get_pore_HI_wetted_wall_surface_area(std::vector<double>& HI_wetted_wall_surface_area) const ;
273 
278  virtual void get_wetted_wall_surface_area(std::vector<double>& wetted_wall_surface_area) const ;
279 
280 
287  virtual void get_pore_knudsen_radius_C1(std::vector<double>& knudsen_radius_C1) const ;
288 
295  virtual void get_pore_knudsen_radius_C3(std::vector<double>& knudsen_radius_C3) const ;
296 
301  virtual void get_knudsen_radius(std::vector<double>& knudsen_radius) const ;
302 
306  virtual void get_diffusivity() const ;
307 
308 
321  virtual void get_critical_radius(std::vector<double>& dst) const;
322 
323  virtual const double get_maximum_cross_sectional_areas() const;
324 
326 
327  protected:
328 
330 
331 
334  virtual boost::shared_ptr<FuelCellShop::MicroScale::BasePSD <dim>> create_replica (const std::string &psd_section_name)
335  {
336  return boost::shared_ptr<FuelCellShop::MicroScale::BasePSD <dim>> (new FuelCellShop::MicroScale::HIPSD<dim> (psd_section_name));
337  }
339 
341 
344  static HIPSD<dim> const* PROTOTYPE;
346 
348  // DATA //
350 
352 
353 
357  std::vector<double> fHI_k;
358 
362  std::vector<double> rHI_k;
363 
367  std::vector<double> sHI_k;
368 
372 
376 
378  double pressure_c;
379 
382  std::vector<double> critical_radius_computed;
383 
390 
392  std::vector<double> saturation_computed;
393 
399 
400 
401 
403  };
404 
405  } // PSD
406 
407 } // FuelCellShop
408 
409 #endif
Definition: system_management.h:77
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_HI.h:334
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_HI.h:156
static HIPSD< dim > const * PROTOTYPE
PROTOTYPE is the pointer is the dynamic pointer pointing to the HIPSD class itself.
Definition: PSD_HI.h:344
VariableNames get_variablename() const
Function to get the VariableNames enumeration corresponding to this struct.
Definition: fcst_variables.h:163
static const std::string concrete_name
Concrete name used for objects of this class.
Definition: PSD_HI.h:206
This structure is used to encapsulate data from constant values and variable solution data that is us...
Definition: fcst_variables.h:86
std::vector< double > rHI_k
The r_k is the characteristic pore size of the distribution k.
Definition: PSD_HI.h:362
std::vector< double > saturation_computed
Saturation_computed by the get_critical_radius function.
Definition: PSD_HI.h:392
void set_saturation()
Member function used to set the saturation at every quadrature point inside the cell.
Definition: PSD_HI.h:187
Definition: system_management.h:75
std::vector< double > critical_radius_computed
Critical_radius_computed by the get_critical_radius function.
Definition: PSD_HI.h:382
bool critical_radius_is_initialized
Check if the critical radius has already been computed by set_critical_radius function.
Definition: PSD_HI.h:389
void set_temperature(const SolutionVariable &T_in)
Member function used to set the temperature [Kelvin] at every quadrature point inside the cell...
Definition: PSD_HI.h:146
bool saturation_is_initialized
Check if the saturation has already been computed by set_saturation function.
Definition: PSD_HI.h:398
std::vector< double > sHI_k
The s_k is the spread of the distribution k.
Definition: PSD_HI.h:367
SolutionVariable T_vector
Temperature at every quadrature point inside the cell.
Definition: PSD_HI.h:371
Hydrophilic Pore Size Distribution.
Definition: PSD_HI.h:105
~HIPSD()
Destructor.
Definition: PSD_HI.h:129
std::vector< double > fHI_k
The f_k is the contribution of the log-normal distribution k to the total PSD.
Definition: PSD_HI.h:357
SolutionVariable Capillary_pressure_vector
Capillary pressure at every quadrature point inside the cell.
Definition: PSD_HI.h:375
Pore Size Distribution.
Definition: PSD_base.h:130
double pressure_c
Constant capillary pressure only for unit_test use.
Definition: PSD_HI.h:378
void set_critical_radius()
Member function used to set the critical radius [nm] at every quadrature point inside the cell...
Definition: PSD_HI.h:174