Loading...
Searching...
No Matches
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"
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
70 mesh.setFluxRequired(alpha1.name());
71}
72
73
74// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
75
77{}
78
79
80// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
81
84{
85 return sigma
94{
95 return Kd
104{
105 return Kdf
108 );
109}
110
111
114{
115 return Vm
116 (
118 );
119}
120
121
123{
124 volScalarField& alpha1 = phase1_;
125 volScalarField& alpha2 = phase2_;
126
127 const dictionary& alphaControls = mesh_.solverDict(alpha1.name());
128
129 label nAlphaSubCycles(alphaControls.get<label>("nAlphaSubCycles"));
130 label nAlphaCorr(alphaControls.get<label>("nAlphaCorr"));
131
132 bool LTS = fv::localEulerDdt::enabled(mesh_);
133
134 word alphaScheme("div(phi," + alpha1.name() + ')');
135 word alpharScheme("div(phir," + alpha1.name() + ')');
136
137 const tmp<surfaceScalarField> tphi1(phase1_.phi());
138 const surfaceScalarField& phi1 = tphi1();
139 const tmp<surfaceScalarField> tphi2(phase2_.phi());
140 const surfaceScalarField& phi2 = tphi2();
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 (
196 mesh_.newIOobject("Sp"),
197 mesh_,
199 );
200
202 (
203 mesh_.newIOobject("Su"),
204 // Divergence term is handled explicitly to be
205 // consistent with the explicit transport solution
206 fvc::div(phi_)*min(alpha1, scalar(1))
207 );
208
209 if (tdgdt.valid())
210 {
211 scalarField& dgdt = tdgdt.ref();
212
213 forAll(dgdt, celli)
214 {
215 if (dgdt[celli] > 0.0)
216 {
217 Sp[celli] -= dgdt[celli]/max(1 - alpha1[celli], 1e-4);
218 Su[celli] += dgdt[celli]/max(1 - alpha1[celli], 1e-4);
219 }
220 else if (dgdt[celli] < 0.0)
221 {
222 Sp[celli] += dgdt[celli]/max(alpha1[celli], 1e-4);
223 }
224 }
225 }
226
228 (
230 (
231 phi_,
232 alpha1,
233 alphaScheme
234 )
235 + fvc::flux
236 (
237 -fvc::flux(-phir, scalar(1) - alpha1, alpharScheme),
238 alpha1,
240 )
241 );
242
243 phase1_.correctInflowOutflow(alphaPhi1);
244
245 if (nAlphaSubCycles > 1)
246 {
247 tmp<volScalarField> trSubDeltaT;
248
249 if (LTS)
250 {
251 trSubDeltaT =
253 }
254
255 for
256 (
258 !(++alphaSubCycle).end();
259 )
260 {
262
264 (
266 alpha1,
267 phi_,
269 (alphaSubCycle.index()*Sp)(),
270 (Su - (alphaSubCycle.index() - 1)*Sp*alpha1)(),
271 UniformField<scalar>(phase1_.alphaMax()),
272 zeroField()
273 );
274
275 if (alphaSubCycle.index() == 1)
276 {
277 phase1_.alphaPhiRef() = alphaPhi10;
278 }
279 else
280 {
281 phase1_.alphaPhiRef() += alphaPhi10;
282 }
283 }
284
285 phase1_.alphaPhiRef() /= nAlphaSubCycles;
286 }
287 else
288 {
290 (
292 alpha1,
293 phi_,
294 alphaPhi1,
295 Sp,
296 Su,
297 UniformField<scalar>(phase1_.alphaMax()),
298 zeroField()
299 );
300
301 phase1_.alphaPhiRef() = alphaPhi1;
302 }
303
304 if (alphaDbyA.valid())
305 {
306 fvScalarMatrix alpha1Eqn
307 (
309 - fvm::laplacian(alphaDbyA(), alpha1, "bounded")
310 );
311
312 alpha1Eqn.relax();
313 alpha1Eqn.solve();
314
315 phase1_.alphaPhiRef() += alpha1Eqn.flux();
316 }
317
318 phase1_.alphaRhoPhiRef() =
319 fvc::interpolate(phase1_.rho())*phase1_.alphaPhi();
320
321 phase2_.alphaPhiRef() = phi_ - phase1_.alphaPhi();
322 phase2_.correctInflowOutflow(phase2_.alphaPhiRef());
323 phase2_.alphaRhoPhiRef() =
324 fvc::interpolate(phase2_.rho())*phase2_.alphaPhi();
325
326 Info<< alpha1.name() << " volume fraction = "
327 << alpha1.weightedAverage(mesh_.V()).value()
328 << " Min(alpha1) = " << min(alpha1).value()
329 << " Max(alpha1) = " << max(alpha1).value()
330 << endl;
331
332 // Ensure the phase-fractions are bounded
333 alpha1.clamp_range(SMALL, 1 - SMALL);
334
335 // Update the phase-fraction of the other phase
336 alpha2 = scalar(1) - alpha1;
337 }
338}
339
340
341// ************************************************************************* //
MULES: Multidimensional universal limiter for explicit solution.
bool found
const volScalarField & alpha1
const auto & alphaPhi1
surfaceScalarField & phi2
phaseModel & phase1
const volScalarField & alpha2
surfaceScalarField & phi1
phaseModel & phase2
dimensioned< Type > weightedAverage(const DimensionedField< scalar, GeoMesh > &weights, const label comm=UPstream::worldComm) const
Return the global weighted average.
DimensionedField< scalar, volMesh > Internal
void clamp_range(const dimensioned< MinMax< Type > > &range)
Clamp field values (in-place) to the specified range.
void correctBoundaryConditions()
Correct boundary field.
const word & name() const noexcept
Return the object name.
Definition IOobjectI.H:205
A class representing the concept of a uniform field which stores only the single value and providing ...
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition fvMatrix.C:1094
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > flux() const
Return the face-flux field from the matrix.
Definition fvMatrix.C:1415
SolverPerformance< Type > solve(const dictionary &)
Solve returning the solution statistics.
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
static tmp< volScalarField > localRSubDeltaT(const fvMesh &mesh, const label nAlphaSubCycles)
Calculate and return the reciprocal of the local sub-cycling.
static bool enabled(const fvMesh &mesh)
Return true if LTS is enabled.
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
An ordered or unorder pair of phase names. Typically specified as follows.
Class to represent a system of phases and model interfacial transfers between them.
Definition phaseSystem.H:72
phaseModelList phaseModels_
Phase models.
virtual const HashPtrTable< surfaceScalarField > & DByAfs() const =0
Return the phase diffusivities divided by the momentum.
phaseSystem(const fvMesh &mesh)
Construct from fvMesh.
bool end() const
Return true if the number of sub-cycles has been reached.
Perform a subCycleTime on a field.
Definition subCycle.H:137
A class for managing temporary objects.
Definition tmp.H:75
bool valid() const noexcept
Identical to good(), or bool operator.
Definition tmp.H:481
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition tmpI.H:235
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.
const volScalarField & dgdt() const
Return the dilatation term.
phaseModel & phase2_
Phase model 2.
const fvMesh & mesh() const
Return the mesh.
tmp< volScalarField > sigma() const
Return the surface tension coefficient.
twoPhaseSystem(const fvMesh &)
Construct from fvMesh.
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:66
A class representing the concept of a field of 0 used to avoid unnecessary manipulations for objects ...
Definition zeroField.H:50
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
alphaPhi10
Definition alphaEqn.H:7
dynamicFvMesh & mesh
bool LTS
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.
auto & name
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:40
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcDiv.C:42
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition fvcDdt.C:40
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
Definition fvcFlux.C:27
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition fvmDdt.C:41
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:40
const dimensionSet dimless
Dimensionless.
fvMatrix< scalar > fvScalarMatrix
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
void Sp(fvMatrix< typename Expr::value_type > &m, const Expr2 &mult, const Expr &expression)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
void Su(fvMatrix< typename Expr::value_type > &m, const Expr &expression)
#define defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.
volScalarField & e
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
const dictionary & alphaControls
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:299