MultiComponentPhaseModel.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2015-2018 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
27 
28 #include "phaseSystem.H"
29 
30 #include "fvmDdt.H"
31 #include "fvmDiv.H"
32 #include "fvmSup.H"
33 #include "fvmLaplacian.H"
34 #include "fvcDdt.H"
35 #include "fvcDiv.H"
36 #include "fvMatrix.H"
37 
38 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
39 
40 template<class BasePhaseModel>
42 (
43  const phaseSystem& fluid,
44  const word& phaseName,
45  const label index
46 )
47 :
48  BasePhaseModel(fluid, phaseName, index),
49  Sct_
50  (
51  "Sct",
52  dimless,
53  fluid.subDict(phaseName)
54  ),
55  residualAlpha_
56  (
57  "residualAlpha",
58  dimless,
59  fluid.mesh().solverDict("Yi")
60  ),
61  inertIndex_(-1)
62 {
63  const word inertSpecie
64  (
65  this->thermo_->lookupOrDefault("inertSpecie", word::null)
66  );
67 
68  if (inertSpecie != word::null)
69  {
70  inertIndex_ = this->thermo_->composition().species()[inertSpecie];
71  }
72 
73  PtrList<volScalarField>& Y = this->thermo_->composition().Y();
74 
75  forAll(Y, i)
76  {
77  if (i != inertIndex_ && this->thermo_->composition().active(i))
78  {
79  const label j = YActive_.size();
80  YActive_.resize(j + 1);
81  YActive_.set(j, &Y[i]);
82  }
83  }
84 }
85 
86 
87 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
88 
89 template<class BasePhaseModel>
91 {}
92 
93 
94 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
95 
96 template<class BasePhaseModel>
98 {
100  (
101  IOobject
102  (
103  IOobject::groupName("Yt", this->name()),
104  this->fluid().mesh().time().timeName(),
105  this->fluid().mesh()
106  ),
107  this->fluid().mesh(),
109  );
110 
111  PtrList<volScalarField>& Yi = YRef();
112 
113  forAll(Yi, i)
114  {
115  if (i != inertIndex_)
116  {
117  Yt += Yi[i];
118  }
119  }
120 
121  if (inertIndex_ != -1)
122  {
123  Yi[inertIndex_] = scalar(1) - Yt;
124  Yi[inertIndex_].max(0);
125  }
126  else
127  {
128  forAll(Yi, i)
129  {
130  Yi[i] /= Yt;
131  Yi[i].max(0);
132  }
133  }
134 
136 }
137 
138 
139 template<class BasePhaseModel>
141 {
142  return false;
143 }
144 
145 
146 template<class BasePhaseModel>
149 {
150  const volScalarField& alpha = *this;
151  const surfaceScalarField alphaRhoPhi(this->alphaRhoPhi());
152  const volScalarField& rho = this->thermo().rho();
153 
154  return
155  (
156  fvm::ddt(alpha, rho, Yi)
157  + fvm::div(alphaRhoPhi, Yi, "div(" + alphaRhoPhi.name() + ",Yi)")
158 
160  (
161  fvc::interpolate(alpha)
162  *fvc::interpolate(this->muEff()/Sct_),
163  Yi
164  )
165  ==
166  alpha*this->R(Yi)
167 
168  + fvc::ddt(residualAlpha_*rho, Yi)
169  - fvm::ddt(residualAlpha_*rho, Yi)
170  );
171 }
172 
173 
174 template<class BasePhaseModel>
177 {
178  return this->thermo_->composition().Y();
179 }
180 
181 
182 template<class BasePhaseModel>
185 {
186  return this->thermo_->composition().Y(name);
187 }
188 
189 
190 template<class BasePhaseModel>
193 {
194  return this->thermo_->composition().Y();
195 }
196 
197 
198 template<class BasePhaseModel>
201 {
202  return YActive_;
203 }
204 
205 
206 template<class BasePhaseModel>
209 {
210  return YActive_;
211 }
212 
213 
214 // ************************************************************************* //
virtual PtrList< volScalarField > & YRef()
Access the species mass fractions.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
multiphaseSystem & fluid
Definition: createFields.H:11
virtual tmp< fvScalarMatrix > YiEqn(volScalarField &Yi)
Return the species fraction equation.
Calculate the matrix for the laplacian of the field.
virtual void correctThermo()
Correct the thermodynamics.
rhoReactionThermo & thermo
Definition: createFields.H:28
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
Calculate the first temporal derivative.
dynamicFvMesh & mesh
virtual const PtrList< volScalarField > & Y() const
Return the species mass fractions.
static word groupName(Name name, const word &group)
mixture correctThermo()
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
Calculate the matrix for the first temporal derivative.
static const word null
An empty word.
Definition: word.H:77
word timeName
Definition: getTimeIndex.H:3
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: UPtrList.H:54
virtual const UPtrList< volScalarField > & YActive() const
Return the active species mass fractions.
Calculate the divergence of the given field.
MultiComponentPhaseModel(const phaseSystem &fluid, const word &phaseName, const label index)
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:46
virtual tmp< volScalarField > rho() const
Density [kg/m^3].
Definition: rhoThermo.C:158
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
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 ~MultiComponentPhaseModel()
Destructor.
Calculate the matrix for the divergence of the given field and flux.
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
#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:46
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:70
void max(const dimensioned< Type > &)
A class for managing temporary objects.
Definition: PtrList.H:53
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
virtual UPtrList< volScalarField > & YActiveRef()
Access the active species mass fractions.
virtual bool pure() const
Return whether the phase is pure (i.e., not multi-component)
const word inertSpecie(thermo.lookup("inertSpecie"))
Calculate the matrix for implicit and explicit sources.
volScalarField Yt(0.0 *Y[0])