NonRandomTwoLiquid.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 "NonRandomTwoLiquid.H"
27 
28 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
29 
30 template<class Thermo, class OtherThermo>
33 (
34  const dictionary& dict,
35  const phasePair& pair
36 )
37 :
38  InterfaceCompositionModel<Thermo, OtherThermo>(dict, pair),
39  gamma1_
40  (
41  IOobject
42  (
43  IOobject::groupName("gamma1", pair.name()),
44  pair.phase1().mesh().time().timeName(),
45  pair.phase1().mesh()
46  ),
47  pair.phase1().mesh(),
49  ),
50  gamma2_
51  (
52  IOobject
53  (
54  IOobject::groupName("gamma2", pair.name()),
55  pair.phase1().mesh().time().timeName(),
56  pair.phase1().mesh()
57  ),
58  pair.phase1().mesh(),
60  ),
61  beta12_("", dimless/dimTemperature, 0),
62  beta21_("", dimless/dimTemperature, 0)
63 {
64  if (this->speciesNames_.size() != 2)
65  {
67  << "NonRandomTwoLiquid model is suitable for two species only."
68  << exit(FatalError);
69  }
70 
71  species1Name_ = this->speciesNames_[0];
72  species2Name_ = this->speciesNames_[1];
73 
74  species1Index_ = this->thermo_.composition().species()[species1Name_];
75  species2Index_ = this->thermo_.composition().species()[species2Name_];
76 
77  alpha12_ = dimensionedScalar
78  (
79  "alpha12",
80  dimless,
81  dict.subDict(species1Name_).lookup("alpha")
82  );
83  alpha21_ = dimensionedScalar
84  (
85  "alpha21",
86  dimless,
87  dict.subDict(species2Name_).lookup("alpha")
88  );
89 
90  beta12_ = dimensionedScalar
91  (
92  "beta12",
94  dict.subDict(species1Name_).lookup("beta")
95  );
96  beta21_ = dimensionedScalar
97  (
98  "beta21",
100  dict.subDict(species2Name_).lookup("beta")
101  );
102 
103  saturationModel12_.reset
104  (
106  (
107  dict.subDict(species1Name_).subDict("interaction"),
108  pair.phase1().mesh()
109  ).ptr()
110  );
111  saturationModel21_.reset
112  (
114  (
115  dict.subDict(species2Name_).subDict("interaction"),
116  pair.phase1().mesh()
117  ).ptr()
118  );
119 
120  speciesModel1_.reset
121  (
123  (
124  dict.subDict(species1Name_),
125  pair
126  ).ptr()
127  );
128 
129  speciesModel2_.reset
130  (
132  (
133  dict.subDict(species2Name_),
134  pair
135  ).ptr()
136  );
137 }
138 
139 
140 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
141 
142 template<class Thermo, class OtherThermo>
145 {}
146 
147 
148 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
149 
150 template<class Thermo, class OtherThermo>
151 void
153 update
154 (
155  const volScalarField& Tf
156 )
157 {
158  volScalarField W(this->thermo_.W());
159 
160  volScalarField X1
161  (
162  this->thermo_.composition().Y(species1Index_)
163  *W
165  (
166  "W",
168  this->thermo_.composition().Wi(species1Index_)
169  )
170  );
171 
172  volScalarField X2
173  (
174  this->thermo_.composition().Y(species2Index_)
175  *W
177  (
178  "W",
180  this->thermo_.composition().Wi(species2Index_)
181  )
182  );
183 
184  volScalarField alpha12(alpha12_ + Tf*beta12_);
185  volScalarField alpha21(alpha21_ + Tf*beta21_);
186 
187  volScalarField tau12(saturationModel12_->lnPSat(Tf));
188  volScalarField tau21(saturationModel21_->lnPSat(Tf));
189 
190  volScalarField G12(exp(- alpha12*tau12));
191  volScalarField G21(exp(- alpha21*tau21));
192 
193  gamma1_ =
194  exp
195  (
196  sqr(X2)
197  *(
198  tau21*sqr(G21)/max(sqr(X1 + X2*G21), small)
199  + tau12*G12/max(sqr(X2 + X1*G12), small)
200  )
201  );
202  gamma2_ =
203  exp
204  (
205  sqr(X1)
206  *(
207  tau12*sqr(G12)/max(sqr(X2 + X1*G12), small)
208  + tau21*G21/max(sqr(X1 + X2*G21), small)
209  )
210  );
211 }
212 
213 
214 template<class Thermo, class OtherThermo>
217 (
218  const word& speciesName,
219  const volScalarField& Tf
220 ) const
221 {
222  if (speciesName == species1Name_)
223  {
224  return
225  this->otherThermo_.composition().Y(speciesName)
226  *speciesModel1_->Yf(speciesName, Tf)
227  *gamma1_;
228  }
229  else if (speciesName == species2Name_)
230  {
231  return
232  this->otherThermo_.composition().Y(speciesName)
233  *speciesModel2_->Yf(speciesName, Tf)
234  *gamma2_;
235  }
236  else
237  {
238  return
239  this->thermo_.composition().Y(speciesName)
240  *(scalar(1) - Yf(species1Name_, Tf) - Yf(species2Name_, Tf));
241  }
242 }
243 
244 
245 template<class Thermo, class OtherThermo>
248 YfPrime
249 (
250  const word& speciesName,
251  const volScalarField& Tf
252 ) const
253 {
254  if (speciesName == species1Name_)
255  {
256  return
257  this->otherThermo_.composition().Y(speciesName)
258  *speciesModel1_->YfPrime(speciesName, Tf)
259  *gamma1_;
260  }
261  else if (speciesName == species2Name_)
262  {
263  return
264  this->otherThermo_.composition().Y(speciesName)
265  *speciesModel2_->YfPrime(speciesName, Tf)
266  *gamma2_;
267  }
268  else
269  {
270  return
271  - this->thermo_.composition().Y(speciesName)
272  *(YfPrime(species1Name_, Tf) + YfPrime(species2Name_, Tf));
273  }
274 }
275 
276 
277 // ************************************************************************* //
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 > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
dimensionedSymmTensor sqr(const dimensionedVector &dv)
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
dynamicFvMesh & mesh
dimensionedScalar exp(const dimensionedScalar &ds)
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:52
NonRandomTwoLiquid(const dictionary &dict, const phasePair &pair)
Construct from components.
word timeName
Definition: getTimeIndex.H:3
phaseModel & phase1
const dimensionSet dimMoles(0, 0, 0, 0, 1, 0, 0)
Definition: dimensionSets.H:53
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const scalarList W(::W(thermo))
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
A class for managing temporary objects.
Definition: PtrList.H:53
virtual void update(const volScalarField &Tf)
Update the composition.