twoPhaseSystem.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2013-2018 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 #include "UniformField.H"
33 
34 #include "fvcDdt.H"
35 #include "fvcDiv.H"
36 #include "fvcSnGrad.H"
37 #include "fvcFlux.H"
38 #include "fvcSup.H"
39 
40 #include "fvmDdt.H"
41 #include "fvmLaplacian.H"
42 #include "fvmSup.H"
43 
44 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
45 
46 namespace Foam
47 {
48  defineTypeNameAndDebug(twoPhaseSystem, 0);
49  defineRunTimeSelectionTable(twoPhaseSystem, dictionary);
50 }
51 
52 
53 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
54 
56 (
57  const fvMesh& mesh
58 )
59 :
60  phaseSystem(mesh),
61  phase1_(phaseModels_[0]),
62  phase2_(phaseModels_[1])
63 {
64  phase2_.volScalarField::operator=(scalar(1) - phase1_);
65 
66  volScalarField& alpha1 = phase1_;
67  mesh.setFluxRequired(alpha1.name());
68 }
69 
70 
71 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
72 
74 {}
75 
76 
77 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
78 
81 {
82  return sigma
83  (
84  phasePairKey(phase1().name(), phase2().name())
85  );
86 }
87 
88 
91 {
92  return Kd
93  (
94  phasePairKey(phase1().name(), phase2().name())
95  );
96 }
97 
98 
101 {
102  return Kdf
103  (
104  phasePairKey(phase1().name(), phase2().name())
105  );
106 }
107 
108 
111 {
112  return Vm
113  (
114  phasePairKey(phase1().name(), phase2().name())
115  );
116 }
117 
118 
120 {
121  const Time& runTime = mesh_.time();
122 
123  volScalarField& alpha1 = phase1_;
124  volScalarField& alpha2 = phase2_;
125 
126  const dictionary& alphaControls = mesh_.solverDict(alpha1.name());
127 
128  label nAlphaSubCycles(readLabel(alphaControls.lookup("nAlphaSubCycles")));
129  label nAlphaCorr(readLabel(alphaControls.lookup("nAlphaCorr")));
130 
131  bool LTS = fv::localEulerDdt::enabled(mesh_);
132 
133  word alphaScheme("div(phi," + alpha1.name() + ')');
134  word alpharScheme("div(phir," + alpha1.name() + ')');
135 
136  const surfaceScalarField& phi1 = phase1_.phi();
137  const surfaceScalarField& phi2 = phase2_.phi();
138 
139  // Construct the dilatation rate source term
140  tmp<volScalarField::Internal> tdgdt;
141 
142  if (phase1_.divU().valid() && phase2_.divU().valid())
143  {
144  tdgdt =
145  (
146  alpha2()
147  *phase1_.divU()()()
148  - alpha1()
149  *phase2_.divU()()()
150  );
151  }
152  else if (phase1_.divU().valid())
153  {
154  tdgdt =
155  (
156  alpha2()
157  *phase1_.divU()()()
158  );
159  }
160  else if (phase2_.divU().valid())
161  {
162  tdgdt =
163  (
164  - alpha1()
165  *phase2_.divU()()()
166  );
167  }
168 
169  alpha1.correctBoundaryConditions();
170 
171  surfaceScalarField phir("phir", phi1 - phi2);
172 
173  tmp<surfaceScalarField> alphaDbyA;
174  if (DByAfs().found(phase1_.name()) && DByAfs().found(phase2_.name()))
175  {
176  surfaceScalarField DbyA
177  (
178  *DByAfs()[phase1_.name()] + *DByAfs()[phase2_.name()]
179  );
180 
181  alphaDbyA =
182  fvc::interpolate(max(alpha1, scalar(0)))
183  *fvc::interpolate(max(alpha2, scalar(0)))
184  *DbyA;
185 
186  phir += DbyA*fvc::snGrad(alpha1, "bounded")*mesh_.magSf();
187  }
188 
189  for (int acorr=0; acorr<nAlphaCorr; acorr++)
190  {
192  (
193  IOobject
194  (
195  "Sp",
196  runTime.timeName(),
197  mesh_
198  ),
199  mesh_,
201  );
202 
204  (
205  IOobject
206  (
207  "Su",
208  runTime.timeName(),
209  mesh_
210  ),
211  // Divergence term is handled explicitly to be
212  // consistent with the explicit transport solution
213  fvc::div(phi_)*min(alpha1, scalar(1))
214  );
215 
216  if (tdgdt.valid())
217  {
218  scalarField& dgdt = tdgdt.ref();
219 
220  forAll(dgdt, celli)
221  {
222  if (dgdt[celli] > 0.0)
223  {
224  Sp[celli] -= dgdt[celli]/max(1 - alpha1[celli], 1e-4);
225  Su[celli] += dgdt[celli]/max(1 - alpha1[celli], 1e-4);
226  }
227  else if (dgdt[celli] < 0.0)
228  {
229  Sp[celli] += dgdt[celli]/max(alpha1[celli], 1e-4);
230  }
231  }
232  }
233 
235  (
236  fvc::flux
237  (
238  phi_,
239  alpha1,
240  alphaScheme
241  )
242  + fvc::flux
243  (
244  -fvc::flux(-phir, scalar(1) - alpha1, alpharScheme),
245  alpha1,
247  )
248  );
249 
250  phase1_.correctInflowOutflow(alphaPhi1);
251 
252  if (nAlphaSubCycles > 1)
253  {
254  tmp<volScalarField> trSubDeltaT;
255 
256  if (LTS)
257  {
258  trSubDeltaT =
260  }
261 
262  for
263  (
264  subCycle<volScalarField> alphaSubCycle(alpha1, nAlphaSubCycles);
265  !(++alphaSubCycle).end();
266  )
267  {
269 
271  (
272  geometricOneField(),
273  alpha1,
274  phi_,
275  alphaPhi10,
276  (alphaSubCycle.index()*Sp)(),
277  (Su - (alphaSubCycle.index() - 1)*Sp*alpha1)(),
278  UniformField<scalar>(phase1_.alphaMax()),
279  zeroField()
280  );
281 
282  if (alphaSubCycle.index() == 1)
283  {
284  phase1_.alphaPhiRef() = alphaPhi10;
285  }
286  else
287  {
288  phase1_.alphaPhiRef() += alphaPhi10;
289  }
290  }
291 
292  phase1_.alphaPhiRef() /= nAlphaSubCycles;
293  }
294  else
295  {
297  (
298  geometricOneField(),
299  alpha1,
300  phi_,
301  alphaPhi1,
302  Sp,
303  Su,
304  UniformField<scalar>(phase1_.alphaMax()),
305  zeroField()
306  );
307 
308  phase1_.alphaPhiRef() = alphaPhi1;
309  }
310 
311  if (alphaDbyA.valid())
312  {
313  fvScalarMatrix alpha1Eqn
314  (
315  fvm::ddt(alpha1) - fvc::ddt(alpha1)
316  - fvm::laplacian(alphaDbyA(), alpha1, "bounded")
317  );
318 
319  alpha1Eqn.relax();
320  alpha1Eqn.solve();
321 
322  phase1_.alphaPhiRef() += alpha1Eqn.flux();
323  }
324 
325  phase1_.alphaRhoPhiRef() =
326  fvc::interpolate(phase1_.rho())*phase1_.alphaPhi();
327 
328  phase2_.alphaPhiRef() = phi_ - phase1_.alphaPhi();
329  phase2_.correctInflowOutflow(phase2_.alphaPhiRef());
330  phase2_.alphaRhoPhiRef() =
331  fvc::interpolate(phase2_.rho())*phase2_.alphaPhi();
332 
333  Info<< alpha1.name() << " volume fraction = "
334  << alpha1.weightedAverage(mesh_.V()).value()
335  << " Min(alpha1) = " << min(alpha1).value()
336  << " Max(alpha1) = " << max(alpha1).value()
337  << endl;
338 
339  // Ensure the phase-fractions are bounded
340  alpha1.maxMin(0, 1);
341 
342  // Update the phase-fraction of the other phase
343  alpha2 = scalar(1) - alpha1;
344  }
345 }
346 
347 
348 // ************************************************************************* //
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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
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.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m^2/K^4].
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.
const volScalarField Kd(fluid.Kd())
label nAlphaCorr(readLabel(alphaControls.lookup("nAlphaCorr")))
tmp< volScalarField > Kd() const
Return the drag coefficient.
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
Calculate the first temporal derivative.
static bool enabled(const fvMesh &mesh)
Return true if LTS is enabled.
Definition: localEulerDdt.C:37
surfaceScalarField phir(fvc::flux(UdmModel.Udm()))
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Calculate the face-flux of the given field.
const volScalarField & alpha1
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
Calculate 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:115
label readLabel(Istream &is)
Definition: label.H:64
phaseModel & phase1
tmp< surfaceScalarField > Kdf() const
Return the face drag coefficient.
Calculate the divergence of the given field.
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
defineTypeNameAndDebug(combustionModel, 0)
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
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.
virtual void solve()
Solve for the phase fractions.
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
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
void correctBoundaryConditions()
Correct boundary field.
const surfaceScalarField Kdf("Kdf", fluid.Kdf())
messageStream Info
MULES: Multidimensional universal limiter for explicit solution.
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:98
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)
bool found
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