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-2021 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  IOobject::NO_READ,
95  IOobject::NO_WRITE
96  ),
97  Y()[i]
98  )
99  );
100  }
101 
102  // Init vol fractions from mass fractions
103  calculateVolumeFractions();
104 }
105 
106 
107 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
108 
109 
110 template<class BasePhaseModel, class phaseThermo>
113 {
114  volScalarField Xtotal(0.0*X_[0]);
115  const volScalarField W(thermo().W());
116 
117  forAll(X_, i)
118  {
119  const dimensionedScalar Wi
120  (
121  "Wi",
123  thermo().composition().W(i)
124  );
125 
126  if (i != inertIndex_)
127  {
128  X_[i] = W*Y()[i]/Wi;
129  Xtotal += X_[i];
131 
132  const volScalarField::Boundary& YBf = Y()[i].boundaryField();
133 
134  forAll(YBf, patchi)
135  {
136  const fvPatchScalarField& YPf = YBf[patchi];
137  if (YPf.fixesValue())
138  {
139  scalarField& xbf = X_[i].boundaryFieldRef()[patchi];
140  const scalarField& ybf = Y()[i].boundaryField()[patchi];
141  const scalarField& Wbf = W.boundaryField()[patchi];
142  forAll(xbf, facei)
143  {
144  xbf[facei] = Wbf[facei]*ybf[facei]/Wi.value();
145  }
146  }
147  }
148  }
149  }
150  X_[inertIndex_] = 1.0 - Xtotal;
151  X_[inertIndex_].correctBoundaryConditions();
152 }
153 
154 
155 template<class BasePhaseModel, class phaseThermo>
158 {
159  volScalarField W(X_[0]*thermo().composition().W(0));
160  for(label i=1; i< species_.size(); i++)
161  {
162  W += X_[i]*thermo().composition().W(i);
163  }
164 
165  forAll(Y(), i)
166  {
167  Y()[i] = X_[i]*thermo().composition().W(i)/W;
168 
169  Info<< Y()[i].name() << " mass fraction = "
170  << " Min(Y) = " << min(Y()[i]).value()
171  << " Max(Y) = " << max(Y()[i]).value()
172  << endl;
173  }
174 }
175 
176 
177 template<class BasePhaseModel, class phaseThermo>
178 const phaseThermo&
180 {
181  return thermoPtr_();
182 }
183 
184 
185 template<class BasePhaseModel, class phaseThermo>
186 phaseThermo&
188 {
189  return thermoPtr_();
190 }
191 
192 
193 template<class BasePhaseModel, class phaseThermo>
195 {
196  return thermo().correct();
197 }
198 
199 
200 template<class BasePhaseModel, class phaseThermo>
202 (
205 )
206 {
207  const volScalarField& alpha1 = *this;
208 
209  const fvMesh& mesh = alpha1.mesh();
210 
211  const dictionary& MULEScontrols = mesh.solverDict(alpha1.name());
212 
213  scalar cAlpha(MULEScontrols.get<scalar>("cYi"));
214 
215  PtrList<surfaceScalarField> phiYiCorrs(species_.size());
216  const surfaceScalarField& phi = this->fluid().phi();
217 
218  surfaceScalarField phic(mag((phi)/mesh.magSf()));
219 
220  phic = min(cAlpha*phic, max(phic));
221 
223 
225  {
226  const volScalarField& alpha2 = iter2()();
227  if (&alpha2 == &alpha1)
228  {
229  continue;
230  }
231 
232  phir += phic*this->fluid().nHatf(alpha1, alpha2);
233  }
234 
235  // Do not compress interface at non-coupled boundary faces
236  // (inlets, outlets etc.)
237  surfaceScalarField::Boundary& phirBf = phir.boundaryFieldRef();
238  forAll(phir.boundaryField(), patchi)
239  {
240  fvsPatchScalarField& phirp = phirBf[patchi];
241 
242  if (!phirp.coupled())
243  {
244  phirp == 0;
245  }
246  }
247 
248  word phirScheme("div(Yiphir," + alpha1.name() + ')');
249 
250  forAll(X_, i)
251  {
252  if (inertIndex_ != i)
253  {
254  volScalarField& Yi = X_[i];
255 
256  phiYiCorrs.set
257  (
258  i,
260  (
261  "phi" + Yi.name() + "Corr",
262  fvc::flux
263  (
264  phi,
265  Yi,
266  "div(phi," + Yi.name() + ')'
267  )
268  )
269  );
270 
271  surfaceScalarField& phiYiCorr = phiYiCorrs[i];
272 
274  (
275  phaseSystem::phaseModelTable, this->fluid().phases(), iter2
276  )
277  {
278  //const volScalarField& alpha2 = iter2()().oldTime();
279  const volScalarField& alpha2 = iter2()();
280 
281  if (&alpha2 == &alpha1)
282  {
283  continue;
284  }
285 
286  phiYiCorr += fvc::flux
287  (
288  -fvc::flux(-phir, alpha2, phirScheme),
289  Yi,
290  phirScheme
291  );
292  }
293 
294  // Ensure that the flux at inflow BCs is preserved
295  forAll(phiYiCorr.boundaryField(), patchi)
296  {
297  fvsPatchScalarField& phiYiCorrp =
298  phiYiCorr.boundaryFieldRef()[patchi];
299 
300  if (!phiYiCorrp.coupled())
301  {
302  const scalarField& phi1p = phi.boundaryField()[patchi];
303  const scalarField& Yip = Yi.boundaryField()[patchi];
304 
305  forAll(phiYiCorrp, facei)
306  {
307  if (phi1p[facei] < 0)
308  {
309  phiYiCorrp[facei] = Yip[facei]*phi1p[facei];
310  }
311  }
312  }
313  }
314 
316  (
317  1.0/mesh.time().deltaT().value(),
319  Yi,
320  phi,
321  phiYiCorr,
322  Sp[i],
323  Su[i],
324  oneField(),
325  zeroField(),
326  true
327  );
328  }
329  }
330 
331  volScalarField Yt(0.0*X_[0]);
332 
333  scalar nYiSubCycles
334  (
335  MULEScontrols.getOrDefault<scalar>("nYiSubCycles", 1)
336  );
337 
338  forAll(X_, i)
339  {
340  if (inertIndex_ != i)
341  {
342  volScalarField& Yi = X_[i];
343 
344  fvScalarMatrix YiEqn
345  (
346  fv::EulerDdtScheme<scalar>(mesh).fvmDdt(Yi)
348  (
349  mesh,
350  phi,
352  ).fvmDiv(phi, Yi)
353  ==
354  Su[i] + fvm::Sp(Sp[i], Yi)
355  );
356 
357  YiEqn.solve();
358 
359  surfaceScalarField& phiYiCorr = phiYiCorrs[i];
360 
361  // Add a bounded upwind U-mean flux
362  phiYiCorr += YiEqn.flux();
363 
364  if (nYiSubCycles > 1)
365  {
366  for
367  (
368  subCycle<volScalarField> YiSubCycle(Yi, nYiSubCycles);
369  !(++YiSubCycle).end();
370  )
371  {
373  (
375  Yi,
376  phi,
377  phiYiCorr,
378  Sp[i],
379  Su[i],
380  oneField(),
381  zeroField()
382  );
383  }
384  }
385  else
386  {
388  (
390  Yi,
391  phi,
392  phiYiCorr,
393  Sp[i],
394  Su[i],
395  oneField(),
396  zeroField()
397  );
398  }
399  Yt += Yi;
400  }
401  }
402 
403  X_[inertIndex_] = scalar(1) - Yt;
404  X_[inertIndex_].max(0.0);
405 
406  calculateMassFractions();
407 }
408 
409 
410 template<class BasePhaseModel, class phaseThermo>
413 {
414  return thermoPtr_->composition().Y();
415 }
416 
417 
418 template<class BasePhaseModel, class phaseThermo>
421 {
422  return thermoPtr_->composition().Y();
423 }
424 
425 
426 template<class BasePhaseModel, class phaseThermo>
427 Foam::label
429 {
430  return inertIndex_;
431 }
432 
433 
434 // ************************************************************************* //
Foam::fvPatchField< scalar >
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
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:65
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
dictName
const word dictName("faMeshDefinition")
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:202
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.
Foam::dimMoles
const dimensionSet dimMoles(0, 0, 0, 0, 1, 0, 0)
Definition: dimensionSets.H:55
fvMatrix.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::dimensioned::value
const Type & value() const
Return const reference to value.
Definition: dimensionedType.C:434
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:107
fluid
twoPhaseSystem & fluid
Definition: setRegionFluidFields.H:3
Foam::fvPatchField::fixesValue
virtual bool fixesValue() const
Return true if this patch field fixes a value.
Definition: fvPatchField.H:332
Foam::MultiComponentPhaseModel::Y
virtual const PtrList< volScalarField > & Y() const
Constant access the species mass fractions.
Definition: MultiComponentPhaseModel.C:412
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:157
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:428
Foam::Field< scalar >
Yt
volScalarField Yt(0.0 *Y[0])
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::MultiComponentPhaseModel::correct
virtual void correct()
Correct phase thermo.
Definition: MultiComponentPhaseModel.C:194
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:59
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:123
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:85
Foam::dimMass
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:51
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:453
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:68
Foam::fvMatrix::solve
SolverPerformance< Type > solve(const dictionary &)
Solve returning the solution statistics.
Definition: fvMatrixSolve.C:319
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:66
Foam::MultiComponentPhaseModel::calculateVolumeFractions
void calculateVolumeFractions()
Transfor mass fraction into volume fractions.
Definition: MultiComponentPhaseModel.C:112
Foam::dictionary::getOrDefault
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:148
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:1534
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:179
Foam::subCycle
Perform a subCycleTime on a field.
Definition: subCycle.H:123