OpenFCST: The open-source Fuel Cell Simulation Toolbox
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
PSD_HO.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__HO__PSD_H
18 #define _FUELCELLSHOP__HO__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 
41 namespace FuelCellShop
42 {
43 
44 
45  namespace MicroScale
46  {
105  template <int dim>
106  class HOPSD : public BasePSD<dim>
107  {
108  public:
110 
112 
121  HOPSD();
122 
126  HOPSD (std::string name);
127 
131  virtual ~HOPSD();
132 
137  void initialize ( ParameterHandler &param );
138 
143  void declare_parameters (ParameterHandler &param) const;
144 
149  inline void set_temperature (const SolutionVariable& T_in)
150  {
151  Assert( T_in.get_variablename() == temperature_of_REV, ExcMessage("Wrong solution variable passed in PSD::set_temperature.") );
152  this->T_vector = T_in;
153  }
154 
159  inline void set_capillary_pressure (const SolutionVariable& C_in)
160  {
161  Assert( C_in.get_variablename() == capillary_pressure, ExcMessage("Wrong solution variable passed in PSD::capillary_pressure.") );
162  this->Capillary_pressure_vector = C_in;
163 
164  critical_radius_is_initialized = false;
165  saturation_is_initialized = false;
166  critical_radius_computed.clear();
167  saturation_computed.clear();
168  }
169 
176  inline void set_critical_radius()
177  {
178  get_critical_radius(critical_radius_computed);
179 
180  critical_radius_is_initialized = true;
181  }
182 
189  inline void set_saturation()
190  {
191  get_saturation(saturation_computed);
192 
193  saturation_is_initialized = true;
194  }
209  static const std::string concrete_name;
211 
213 
222  virtual void get_saturation(std::vector<double>& S) const ;
223 
230  virtual void get_global_saturated_permeability(double& saturated_permeability) const ;
231 
239  virtual void get_pore_HO_liquid_saturated_permeability(std::vector<double>& saturated_HO_permeability) const ;
240 
245  virtual void get_relative_liquid_permeability(std::vector<double>& liquid_permeability) const ;
246 
253  virtual void get_pore_HO_gas_saturated_permeability(std::vector<double>& saturated_HO_permeability) const ;
254 
259  virtual void get_relative_gas_permeability(std::vector<double>& gas_permeability) const ;
260 
267  virtual void get_liquid_gas_interfacial_surface(std::vector<double>& HO_liquid_gas_interfacial_surface) const ;
268 
275  virtual void get_pore_HO_wetted_wall_surface_area(std::vector<double>& HO_wetted_wall_surface_area) const ;
276 
281  virtual void get_wetted_wall_surface_area(std::vector<double>& wetted_wall_surface_area) const ;
282 
289  virtual void get_pore_knudsen_radius_C2(std::vector<double>& knudsen_radius_C2) const ;
290 
297  virtual void get_pore_knudsen_radius_C4(std::vector<double>& knudsen_radius_C4) const ;
298 
303  virtual void get_knudsen_radius(std::vector<double>& knudsen_radius) const ;
304 
308  virtual void get_diffusivity() const ;
309 
310 
323  virtual void get_critical_radius(std::vector<double>& dst) const;
324 
325  virtual const double get_maximum_cross_sectional_areas() const;
326 
327 
329 
330  protected:
331 
333 
334 
339  virtual boost::shared_ptr<FuelCellShop::MicroScale::BasePSD <dim>> create_replica (const std::string &psd_section_name)
340  {
341  return boost::shared_ptr<FuelCellShop::MicroScale::BasePSD <dim>> (new FuelCellShop::MicroScale::HOPSD<dim> (psd_section_name));
342  }
344 
346 
347 
350  static HOPSD<dim> const* PROTOTYPE;
352 
354  // DATA //
356 
358 
359 
363  std::vector<double> fHO_k;
364 
368  std::vector<double> rHO_k;
369 
373  std::vector<double> sHO_k;
374 
378 
382 
385  double pressure_c;
386 
389  std::vector<double> critical_radius_computed;
390 
397 
399  std::vector<double> saturation_computed;
400 
406 
408  };
409 
410  } // PSD
411 
412 } // FuelCellShop
413 
414 #endif
static const std::string concrete_name
Concrete name used for objects of this class.
Definition: PSD_HO.h:209
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_HO.h:339
std::vector< double > rHO_k
The r_k is the characteristic pore size of the distribution k.
Definition: PSD_HO.h:368
Definition: system_management.h:77
Hydrophobic Pore Size Distribution.
Definition: PSD_HO.h:106
std::vector< double > critical_radius_computed
Critical_radius_computed by the get_critical_radius function.
Definition: PSD_HO.h:389
std::vector< double > sHO_k
The s_k is the spread of the distribution k.
Definition: PSD_HO.h:373
VariableNames get_variablename() const
Function to get the VariableNames enumeration corresponding to this struct.
Definition: fcst_variables.h:163
void set_saturation()
Member function used to set the saturation at every quadrature point inside the cell.
Definition: PSD_HO.h:189
This structure is used to encapsulate data from constant values and variable solution data that is us...
Definition: fcst_variables.h:86
bool critical_radius_is_initialized
Check if the critical radius has already been computed by set_critical_radius function.
Definition: PSD_HO.h:396
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_HO.h:159
static HOPSD< dim > const * PROTOTYPE
PROTOTYPE is the pointer is the dynamic pointer pointing to the HOPSD class itself.
Definition: PSD_HO.h:350
void set_temperature(const SolutionVariable &T_in)
Member function used to set the temperature [Kelvin] at every quadrature point inside the cell...
Definition: PSD_HO.h:149
SolutionVariable Capillary_pressure_vector
Capillary pressure at every quadrature point inside the cell.
Definition: PSD_HO.h:381
Definition: system_management.h:75
SolutionVariable T_vector
Temperature at every quadrature point inside the cell.
Definition: PSD_HO.h:377
std::vector< double > saturation_computed
Saturation_computed by the get_critical_radius function.
Definition: PSD_HO.h:399
void set_critical_radius()
Member function used to set the critical radius [nm] at every quadrature point inside the cell...
Definition: PSD_HO.h:176
bool saturation_is_initialized
Check if the saturation has already been computed by set_saturation function.
Definition: PSD_HO.h:405
Pore Size Distribution.
Definition: PSD_base.h:130
std::vector< double > fHO_k
The f_k is the contribution of the log-normal distribution k to the total PSD.
Definition: PSD_HO.h:363
double pressure_c
Constant capillary pressure only for unit_test use.
Definition: PSD_HO.h:385