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-2020 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 "phaseSystem.H"
31 #include "multiphaseSystem.H"
32 #include "fvmDdt.H"
33 #include "fvmDiv.H"
34 #include "fvmSup.H"
35 #include "fvmLaplacian.H"
36 #include "fvcDdt.H"
37 #include "fvcDiv.H"
38 #include "fvcDDt.H"
39 #include "fvMatrix.H"
40 #include "fvcFlux.H"
41 #include "CMULES.H"
42 #include "subCycle.H"
43 
44 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
45 
46 template<class BasePhaseModel, class phaseThermo>
49 (
50  const phaseSystem& fluid,
51  const word& phaseName
52 )
53 :
54  BasePhaseModel(fluid, phaseName),
55  species_(),
56  inertIndex_(-1)
57 {
58  thermoPtr_.set
59  (
61  (
62  fluid.mesh(),
63  phaseName,
64  basicThermo::phasePropertyName(basicThermo::dictName, phaseName)
65  ).ptr()
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_[thermoPtr_().template get<word>("inertSpecie")];
78 
79  X_.setSize(thermoPtr_->composition().species().size());
80 
81  // Initiate X's using Y's to set BC's
82  forAll(species_, i)
83  {
84  X_.set
85  (
86  i,
87  new volScalarField
88  (
89  IOobject
90  (
91  IOobject::groupName("X" + species_[i], phaseName),
92  fluid.mesh().time().timeName(),
93  fluid.mesh()
94  ),
95  Y()[i]
96  )
97  );
98  }
99 
100  // Init vol fractions from mass fractions
101  calculateVolumeFractions();
102 }
103 
104 
105 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
106 
107 
108 template<class BasePhaseModel, class phaseThermo>
111 {
112  volScalarField Xtotal(0.0*X_[0]);
113  const volScalarField W(thermo().W());
114 
115  forAll(X_, i)
116  {
117  const dimensionedScalar Wi
118  (
119  "W",
121  thermo().composition().W(i)
122  );
123 
124  if (i != inertIndex_)
125  {
126  X_[i] = W*Y()[i]/Wi;
127  Xtotal += X_[i];
129  }
130  }
131  X_[inertIndex_] = 1.0 - Xtotal;
132  X_[inertIndex_].correctBoundaryConditions();
133 }
134 
135 
136 template<class BasePhaseModel, class phaseThermo>
139 {
140  volScalarField W(X_[0]*thermo().composition().W(0));
141  for(label i=1; i< species_.size(); i++)
142  {
143  W += X_[i]*thermo().composition().W(i);
144  }
145 
146  forAll(Y(), i)
147  {
148  Y()[i] = X_[i]*thermo().composition().W(i)/W;
149 
150  Info<< Y()[i].name() << " mass fraction = "
151  << " Min(Y) = " << min(Y()[i]).value()
152  << " Max(Y) = " << max(Y()[i]).value()
153  << endl;
154  }
155 }
156 
157 
158 template<class BasePhaseModel, class phaseThermo>
159 const phaseThermo&
161 {
162  return thermoPtr_();
163 }
164 
165 
166 template<class BasePhaseModel, class phaseThermo>
167 phaseThermo&
169 {
170  return thermoPtr_();
171 }
172 
173 
174 template<class BasePhaseModel, class phaseThermo>
176 {
177  return thermo().correct();
178 }
179 
180 
181 template<class BasePhaseModel, class phaseThermo>
183 (
186 )
187 {
188  const volScalarField& alpha1 = *this;
189 
190  const fvMesh& mesh = alpha1.mesh();
191 
192  const dictionary& MULEScontrols = mesh.solverDict(alpha1.name());
193 
194  scalar cAlpha(MULEScontrols.get<scalar>("cYi"));
195 
196  PtrList<surfaceScalarField> phiYiCorrs(species_.size());
197  const surfaceScalarField& phi = this->fluid().phi();
198 
199  surfaceScalarField phic(mag((phi)/mesh.magSf()));
200 
201  phic = min(cAlpha*phic, max(phic));
202 
204 
206  {
207  const volScalarField& alpha2 = iter2()();
208  if (&alpha2 == &alpha1)
209  {
210  continue;
211  }
212 
213  phir += phic*this->fluid().nHatf(alpha1, alpha2);
214  }
215 
216  // Do not compress interface at non-coupled boundary faces
217  // (inlets, outlets etc.)
218  surfaceScalarField::Boundary& phirBf = phir.boundaryFieldRef();
219  forAll(phir.boundaryField(), patchi)
220  {
221  fvsPatchScalarField& phirp = phirBf[patchi];
222 
223  if (!phirp.coupled())
224  {
225  phirp == 0;
226  }
227  }
228 
229  word phirScheme("div(Yiphir," + alpha1.name() + ')');
230 
231  forAll(X_, i)
232  {
233  if (inertIndex_ != i)
234  {
235  volScalarField& Yi = X_[i];
236 
237  phiYiCorrs.set
238  (
239  i,
241  (
242  "phi" + Yi.name() + "Corr",
243  fvc::flux
244  (
245  phi,
246  Yi,
247  "div(phi," + Yi.name() + ')'
248  )
249  )
250  );
251 
252  surfaceScalarField& phiYiCorr = phiYiCorrs[i];
253 
255  (
256  phaseSystem::phaseModelTable, this->fluid().phases(), iter2
257  )
258  {
259  //const volScalarField& alpha2 = iter2()().oldTime();
260  const volScalarField& alpha2 = iter2()();
261 
262  if (&alpha2 == &alpha1)
263  {
264  continue;
265  }
266 
267  phiYiCorr += fvc::flux
268  (
269  -fvc::flux(-phir, alpha2, phirScheme),
270  Yi,
271  phirScheme
272  );
273  }
274 
275  // Ensure that the flux at inflow BCs is preserved
276  forAll(phiYiCorr.boundaryField(), patchi)
277  {
278  fvsPatchScalarField& phiYiCorrp =
279  phiYiCorr.boundaryFieldRef()[patchi];
280 
281  if (!phiYiCorrp.coupled())
282  {
283  const scalarField& phi1p = phi.boundaryField()[patchi];
284  const scalarField& Yip = Yi.boundaryField()[patchi];
285 
286  forAll(phiYiCorrp, facei)
287  {
288  if (phi1p[facei] < 0)
289  {
290  phiYiCorrp[facei] = Yip[facei]*phi1p[facei];
291  }
292  }
293  }
294  }
295 
297  (
298  1.0/mesh.time().deltaT().value(),
300  Yi,
301  phi,
302  phiYiCorr,
303  Sp[i],
304  Su[i],
305  oneField(),
306  zeroField(),
307  true
308  );
309  }
310  }
311 
312  volScalarField Yt(0.0*X_[0]);
313 
314  scalar nYiSubCycles
315  (
316  MULEScontrols.getOrDefault<scalar>("nYiSubCycles", 1)
317  );
318 
319  forAll(X_, i)
320  {
321  if (inertIndex_ != i)
322  {
323  volScalarField& Yi = X_[i];
324 
325  fvScalarMatrix YiEqn
326  (
327  fv::EulerDdtScheme<scalar>(mesh).fvmDdt(Yi)
329  (
330  mesh,
331  phi,
333  ).fvmDiv(phi, Yi)
334  ==
335  Su[i] + fvm::Sp(Sp[i], Yi)
336  );
337 
338  YiEqn.solve();
339 
340  surfaceScalarField& phiYiCorr = phiYiCorrs[i];
341 
342  // Add a bounded upwind U-mean flux
343  phiYiCorr += YiEqn.flux();
344 
345  if (nYiSubCycles > 1)
346  {
347  for
348  (
349  subCycle<volScalarField> YiSubCycle(Yi, nYiSubCycles);
350  !(++YiSubCycle).end();
351  )
352  {
354  (
356  Yi,
357  phi,
358  phiYiCorr,
359  Sp[i],
360  Su[i],
361  oneField(),
362  zeroField()
363  );
364  }
365  }
366  else
367  {
369  (
371  Yi,
372  phi,
373  phiYiCorr,
374  Sp[i],
375  Su[i],
376  oneField(),
377  zeroField()
378  );
379  }
380  Yt += Yi;
381  }
382  }
383 
384  X_[inertIndex_] = scalar(1) - Yt;
385  X_[inertIndex_].max(0.0);
386 
387  calculateMassFractions();
388 }
389 
390 
391 template<class BasePhaseModel, class phaseThermo>
394 {
395  return thermoPtr_->composition().Y();
396 }
397 
398 
399 template<class BasePhaseModel, class phaseThermo>
402 {
403  return thermoPtr_->composition().Y();
404 }
405 
406 
407 template<class BasePhaseModel, class phaseThermo>
408 Foam::label
410 {
411  return inertIndex_;
412 }
413 
414 
415 // ************************************************************************* //
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::fvc::flux
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
Foam::MULES::explicitSolve
void explicitSolve(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiPsi, const SpType &Sp, const SuType &Su)
Definition: MULESTemplates.C:41
CMULES.H
CMULES: Multidimensional universal limiter for explicit corrected implicit solution.
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::MULES::limit
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)
Definition: MULESTemplates.C:581
Foam::upwind
Upwind differencing scheme class.
Definition: upwind.H:56
alpha2
const volScalarField & alpha2
Definition: setRegionFluidFields.H:9
subCycle.H
thermo
psiReactionThermo & thermo
Definition: createFields.H:28
Foam::vtk::Tools::zeroField
vtkSmartPointer< vtkFloatArray > zeroField(const word &name, const label size)
Create named field initialized to zero.
Definition: foamVtkToolsTemplates.C:327
Foam::fvsPatchField
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:68
Sp
zeroField Sp
Definition: alphaSuSp.H:2
Foam::MultiComponentPhaseModel::solveYi
virtual void solveYi(PtrList< volScalarField::Internal > &, PtrList< volScalarField::Internal > &)
Solve species fraction equation.
Definition: MultiComponentPhaseModel.C:183
Foam::geometricOneField
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
Definition: geometricOneField.H:55
fvcDiv.H
Calculate the divergence of the given field.
dictName
const word dictName("blockMeshDict")
Foam::dimMoles
const dimensionSet dimMoles(0, 0, 0, 0, 1, 0, 0)
Definition: dimensionSets.H:56
fvMatrix.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
alpha1
const volScalarField & alpha1
Definition: setRegionFluidFields.H:8
Foam::dictionary::get
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:81
fluid
twoPhaseSystem & fluid
Definition: setRegionFluidFields.H:3
Foam::MultiComponentPhaseModel::Y
virtual const PtrList< volScalarField > & Y() const
Constant access the species mass fractions.
Definition: MultiComponentPhaseModel.C:393
Foam::fv::gaussConvectionScheme
Basic second-order convection using face-gradients and Gauss' theorem.
Definition: gaussConvectionScheme.H:59
composition
basicSpecieMixture & composition
Definition: createFieldRefs.H:6
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
forAllConstIter
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object.
Definition: stdFoam.H:344
fvmDiv.H
Calculate the matrix for the divergence of the given field and flux.
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
phic
surfaceScalarField phic(mixture.cAlpha() *mag(alphaPhic/mesh.magSf()))
Foam::MultiComponentPhaseModel::calculateMassFractions
void calculateMassFractions()
Transfor volume fraction into mass fractions.
Definition: MultiComponentPhaseModel.C:138
Su
zeroField Su
Definition: alphaSuSp.H:1
Foam::MultiComponentPhaseModel::MultiComponentPhaseModel
MultiComponentPhaseModel(const phaseSystem &fluid, const word &phaseName)
Definition: MultiComponentPhaseModel.C:49
phir
surfaceScalarField phir(fvc::flux(UdmModel.Udm()))
Foam::MultiComponentPhaseModel::inertIndex
label inertIndex() const
Return inert species index.
Definition: MultiComponentPhaseModel.C:409
Foam::Field< scalar >
Yt
volScalarField Yt(0.0 *Y[0])
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::MultiComponentPhaseModel::correct
virtual void correct()
Correct phase thermo.
Definition: MultiComponentPhaseModel.C:175
fvmLaplacian.H
Calculate the matrix for the laplacian of the field.
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:62
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::dimensioned< scalar >
stdFoam::end
constexpr auto end(C &c) -> decltype(c.end())
Return iterator to the end of the container c.
Definition: stdFoam.H:121
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:83
Foam::dimMass
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
Y
PtrList< volScalarField > & Y
Definition: createFieldRefs.H:7
fvmSup.H
Calculate the matrix for implicit and explicit sources.
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::GeometricField::correctBoundaryConditions
void correctBoundaryConditions()
Correct boundary field.
Definition: GeometricField.C:940
Foam::HashTable< autoPtr< phaseModel > >
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:381
Foam::GeometricField::boundaryFieldRef
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
Definition: GeometricField.C:783
fvcDDt.H
Calculate the substantive (total) derivative.
MultiComponentPhaseModel.H
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
fvcFlux.H
Calculate the face-flux of the given field.
Foam::oneField
A class representing the concept of a field of 1 used to avoid unnecessary manipulations for objects ...
Definition: oneField.H:53
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:76
Foam::fvMatrix::solve
SolverPerformance< Type > solve(const dictionary &)
Solve returning the solution statistics.
Definition: fvMatrixSolve.C:301
fvcDdt.H
Calculate the first temporal derivative.
Foam::phaseSystem
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:63
Foam::MultiComponentPhaseModel::calculateVolumeFractions
void calculateVolumeFractions()
Transfor mass fraction into volume fractions.
Definition: MultiComponentPhaseModel.C:110
Foam::dictionary::getOrDefault
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:122
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::fvsPatchField::coupled
virtual bool coupled() const
Return true if this patch field is coupled.
Definition: fvsPatchField.H:310
Foam::fv::EulerDdtScheme
Basic first-order Euler implicit/explicit ddt using only the current and previous time-step values.
Definition: EulerDdtScheme.H:60
Foam::fvMatrix::flux
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > flux() const
Return the face-flux field from the matrix.
Definition: fvMatrix.C:909
Foam::GeometricField::boundaryField
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Definition: GeometricFieldI.H:62
fvmDdt.H
Calculate the matrix for the first temporal derivative.
phases
multiphaseSystem::phaseModelList & phases
Definition: createFields.H:12
Foam::MultiComponentPhaseModel::thermo
virtual const phaseThermo & thermo() const
Access to thermo.
Definition: MultiComponentPhaseModel.C:160
Foam::subCycle
Perform a subCycleTime on a field.
Definition: subCycle.H:123