Saturated.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2015 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "Saturated.H"
27 
28 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
29 
30 template<class Thermo, class OtherThermo>
33 wRatioByP() const
34 {
35  return
36  this->thermo_.composition().W(saturatedIndex_)
37  /this->thermo_.composition().W()
38  /this->thermo_.p();
39 }
40 
41 
42 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
43 
44 template<class Thermo, class OtherThermo>
46 (
47  const dictionary& dict,
48  const phasePair& pair
49 )
50 :
51  InterfaceCompositionModel<Thermo, OtherThermo>(dict, pair),
52  saturatedName_(this->speciesNames_[0]),
53  saturatedIndex_
54  (
55  this->thermo_.composition().species()[saturatedName_]
56  ),
57  saturationModel_
58  (
59  saturationModel::New
60  (
61  dict.subDict("saturationPressure")
62  )
63  )
64 {
65  if (this->speciesNames_.size() != 1)
66  {
68  << "Saturated model is suitable for one species only."
69  << exit(FatalError);
70  }
71 }
72 
73 
74 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
75 
76 template<class Thermo, class OtherThermo>
78 {}
79 
80 
81 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
82 
83 template<class Thermo, class OtherThermo>
84 void
86 (
87  const volScalarField& Tf
88 )
89 {
90  // do nothing
91 }
92 
93 
94 template<class Thermo, class OtherThermo>
97 (
98  const word& speciesName,
99  const volScalarField& Tf
100 ) const
101 {
102  if (saturatedName_ == speciesName)
103  {
104  return wRatioByP()*saturationModel_->pSat(Tf);
105  }
106  else
107  {
108  const label speciesIndex
109  (
110  this->thermo_.composition().species()[speciesName]
111  );
112 
113  return
114  this->thermo_.Y()[speciesIndex]
115  *(scalar(1) - wRatioByP()*saturationModel_->pSat(Tf))
116  /max(scalar(1) - this->thermo_.Y()[saturatedIndex_], SMALL);
117  }
118 }
119 
120 
121 template<class Thermo, class OtherThermo>
124 (
125  const word& speciesName,
126  const volScalarField& Tf
127 ) const
128 {
129  if (saturatedName_ == speciesName)
130  {
131  return wRatioByP()*saturationModel_->pSatPrime(Tf);
132  }
133  else
134  {
135  const label speciesIndex
136  (
137  this->thermo_.composition().species()[speciesName]
138  );
139 
140  return
141  - this->thermo_.Y()[speciesIndex]
142  *wRatioByP()*saturationModel_->pSatPrime(Tf)
143  /max(scalar(1) - this->thermo_.Y()[saturatedIndex_], SMALL);
144  }
145 }
146 
147 
148 // ************************************************************************* //
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
virtual void update(const volScalarField &Tf)
Update the composition.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
basicMultiComponentMixture & composition
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
virtual tmp< volScalarField > Yf(const word &speciesName, const volScalarField &Tf) const
The interface species fraction.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
virtual tmp< volScalarField > YfPrime(const word &speciesName, const volScalarField &Tf) const
The interface species fraction derivative w.r.t. temperature.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
Saturated(const dictionary &dict, const phasePair &pair)
Construct from components.
const hashedWordList speciesNames_
Names of the transferring species.
tmp< volScalarField > wRatioByP() const
Constant of propotionality between partial pressure and mass.
A class for managing temporary objects.
Definition: PtrList.H:54