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