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-2023 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 {
37  (
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().name(),
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().name(),
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_ = thermo().species()[species1Name_];
91  species2Index_ = thermo().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  "interaction",
124  dict.subDict(species1Name_)
125  ).ptr()
126  );
127  saturationModel21_.reset
128  (
130  (
131  "interaction",
132  dict.subDict(species2Name_)
133  ).ptr()
134  );
135 
136  speciesModel1_.reset
137  (
139  (
140  dict.subDict(species1Name_),
141  interface,
142  false
143  ).ptr()
144  );
145  speciesModel2_.reset
146  (
148  (
149  dict.subDict(species2Name_),
150  interface,
151  false
152  ).ptr()
153  );
154 }
155 
156 
157 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
158 
160 {}
161 
162 
163 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
164 
166 (
167  const volScalarField& Tf
168 )
169 {
170  const volScalarField W(thermo().W());
171 
172  const volScalarField X1
173  (
174  thermo().Y(species1Index_)*W/thermo().Wi(species1Index_)
175  );
176 
177  const volScalarField X2
178  (
179  thermo().Y(species2Index_)*W/thermo().Wi(species2Index_)
180  );
181 
182  const volScalarField alpha12(alpha12_ + Tf*beta12_);
183  const volScalarField alpha21(alpha21_ + Tf*beta21_);
184 
185  const volScalarField tau12(saturationModel12_->lnPSat(Tf));
186  const volScalarField tau21(saturationModel21_->lnPSat(Tf));
187 
188  const volScalarField G12(exp(- alpha12*tau12));
189  const volScalarField G21(exp(- alpha21*tau21));
190 
191  gamma1_ =
192  exp
193  (
194  sqr(X2)
195  *(
196  tau21*sqr(G21)/max(sqr(X1 + X2*G21), small)
197  + tau12*G12/max(sqr(X2 + X1*G12), small)
198  )
199  );
200  gamma2_ =
201  exp
202  (
203  sqr(X1)
204  *(
205  tau12*sqr(G12)/max(sqr(X2 + X1*G12), small)
206  + tau21*G21/max(sqr(X1 + X2*G21), small)
207  )
208  );
209 }
210 
211 
214 (
215  const word& speciesName,
216  const volScalarField& Tf
217 ) const
218 {
219  if (speciesName == species1Name_)
220  {
221  return
222  otherMulticomponentThermo().Y(speciesName)
223  *speciesModel1_->Yf(speciesName, Tf)
224  *gamma1_;
225  }
226  else if (speciesName == species2Name_)
227  {
228  return
229  otherMulticomponentThermo().Y(speciesName)
230  *speciesModel2_->Yf(speciesName, Tf)
231  *gamma2_;
232  }
233  else
234  {
235  return
236  thermo().Y(speciesName)
237  *(scalar(1) - Yf(species1Name_, Tf) - Yf(species2Name_, Tf));
238  }
239 }
240 
241 
244 (
245  const word& speciesName,
246  const volScalarField& Tf
247 ) const
248 {
249  if (speciesName == species1Name_)
250  {
251  return
252  otherMulticomponentThermo().Y(speciesName)
253  *speciesModel1_->YfPrime(speciesName, Tf)
254  *gamma1_;
255  }
256  else if (speciesName == species2Name_)
257  {
258  return
259  otherMulticomponentThermo().Y(speciesName)
260  *speciesModel2_->YfPrime(speciesName, Tf)
261  *gamma2_;
262  }
263  else
264  {
265  return
266  - thermo().Y(speciesName)
267  *(YfPrime(species1Name_, Tf) + YfPrime(species2Name_, Tf));
268  }
269 }
270 
271 
272 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Generic base class for interface composition models. These models describe the composition in phase 1...
const sidedPhaseInterface & interface() const
Return the interface.
static autoPtr< interfaceCompositionModel > New(const dictionary &dict, const phaseInterface &interface, const bool outer=true)
const rhoFluidMulticomponentThermo & thermo() const
Return the thermo.
const hashedWordList & species() const
Return the transferring species names.
Non ideal law for the mixing of two species. A separate composition model is given for each species....
virtual void update(const volScalarField &Tf)
Update the composition.
nonRandomTwoLiquid(const dictionary &dict, const phaseInterface &interface)
Construct from a dictionary and an interface.
virtual tmp< volScalarField > Yf(const word &speciesName, const volScalarField &Tf) const
The interface species fraction.
virtual tmp< volScalarField > YfPrime(const word &speciesName, const volScalarField &Tf) const
The interface species fraction derivative w.r.t. temperature.
virtual const speciesTable & species() const =0
The table of species.
Class to represent an interface between phases. Derivations can further specify the configuration of ...
static autoPtr< saturationPressureModel > New(const dictionary &dict)
Select with dictionary.
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
addToRunTimeSelectionTable(interfaceCompositionModel, Henry, dictionary)
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
dimensionedScalar exp(const dimensionedScalar &ds)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
const dimensionSet dimless
const dimensionSet dimTemperature
scalarList W(const fluidMulticomponentThermo &thermo)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
error FatalError
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dictionary dict
PtrList< volScalarField > & Y
fluidMulticomponentThermo & thermo
Definition: createFields.H:31