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 -------------------------------------------------------------------------------
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 
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 "fvMatrix.H"
39 
40 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
41 
42 template<class BasePhaseModel>
44 (
45  const phaseSystem& fluid,
46  const word& phaseName,
47  const label index
48 )
49 :
50  BasePhaseModel(fluid, phaseName, index),
51  Sct_
52  (
53  "Sct",
54  dimless,
55  fluid.subDict(phaseName)
56  ),
57  residualAlpha_
58  (
59  "residualAlpha",
60  dimless,
61  fluid.mesh().solverDict("Yi")
62  ),
63  inertIndex_(-1)
64 {
65  const word inertSpecie
66  (
67  this->thermo_->lookupOrDefault("inertSpecie", word::null)
68  );
69 
70  if (inertSpecie != word::null)
71  {
72  inertIndex_ = this->thermo_->composition().species()[inertSpecie];
73  }
74 
75  PtrList<volScalarField>& Y = this->thermo_->composition().Y();
76 
77  forAll(Y, i)
78  {
79  if (i != inertIndex_ && this->thermo_->composition().active(i))
80  {
81  const label j = YActive_.size();
82  YActive_.resize(j + 1);
83  YActive_.set(j, &Y[i]);
84  }
85  }
86 }
87 
88 
89 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
90 
91 template<class BasePhaseModel>
93 {}
94 
95 
96 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
97 
98 template<class BasePhaseModel>
100 {
102  (
103  IOobject
104  (
105  IOobject::groupName("Yt", this->name()),
106  this->fluid().mesh().time().timeName(),
107  this->fluid().mesh()
108  ),
109  this->fluid().mesh(),
111  );
112 
113  PtrList<volScalarField>& Yi = YRef();
114 
115  forAll(Yi, i)
116  {
117  if (i != inertIndex_)
118  {
119  Yt += Yi[i];
120  }
121  }
122 
123  if (inertIndex_ != -1)
124  {
125  Yi[inertIndex_] = scalar(1) - Yt;
126  Yi[inertIndex_].max(0);
127  }
128  else
129  {
130  forAll(Yi, i)
131  {
132  Yi[i] /= Yt;
133  Yi[i].max(0);
134  }
135  }
136 
138 }
139 
140 
141 template<class BasePhaseModel>
143 {
144  return false;
145 }
146 
147 
148 template<class BasePhaseModel>
151 {
152  const volScalarField& alpha = *this;
153  const surfaceScalarField alphaRhoPhi(this->alphaRhoPhi());
154  const volScalarField& rho = this->thermo().rho();
155 
156  return
157  (
158  fvm::ddt(alpha, rho, Yi)
159  + fvm::div(alphaRhoPhi, Yi, "div(" + alphaRhoPhi.name() + ",Yi)")
160 
162  (
164  *fvc::interpolate(this->muEff()/Sct_),
165  Yi
166  )
167  ==
168  alpha*this->R(Yi)
169 
170  + fvc::ddt(residualAlpha_*rho, Yi)
171  - fvm::ddt(residualAlpha_*rho, Yi)
172  );
173 }
174 
175 
176 template<class BasePhaseModel>
179 {
180  return this->thermo_->composition().Y();
181 }
182 
183 
184 template<class BasePhaseModel>
187 {
188  return this->thermo_->composition().Y(name);
189 }
190 
191 
192 template<class BasePhaseModel>
195 {
196  return this->thermo_->composition().Y();
197 }
198 
199 
200 template<class BasePhaseModel>
203 {
204  return YActive_;
205 }
206 
207 
208 template<class BasePhaseModel>
211 {
212  return YActive_;
213 }
214 
215 
216 // ************************************************************************* //
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Dimensionless.
Definition: dimensionSets.H:50
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
thermo
psiReactionThermo & thermo
Definition: createFields.H:28
fvcDiv.H
Calculate the divergence of the given field.
fvMatrix.H
fluid
twoPhaseSystem & fluid
Definition: setRegionFluidFields.H:3
Foam::fac::div
tmp< GeometricField< Type, faPatchField, areaMesh > > div(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Definition: facDiv.C:50
rho
rho
Definition: readInitialConditions.H:96
Foam::MultiComponentPhaseModel::Y
virtual const PtrList< volScalarField > & Y() const
Return the species mass fractions.
Definition: MultiComponentPhaseModel.C:178
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:290
Foam::MultiComponentPhaseModel::YRef
virtual PtrList< volScalarField > & YRef()
Access the species mass fractions.
Definition: MultiComponentPhaseModel.C:194
R
#define R(A, B, C, D, E, F, K, M)
Foam::label
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:62
MultiComponentPhaseModel.H
Yt
volScalarField Yt(0.0 *Y[0])
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
fvmLaplacian.H
Calculate the matrix for the laplacian of the field.
Foam::MultiComponentPhaseModel::pure
virtual bool pure() const
Return whether the phase is pure (i.e., not multi-component)
Definition: MultiComponentPhaseModel.C:142
Foam::UPtrList
A list of pointers to objects of type <T>, without allocation/deallocation management of the pointers...
Definition: UPtrList.H:63
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:43
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:65
timeName
word timeName
Definition: getTimeIndex.H:3
inertSpecie
const word inertSpecie(thermo.get< word >("inertSpecie"))
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
correctThermo
mixture correctThermo()
Y
PtrList< volScalarField > & Y
Definition: createFieldRefs.H:7
Foam::MultiComponentPhaseModel::~MultiComponentPhaseModel
virtual ~MultiComponentPhaseModel()
Destructor.
Definition: MultiComponentPhaseModel.C:92
fvmSup.H
Calculate the matrix for implicit and explicit sources.
Foam::MultiComponentPhaseModel::YActive
virtual const UPtrList< volScalarField > & YActive() const
Return the active species mass fractions.
Definition: MultiComponentPhaseModel.C:202
Foam::MultiComponentPhaseModel::YiEqn
virtual tmp< fvScalarMatrix > YiEqn(volScalarField &Yi)
Return the species fraction equation.
Definition: MultiComponentPhaseModel.C:150
Foam::fac::ddt
tmp< GeometricField< Type, faPatchField, areaMesh > > ddt(const dimensioned< Type > dt, const faMesh &mesh)
Definition: facDdt.C:47
Foam::MultiComponentPhaseModel::correctThermo
virtual void correctThermo()
Correct the thermodynamics.
Definition: MultiComponentPhaseModel.C:99
Foam::MultiComponentPhaseModel::YActiveRef
virtual UPtrList< volScalarField > & YActiveRef()
Access the active species mass fractions.
Definition: MultiComponentPhaseModel.C:210
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:69
Foam::fac::laplacian
tmp< GeometricField< Type, faPatchField, areaMesh > > laplacian(const GeometricField< Type, faPatchField, areaMesh > &vf, const word &name)
Definition: facLaplacian.C:47
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::MultiComponentPhaseModel::MultiComponentPhaseModel
MultiComponentPhaseModel(const phaseSystem &fluid, const word &phaseName, const label index)
Definition: MultiComponentPhaseModel.C:44
Foam::fac::interpolate
static tmp< GeometricField< Type, faePatchField, edgeMesh > > interpolate(const GeometricField< Type, faPatchField, areaMesh > &tvf, const edgeScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
fvmDdt.H
Calulate the matrix for the first temporal derivative.