MultiComponentPhaseModel.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 
29 
30 #include "multiphaseInterSystem.H"
31 #include "fvmDdt.H"
32 #include "fvmDiv.H"
33 #include "fvmSup.H"
34 #include "fvmLaplacian.H"
35 #include "fvcDdt.H"
36 #include "fvcDiv.H"
37 #include "fvcDDt.H"
38 #include "fvMatrix.H"
39 #include "fvcFlux.H"
40 #include "CMULES.H"
41 #include "subCycle.H"
42 
43 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44 
45 template<class BasePhaseModel, class phaseThermo>
48 (
50  const word& phaseName
51 )
52 :
53  BasePhaseModel(fluid, phaseName),
54  species_(),
55  inertIndex_(-1),
56  addDiffusion_(false)
57 {
59  (
61  (
62  fluid.mesh(),
63  phaseName,
65  )
66  );
67 
68  if (thermoPtr_->composition().species().empty())
69  {
71  << " The selected thermo is pure. Use a multicomponent thermo."
72  << exit(FatalError);
73  }
74 
75  species_ = thermoPtr_->composition().species();
76 
77  inertIndex_ = species_.find(thermoPtr_().template get<word>("inertSpecie"));
78 
80  thermoPtr_().template getOrDefault<bool>("addDiffusion", false);
81 
82  Sct_ = thermoPtr_().template getOrDefault<scalar>("Sct", 1.0);
83 
84  X_.setSize(thermoPtr_->composition().species().size());
85 
86  // Initiate X's using Y's to set BC's
87  forAll(species_, i)
88  {
89  X_.set
90  (
91  i,
92  new volScalarField
93  (
94  IOobject
95  (
96  IOobject::groupName("X" + species_[i], phaseName),
97  fluid.mesh().time().timeName(),
98  fluid.mesh(),
101  ),
102  Y()[i]
103  )
104  );
105  }
106 
107  // Init vol fractions from mass fractions
109 }
111 
112 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
113 
114 
115 template<class BasePhaseModel, class phaseThermo>
118 {
119  volScalarField Xtotal(0.0*X_[0]);
120  const volScalarField W(thermo().W());
121 
122  forAll(X_, i)
123  {
124  const dimensionedScalar Wi
125  (
126  "Wi",
128  thermo().composition().W(i)
129  );
130 
131  if (i != inertIndex_)
132  {
133  X_[i] = W*Y()[i]/Wi;
134  Xtotal += X_[i];
136 
137  const volScalarField::Boundary& YBf = Y()[i].boundaryField();
138 
139  forAll(YBf, patchi)
140  {
141  const fvPatchScalarField& YPf = YBf[patchi];
142  if (YPf.fixesValue())
143  {
144  scalarField& xbf = X_[i].boundaryFieldRef()[patchi];
145  const scalarField& ybf = Y()[i].boundaryField()[patchi];
146  const scalarField& Wbf = W.boundaryField()[patchi];
147  forAll(xbf, facei)
148  {
149  xbf[facei] = Wbf[facei]*ybf[facei]/Wi.value();
150  }
151  }
152  }
153  }
154  }
155  X_[inertIndex_] = 1.0 - Xtotal;
156  X_[inertIndex_].correctBoundaryConditions();
157 }
158 
159 
160 template<class BasePhaseModel, class phaseThermo>
163 {
164  volScalarField W(X_[0]*thermo().composition().W(0));
165  for(label i=1; i< species_.size(); i++)
166  {
167  W += X_[i]*thermo().composition().W(i);
168  }
169 
170  forAll(Y(), i)
171  {
172  Y()[i] = X_[i]*thermo().composition().W(i)/W;
173 
174  Info<< Y()[i].name() << " mass fraction = "
175  << " Min(Y) = " << min(Y()[i]).value()
176  << " Max(Y) = " << max(Y()[i]).value()
177  << endl;
178  }
179 }
180 
181 
182 template<class BasePhaseModel, class phaseThermo>
183 const phaseThermo&
185 {
186  return thermoPtr_();
187 }
188 
189 
190 template<class BasePhaseModel, class phaseThermo>
191 phaseThermo&
193 {
194  return thermoPtr_();
195 }
196 
197 
198 template<class BasePhaseModel, class phaseThermo>
200 {
201  return thermo().correct();
202 }
203 
204 
205 template<class BasePhaseModel, class phaseThermo>
207 (
210 )
211 {
212  const volScalarField& alpha1 = *this;
213 
214  const fvMesh& mesh = alpha1.mesh();
215 
216  const dictionary& MULEScontrols = mesh.solverDict(alpha1.name());
217 
218  scalar cAlpha(MULEScontrols.get<scalar>("cYi"));
219 
220  PtrList<surfaceScalarField> phiYiCorrs(species_.size());
221  const surfaceScalarField& phi = this->fluid().phi();
222 
224 
225  phic = min(cAlpha*phic, max(phic));
226 
228 
230  (
232  )
233  {
234  const volScalarField& alpha2 = iter2()();
235  if (&alpha2 == &alpha1)
236  {
237  continue;
238  }
239 
240  phir += phic*this->fluid().nHatf(alpha1, alpha2);
241  }
242 
243  // Do not compress interface at non-coupled boundary faces
244  // (inlets, outlets etc.)
246  forAll(phir.boundaryField(), patchi)
247  {
248  fvsPatchScalarField& phirp = phirBf[patchi];
249 
250  if (!phirp.coupled())
251  {
252  phirp == 0;
253  }
254  }
255 
256  word phirScheme("div(Yiphir," + alpha1.name() + ')');
257 
258  forAll(X_, i)
259  {
260  if (inertIndex_ != i)
261  {
262  volScalarField& Yi = X_[i];
263 
264  phiYiCorrs.set
265  (
266  i,
268  (
269  "phi" + Yi.name() + "Corr",
270  fvc::flux
271  (
272  phi,
273  Yi,
274  "div(phi," + Yi.name() + ')'
275  )
276  )
277  );
278 
279  surfaceScalarField& phiYiCorr = phiYiCorrs[i];
280 
282  (
284  this->fluid().phases(),
285  iter2
286  )
287  {
288  //const volScalarField& alpha2 = iter2()().oldTime();
289  const volScalarField& alpha2 = iter2()();
290 
291  if (&alpha2 == &alpha1)
292  {
293  continue;
294  }
295 
296  phiYiCorr += fvc::flux
297  (
298  -fvc::flux(-phir, alpha2, phirScheme),
299  Yi,
300  phirScheme
301  );
302  }
303 
304  // Ensure that the flux at inflow BCs is preserved
305  forAll(phiYiCorr.boundaryField(), patchi)
306  {
307  fvsPatchScalarField& phiYiCorrp =
308  phiYiCorr.boundaryFieldRef()[patchi];
309 
310  if (!phiYiCorrp.coupled())
311  {
312  const scalarField& phi1p = phi.boundaryField()[patchi];
313  const scalarField& Yip = Yi.boundaryField()[patchi];
314 
315  forAll(phiYiCorrp, facei)
316  {
317  if (phi1p[facei] < 0)
318  {
319  phiYiCorrp[facei] = Yip[facei]*phi1p[facei];
320  }
321  }
322  }
323  }
324 
326  (
327  1.0/mesh.time().deltaTValue(),
328  geometricOneField(),
329  Yi,
330  phi,
331  phiYiCorr,
332  Sp[i],
333  Su[i],
334  oneField(),
335  zeroField(),
336  true
337  );
338  }
339  }
340 
341  volScalarField Yt(0.0*X_[0]);
342 
343  scalar nYiSubCycles
344  (
345  MULEScontrols.getOrDefault<scalar>("nYiSubCycles", 1)
346  );
347 
348  forAll(X_, i)
349  {
350  if (inertIndex_ != i)
351  {
352  volScalarField& Yi = X_[i];
353 
354  fvScalarMatrix YiEqn
355  (
356  fv::EulerDdtScheme<scalar>(mesh).fvmDdt(Yi)
357  + fv::gaussConvectionScheme<scalar>
358  (
359  mesh,
360  phi,
361  upwind<scalar>(mesh, phi)
362  ).fvmDiv(phi, Yi)
363  ==
364  Su[i] + fvm::Sp(Sp[i], Yi)
365  );
366 
367  YiEqn.solve();
368 
369  surfaceScalarField& phiYiCorr = phiYiCorrs[i];
370 
371  // Add a bounded upwind U-mean flux
372  phiYiCorr += YiEqn.flux();
373 
374  if (nYiSubCycles > 1)
375  {
376  for
377  (
378  subCycle<volScalarField> YiSubCycle(Yi, nYiSubCycles);
379  !(++YiSubCycle).end();
380  )
381  {
383  (
384  geometricOneField(),
385  Yi,
386  phi,
387  phiYiCorr,
388  Sp[i],
389  Su[i],
390  oneField(),
391  zeroField()
392  );
393  }
394  }
395  else
396  {
398  (
399  geometricOneField(),
400  Yi,
401  phi,
402  phiYiCorr,
403  Sp[i],
404  Su[i],
405  oneField(),
406  zeroField()
407  );
408  }
409 
410  if (addDiffusion_)
411  {
412  const volScalarField& alpha = *this;
413  fvScalarMatrix YiDiffEqn
414  (
415  fvm::ddt(Yi) - fvc::ddt(Yi)
417  (
418  alpha*this->fluid().turbulence()->nut()/Sct_,
419  Yi
420  )
421  );
422 
423  YiDiffEqn.solve("diffusion" + Yi.name());
424  }
425 
426  Yt += Yi;
427  }
428  }
429 
430  X_[inertIndex_] = scalar(1) - Yt;
431  X_[inertIndex_].clamp_min(0);
433  calculateMassFractions();
434 }
435 
436 
437 template<class BasePhaseModel, class phaseThermo>
440 {
441  return thermoPtr_->composition().Y();
442 }
443 
444 
445 template<class BasePhaseModel, class phaseThermo>
448 {
449  return thermoPtr_->composition().Y();
450 }
451 
452 
453 template<class BasePhaseModel, class phaseThermo>
454 Foam::label
456 {
457  return inertIndex_;
458 }
459 
460 
461 // ************************************************************************* //
word dictName() const
The local dictionary name (final part of scoped name)
Definition: dictionaryI.H:53
void clamp_min(const Type &lower)
Impose lower (floor) clamp on the field values (in-place)
twoPhaseSystem & fluid
fvsPatchField< scalar > fvsPatchScalarField
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)
label inertIndex_
Inert species index.
const Type & value() const noexcept
Return const reference to value.
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:37
virtual void solveYi(PtrList< volScalarField::Internal > &, PtrList< volScalarField::Internal > &)
Solve species fraction equation.
multiphaseSystem::phaseModelList & phases
Definition: createFields.H:12
static word phasePropertyName(const word &name, const word &phaseName)
The phase property name as property.phase.
Definition: basicThermo.H:317
scalar deltaTValue() const noexcept
Return time step value.
Definition: TimeStateI.H:49
zeroField Su
Definition: alphaSuSp.H:1
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
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:598
basicSpecieMixture & composition
bool empty() const noexcept
Deprecated(2020-07) True if the managed pointer is null.
Definition: autoPtr.H:424
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
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:489
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
Calculate the matrix for the laplacian of the field.
vtkSmartPointer< vtkFloatArray > zeroField(const word &name, const label size)
Create named field initialized to zero.
void calculateVolumeFractions()
Transfor mass fraction into volume fractions.
const volScalarField & alpha2
GeometricBoundaryField< scalar, fvsPatchField, surfaceMesh > Boundary
Type of boundary fields.
Ignore writing from objectRegistry::writeObject()
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.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
Calculate the substantive (total) derivative.
const dictionary & solverDict(const word &name) const
The solver controls dictionary for the given field. Same as solversDict().subDict(...)
Definition: solution.C:428
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object.
Definition: stdFoam.H:478
PtrList< volScalarField > X_
Ptr list of volumetric fractions for species.
virtual void correct()
Correct phase thermo.
virtual const PtrList< volScalarField > & Y() const
Constant access the species mass fractions.
#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
void reset(T *p=nullptr) noexcept
Delete managed object and set to new given pointer.
Definition: autoPtrI.H:37
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
Calculate the first temporal derivative.
psiReactionThermo & thermo
Definition: createFields.H:28
MultiComponentPhaseModel(const multiphaseInterSystem &fluid, const word &phaseName)
HashTable< autoPtr< multiphaseInter::phaseModel > > phaseModelTable
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
autoPtr< phaseThermo > thermoPtr_
Thermophysical model.
dynamicFvMesh & mesh
surfaceScalarField phir(fvc::flux(UdmModel.Udm()))
void calculateMassFractions()
Transfor volume fraction into mass fractions.
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.
bool addDiffusion_
Add diffusion transport on Yi&#39;s Eq.
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:41
Calculate the matrix for the first temporal derivative.
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.
const Mesh & mesh() const noexcept
Return mesh.
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
constexpr auto end(C &c) -> decltype(c.end())
Return iterator to the end of the container c.
Definition: stdFoam.H:201
const dimensionSet dimMoles(0, 0, 0, 0, 1, 0, 0)
Definition: dimensionSets.H:53
Info<< "Reading field U\"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);volScalarField rho(IOobject("rho", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), thermo.rho());volVectorField rhoU(IOobject("rhoU", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *U);volScalarField rhoE(IOobject("rhoE", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *(e+0.5 *magSqr(U)));surfaceScalarField pos(IOobject("pos", runTime.timeName(), mesh), mesh, dimensionedScalar("pos", dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar("neg", dimless, -1.0));surfaceScalarField phi("phi", fvc::flux(rhoU));Info<< "Creating turbulence model\"<< endl;autoPtr< compressible::turbulenceModel > turbulence(compressible::turbulenceModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
surfaceScalarField phic(mixture.cAlpha() *mag(alphaPhic/mesh.magSf()))
Calculate the matrix for the divergence of the given field and flux.
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:41
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
PtrList< volScalarField > & Y
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
virtual bool fixesValue() const
True if the patch field fixes a value.
Definition: fvPatchField.H:245
Nothing to be read.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
void correctBoundaryConditions()
Correct boundary field.
messageStream Info
Information stream (stdout output on master, null elsewhere)
virtual basicSpecieMixture & composition()=0
Return the composition of the multi-component mixture.
const surfaceScalarField & phi() const
Return the mixture flux.
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
label find(const word &val) const
Find index of the value (searches the hash).
virtual void correct()=0
Update properties.
virtual scalar W(const label speciei) const =0
Molecular weight of the given specie [kg/kmol].
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)
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
virtual const phaseThermo & thermo() const
Access to thermo.
label inertIndex() const
Return inert species index.
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
hashedWordList species_
Species table.
CMULES: Multidimensional universal limiter for explicit corrected implicit solution.
Calculate the finiteVolume matrix for implicit and explicit sources.
scalar nut
volScalarField Yt(0.0 *Y[0])
zeroField Sp
Definition: alphaSuSp.H:2
const volScalarField & alpha1