OpenFCST: The open-source Fuel Cell Simulation Toolbox
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
GasMixture.h
Go to the documentation of this file.
1 // ----------------------------------------------------------------------------
2 //
3 // FCST: Fuel Cell Simulation Toolbox
4 //
5 // Copyright (C) 2006-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: GasMixture.h
11 // - Description: This class describes properties of gas mixtures
12 // - Developers: Valentin N. Zingan, University of Alberta
13 //
14 // ----------------------------------------------------------------------------
15 
16 #ifndef _FCST_FUELCELLSHOP_MATERIAL_GASMIXTURE_H_
17 #define _FCST_FUELCELLSHOP_MATERIAL_GASMIXTURE_H_
18 
19 #include <materials/PureGas.h>
20 
21 namespace FuelCellShop
22 {
23  namespace Material
24  {
25 
102  class GasMixture : public BaseMaterial
103  {
104  public:
105 
107 
108 
111  GasMixture();
115  GasMixture(const std::string& name);
116 
120  virtual ~GasMixture();
121 
125  virtual void declare_parameters(ParameterHandler& param) const;
126 
130  virtual void initialize(ParameterHandler& param);
131 
136  void set_gases(const std::vector< PureGas* >& rgases)
137  {
138  gases.clear();
139  gases = rgases;
140  }
141 
147  void set_total_pressure(const double& rtotal_pressure)
148  {
149  total_pressure = rtotal_pressure;
150  }
151 
157  void set_temperature(const double& rtemperature)
158  {
159  temperature = rtemperature;
160  }
161 
163 
165 
166 
169  inline unsigned int n_gases() const
170  {
171  return gases.size();
172  }
173 
177  inline PureGas* get_gas(unsigned int& ind) const
178  {
179  return gases[ind];
180  }
181 
186  const std::vector< PureGas* >& get_gases() const
187  {
188  return gases;
189  }
190 
195  const double& get_total_pressure() const
196  {
197  return total_pressure;
198  }
199 
204  const double& get_temperature() const
205  {
206  return temperature;
207  }
208 
213  virtual void print_material_properties() const;
214 
216 
218 
219 
226 
236  void get_ChapmanEnskog_isobaric_diffusion_coefficient(std::vector<double>& diffusion_coefficient) const;
237 
245  const double get_ChapmanEnskog_isobaric_diffusion_coefficient(const double& temperature) const;
246 
257  void get_ChapmanEnskog_isobaric_diffusion_coefficient(const std::vector<double>& temperature,
258  std::vector<double>& diffusion_coefficient) const;
259 
261 
263 
264 
274 
287  std::vector<double>& dst) const;
288 
290 
292 
293 
299  const double get_ChapmanEnskog_diffusion_coefficient() const;
300 
310  void get_ChapmanEnskog_diffusion_coefficient(std::vector<double>& diffusion_coefficient) const;
311 
320 
332  std::vector<double>& diffusion_coefficient) const;
333 
342 
354  std::vector<double>& diffusion_coefficient) const;
355 
364  const double get_ChapmanEnskog_diffusion_coefficient(const double& total_pressure,
365  const double& temperature) const;
366 
378  void get_ChapmanEnskog_diffusion_coefficient(const std::vector<double>& total_pressure,
379  const std::vector<double>& temperature,
380  std::vector<double>& diffusion_coefficient) const;
381 
383 
385 
386 
395  const double get_DChapmanEnskog_diffusion_coefficient_Dpressure(const double& total_pressure) const;
396 
409  std::vector<double>& dst) const;
410 
421  const double& temperature) const;
422 
436  const std::vector<double>& temperature,
437  std::vector<double>& dst) const;
438 
447  const double get_DChapmanEnskog_diffusion_coefficient_Dtemperature(const double& temperature) const;
448 
461  std::vector<double>& dst) const;
462 
473  const double& temperature) const;
474 
488  const std::vector<double>& temperature,
489  std::vector<double>& dst) const;
490 
492 
494 
495 
501  const Table< 2, double > get_ChapmanEnskog_isobaric_diffusion_coefficients() const;
502 
512  void get_ChapmanEnskog_isobaric_diffusion_coefficients(std::vector< Table< 2, double > >& diffusion_coefficients) const;
513 
521  const Table< 2, double > get_ChapmanEnskog_isobaric_diffusion_coefficients(const double& temperature) const;
522 
533  void get_ChapmanEnskog_isobaric_diffusion_coefficients(const std::vector<double>& temperature,
534  std::vector< Table< 2, double > >& diffusion_coefficients) const;
535 
537 
539 
540 
549  const Table< 2, double > get_DChapmanEnskog_isobaric_diffusion_coefficients_Dtemperature(const double& temperature) const;
550 
563  std::vector< Table< 2, double > >& dst) const;
564 
566 
568 
569 
575  const Table< 2, double > get_ChapmanEnskog_diffusion_coefficients() const;
576 
586  void get_ChapmanEnskog_diffusion_coefficients(std::vector< Table< 2, double > >& diffusion_coefficients) const;
587 
595  const Table< 2, double > get_ChapmanEnskog_diffusion_coefficients_at_constant_pressure(const double& temperature) const;
596 
608  std::vector< Table< 2, double > >& diffusion_coefficients) const;
609 
617  const Table< 2, double > get_ChapmanEnskog_diffusion_coefficients_at_constant_temperature(const double& total_pressure) const;
618 
630  std::vector< Table< 2, double > >& diffusion_coefficients) const;
631 
640  const Table< 2, double > get_ChapmanEnskog_diffusion_coefficients(const double& total_pressure,
641  const double& temperature) const;
642 
654  void get_ChapmanEnskog_diffusion_coefficients(const std::vector<double>& total_pressure,
655  const std::vector<double>& temperature,
656  std::vector< Table< 2, double > >& diffusion_coefficients) const;
657 
659 
661 
662 
671  const Table< 2, double > get_DChapmanEnskog_diffusion_coefficients_Dpressure(const double& total_pressure) const;
672 
685  std::vector< Table< 2, double > >& dst) const;
686 
696  const Table< 2, double > get_DChapmanEnskog_diffusion_coefficients_Dpressure(const double& total_pressure,
697  const double& temperature) const;
698 
712  const std::vector<double>& temperature,
713  std::vector< Table< 2, double > >& dst) const;
714 
723  const Table< 2, double > get_DChapmanEnskog_diffusion_coefficients_Dtemperature(const double& temperature) const;
724 
737  std::vector< Table< 2, double > >& dst) const;
738 
748  const Table< 2, double > get_DChapmanEnskog_diffusion_coefficients_Dtemperature(const double& total_pressure,
749  const double& temperature) const;
750 
764  const std::vector<double>& temperature,
765  std::vector< Table< 2, double > >& dst) const;
767 
769 
770 
791  const std::vector< std::vector<double> > get_isothermal_nonisobaric_partial_viscosity(const double& tempOfMixture,
792  const std::vector< std::vector<double> >& density,
793  const std::vector<double>& molarMass,
794  const std::vector<double>& dynamicViscosity,
795  const std::vector<double>& collisionDiameter,
796  const std::vector<double>& porosity,
797  std::vector< std::vector< std::vector<double> > >& paramMatrix,
798  std::vector< FullMatrix<double> >& PInv) const;
799 
820  const std::vector< std::vector<double> > get_nonisothermal_nonisobaric_partial_viscosity(const std::vector<double>& temperature,
821  const std::vector< std::vector<double> >& density,
822  const std::vector<double>& molarMass,
823  const std::vector< std::vector<double> >& dynamicViscosity,
824  const std::vector<double>& collisionDiameter,
825  const std::vector<double>& porosity,
826  std::vector< std::vector< std::vector<double> > >& paramMatrix,
827  std::vector< FullMatrix<double> >& PInv) const;
828 
849  std::vector<double> get_isothermal_nonisobaric_Wilke_partial_viscosity(const std::vector<double>& density,
850  const std::vector<double>& dynamicViscosity,
851  const std::vector<double>& molarMass,
852  const double& porosity) const;
853 
879  std::vector<double> get_isothermal_nonisobaric_Wilke_partial_viscosity(const std::vector<double>& density,
880  const std::vector<double>& dynamicViscosity,
881  const std::vector<double>& molarMass,
882  const double& porosity,
883  std::vector< std::vector<double> >& xi) const;
884 
889  std::vector<double> get_nonisothermal_nonisobaric_Wilke_partial_viscosity(const std::vector<double>& density,
890  const std::vector<double>& dynamicViscosity,
891  const std::vector<double>& molarMass,
892  const double& porosity,
893  std::vector< std::vector<double> >& xi) const;
894 
948  std::vector<double> get_isothermal_nonisobaric_OmegaKG_partial_viscosity(const std::vector<double>& density,
949  const std::vector<double>& collisionDiameter,
950  const std::vector<double>& molarMass,
951  const double& temperature,
952  const double& porosity) const;
953 
1013  std::vector<double> get_isothermal_nonisobaric_OmegaKG_partial_viscosity(const std::vector<double>& density,
1014  const std::vector<double>& collisionDiameter,
1015  const std::vector<double>& molarMass,
1016  const double& temperature,
1017  const double& porosity,
1018  std::vector< std::vector<double> >& omegaIntegralTable,
1019  FullMatrix<double>& PInv) const;
1024  std::vector<double> get_nonisothermal_nonisobaric_OmegaKG_partial_viscosity(const std::vector<double>& density,
1025  const std::vector<double>& collisionDiameter,
1026  const std::vector<double>& molarMass,
1027  const double& temperature,
1028  const double& porosity,
1029  std::vector< std::vector<double> >& omegaIntegralTable,
1030  FullMatrix<double>& PInv) const;
1031 
1050  void get_isothermal_nonisobaric_delta_partial_viscosity(const double& tempOfMixture,
1051  const std::vector< std::vector< std::vector<double> > >& paramMatrix,
1052  const std::vector< FullMatrix<double> >& PInv,
1053  const std::vector<double>& porosity,
1054  const std::vector<double>& molarMass,
1055  const std::vector<double>& dynamicViscosity,
1056  const std::vector<double>& collisionDiameter,
1057  const std::vector< std::vector< std::vector<double> > >& deltaDensity,
1058  const std::vector< std::vector<double> >& density,
1059  std::vector< std::vector< std::vector<double> > >& deltaPartialViscosity) const;
1060 
1079  const std::vector< std::vector< std::vector<double> > >& paramMatrix,
1080  const std::vector< FullMatrix<double> >& PInv,
1081  const std::vector<double>& porosity,
1082  const std::vector<double>& molarMass,
1083  const std::vector< std::vector<double> >& dynamicViscosity,
1084  const std::vector<double>& collisionDiameter,
1085  const std::vector< std::vector< std::vector<double> > >& deltaDensity,
1086  const std::vector< std::vector<double> >& density,
1087  std::vector< std::vector< std::vector<double> > >& deltaPartialViscosity) const;
1088 
1114  void get_Wilke_delta_partial_viscosity_wrt_density(const std::vector< std::vector< std::vector<double> > >& xi,
1115  const std::vector<double>& porosity,
1116  const std::vector<double>& molarMass,
1117  const std::vector<double>& dynamicViscosity,
1118  const std::vector< std::vector< std::vector<double> > >& deltaDensity,
1119  const std::vector< std::vector<double> >& density,
1120  std::vector< std::vector< std::vector<double> > >& deltaPartialViscosity) const;
1121 
1157  void get_Wilke_delta_partial_viscosity_wrt_temperature(const std::vector< std::vector< std::vector<double> > >& xi,
1158  const std::vector<double>& porosity,
1159  const std::vector<double>& molarMass,
1160  const std::vector<double>& dynamicViscosity,
1161  const std::vector< std::vector< std::vector<double> > >& deltaTemperature,
1162  const std::vector< std::vector<double> >& density,
1163  std::vector< std::vector< std::vector<double> > >& deltaPartialViscosity) const;
1164 
1198  void get_OmegaKG_delta_partial_viscosity_wrt_density(const double& tempOfMixture,
1199  const std::vector< std::vector< std::vector<double> > >& omegaIntegralTable,
1200  const std::vector< FullMatrix<double> >& PInv,
1201  const std::vector<double>& porosity,
1202  const std::vector<double>& molarMass,
1203  const std::vector<double>& collisionDiameter,
1204  const std::vector< std::vector< std::vector<double> > >& deltaDensity,
1205  const std::vector< std::vector<double> >& density,
1206  std::vector< std::vector< std::vector<double> > >& deltaPartialViscosity) const;
1207 
1234  void get_OmegaKG_delta_partial_viscosity_wrt_temperature(const std::vector<double>& temperature,
1235  const std::vector< std::vector< std::vector<double> > >& omegaIntegralTable,
1236  const std::vector< FullMatrix<double> >& PInv,
1237  const std::vector<double>& porosity,
1238  const std::vector<double>& molarMass,
1239  const std::vector<double>& collisionDiameter,
1240  const std::vector< std::vector< std::vector<double> > >& deltaDensity,
1241  const std::vector< std::vector<double> >& density,
1242  std::vector< std::vector< std::vector<double> > >& deltaPartialViscosity) const;
1243 
1249  void get_Null_delta_viscosity(std::vector< std::vector< std::vector<double> > >& deltaViscosity) const;
1250 
1255  void get_delta_bulk_viscosity(const std::vector< PureGas* >& gases,
1256  const std::vector< std::vector< std::vector<double> > >& deltaPartialViscosity,
1257  std::vector< std::vector< std::vector<double> > >& deltaBulkViscosity) const;
1258 
1260 
1261  protected:
1262 
1264 
1265 
1278  const double get_binary_collision_integral(const unsigned int& N1 = 0,
1279  const unsigned int& N2 = 1) const;
1280 
1296  void get_binary_collision_integral(std::vector<double>& binary_collision_integral,
1297  const unsigned int& N1 = 0,
1298  const unsigned int& N2 = 1) const;
1299 
1313  const double get_binary_collision_integral(const double& temperature,
1314  const unsigned int& N1 = 0,
1315  const unsigned int& N2 = 1) const;
1316 
1333  void get_binary_collision_integral(const std::vector<double>& temperature,
1334  std::vector<double>& binary_collision_integral,
1335  const unsigned int& N1 = 0,
1336  const unsigned int& N2 = 1) const;
1337 
1347  const double get_Dbinary_collision_integral_Dtemperature(const double& temperature,
1348  const unsigned int& N1 = 0,
1349  const unsigned int& N2 = 1) const;
1350 
1362  void get_Dbinary_collision_integral_Dtemperature(const std::vector<double>& temperature,
1363  std::vector<double>& dst,
1364  const unsigned int& N1 = 0,
1365  const unsigned int& N2 = 1) const;
1366 
1368 
1370 
1371 
1391  const double get_omega_star_11_integral(const double& temperature,
1392  const unsigned int& N1 = 0,
1393  const unsigned int& N2 = 1) const;
1394 
1416  void get_omega_star_11_integral(const double & temperature,
1417  std::vector<double>& omega_integral,
1418  const unsigned int& N1 = 0,
1419  const unsigned int& N2 = 1) const;
1420 
1440  const double get_omega_star_22_integral(const double& temperature,
1441  const unsigned int& N1 = 0,
1442  const unsigned int& N2 = 1) const;
1443 
1465  void get_omega_star_22_integral(const double& temperature,
1466  std::vector<double>& omega_integral,
1467  const unsigned int& N1 = 0,
1468  const unsigned int& N2 = 1) const;
1469 
1471 
1473  // DATA //
1475 
1477 
1478 
1482  std::vector< PureGas* > gases;
1483 
1487  bool tempIsoTherm = false;
1488 
1492  bool pressIsoBaric = false;
1493 
1498  double total_pressure = 1.0;
1499 
1504  double temperature = 1.0;
1505 
1507 
1509 
1510 
1517 
1519 
1520  };
1521 
1522  } // Material
1523 
1524 } // FuelCellShop
1525 
1526 #endif
std::string mixture_viscosity_mode
Definition: GasMixture.h:1516
const double & get_temperature() const
This function returns temperature.
Definition: GasMixture.h:204
void get_nonisothermal_nonisobaric_delta_partial_viscosity(const std::vector< double > &temperature, const std::vector< std::vector< std::vector< double > > > &paramMatrix, const std::vector< FullMatrix< double > > &PInv, const std::vector< double > &porosity, const std::vector< double > &molarMass, const std::vector< std::vector< double > > &dynamicViscosity, const std::vector< double > &collisionDiameter, const std::vector< std::vector< std::vector< double > > > &deltaDensity, const std::vector< std::vector< double > > &density, std::vector< std::vector< std::vector< double > > > &deltaPartialViscosity) const
This function calculates the variation in partial viscosity based on value stored in mixture_viscosit...
void get_isothermal_nonisobaric_delta_partial_viscosity(const double &tempOfMixture, const std::vector< std::vector< std::vector< double > > > &paramMatrix, const std::vector< FullMatrix< double > > &PInv, const std::vector< double > &porosity, const std::vector< double > &molarMass, const std::vector< double > &dynamicViscosity, const std::vector< double > &collisionDiameter, const std::vector< std::vector< std::vector< double > > > &deltaDensity, const std::vector< std::vector< double > > &density, std::vector< std::vector< std::vector< double > > > &deltaPartialViscosity) const
This function calculates the variation in partial viscosity based on value stored in mixture_viscosit...
double temperature
Temperature of the whole gas mixture, .
Definition: GasMixture.h:1504
const double get_DChapmanEnskog_isobaric_diffusion_coefficient_Dtemperature(const double &temperature) const
This function returns the first derivative of the Maxwell-Stefan isobaric diffusion coefficient of g...
virtual void print_material_properties() const
This function prints out the material properties.
const double get_binary_collision_integral(const unsigned int &N1=0, const unsigned int &N2=1) const
This function returns binary collision integral at a constant temperature.
const std::string name
Name of the layer.
Definition: base_material.h:155
const std::vector< PureGas * > & get_gases() const
This function returns gases.
Definition: GasMixture.h:186
const double get_DChapmanEnskog_diffusion_coefficient_Dpressure(const double &total_pressure) const
This function returns the first derivative of the Maxwell-Stefan diffusion coefficient of gas in ga...
const double get_omega_star_11_integral(const double &temperature, const unsigned int &N1=0, const unsigned int &N2=1) const
This function returns integral at a variable temperature.
std::vector< double > get_isothermal_nonisobaric_OmegaKG_partial_viscosity(const std::vector< double > &density, const std::vector< double > &collisionDiameter, const std::vector< double > &molarMass, const double &temperature, const double &porosity) const
This is an overloaded function that calculates the Kerkhof and Geboers partial viscosity model based ...
const Table< 2, double > get_DChapmanEnskog_diffusion_coefficients_Dtemperature(const double &temperature) const
This function returns the first derivative of the Maxwell-Stefan diffusion coefficients of gas in g...
const double get_DChapmanEnskog_diffusion_coefficient_Dtemperature(const double &temperature) const
This function returns the first derivative of the Maxwell-Stefan diffusion coefficient of gas in ga...
std::vector< double > get_nonisothermal_nonisobaric_OmegaKG_partial_viscosity(const std::vector< double > &density, const std::vector< double > &collisionDiameter, const std::vector< double > &molarMass, const double &temperature, const double &porosity, std::vector< std::vector< double > > &omegaIntegralTable, FullMatrix< double > &PInv) const
This function calculates the partial viscosity using the OmegaKG model based on value stored in mixtu...
virtual void initialize(ParameterHandler &param)
Initialize parameters.
const std::vector< std::vector< double > > get_isothermal_nonisobaric_partial_viscosity(const double &tempOfMixture, const std::vector< std::vector< double > > &density, const std::vector< double > &molarMass, const std::vector< double > &dynamicViscosity, const std::vector< double > &collisionDiameter, const std::vector< double > &porosity, std::vector< std::vector< std::vector< double > > > &paramMatrix, std::vector< FullMatrix< double > > &PInv) const
This function calculates the partial viscosity based on value stored in mixture_viscosity_mode and AS...
const std::vector< std::vector< double > > get_nonisothermal_nonisobaric_partial_viscosity(const std::vector< double > &temperature, const std::vector< std::vector< double > > &density, const std::vector< double > &molarMass, const std::vector< std::vector< double > > &dynamicViscosity, const std::vector< double > &collisionDiameter, const std::vector< double > &porosity, std::vector< std::vector< std::vector< double > > > &paramMatrix, std::vector< FullMatrix< double > > &PInv) const
This function calculates the partial viscosity based on value stored in mixture_viscosity_mode and AS...
void get_OmegaKG_delta_partial_viscosity_wrt_temperature(const std::vector< double > &temperature, const std::vector< std::vector< std::vector< double > > > &omegaIntegralTable, const std::vector< FullMatrix< double > > &PInv, const std::vector< double > &porosity, const std::vector< double > &molarMass, const std::vector< double > &collisionDiameter, const std::vector< std::vector< std::vector< double > > > &deltaDensity, const std::vector< std::vector< double > > &density, std::vector< std::vector< std::vector< double > > > &deltaPartialViscosity) const
Calculates the variation in Kerkhof and Geboers partial viscosity model using Omega integrals w...
void get_Null_delta_viscosity(std::vector< std::vector< std::vector< double > > > &deltaViscosity) const
Returns a std::vector&lt; std::vector&lt; std::vector&lt;double&gt; &gt; &gt; of 0s.
const Table< 2, double > get_ChapmanEnskog_diffusion_coefficients_at_constant_temperature(const double &total_pressure) const
This function returns Maxwell-Stefan diffusion coefficients of gas in gas written in the Chapman En...
std::vector< double > get_isothermal_nonisobaric_Wilke_partial_viscosity(const std::vector< double > &density, const std::vector< double > &dynamicViscosity, const std::vector< double > &molarMass, const double &porosity) const
This is an overloaded function that calculates the Wilke partial viscosity model and ASSUMES ISOTHERM...
const double get_omega_star_22_integral(const double &temperature, const unsigned int &N1=0, const unsigned int &N2=1) const
This function returns integral at a variable temperature.
std::vector< PureGas * > gases
This std::vector contains all pure gases which form the whole gas mixture of a problem at hand...
Definition: GasMixture.h:1482
double total_pressure
Total pressure of the whole mixture, .
Definition: GasMixture.h:1498
const Table< 2, double > get_ChapmanEnskog_diffusion_coefficients() const
This function returns Maxwell-Stefan diffusion coefficients of gas in gas written in the Chapman En...
virtual void declare_parameters(ParameterHandler &param) const
Declare parameters.
const double get_ChapmanEnskog_diffusion_coefficient_at_constant_temperature(const double &total_pressure) const
This function returns Maxwell-Stefan diffusion coefficient of gas in gas (or vice-versa) written in...
const Table< 2, double > get_DChapmanEnskog_isobaric_diffusion_coefficients_Dtemperature(const double &temperature) const
This function returns the first derivative of the Maxwell-Stefan isobaric diffusion coefficients of ...
const double get_ChapmanEnskog_isobaric_diffusion_coefficient() const
This function returns Maxwell-Stefan isobaric diffusion coefficient of gas in gas (or vice-versa) w...
std::vector< double > get_nonisothermal_nonisobaric_Wilke_partial_viscosity(const std::vector< double > &density, const std::vector< double > &dynamicViscosity, const std::vector< double > &molarMass, const double &porosity, std::vector< std::vector< double > > &xi) const
This function calculates the partial viscosity using the Wilke model based on value stored in mixture...
virtual ~GasMixture()
Destructor.
const double get_ChapmanEnskog_diffusion_coefficient() const
This function returns Maxwell-Stefan diffusion coefficient of gas in gas (or vice-versa) written in...
const Table< 2, double > get_DChapmanEnskog_diffusion_coefficients_Dpressure(const double &total_pressure) const
This function returns the first derivative of the Maxwell-Stefan diffusion coefficients of gas in g...
void get_Wilke_delta_partial_viscosity_wrt_temperature(const std::vector< std::vector< std::vector< double > > > &xi, const std::vector< double > &porosity, const std::vector< double > &molarMass, const std::vector< double > &dynamicViscosity, const std::vector< std::vector< std::vector< double > > > &deltaTemperature, const std::vector< std::vector< double > > &density, std::vector< std::vector< std::vector< double > > > &deltaPartialViscosity) const
Calculates the variation in Wilke&#39;s partial viscosity model w.r.t.
void set_total_pressure(const double &rtotal_pressure)
This function takes a pressure (Pa) and if Isobaric fluid flow is set to true in the data file then t...
Definition: GasMixture.h:147
This class describes properties of gas mixtures.
Definition: GasMixture.h:102
This class is a base class for all pure gases used in OpenFCST.
Definition: PureGas.h:88
unsigned int n_gases() const
Function returning the number of gases in the mixture.
Definition: GasMixture.h:169
bool pressIsoBaric
This bool tells GasMisture if it should use isobaric assumption.
Definition: GasMixture.h:1492
const double & get_total_pressure() const
This function returns total_pressure.
Definition: GasMixture.h:195
bool tempIsoTherm
This bool tells GasMisture if it should use isothermal assumption.
Definition: GasMixture.h:1487
const double get_ChapmanEnskog_diffusion_coefficient_at_constant_pressure(const double &temperature) const
This function returns Maxwell-Stefan diffusion coefficient of gas in gas (or vice-versa) written in...
const Table< 2, double > get_ChapmanEnskog_isobaric_diffusion_coefficients() const
This function returns Maxwell-Stefan isobaric diffusion coefficients of gas in gas written in the C...
const double get_Dbinary_collision_integral_Dtemperature(const double &temperature, const unsigned int &N1=0, const unsigned int &N2=1) const
This function returns the first derivative of the binary collision integral at a variable temperatur...
PureGas * get_gas(unsigned int &ind) const
Return gas stored in index ind.
Definition: GasMixture.h:177
void get_OmegaKG_delta_partial_viscosity_wrt_density(const double &tempOfMixture, const std::vector< std::vector< std::vector< double > > > &omegaIntegralTable, const std::vector< FullMatrix< double > > &PInv, const std::vector< double > &porosity, const std::vector< double > &molarMass, const std::vector< double > &collisionDiameter, const std::vector< std::vector< std::vector< double > > > &deltaDensity, const std::vector< std::vector< double > > &density, std::vector< std::vector< std::vector< double > > > &deltaPartialViscosity) const
Calculates the variation in Kerkhof and Geboers partial viscosity model using Omega integrals w...
void set_gases(const std::vector< PureGas * > &rgases)
This function sets gases.
Definition: GasMixture.h:136
void get_delta_bulk_viscosity(const std::vector< PureGas * > &gases, const std::vector< std::vector< std::vector< double > > > &deltaPartialViscosity, std::vector< std::vector< std::vector< double > > > &deltaBulkViscosity) const
Loops through all variations of partial viscosity and calculates the varation in bulk viscosity accor...
void set_temperature(const double &rtemperature)
This function takes a temperature (K) and if Isothermal fluid flow is set to true in the data file th...
Definition: GasMixture.h:157
Virtual class used to provide the interface for all material classes.
Definition: base_material.h:54
void get_Wilke_delta_partial_viscosity_wrt_density(const std::vector< std::vector< std::vector< double > > > &xi, const std::vector< double > &porosity, const std::vector< double > &molarMass, const std::vector< double > &dynamicViscosity, const std::vector< std::vector< std::vector< double > > > &deltaDensity, const std::vector< std::vector< double > > &density, std::vector< std::vector< std::vector< double > > > &deltaPartialViscosity) const
Calculates the variation in Wilke&#39;s partial viscosity model w.r.t.
const Table< 2, double > get_ChapmanEnskog_diffusion_coefficients_at_constant_pressure(const double &temperature) const
This function returns Maxwell-Stefan diffusion coefficients of gas in gas written in the Chapman En...