twoPhaseSystem.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) 2013-2018 OpenFOAM Foundation
9 Copyright (C) 2020 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "twoPhaseSystem.H"
30#include "dragModel.H"
31#include "virtualMassModel.H"
32
33#include "MULES.H"
34#include "subCycle.H"
35#include "UniformField.H"
36
37#include "fvcDdt.H"
38#include "fvcDiv.H"
39#include "fvcSnGrad.H"
40#include "fvcFlux.H"
41#include "fvcSup.H"
42
43#include "fvmDdt.H"
44#include "fvmLaplacian.H"
45#include "fvmSup.H"
46
47// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
48
49namespace Foam
50{
53}
54
55
56// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
57
59(
60 const fvMesh& mesh
61)
62:
64 phase1_(phaseModels_[0]),
65 phase2_(phaseModels_[1])
66{
67 phase2_.volScalarField::operator=(scalar(1) - phase1_);
68
71}
72
73
74// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
75
77{}
78
79
80// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
81
84{
85 return sigma
86 (
88 );
89}
90
91
94{
95 return Kd
96 (
98 );
99}
100
101
104{
105 return Kdf
106 (
108 );
109}
110
111
114{
115 return Vm
116 (
118 );
119}
120
121
123{
124 const Time& runTime = mesh_.time();
125
126 volScalarField& alpha1 = phase1_;
127 volScalarField& alpha2 = phase2_;
128
129 const dictionary& alphaControls = mesh_.solverDict(alpha1.name());
130
131 label nAlphaSubCycles(alphaControls.get<label>("nAlphaSubCycles"));
132 label nAlphaCorr(alphaControls.get<label>("nAlphaCorr"));
133
134 bool LTS = fv::localEulerDdt::enabled(mesh_);
135
136 word alphaScheme("div(phi," + alpha1.name() + ')');
137 word alpharScheme("div(phir," + alpha1.name() + ')');
138
139 const surfaceScalarField& phi1 = phase1_.phi();
140 const surfaceScalarField& phi2 = phase2_.phi();
141
142 // Construct the dilatation rate source term
144
145 if (phase1_.divU().valid() && phase2_.divU().valid())
146 {
147 tdgdt =
148 (
149 alpha2()
150 *phase1_.divU()()()
151 - alpha1()
152 *phase2_.divU()()()
153 );
154 }
155 else if (phase1_.divU().valid())
156 {
157 tdgdt =
158 (
159 alpha2()
160 *phase1_.divU()()()
161 );
162 }
163 else if (phase2_.divU().valid())
164 {
165 tdgdt =
166 (
167 - alpha1()
168 *phase2_.divU()()()
169 );
170 }
171
173
174 surfaceScalarField phir("phir", phi1 - phi2);
175
176 tmp<surfaceScalarField> alphaDbyA;
177 if (DByAfs().found(phase1_.name()) && DByAfs().found(phase2_.name()))
178 {
180 (
181 *DByAfs()[phase1_.name()] + *DByAfs()[phase2_.name()]
182 );
183
184 alphaDbyA =
185 fvc::interpolate(max(alpha1, scalar(0)))
186 *fvc::interpolate(max(alpha2, scalar(0)))
187 *DbyA;
188
189 phir += DbyA*fvc::snGrad(alpha1, "bounded")*mesh_.magSf();
190 }
191
192 for (int acorr=0; acorr<nAlphaCorr; acorr++)
193 {
195 (
197 (
198 "Sp",
200 mesh_
201 ),
202 mesh_,
204 );
205
207 (
209 (
210 "Su",
212 mesh_
213 ),
214 // Divergence term is handled explicitly to be
215 // consistent with the explicit transport solution
216 fvc::div(phi_)*min(alpha1, scalar(1))
217 );
218
219 if (tdgdt.valid())
220 {
221 scalarField& dgdt = tdgdt.ref();
222
223 forAll(dgdt, celli)
224 {
225 if (dgdt[celli] > 0.0)
226 {
227 Sp[celli] -= dgdt[celli]/max(1 - alpha1[celli], 1e-4);
228 Su[celli] += dgdt[celli]/max(1 - alpha1[celli], 1e-4);
229 }
230 else if (dgdt[celli] < 0.0)
231 {
232 Sp[celli] += dgdt[celli]/max(alpha1[celli], 1e-4);
233 }
234 }
235 }
236
238 (
240 (
241 phi_,
242 alpha1,
243 alphaScheme
244 )
245 + fvc::flux
246 (
247 -fvc::flux(-phir, scalar(1) - alpha1, alpharScheme),
248 alpha1,
250 )
251 );
252
253 phase1_.correctInflowOutflow(alphaPhi1);
254
255 if (nAlphaSubCycles > 1)
256 {
257 tmp<volScalarField> trSubDeltaT;
258
259 if (LTS)
260 {
261 trSubDeltaT =
263 }
264
265 for
266 (
268 !(++alphaSubCycle).end();
269 )
270 {
272
274 (
276 alpha1,
277 phi_,
279 (alphaSubCycle.index()*Sp)(),
280 (Su - (alphaSubCycle.index() - 1)*Sp*alpha1)(),
281 UniformField<scalar>(phase1_.alphaMax()),
282 zeroField()
283 );
284
285 if (alphaSubCycle.index() == 1)
286 {
287 phase1_.alphaPhiRef() = alphaPhi10;
288 }
289 else
290 {
291 phase1_.alphaPhiRef() += alphaPhi10;
292 }
293 }
294
295 phase1_.alphaPhiRef() /= nAlphaSubCycles;
296 }
297 else
298 {
300 (
302 alpha1,
303 phi_,
304 alphaPhi1,
305 Sp,
306 Su,
307 UniformField<scalar>(phase1_.alphaMax()),
308 zeroField()
309 );
310
311 phase1_.alphaPhiRef() = alphaPhi1;
312 }
313
314 if (alphaDbyA.valid())
315 {
316 fvScalarMatrix alpha1Eqn
317 (
319 - fvm::laplacian(alphaDbyA(), alpha1, "bounded")
320 );
321
322 alpha1Eqn.relax();
323 alpha1Eqn.solve();
324
325 phase1_.alphaPhiRef() += alpha1Eqn.flux();
326 }
327
328 phase1_.alphaRhoPhiRef() =
329 fvc::interpolate(phase1_.rho())*phase1_.alphaPhi();
330
331 phase2_.alphaPhiRef() = phi_ - phase1_.alphaPhi();
332 phase2_.correctInflowOutflow(phase2_.alphaPhiRef());
333 phase2_.alphaRhoPhiRef() =
334 fvc::interpolate(phase2_.rho())*phase2_.alphaPhi();
335
336 Info<< alpha1.name() << " volume fraction = "
337 << alpha1.weightedAverage(mesh_.V()).value()
338 << " Min(alpha1) = " << min(alpha1).value()
339 << " Max(alpha1) = " << max(alpha1).value()
340 << endl;
341
342 // Ensure the phase-fractions are bounded
343 alpha1.clip(SMALL, 1 - SMALL);
344
345 // Update the phase-fraction of the other phase
346 alpha2 = scalar(1) - alpha1;
347 }
348}
349
350
351// ************************************************************************* //
MULES: Multidimensional universal limiter for explicit solution.
bool found
const volScalarField & alpha1
surfaceScalarField & phi2
phaseModel & phase1
const volScalarField & alpha2
const surfaceScalarField & alphaPhi1
surfaceScalarField & phi1
phaseModel & phase2
dimensioned< Type > weightedAverage(const DimensionedField< scalar, GeoMesh > &weightField) const
Calculate and return weighted average.
void clip(const dimensioned< MinMax< Type > > &range)
Clip the field to be bounded within the specified range.
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
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:80
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
A class representing the concept of a uniform field which stores only the single value and providing ...
Definition: UniformField.H:54
virtual bool enabled() const
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
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:121
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:1092
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
static tmp< volScalarField > localRSubDeltaT(const fvMesh &mesh, const label nAlphaSubCycles)
Calculate and return the reciprocal of the local sub-cycling.
Definition: localEulerDdt.C:72
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
const Time & time() const noexcept
Return time registry.
An ordered or unorder pair of phase names. Typically specified as follows.
Definition: phasePairKey.H:68
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:76
void setFluxRequired(const word &name) const
Get flux-required for given name, or default.
Perform a subCycleTime on a field.
Definition: subCycle.H:127
A class for managing temporary objects.
Definition: tmp.H:65
bool valid() const noexcept
Identical to good(), or bool operator.
Definition: tmp.H:292
T & ref() const
Definition: tmpI.H:227
Class which solves the volume fraction equations for two phases.
tmp< volScalarField > Kd() const
Return the drag coefficient.
virtual ~twoPhaseSystem()
Destructor.
phaseModel & phase1_
Phase model 1.
tmp< surfaceScalarField > Kdf() const
Return the face drag coefficient.
phaseModel & phase2_
Phase model 2.
const fvMesh & mesh() const
Return the mesh.
tmp< volScalarField > sigma() const
Return the surface tension coefficient.
virtual void solve()
Solve for the phase fractions.
tmp< volScalarField > Vm() const
Return the virtual mass coefficient.
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
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
alphaPhi10
Definition: alphaEqn.H:7
dynamicFvMesh & mesh
engineTime & runTime
bool LTS
Definition: createRDeltaT.H:1
surfaceScalarField phir(fvc::flux(UdmModel.Udm()))
word alpharScheme("div(phirb,alpha)")
Calculate the first temporal derivative.
Calculate the divergence of the given field.
Calculate the face-flux of the given field.
Calculate the snGrad of the given volField.
Calculate the field for explicit evaluation of implicit and explicit sources.
Calculate the matrix for the first temporal derivative.
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
const volScalarField Kd(fluid.Kd())
const surfaceScalarField Kdf("Kdf", fluid.Kdf())
void explicitSolve(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiPsi, const SpType &Sp, const SuType &Su)
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:47
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:49
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
Namespace for OpenFOAM.
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 dimless
Dimensionless.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
#define defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.
volScalarField & e
Definition: createFields.H:11
const dictionary & alphaControls
Definition: alphaControls.H:1
label nAlphaSubCycles(alphaControls.get< label >("nAlphaSubCycles"))
label nAlphaCorr(alphaControls.get< label >("nAlphaCorr"))
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333