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 -------------------------------------------------------------------------------
11 License
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 
49 namespace Foam
50 {
51  defineTypeNameAndDebug(twoPhaseSystem, 0);
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
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  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  {
179  surfaceScalarField DbyA
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  (
229  fvc::flux
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  (
257  subCycle<volScalarField> alphaSubCycle(alpha1, nAlphaSubCycles);
258  !(++alphaSubCycle).end();
259  )
260  {
262 
264  (
265  geometricOneField(),
266  alpha1,
267  phi_,
268  alphaPhi10,
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  (
291  geometricOneField(),
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 // ************************************************************************* //
virtual ~twoPhaseSystem()
Destructor.
phaseModel & phase1_
Phase model 1.
void clamp_range(const dimensioned< MinMax< Type >> &range)
Clamp field values (in-place) to the specified range.
bool valid() const noexcept
Identical to good(), or bool operator.
Definition: tmp.H:472
tmp< volScalarField > Kd() const
Return the drag coefficient.
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:37
zeroField Su
Definition: alphaSuSp.H:1
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:195
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition: tmpI.H:235
tmp< volScalarField > sigma() const
Return the surface tension coefficient.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:42
const volScalarField Kd(fluid.Kd())
static tmp< volScalarField > localRSubDeltaT(const fvMesh &mesh, const label nAlphaSubCycles)
Calculate and return the reciprocal of the local sub-cycling.
Definition: localEulerDdt.C:65
word alpharScheme("div(phirb,alpha)")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
Calculate the matrix for the laplacian of the field.
const fvMesh & mesh() const
Return the mesh.
vtkSmartPointer< vtkFloatArray > zeroField(const word &name, const label size)
Create named field initialized to zero.
Calculate the snGrad of the given volField.
const volScalarField & alpha2
const dimensionSet dimless
Dimensionless.
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
surfaceScalarField & phi1
phaseModel & phase2
DimensionedField< scalar, volMesh > Internal
The internal field type from which this GeometricField is derived.
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
void setFluxRequired(const word &name) const
Set flux-required for given name (mutable)
Class which solves the volume fraction equations for two phases.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:40
virtual void solve()
Solve for the phase fractions.
Calculate the first temporal derivative.
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
phaseModel & phase2_
Phase model 2.
dynamicFvMesh & mesh
static bool enabled(const fvMesh &mesh)
Return true if LTS is enabled.
Definition: localEulerDdt.C:32
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
const dictionary & alphaControls
Definition: alphaControls.H:1
surfaceScalarField phir(fvc::flux(UdmModel.Udm()))
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:69
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Calculate the face-flux of the given field.
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:41
Calculate the matrix for the first temporal derivative.
tmp< volScalarField > Vm() const
Return the virtual mass coefficient.
twoPhaseSystem(const fvMesh &)
Construct from fvMesh.
Calculate the field for explicit evaluation of implicit and explicit sources.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
surfaceScalarField & phi2
Calculate the divergence of the given field.
An ordered or unorder pair of phase names. Typically specified as follows.
Definition: phasePairKey.H:61
const surfaceScalarField Kdf("Kdf", fluid.Kdf())
phaseModel & phase1
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
bool LTS
Definition: createRDeltaT.H:1
constexpr auto end(C &c) -> decltype(c.end())
Return iterator to the end of the container c.
Definition: stdFoam.H:201
tmp< surfaceScalarField > Kdf() const
Return the face drag coefficient.
defineTypeNameAndDebug(combustionModel, 0)
alphaPhi10
Definition: alphaEqn.H:7
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.
label nAlphaSubCycles(alphaControls.get< label >("nAlphaSubCycles"))
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:41
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
void correctBoundaryConditions()
Correct boundary field.
messageStream Info
Information stream (stdout output on master, null elsewhere)
MULES: Multidimensional universal limiter for explicit solution.
dimensioned< Type > weightedAverage(const DimensionedField< scalar, GeoMesh > &weightField) const
Calculate and return weighted average.
A class for managing temporary objects.
Definition: HashPtrTable.H:50
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
const surfaceScalarField & alphaPhi1
void explicitSolve(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiPsi, const SpType &Sp, const SuType &Su)
label nAlphaCorr(alphaControls.get< label >("nAlphaCorr"))
bool found
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:40
Calculate the finiteVolume matrix for implicit and explicit sources.
Namespace for OpenFOAM.
zeroField Sp
Definition: alphaSuSp.H:2
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127
const volScalarField & alpha1