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) 2015-2018 OpenFOAM Foundation
9  Copyright (C) 2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
30 
31 #include "phaseSystem.H"
32 
33 #include "fvmDdt.H"
34 #include "fvmDiv.H"
35 #include "fvmSup.H"
36 #include "fvmLaplacian.H"
37 #include "fvcDdt.H"
38 #include "fvcDiv.H"
39 #include "fvMatrix.H"
40 
41 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
42 
43 template<class BasePhaseModel>
45 (
46  const phaseSystem& fluid,
47  const word& phaseName,
48  const label index
49 )
50 :
51  BasePhaseModel(fluid, phaseName, index),
52  Sct_
53  (
54  "Sct",
55  dimless,
56  fluid.subDict(phaseName)
57  ),
58  residualAlpha_
59  (
60  "residualAlpha",
61  dimless,
62  fluid.mesh().solverDict("Yi")
63  ),
64  inertIndex_(-1)
65 {
67  if
68  (
69  this->thermo_->readIfPresent("inertSpecie", inertSpecie)
70  && !inertSpecie.empty()
71  )
72  {
73  inertIndex_ = this->thermo_->composition().species().find(inertSpecie);
74  }
75 
76  PtrList<volScalarField>& Y = this->thermo_->composition().Y();
77 
78  forAll(Y, i)
79  {
80  if (i != inertIndex_ && this->thermo_->composition().active(i))
81  {
82  const label j = YActive_.size();
83  YActive_.resize(j + 1);
84  YActive_.set(j, &Y[i]);
85  }
86  }
87 }
88 
89 
90 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
91 
92 template<class BasePhaseModel>
94 {}
95 
96 
97 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
98 
99 template<class BasePhaseModel>
101 {
103  (
104  IOobject
105  (
106  IOobject::groupName("Yt", this->name()),
107  this->fluid().mesh().time().timeName(),
108  this->fluid().mesh()
109  ),
110  this->fluid().mesh(),
112  );
113 
114  PtrList<volScalarField>& Yi = YRef();
115 
116  forAll(Yi, i)
117  {
118  if (i != inertIndex_)
119  {
120  Yt += Yi[i];
121  }
122  }
123 
124  if (inertIndex_ != -1)
125  {
126  Yi[inertIndex_] = scalar(1) - Yt;
127  Yi[inertIndex_].clamp_min(0);
128  }
129  else
130  {
131  forAll(Yi, i)
132  {
133  Yi[i] /= Yt;
134  Yi[i].clamp_min(0);
135  }
136  }
137 
139 }
140 
141 
142 template<class BasePhaseModel>
144 {
145  return false;
146 }
147 
148 
149 template<class BasePhaseModel>
152 {
153  const volScalarField& alpha = *this;
154  const surfaceScalarField alphaRhoPhi(this->alphaRhoPhi());
155  const volScalarField& rho = this->thermo().rho();
156 
157  return
158  (
159  fvm::ddt(alpha, rho, Yi)
160  + fvm::div(alphaRhoPhi, Yi, "div(" + alphaRhoPhi.name() + ",Yi)")
161 
163  (
165  *fvc::interpolate(this->muEff()/Sct_),
166  Yi
167  )
168  ==
169  alpha*this->R(Yi)
170 
171  + fvc::ddt(residualAlpha_*rho, Yi)
172  - fvm::ddt(residualAlpha_*rho, Yi)
173  );
174 }
175 
176 
177 template<class BasePhaseModel>
180 {
181  return this->thermo_->composition().Y();
182 }
183 
184 
185 template<class BasePhaseModel>
188 {
189  return this->thermo_->composition().Y(name);
190 }
191 
192 
193 template<class BasePhaseModel>
196 {
197  return this->thermo_->composition().Y();
198 }
199 
200 
201 template<class BasePhaseModel>
204 {
205  return YActive_;
206 }
207 
208 
209 template<class BasePhaseModel>
212 {
213  return YActive_;
214 }
215 
216 
217 // ************************************************************************* //
twoPhaseSystem & fluid
virtual UPtrList< volScalarField > & YActiveRef()
Access the active species mass fractions.
label inertIndex_
Inert species index.
virtual tmp< fvScalarMatrix > YiEqn(volScalarField &Yi)
Return the species fraction equation.
const word inertSpecie(thermo.get< word >("inertSpecie"))
virtual ~MultiComponentPhaseModel()=default
Destructor.
Calculate the matrix for the laplacian of the field.
const dimensionSet dimless
Dimensionless.
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
Calculate the first temporal derivative.
psiReactionThermo & thermo
Definition: createFields.H:28
MultiComponentPhaseModel(const multiphaseInterSystem &fluid, const word &phaseName)
word timeName
Definition: getTimeIndex.H:3
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:69
A class for handling words, derived from Foam::string.
Definition: word.H:63
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:41
Calculate the matrix for the first temporal derivative.
virtual bool pure() const
Return whether the phase is pure (i.e., not multi-component)
A list of pointers to objects of type <T>, without allocation/deallocation management of the pointers...
Definition: HashTable.H:106
virtual tmp< volScalarField > rho() const
Density [kg/m^3] - uses current value of pressure.
Definition: psiThermo.C:139
Calculate the divergence of the given field.
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:41
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 correctThermo()
Correct the thermodynamics.
virtual const UPtrList< volScalarField > & YActive() const
Return the active species mass fractions.
Calculate the matrix for the divergence of the given field and flux.
virtual PtrList< volScalarField > & YRef()
Access the species mass fractions.
#define R(A, B, C, D, E, F, K, M)
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:41
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
A class for managing temporary objects.
Definition: HashPtrTable.H:50
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
UPtrList< volScalarField > YActive_
Pointer list to active species.
Calculate the finiteVolume matrix for implicit and explicit sources.
mixture correctThermo()
volScalarField Yt(0.0 *Y[0])