MultiComponentPhaseModel.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2017-2022 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
29
31#include "fvmDdt.H"
32#include "fvmDiv.H"
33#include "fvmSup.H"
34#include "fvmLaplacian.H"
35#include "fvcDdt.H"
36#include "fvcDiv.H"
37#include "fvcDDt.H"
38#include "fvMatrix.H"
39#include "fvcFlux.H"
40#include "CMULES.H"
41#include "subCycle.H"
42
43// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44
45template<class BasePhaseModel, class phaseThermo>
48(
50 const word& phaseName
51)
52:
53 BasePhaseModel(fluid, phaseName),
54 species_(),
55 inertIndex_(-1),
56 addDiffusion_(false)
57{
59 (
60 phaseThermo::New
61 (
62 fluid.mesh(),
63 phaseName,
65 )
66 );
67
68 if (thermoPtr_->composition().species().empty())
69 {
71 << " The selected thermo is pure. Use a multicomponent thermo."
72 << exit(FatalError);
73 }
74
75 species_ = thermoPtr_->composition().species();
76
77 inertIndex_ = species_.find(thermoPtr_().template get<word>("inertSpecie"));
78
80 thermoPtr_().template getOrDefault<bool>("addDiffusion", false);
81
82 Sct_ = thermoPtr_().template getOrDefault<scalar>("Sct", 1.0);
83
84 X_.setSize(thermoPtr_->composition().species().size());
85
86 // Initiate X's using Y's to set BC's
88 {
89 X_.set
90 (
91 i,
93 (
95 (
96 IOobject::groupName("X" + species_[i], phaseName),
97 fluid.mesh().time().timeName(),
98 fluid.mesh(),
101 ),
102 Y()[i]
103 )
104 );
105 }
106
107 // Init vol fractions from mass fractions
109}
110
111
112// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
113
114
115template<class BasePhaseModel, class phaseThermo>
118{
119 volScalarField Xtotal(0.0*X_[0]);
120 const volScalarField W(thermo().W());
121
122 forAll(X_, i)
123 {
124 const dimensionedScalar Wi
125 (
126 "Wi",
128 thermo().composition().W(i)
129 );
130
131 if (i != inertIndex_)
132 {
133 X_[i] = W*Y()[i]/Wi;
134 Xtotal += X_[i];
136
137 const volScalarField::Boundary& YBf = Y()[i].boundaryField();
138
139 forAll(YBf, patchi)
140 {
141 const fvPatchScalarField& YPf = YBf[patchi];
142 if (YPf.fixesValue())
143 {
144 scalarField& xbf = X_[i].boundaryFieldRef()[patchi];
145 const scalarField& ybf = Y()[i].boundaryField()[patchi];
146 const scalarField& Wbf = W.boundaryField()[patchi];
147 forAll(xbf, facei)
148 {
149 xbf[facei] = Wbf[facei]*ybf[facei]/Wi.value();
150 }
151 }
152 }
153 }
154 }
155 X_[inertIndex_] = 1.0 - Xtotal;
156 X_[inertIndex_].correctBoundaryConditions();
157}
158
159
160template<class BasePhaseModel, class phaseThermo>
163{
164 volScalarField W(X_[0]*thermo().composition().W(0));
165 for(label i=1; i< species_.size(); i++)
166 {
167 W += X_[i]*thermo().composition().W(i);
168 }
169
170 forAll(Y(), i)
171 {
172 Y()[i] = X_[i]*thermo().composition().W(i)/W;
173
174 Info<< Y()[i].name() << " mass fraction = "
175 << " Min(Y) = " << min(Y()[i]).value()
176 << " Max(Y) = " << max(Y()[i]).value()
177 << endl;
178 }
179}
180
181
182template<class BasePhaseModel, class phaseThermo>
183const phaseThermo&
185{
186 return thermoPtr_();
187}
188
189
190template<class BasePhaseModel, class phaseThermo>
191phaseThermo&
193{
194 return thermoPtr_();
195}
196
197
198template<class BasePhaseModel, class phaseThermo>
200{
201 return thermo().correct();
202}
203
204
205template<class BasePhaseModel, class phaseThermo>
207(
210)
211{
212 const volScalarField& alpha1 = *this;
213
214 const fvMesh& mesh = alpha1.mesh();
215
216 const dictionary& MULEScontrols = mesh.solverDict(alpha1.name());
217
218 scalar cAlpha(MULEScontrols.get<scalar>("cYi"));
219
220 PtrList<surfaceScalarField> phiYiCorrs(species_.size());
221 const surfaceScalarField& phi = this->fluid().phi();
222
224
225 phic = min(cAlpha*phic, max(phic));
226
228
230 (
232 )
233 {
234 const volScalarField& alpha2 = iter2()();
235 if (&alpha2 == &alpha1)
236 {
237 continue;
238 }
239
240 phir += phic*this->fluid().nHatf(alpha1, alpha2);
241 }
242
243 // Do not compress interface at non-coupled boundary faces
244 // (inlets, outlets etc.)
246 forAll(phir.boundaryField(), patchi)
247 {
248 fvsPatchScalarField& phirp = phirBf[patchi];
249
250 if (!phirp.coupled())
251 {
252 phirp == 0;
253 }
254 }
255
256 word phirScheme("div(Yiphir," + alpha1.name() + ')');
257
258 forAll(X_, i)
259 {
260 if (inertIndex_ != i)
261 {
262 volScalarField& Yi = X_[i];
263
264 phiYiCorrs.set
265 (
266 i,
268 (
269 "phi" + Yi.name() + "Corr",
271 (
272 phi,
273 Yi,
274 "div(phi," + Yi.name() + ')'
275 )
276 )
277 );
278
279 surfaceScalarField& phiYiCorr = phiYiCorrs[i];
280
282 (
284 this->fluid().phases(),
285 iter2
286 )
287 {
288 //const volScalarField& alpha2 = iter2()().oldTime();
289 const volScalarField& alpha2 = iter2()();
290
291 if (&alpha2 == &alpha1)
292 {
293 continue;
294 }
295
296 phiYiCorr += fvc::flux
297 (
298 -fvc::flux(-phir, alpha2, phirScheme),
299 Yi,
300 phirScheme
301 );
302 }
303
304 // Ensure that the flux at inflow BCs is preserved
305 forAll(phiYiCorr.boundaryField(), patchi)
306 {
307 fvsPatchScalarField& phiYiCorrp =
308 phiYiCorr.boundaryFieldRef()[patchi];
309
310 if (!phiYiCorrp.coupled())
311 {
312 const scalarField& phi1p = phi.boundaryField()[patchi];
313 const scalarField& Yip = Yi.boundaryField()[patchi];
314
315 forAll(phiYiCorrp, facei)
316 {
317 if (phi1p[facei] < 0)
318 {
319 phiYiCorrp[facei] = Yip[facei]*phi1p[facei];
320 }
321 }
322 }
323 }
324
326 (
327 1.0/mesh.time().deltaT().value(),
329 Yi,
330 phi,
331 phiYiCorr,
332 Sp[i],
333 Su[i],
334 oneField(),
335 zeroField(),
336 true
337 );
338 }
339 }
340
341 volScalarField Yt(0.0*X_[0]);
342
343 scalar nYiSubCycles
344 (
345 MULEScontrols.getOrDefault<scalar>("nYiSubCycles", 1)
346 );
347
348 forAll(X_, i)
349 {
350 if (inertIndex_ != i)
351 {
352 volScalarField& Yi = X_[i];
353
354 fvScalarMatrix YiEqn
355 (
358 (
359 mesh,
360 phi,
362 ).fvmDiv(phi, Yi)
363 ==
364 Su[i] + fvm::Sp(Sp[i], Yi)
365 );
366
367 YiEqn.solve();
368
369 surfaceScalarField& phiYiCorr = phiYiCorrs[i];
370
371 // Add a bounded upwind U-mean flux
372 phiYiCorr += YiEqn.flux();
373
374 if (nYiSubCycles > 1)
375 {
376 for
377 (
378 subCycle<volScalarField> YiSubCycle(Yi, nYiSubCycles);
379 !(++YiSubCycle).end();
380 )
381 {
383 (
385 Yi,
386 phi,
387 phiYiCorr,
388 Sp[i],
389 Su[i],
390 oneField(),
391 zeroField()
392 );
393 }
394 }
395 else
396 {
398 (
400 Yi,
401 phi,
402 phiYiCorr,
403 Sp[i],
404 Su[i],
405 oneField(),
406 zeroField()
407 );
408 }
409
410 if (addDiffusion_)
411 {
412 const volScalarField& alpha = *this;
413 fvScalarMatrix YiDiffEqn
414 (
415 fvm::ddt(Yi) - fvc::ddt(Yi)
417 (
418 alpha*this->fluid().turbulence()->nut()/Sct_,
419 Yi
420 )
421 );
422
423 YiDiffEqn.solve(mesh.solver("diffusion" + Yi.name()));
424 }
425
426 Yt += Yi;
427 }
428 }
429
430 X_[inertIndex_] = scalar(1) - Yt;
431 X_[inertIndex_].max(0.0);
432
433 calculateMassFractions();
434}
435
436
437template<class BasePhaseModel, class phaseThermo>
440{
441 return thermoPtr_->composition().Y();
442}
443
444
445template<class BasePhaseModel, class phaseThermo>
448{
449 return thermoPtr_->composition().Y();
450}
451
452
453template<class BasePhaseModel, class phaseThermo>
454Foam::label
456{
457 return inertIndex_;
458}
459
460
461// ************************************************************************* //
CMULES: Multidimensional universal limiter for explicit corrected implicit solution.
surfaceScalarField phic(mixture.cAlpha() *mag(alphaPhic/mesh.magSf()))
volScalarField Yt(0.0 *Y[0])
surfaceScalarField & phi
const volScalarField & alpha1
twoPhaseSystem & fluid
const volScalarField & alpha2
const Mesh & mesh() const
Return mesh.
void max(const dimensioned< Type > &dt)
Use the maximum of the field and specified value.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
void correctBoundaryConditions()
Correct boundary field.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:65
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
Class which represents a phase with multiple species. Returns the species' mass fractions,...
void calculateMassFractions()
Transfor volume fraction into mass fractions.
virtual void solveYi(PtrList< volScalarField::Internal > &, PtrList< volScalarField::Internal > &)
Solve species fraction equation.
virtual void correct()
Correct phase thermo.
label inertIndex() const
Return inert species index.
autoPtr< phaseThermo > thermoPtr_
Thermophysical model.
PtrList< volScalarField > X_
Ptr list of volumetric fractions for species.
virtual const phaseThermo & thermo() const
Access to thermo.
void calculateVolumeFractions()
Transfor mass fraction into volume fractions.
label inertIndex_
Inert species index.
bool addDiffusion_
Add diffusion transport on Yi's Eq.
hashedWordList species_
Species table.
virtual const PtrList< volScalarField > & Y() const
Constant access the species mass fractions.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
const T * set(const label i) const
Definition: PtrList.H:138
dimensionedScalar deltaT() const
Return time step.
Definition: TimeStateI.H:55
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
bool empty() const noexcept
Deprecated(2020-07) True if the managed pointer is null.
Definition: autoPtr.H:278
void reset(autoPtr< T > &&other) noexcept
Delete managed object and set to new given pointer.
Definition: autoPtrI.H:117
static word phasePropertyName(const word &name, const word &phaseName)
Definition: basicThermo.H:259
static const word dictName
Definition: basicThermo.H:256
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
const Type & value() const
Return const reference to value.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:121
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > flux() const
Return the face-flux field from the matrix.
Definition: fvMatrix.C:1443
SolverPerformance< Type > solve(const dictionary &)
Solve returning the solution statistics.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:290
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
virtual bool fixesValue() const
Return true if this patch field fixes a value.
Definition: fvPatchField.H:337
Basic first-order Euler implicit/explicit ddt using only the current and previous time-step values.
Basic second-order convection using face-gradients and Gauss' theorem.
virtual bool coupled() const
Return true if this patch field is coupled.
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
label find(const word &val) const
Find index of the value (searches the hash).
A class representing the concept of a field of 1 used to avoid unnecessary manipulations for objects ...
Definition: oneField.H:56
const dictionary & solver(const word &name) const
Return the solver controls dictionary for the given field.
Definition: solution.C:373
const dictionary & solverDict(const word &name) const
Return the solver controls dictionary for the given field.
Definition: solution.C:366
Perform a subCycleTime on a field.
Definition: subCycle.H:127
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
const surfaceScalarField & phi() const
Return the mixture flux.
const fvMesh & mesh() const
Return the mesh.
Upwind differencing scheme class.
Definition: upwind.H:59
A class for handling words, derived from Foam::string.
Definition: word.H:68
A class representing the concept of a field of 0 used to avoid unnecessary manipulations for objects ...
Definition: zeroField.H:56
basicSpecieMixture & composition
PtrList< volScalarField > & Y
dynamicFvMesh & mesh
surfaceScalarField phir(fvc::flux(UdmModel.Udm()))
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
scalar nut
compressible::turbulenceModel & turbulence
Calculate the substantive (total) derivative.
Calculate the first temporal derivative.
Calculate the divergence of the given field.
Calculate the face-flux of the given field.
Calculate the matrix for the first temporal derivative.
Calculate the matrix for the divergence of the given field and flux.
Calculate the matrix for the laplacian of the field.
Calculate the finiteVolume matrix for implicit and explicit sources.
zeroField Su
Definition: alphaSuSp.H:1
zeroField Sp
Definition: alphaSuSp.H:2
void explicitSolve(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiPsi, const SpType &Sp, const SuType &Su)
void limit(const RdeltaTType &rDeltaT, const RhoType &rho, const volScalarField &psi, const surfaceScalarField &phi, surfaceScalarField &phiPsi, const SpType &Sp, const SuType &Su, const PsiMaxType &psiMax, const PsiMinType &psiMin, const bool returnCorr)
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:47
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
Definition: fvcFlux.C:34
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:48
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:48
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
const dimensionSet dimMoles(0, 0, 0, 0, 1, 0, 0)
Definition: dimensionSets.H:55
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
error FatalError
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:51
volScalarField & alpha
multiphaseSystem::phaseModelList & phases
Definition: createFields.H:12
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object.
Definition: stdFoam.H:381