multiphaseSystem.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) 2017-2022 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "multiphaseSystem.H"
29 
31 #include "Time.H"
32 #include "subCycle.H"
33 #include "fvcMeshPhi.H"
34 
35 #include "surfaceInterpolate.H"
36 #include "fvcGrad.H"
37 #include "fvcSnGrad.H"
38 #include "fvcDiv.H"
39 #include "fvcDdt.H"
40 #include "fvcFlux.H"
41 #include "fvmDdt.H"
42 #include "fvcAverage.H"
43 #include "fvMatrix.H"
44 #include "fvmSup.H"
45 #include "CMULES.H"
46 
47 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
48 
49 namespace Foam
50 {
51 namespace multiphaseInter
52 {
53  defineTypeNameAndDebug(multiphaseSystem, 0);
55 }
56 }
57 
58 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
59 
61 (
62  const fvMesh& mesh
63 )
64 :
66  cAlphas_(),
67  ddtAlphaMax_(0.0),
68  limitedPhiAlphas_(phaseModels_.size()),
69  Su_(phaseModels_.size()),
70  Sp_(phaseModels_.size())
71 {
72  label phasei = 0;
75  {
76  multiphaseInter::phaseModel& pm = iter()();
77  phases_.set(phasei++, &pm);
78  }
79 
80  mesh.solverDict("alpha").readEntry("cAlphas", cAlphas_);
81 
82  // Initiate Su and Sp
84  {
85  const multiphaseInter::phaseModel& pm = iter()();
86 
87  Su_.emplace
88  (
89  pm.name(),
90  // volScalarField::Internal
91  mesh_.newIOobject("Su" + pm.name()),
92  mesh_,
94  );
95 
96  Sp_.emplace
97  (
98  pm.name(),
99  // volScalarField::Internal
100  mesh_.newIOobject("Sp" + pm.name()),
101  mesh_,
103  );
104  }
105 }
106 
107 
108 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
111 {
112  this->alphaTransfer(Su_, Sp_);
113 }
114 
115 
117 {
118  const dictionary& alphaControls = mesh_.solverDict("alpha");
119  label nAlphaSubCycles(alphaControls.get<label>("nAlphaSubCycles"));
120 
121  volScalarField& alpha = phases_.first();
122 
123  if (nAlphaSubCycles > 1)
124  {
125  surfaceScalarField rhoPhiSum
126  (
127  mesh_.newIOobject("rhoPhiSum"),
128  mesh_,
129  dimensionedScalar(rhoPhi_.dimensions(), Zero)
130  );
131 
132  dimensionedScalar totalDeltaT = mesh_.time().deltaT();
133 
134  for
135  (
137  !(++alphaSubCycle).end();
138  )
139  {
140  solveAlphas();
141  rhoPhiSum += (mesh_.time().deltaT()/totalDeltaT)*rhoPhi_;
142  }
143 
144  rhoPhi_ = rhoPhiSum;
145  }
146  else
147  {
148  solveAlphas();
149  }
150 
151 }
152 
154 {
155 
156  const dictionary& alphaControls = mesh_.solverDict("alpha");
157  alphaControls.readEntry("cAlphas", cAlphas_);
158  label nAlphaCorr(alphaControls.get<label>("nAlphaCorr"));
159 
160  PtrList<surfaceScalarField> phiAlphaCorrs(phases_.size());
161 
162  const surfaceScalarField& phi = this->phi();
163 
164  surfaceScalarField phic(mag((phi)/mesh_.magSf()));
165 
166  // Do not compress interface at non-coupled boundary faces
167  // (inlets, outlets etc.)
168  surfaceScalarField::Boundary& phicBf = phic.boundaryFieldRef();
169  forAll(phic.boundaryField(), patchi)
170  {
171  fvsPatchScalarField& phicp = phicBf[patchi];
172 
173  if (!phicp.coupled())
174  {
175  phicp == 0;
176  }
177  }
178 
179  for (int acorr=0; acorr<nAlphaCorr; acorr++)
180  {
181  label phasei = 0;
182  for (multiphaseInter::phaseModel& phase1 : phases_)
183  {
184  const volScalarField& alpha1 = phase1;
185 
186  phiAlphaCorrs.set
187  (
188  phasei,
190  (
191  "phi" + alpha1.name() + "Corr",
192  fvc::flux
193  (
194  phi,
195  alpha1,
196  "div(phi," + alpha1.name() + ')'
197  )
198  )
199  );
200 
201  surfaceScalarField& phiAlphaCorr = phiAlphaCorrs[phasei];
202 
203  for (multiphaseInter::phaseModel& phase2 : phases_)
204  {
205  const volScalarField& alpha2 = phase2;
206 
207  if (&phase2 == &phase1) continue;
208 
209  const phasePairKey key12(phase1.name(), phase2.name());
210 
211  if (!cAlphas_.found(key12))
212  {
214  << "Phase compression factor (cAlpha) not found for : "
215  << key12
216  << exit(FatalError);
217  }
218  scalar cAlpha = cAlphas_.find(key12)();
219 
220  phic = min(cAlpha*phic, max(phic));
221 
223 
224  word phirScheme
225  (
226  "div(phir," + alpha2.name() + ',' + alpha1.name() + ')'
227  );
228 
229  phiAlphaCorr += fvc::flux
230  (
231  -fvc::flux(-phir, alpha2, phirScheme),
232  alpha1,
233  phirScheme
234  );
235  }
236 
237  // Ensure that the flux at inflow BCs is preserved
238  forAll(phiAlphaCorr.boundaryField(), patchi)
239  {
240  fvsPatchScalarField& phiAlphaCorrp =
241  phiAlphaCorr.boundaryFieldRef()[patchi];
242 
243  if (!phiAlphaCorrp.coupled())
244  {
245  const scalarField& phi1p = phi.boundaryField()[patchi];
246  const scalarField& alpha1p =
247  alpha1.boundaryField()[patchi];
248 
249  forAll(phiAlphaCorrp, facei)
250  {
251  if (phi1p[facei] < 0)
252  {
253  phiAlphaCorrp[facei] = alpha1p[facei]*phi1p[facei];
254  }
255  }
256  }
257  }
258 
259  ++phasei;
260  }
261 
262  // Set Su and Sp to zero
263  for (const multiphaseInter::phaseModel& phase : phases_)
264  {
265  Su_[phase.name()] = dimensionedScalar("Su", dimless/dimTime, Zero);
266  Sp_[phase.name()] = dimensionedScalar("Sp", dimless/dimTime, Zero);
267 
268  // Add alpha*div(U)
269  //const volScalarField& alpha = phase;
270  //Su_[phase.name()] += fvc::div(phi)*clamp(alpha, zero_one{});
271  }
272 
273  // Fill Su and Sp
274  calculateSuSp();
275 
276  // Limit phiAlphaCorr on each phase
277  phasei = 0;
278  for (multiphaseInter::phaseModel& phase : phases_)
279  {
280  volScalarField& alpha1 = phase;
281 
282  surfaceScalarField& phiAlphaCorr = phiAlphaCorrs[phasei];
283 
284  volScalarField::Internal& Su = Su_[phase.name()];
285  volScalarField::Internal& Sp = Sp_[phase.name()];
286 
288  (
289  1.0/mesh_.time().deltaTValue(),
290  geometricOneField(),
291  alpha1,
292  phi,
293  phiAlphaCorr,
294  Sp,
295  Su,
296  oneField(),
297  zeroField(),
298  true
299  );
300  ++phasei;
301  }
302 
303  MULES::limitSum(phiAlphaCorrs);
304 
305  volScalarField sumAlpha
306  (
307  mesh_.newIOobject("sumAlpha"),
308  mesh_,
310  );
311 
312  phasei = 0;
313  for (multiphaseInter::phaseModel& phase : phases_)
314  {
315  volScalarField& alpha1 = phase;
316 
317  const volScalarField::Internal& Su = Su_[phase.name()];
318 
319  const volScalarField::Internal& Sp = Sp_[phase.name()];
320 
321  surfaceScalarField& phiAlpha = phiAlphaCorrs[phasei];
322 
323  // Add a bounded upwind U-mean flux
324  //phiAlpha += upwind<scalar>(mesh_, phi).flux(alpha1);
325  fvScalarMatrix alpha1Eqn
326  (
327  fv::EulerDdtScheme<scalar>(mesh_).fvmDdt(alpha1)
328  + fv::gaussConvectionScheme<scalar>
329  (
330  mesh_,
331  phi,
332  upwind<scalar>(mesh_, phi)
333  ).fvmDiv(phi, alpha1)
334  ==
335  Su + fvm::Sp(Sp, alpha1)
336  );
337 
338  alpha1Eqn.boundaryManipulate(alpha1.boundaryFieldRef());
339 
340  alpha1Eqn.solve();
341 
342  phiAlpha += alpha1Eqn.flux();
343 
345  (
346  geometricOneField(),
347  alpha1,
348  phi,
349  phiAlpha,
350  Sp,
351  Su,
352  oneField(),
353  zeroField()
354  );
355 
356  phase.alphaPhi() = phiAlpha;
357 
358  ++phasei;
359  }
360 
361  if (acorr == nAlphaCorr - 1)
362  {
363  volScalarField sumAlpha
364  (
365  mesh_.newIOobject("sumAlpha"),
366  mesh_,
368  );
369 
370  // Reset rhoPhi
371  rhoPhi_ = dimensionedScalar("rhoPhi", dimMass/dimTime, Zero);
372 
373  for (multiphaseInter::phaseModel& phase : phases_)
374  {
375  volScalarField& alpha1 = phase;
376  sumAlpha += alpha1;
377 
378  // Update rhoPhi
379  rhoPhi_ += fvc::interpolate(phase.rho()) * phase.alphaPhi();
380  }
381 
382  Info<< "Phase-sum volume fraction, min, max = "
383  << sumAlpha.weightedAverage(mesh_.V()).value()
384  << ' ' << min(sumAlpha).value()
385  << ' ' << max(sumAlpha).value()
386  << endl;
387 
388  volScalarField sumCorr(1.0 - sumAlpha);
389 
390  for (multiphaseInter::phaseModel& phase : phases_)
391  {
392  volScalarField& alpha = phase;
393  alpha += alpha*sumCorr;
394 
395  Info<< alpha.name() << " volume fraction = "
396  << alpha.weightedAverage(mesh_.V()).value()
397  << " Min(alpha) = " << min(alpha).value()
398  << " Max(alpha) = " << max(alpha).value()
399  << endl;
400  }
401  }
402  }
403 }
404 
405 
408 {
409  return phases_;
410 }
411 
412 
415 {
416  return phases_;
417 }
418 
419 
422 {
423  return phases_[i];
424 }
425 
426 
429 {
430  return phases_[i];
431 }
432 
433 
436 {
437  return ddtAlphaMax_;
438 }
439 
440 
441 Foam::scalar
443 {
444  auto iter = phaseModels_.cbegin();
445 
446  scalar maxVal = max(iter()->diffNo()).value();
447 
448  for (++iter; iter != phaseModels_.cend(); ++iter)
449  {
450  maxVal = max(maxVal, max(iter()->diffNo()).value());
451  }
452 
453  return maxVal * mesh_.time().deltaTValue();
454 }
455 
456 
459 {
460  return limitedPhiAlphas_;
461 }
462 
463 
466 {
467  return Su_;
468 }
469 
470 
473 {
474  return Sp_;
475 }
476 
477 
479 {
480  return true;
481 }
482 
483 
484 // ************************************************************************* //
const fvMesh & mesh() const
Return mesh.
fvsPatchField< scalar > fvsPatchScalarField
UPtrList< phaseModel > phases_
Unallocated phase list.
void limit(const RdeltaTType &rDeltaT, const RhoType &rho, const volScalarField &psi, const surfaceScalarField &phi, surfaceScalarField &phiPsi, const SpType &Sp, const SuType &Su, const PsiMaxType &psiMax, const PsiMinType &psiMin, const bool returnCorr)
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:37
bool emplace(const Key &key, Args &&... args)
Emplace insert a new entry, not overwriting existing entries.
Definition: HashTableI.H:129
zeroField Su
Definition: alphaSuSp.H:1
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
label phasei
Definition: pEqn.H:27
virtual bool read()
Read thermophysical properties dictionary.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:608
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
std::enable_if< std::is_same< bool, TypeT >::value, bool >::type set(const label i, bool val=true)
A bitSet::set() method for a list of bool.
Definition: List.H:493
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
vtkSmartPointer< vtkFloatArray > zeroField(const word &name, const label size)
Create named field initialized to zero.
Calculate the snGrad of the given volField.
const word & name() const
The name of this phase.
Definition: phaseModel.H:122
const volScalarField & alpha2
const dimensionSet dimless
Dimensionless.
const dictionary & solverDict(const word &name) const
The solver controls dictionary for the given field. Same as solversDict().subDict(...)
Definition: solution.C:467
phaseModel & phase2
DimensionedField< scalar, volMesh > Internal
The internal field type from which this GeometricField is derived.
Area-weighted average a surfaceField creating a volField.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:72
Calculate the first temporal derivative.
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, IOobjectOption::readOption readOpt=IOobjectOption::MUST_READ) const
Find entry and assign to T val. FatalIOError if it is found and the number of tokens is incorrect...
label size() const noexcept
The number of elements in table.
Definition: HashTable.H:358
const fvMesh & mesh_
Reference to the mesh.
virtual void solve()
Solve for the phase fractions.
dynamicFvMesh & mesh
const dictionary & alphaControls
Definition: alphaControls.H:1
void limitSum(UPtrList< scalarField > &phiPsiCorrs)
Definition: MULES.C:27
multiphaseSystem(const fvMesh &)
Construct from fvMesh.
const phaseModel & phase(const label i) const
Constant access phase model i.
Calculate the gradient of the given field.
surfaceScalarField phir(fvc::flux(UdmModel.Udm()))
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Calculate the face-flux of the given field.
Calculate the matrix for the first temporal derivative.
#define forAllIters(container, iter)
Iterate across all elements in the container object.
Definition: stdFoam.H:336
dimensionedScalar ddtAlphaMax() const
Access to ddtAlphaMax.
Perform a subCycleTime on a field.
Definition: subCycle.H:133
const compressionFluxTable & limitedPhiAlphas() const
Access to compression fluxes for phaes.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
Calculate the divergence of the given field.
defineTypeNameAndDebug(interfaceCompositionModel, 0)
phaseModel & phase1
constexpr auto end(C &c) -> decltype(c.end())
Return iterator to the end of the container c.
Definition: stdFoam.H:201
scalarTable cAlphas_
Table for compression factors between phases.
virtual bool coupled() const
True if the patch field is coupled.
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.
Calculate the mesh motion flux and convert fluxes from absolute to relative and back.
surfaceScalarField::Boundary & phicBf
Definition: alphaEqn.H:75
const word & name() const noexcept
Return const reference to name.
surfaceScalarField phic(mixture.cAlpha() *mag(alphaPhic/mesh.magSf()))
SuSpTable Su_
Su phase source terms.
phaseModelTable phaseModels_
Phase models.
label nAlphaSubCycles(alphaControls.get< label >("nAlphaSubCycles"))
SuSpTable Sp_
Sp phase source terms.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
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
scalar maxDiffNo() const
Maximum diffusion number.
messageStream Info
Information stream (stdout output on master, null elsewhere)
defineRunTimeSelectionTable(interfaceCompositionModel, dictionary)
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
IOobject newIOobject(const word &name, IOobjectOption ioOpt) const
Create an IOobject at the current time instance (timeName) with the specified options.
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
void setSize(const label n)
Alias for resize()
Definition: UPtrList.H:840
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
void explicitSolve(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiPsi, const SpType &Sp, const SuType &Su)
const UPtrList< phaseModel > & phases() const
Return phases.
label nAlphaCorr(alphaControls.get< label >("nAlphaCorr"))
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
void calculateSuSp()
Calculate Sp and Su.
CMULES: Multidimensional universal limiter for explicit corrected implicit solution.
Calculate the finiteVolume matrix for implicit and explicit sources.
Namespace for OpenFOAM.
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie...
Definition: UPtrList.H:366
zeroField Sp
Definition: alphaSuSp.H:2
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127
const volScalarField & alpha1