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-2020 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 #include "phasePair.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace interfaceCompositionModels
35 {
36  defineTypeNameAndDebug(nonRandomTwoLiquid, 0);
38  (
39  interfaceCompositionModel,
40  nonRandomTwoLiquid,
41  dictionary
42  );
43 }
44 }
45 
46 
47 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48 
50 (
51  const dictionary& dict,
52  const phasePair& pair
53 )
54 :
55  interfaceCompositionModel(dict, pair),
56  gamma1_
57  (
58  IOobject
59  (
60  IOobject::groupName("gamma1", pair.name()),
61  pair.phase1().mesh().time().timeName(),
62  pair.phase1().mesh()
63  ),
64  pair.phase1().mesh(),
66  ),
67  gamma2_
68  (
69  IOobject
70  (
71  IOobject::groupName("gamma2", pair.name()),
72  pair.phase1().mesh().time().timeName(),
73  pair.phase1().mesh()
74  ),
75  pair.phase1().mesh(),
77  ),
78  beta12_("", dimless/dimTemperature, 0),
79  beta21_("", dimless/dimTemperature, 0)
80 {
81  if (species().size() != 2)
82  {
84  << "nonRandomTwoLiquid model is suitable for two species only."
85  << exit(FatalError);
86  }
87 
88  species1Name_ = species()[0];
89  species2Name_ = species()[1];
90 
91  species1Index_ = composition().species()[species1Name_];
92  species2Index_ = composition().species()[species2Name_];
93 
94  alpha12_ = dimensionedScalar
95  (
96  "alpha12",
97  dimless,
98  dict.subDict(species1Name_).lookup("alpha")
99  );
100  alpha21_ = dimensionedScalar
101  (
102  "alpha21",
103  dimless,
104  dict.subDict(species2Name_).lookup("alpha")
105  );
106 
107  beta12_ = dimensionedScalar
108  (
109  "beta12",
111  dict.subDict(species1Name_).lookup("beta")
112  );
113  beta21_ = dimensionedScalar
114  (
115  "beta21",
117  dict.subDict(species2Name_).lookup("beta")
118  );
119 
120  saturationModel12_.reset
121  (
123  (
124  dict.subDict(species1Name_).subDict("interaction"),
125  pair
126  ).ptr()
127  );
128  saturationModel21_.reset
129  (
131  (
132  dict.subDict(species2Name_).subDict("interaction"),
133  pair
134  ).ptr()
135  );
136 
137  speciesModel1_.reset
138  (
140  (
141  dict.subDict(species1Name_),
142  pair
143  ).ptr()
144  );
145  speciesModel2_.reset
146  (
148  (
149  dict.subDict(species2Name_),
150  pair
151  ).ptr()
152  );
153 }
154 
155 
156 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
157 
159 {}
160 
161 
162 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
163 
165 (
166  const volScalarField& Tf
167 )
168 {
169  const volScalarField W(thermo().W());
170 
171  const volScalarField X1
172  (
173  composition().Y(species1Index_)
174  *W
176  (
177  "W",
179  composition().Wi(species1Index_)
180  )
181  );
182 
183  const volScalarField X2
184  (
185  composition().Y(species2Index_)
186  *W
188  (
189  "W",
191  composition().Wi(species2Index_)
192  )
193  );
194 
195  const volScalarField alpha12(alpha12_ + Tf*beta12_);
196  const volScalarField alpha21(alpha21_ + Tf*beta21_);
197 
198  const volScalarField tau12(saturationModel12_->lnPSat(Tf));
199  const volScalarField tau21(saturationModel21_->lnPSat(Tf));
200 
201  const volScalarField G12(exp(- alpha12*tau12));
202  const volScalarField G21(exp(- alpha21*tau21));
203 
204  gamma1_ =
205  exp
206  (
207  sqr(X2)
208  *(
209  tau21*sqr(G21)/max(sqr(X1 + X2*G21), small)
210  + tau12*G12/max(sqr(X2 + X1*G12), small)
211  )
212  );
213  gamma2_ =
214  exp
215  (
216  sqr(X1)
217  *(
218  tau12*sqr(G12)/max(sqr(X2 + X1*G12), small)
219  + tau21*G21/max(sqr(X1 + X2*G21), small)
220  )
221  );
222 }
223 
224 
227 (
228  const word& speciesName,
229  const volScalarField& Tf
230 ) const
231 {
232  if (speciesName == species1Name_)
233  {
234  return
235  otherComposition().Y(speciesName)
236  *speciesModel1_->Yf(speciesName, Tf)
237  *gamma1_;
238  }
239  else if (speciesName == species2Name_)
240  {
241  return
242  otherComposition().Y(speciesName)
243  *speciesModel2_->Yf(speciesName, Tf)
244  *gamma2_;
245  }
246  else
247  {
248  return
249  composition().Y(speciesName)
250  *(scalar(1) - Yf(species1Name_, Tf) - Yf(species2Name_, Tf));
251  }
252 }
253 
254 
257 (
258  const word& speciesName,
259  const volScalarField& Tf
260 ) const
261 {
262  if (speciesName == species1Name_)
263  {
264  return
265  otherComposition().Y(speciesName)
266  *speciesModel1_->YfPrime(speciesName, Tf)
267  *gamma1_;
268  }
269  else if (speciesName == species2Name_)
270  {
271  return
272  otherComposition().Y(speciesName)
273  *speciesModel2_->YfPrime(speciesName, Tf)
274  *gamma2_;
275  }
276  else
277  {
278  return
279  - composition().Y(speciesName)
280  *(YfPrime(species1Name_, Tf) + YfPrime(species2Name_, Tf));
281  }
282 }
283 
284 
285 // ************************************************************************* //
const rhoReactionThermo & thermo() const
Return the thermo.
const hashedWordList & species() const
Return the transferring species names.
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)
static autoPtr< saturationModel > New(const dictionary &dict, const phasePair &pair)
Select null constructed.
static autoPtr< interfaceCompositionModel > New(const dictionary &dict, const phasePair &pair)
const basicSpecieMixture & otherComposition() const
Return the other composition.
virtual tmp< volScalarField > Yf(const word &speciesName, const volScalarField &Tf) const
The interface species fraction.
Macros for easy insertion into run-time selection tables.
const phasePair & pair() const
Return the phase pair.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
dynamicFvMesh & mesh
dimensionedScalar exp(const dimensionedScalar &ds)
virtual tmp< volScalarField > Yf(const word &speciesName, const volScalarField &Tf) const =0
Interface mass fraction.
virtual tmp< volScalarField > YfPrime(const word &speciesName, const volScalarField &Tf) const
The interface species fraction derivative w.r.t. temperature.
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:52
word timeName
Definition: getTimeIndex.H:3
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
const dimensionSet dimMoles(0, 0, 0, 0, 1, 0, 0)
Definition: dimensionSets.H:53
defineTypeNameAndDebug(combustionModel, 0)
nonRandomTwoLiquid(const dictionary &dict, const phasePair &pair)
Construct from components.
virtual void update(const volScalarField &Tf)
Update the composition.
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
const basicSpecieMixture & composition() const
Return the composition.
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))
PtrList< volScalarField > & Y
virtual tmp< volScalarField > YfPrime(const word &speciesName, const volScalarField &Tf) const =0
The interface mass fraction derivative w.r.t. temperature.
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
A class for managing temporary objects.
Definition: PtrList.H:53
PtrList< volScalarField > & Y()
Return the mass-fraction fields.
Namespace for OpenFOAM.
const speciesTable & species() const
Return the table of species.