Saturated.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2015-2018 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  const dimensionedScalar Wi
36  (
37  "W",
39  this->thermo_.composition().Wi(saturatedIndex_)
40  );
41 
42  return Wi/this->thermo_.W()/this->thermo_.p();
43 }
44 
45 
46 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
47 
48 template<class Thermo, class OtherThermo>
50 (
51  const dictionary& dict,
52  const phasePair& pair
53 )
54 :
55  InterfaceCompositionModel<Thermo, OtherThermo>(dict, pair),
56  saturatedName_(this->speciesNames_[0]),
57  saturatedIndex_
58  (
59  this->thermo_.composition().species()[saturatedName_]
60  ),
61  saturationModel_
62  (
63  saturationModel::New
64  (
65  dict.subDict("saturationPressure"),
66  pair.phase1().mesh()
67  )
68  )
69 {
70  if (this->speciesNames_.size() != 1)
71  {
73  << "Saturated model is suitable for one species only."
74  << exit(FatalError);
75  }
76 }
77 
78 
79 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
80 
81 template<class Thermo, class OtherThermo>
83 {}
84 
85 
86 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
87 
88 template<class Thermo, class OtherThermo>
89 void
91 (
92  const volScalarField& Tf
93 )
94 {}
95 
96 
97 template<class Thermo, class OtherThermo>
100 (
101  const word& speciesName,
102  const volScalarField& Tf
103 ) const
104 {
105  if (saturatedName_ == speciesName)
106  {
107  return wRatioByP()*saturationModel_->pSat(Tf);
108  }
109  else
110  {
111  const label speciesIndex
112  (
113  this->thermo_.composition().species()[speciesName]
114  );
115 
116  return
117  this->thermo_.Y()[speciesIndex]
118  *(scalar(1) - wRatioByP()*saturationModel_->pSat(Tf))
119  /max(scalar(1) - this->thermo_.Y()[saturatedIndex_], small);
120  }
121 }
122 
123 
124 template<class Thermo, class OtherThermo>
127 (
128  const word& speciesName,
129  const volScalarField& Tf
130 ) const
131 {
132  if (saturatedName_ == speciesName)
133  {
134  return wRatioByP()*saturationModel_->pSatPrime(Tf);
135  }
136  else
137  {
138  const label speciesIndex
139  (
140  this->thermo_.composition().species()[speciesName]
141  );
142 
143  return
144  - this->thermo_.Y()[speciesIndex]
145  *wRatioByP()*saturationModel_->pSatPrime(Tf)
146  /max(scalar(1) - this->thermo_.Y()[saturatedIndex_], small);
147  }
148 }
149 
150 
151 // ************************************************************************* //
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 tmp< volScalarField > Yf(const word &speciesName, const volScalarField &Tf) const
The interface species fraction.
virtual void update(const volScalarField &Tf)
Update the composition.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
virtual tmp< volScalarField > YfPrime(const word &speciesName, const volScalarField &Tf) const
The interface species fraction derivative w.r.t. temperature.
basicSpecieMixture & composition
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:163
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
dynamicFvMesh & mesh
Saturated(const dictionary &dict, const phasePair &pair)
Construct from components.
tmp< volScalarField > wRatioByP() const
Constant of proportionality between partial pressure and mass.
phaseModel & phase1
const dimensionSet dimMoles(0, 0, 0, 0, 1, 0, 0)
Definition: dimensionSets.H:53
const hashedWordList speciesNames_
Names of the transferring species.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
A class for managing temporary objects.
Definition: PtrList.H:53