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 | Copyright (C) 2013-2017 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 "twoPhaseSystem.H"
27 #include "dragModel.H"
28 #include "virtualMassModel.H"
29 
30 #include "MULES.H"
31 #include "subCycle.H"
32 
33 #include "fvcDdt.H"
34 #include "fvcDiv.H"
35 #include "fvcSnGrad.H"
36 #include "fvcFlux.H"
37 #include "fvcSup.H"
38 
39 #include "fvmDdt.H"
40 #include "fvmLaplacian.H"
41 #include "fvmSup.H"
42 
43 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
44 
45 namespace Foam
46 {
47  defineTypeNameAndDebug(twoPhaseSystem, 0);
48  defineRunTimeSelectionTable(twoPhaseSystem, dictionary);
49 }
50 
51 
52 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
53 
55 (
56  const fvMesh& mesh
57 )
58 :
59  phaseSystem(mesh),
60  phase1_(phaseModels_[0]),
61  phase2_(phaseModels_[1])
62 {
63  phase2_.volScalarField::operator=(scalar(1) - phase1_);
64 
65  volScalarField& alpha1 = phase1_;
66  mesh.setFluxRequired(alpha1.name());
67 }
68 
69 
70 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
71 
73 {}
74 
75 
76 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
77 
80 {
81  return sigma
82  (
83  phasePairKey(phase1().name(), phase2().name())
84  );
85 }
86 
87 
88 const Foam::dragModel& Foam::twoPhaseSystem::drag(const phaseModel& phase) const
89 {
90  return lookupSubModel<dragModel>(phase, otherPhase(phase));
91 }
92 
93 
96 {
97  return Kd
98  (
99  phasePairKey(phase1().name(), phase2().name())
100  );
101 }
102 
103 
106 {
107  return Kdf
108  (
109  phasePairKey(phase1().name(), phase2().name())
110  );
111 }
112 
113 
115 Foam::twoPhaseSystem::virtualMass(const phaseModel& phase) const
116 {
117  return lookupSubModel<virtualMassModel>(phase, otherPhase(phase));
118 }
119 
120 
123 {
124  return Vm
125  (
126  phasePairKey(phase1().name(), phase2().name())
127  );
128 }
129 
130 
133 {
134  return Vmf
135  (
136  phasePairKey(phase1().name(), phase2().name())
137  );
138 }
139 
140 
143 {
144  return F
145  (
146  phasePairKey(phase1().name(), phase2().name())
147  );
148 }
149 
150 
153 {
154  return Ff
155  (
156  phasePairKey(phase1().name(), phase2().name())
157  );
158 }
159 
160 
163 {
164  return D
165  (
166  phasePairKey(phase1().name(), phase2().name())
167  );
168 }
169 
170 
172 {
173  return transfersMass(phase1());
174 }
175 
176 
179 {
180  return dmdt
181  (
182  phasePairKey(phase1().name(), phase2().name())
183  );
184 }
185 
186 
188 {
189  const Time& runTime = mesh_.time();
190 
191  volScalarField& alpha1 = phase1_;
192  volScalarField& alpha2 = phase2_;
193 
194  const dictionary& alphaControls = mesh_.solverDict(alpha1.name());
195 
196  label nAlphaSubCycles(readLabel(alphaControls.lookup("nAlphaSubCycles")));
197  label nAlphaCorr(readLabel(alphaControls.lookup("nAlphaCorr")));
198 
199  bool LTS = fv::localEulerDdt::enabled(mesh_);
200 
201  word alphaScheme("div(phi," + alpha1.name() + ')');
202  word alpharScheme("div(phir," + alpha1.name() + ')');
203 
204  const surfaceScalarField& phi1 = phase1_.phi();
205  const surfaceScalarField& phi2 = phase2_.phi();
206 
207  // Construct the dilatation rate source term
208  tmp<volScalarField::Internal> tdgdt;
209 
210  if (phase1_.divU().valid() && phase2_.divU().valid())
211  {
212  tdgdt =
213  (
214  alpha2()
215  *phase1_.divU()()()
216  - alpha1()
217  *phase2_.divU()()()
218  );
219  }
220  else if (phase1_.divU().valid())
221  {
222  tdgdt =
223  (
224  alpha2()
225  *phase1_.divU()()()
226  );
227  }
228  else if (phase2_.divU().valid())
229  {
230  tdgdt =
231  (
232  - alpha1()
233  *phase2_.divU()()()
234  );
235  }
236 
237  alpha1.correctBoundaryConditions();
238 
239  surfaceScalarField phir("phir", phi1 - phi2);
240 
241  tmp<surfaceScalarField> alphaDbyA;
242 
243  if (notNull(phase1_.DbyA()) && notNull(phase2_.DbyA()))
244  {
245  surfaceScalarField DbyA(phase1_.DbyA() + phase2_.DbyA());
246 
247  alphaDbyA =
248  fvc::interpolate(max(alpha1, scalar(0)))
249  *fvc::interpolate(max(alpha2, scalar(0)))
250  *DbyA;
251 
252  phir += DbyA*fvc::snGrad(alpha1, "bounded")*mesh_.magSf();
253  }
254 
255  for (int acorr=0; acorr<nAlphaCorr; acorr++)
256  {
258  (
259  IOobject
260  (
261  "Sp",
262  runTime.timeName(),
263  mesh_
264  ),
265  mesh_,
266  dimensionedScalar("Sp", dimless/dimTime, 0.0)
267  );
268 
270  (
271  IOobject
272  (
273  "Su",
274  runTime.timeName(),
275  mesh_
276  ),
277  // Divergence term is handled explicitly to be
278  // consistent with the explicit transport solution
279  fvc::div(phi_)*min(alpha1, scalar(1))
280  );
281 
282  if (tdgdt.valid())
283  {
284  scalarField& dgdt = tdgdt.ref();
285 
286  forAll(dgdt, celli)
287  {
288  if (dgdt[celli] > 0.0)
289  {
290  Sp[celli] -= dgdt[celli]/max(1.0 - alpha1[celli], 1e-4);
291  Su[celli] += dgdt[celli]/max(1.0 - alpha1[celli], 1e-4);
292  }
293  else if (dgdt[celli] < 0.0)
294  {
295  Sp[celli] += dgdt[celli]/max(alpha1[celli], 1e-4);
296  }
297  }
298  }
299 
301  (
302  fvc::flux
303  (
304  phi_,
305  alpha1,
306  alphaScheme
307  )
308  + fvc::flux
309  (
310  -fvc::flux(-phir, scalar(1) - alpha1, alpharScheme),
311  alpha1,
313  )
314  );
315 
316  phase1_.correctInflowOutflow(alphaPhi1);
317 
318  if (nAlphaSubCycles > 1)
319  {
320  tmp<volScalarField> trSubDeltaT;
321 
322  if (LTS)
323  {
324  trSubDeltaT =
326  }
327 
328  for
329  (
330  subCycle<volScalarField> alphaSubCycle(alpha1, nAlphaSubCycles);
331  !(++alphaSubCycle).end();
332  )
333  {
335 
337  (
338  geometricOneField(),
339  alpha1,
340  phi_,
341  alphaPhi10,
342  (alphaSubCycle.index()*Sp)(),
343  (Su - (alphaSubCycle.index() - 1)*Sp*alpha1)(),
344  phase1_.alphaMax(),
345  0
346  );
347 
348  if (alphaSubCycle.index() == 1)
349  {
350  phase1_.alphaPhi() = alphaPhi10;
351  }
352  else
353  {
354  phase1_.alphaPhi() += alphaPhi10;
355  }
356  }
357 
358  phase1_.alphaPhi() /= nAlphaSubCycles;
359  }
360  else
361  {
363  (
364  geometricOneField(),
365  alpha1,
366  phi_,
367  alphaPhi1,
368  Sp,
369  Su,
370  phase1_.alphaMax(),
371  0
372  );
373 
374  phase1_.alphaPhi() = alphaPhi1;
375  }
376 
377  if (alphaDbyA.valid())
378  {
379  fvScalarMatrix alpha1Eqn
380  (
381  fvm::ddt(alpha1) - fvc::ddt(alpha1)
382  - fvm::laplacian(alphaDbyA(), alpha1, "bounded")
383  );
384 
385  alpha1Eqn.relax();
386  alpha1Eqn.solve();
387 
388  phase1_.alphaPhi() += alpha1Eqn.flux();
389  }
390 
391  phase1_.alphaRhoPhi() =
392  fvc::interpolate(phase1_.rho())*phase1_.alphaPhi();
393 
394  phase2_.alphaPhi() = phi_ - phase1_.alphaPhi();
395  phase2_.correctInflowOutflow(phase2_.alphaPhi());
396  phase2_.alphaRhoPhi() =
397  fvc::interpolate(phase2_.rho())*phase2_.alphaPhi();
398 
399  Info<< alpha1.name() << " volume fraction = "
400  << alpha1.weightedAverage(mesh_.V()).value()
401  << " Min(alpha1) = " << min(alpha1).value()
402  << " Max(alpha1) = " << max(alpha1).value()
403  << endl;
404 
405  // Ensure the phase-fractions are bounded
406  alpha1.maxMin(0, 1);
407 
408  // Update the phase-fraction of the other phase
409  alpha2 = scalar(1) - alpha1;
410  }
411 }
412 
413 
414 // ************************************************************************* //
tmp< volVectorField > F() const
Return the combined force (lift + wall-lubrication)
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
zeroField Su
Definition: alphaSuSp.H:1
const double e
Elementary charge.
Definition: doubleFloat.H:78
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
static tmp< volScalarField > localRSubDeltaT(const fvMesh &mesh, const label nAlphaSubCycles)
Calculate and return the reciprocal of the local sub-cycling.
Definition: localEulerDdt.C:70
word alpharScheme("div(phirb,alpha)")
tmp< volScalarField > sigma() const
Return the surface tension coefficient.
const phaseModel & otherPhase() const
Definition: IATEsource.H:148
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m2/K4].
Calculate the matrix for the laplacian of the field.
Calculate the snGrad of the given volField.
const surfaceScalarField & alphaPhi1
surfaceScalarField alphaPhi10(alphaPhi10Header, phi *fvc::interpolate(alpha1))
DimensionedField< scalar, volMesh > Internal
Type of the internal field from which this GeometricField is derived.
label nAlphaCorr(readLabel(alphaControls.lookup("nAlphaCorr")))
tmp< volScalarField > Kd() const
Return the drag coefficient.
surfaceScalarField Ff(fluid.Ff())
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
surfaceScalarField Kdf("Kdf", fluid.Kdf())
Calculate the first temporal derivative.
bool transfersMass() const
Return true if there is mass transfer.
volScalarField Kd(fluid.Kd())
static bool enabled(const fvMesh &mesh)
Return true if LTS is enabled.
Definition: localEulerDdt.C:37
surfaceScalarField phir(fvc::flux(UdmModel.Udm()))
tmp< volScalarField > dmdt() const
Return the interfacial mass flow rate.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Calculate the face-flux of the given field.
const phaseModel & phase() const
Definition: IATEsource.H:138
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
Calulate the matrix for the first temporal derivative.
twoPhaseSystem(const fvMesh &)
Construct from fvMesh.
Calculate the field for explicit evaluation of implicit and explicit sources.
alpha2
Definition: alphaEqn.H:112
label readLabel(Istream &is)
Definition: label.H:64
phaseModel & phase1
tmp< surfaceScalarField > Kdf() const
Return the face drag coefficient.
tmp< volScalarField > D() const
Return the turbulent diffusivity.
Calculate the divergence of the given field.
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
defineTypeNameAndDebug(combustionModel, 0)
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
volScalarField & alpha1
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
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.
bool notNull(const T &t)
Return true if t is not a reference to the nullObject of type T.
Definition: nullObjectI.H:46
virtual void solve()
Solve for the phase fractions.
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
const dimensionedScalar F
Faraday constant: default SI units: [C/mol].
virtual ~twoPhaseSystem()
Destructor.
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:46
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
tmp< surfaceScalarField > Vmf() const
Return the face virtual mass coefficient.
tmp< surfaceScalarField > Ff() const
Return the combined face-force (lift + wall-lubrication)
void correctBoundaryConditions()
Correct boundary field.
messageStream Info
const dragModel & drag(const phaseModel &phase) const
Return the drag model for the given phase.
const virtualMassModel & virtualMass(const phaseModel &phase) const
Return the virtual mass model for the given phase.
MULES: Multidimensional universal limiter for explicit solution.
surfaceScalarField Vmf("Vmf", fluid.Vmf())
phaseModel & phase2
A class for managing temporary objects.
Definition: PtrList.H:53
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
Definition: fvcFlux.C:32
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
void explicitSolve(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiPsi, const SpType &Sp, const SuType &Su)
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:45
Calculate the matrix for implicit and explicit sources.
Namespace for OpenFOAM.
tmp< volScalarField > Vm() const
Return the virtual mass coefficient.
label nAlphaSubCycles(readLabel(alphaControls.lookup("nAlphaSubCycles")))
zeroField Sp
Definition: alphaSuSp.H:2