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 | Copyright (C) 2015 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  Sc_
50  (
51  "Sc",
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 
74 
75 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
76 
77 template<class BasePhaseModel>
79 {}
80 
81 
82 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
83 
84 template<class BasePhaseModel>
86 {
88  (
89  IOobject
90  (
91  IOobject::groupName("Yt", this->name()),
92  this->fluid().mesh().time().timeName(),
93  this->fluid().mesh()
94  ),
95  this->fluid().mesh(),
96  dimensionedScalar("zero", dimless, 0)
97  );
98 
99  PtrList<volScalarField>& Yi = Y();
100 
101  forAll(Yi, i)
102  {
103  if (i != inertIndex_)
104  {
105  Yt += Yi[i];
106  }
107  }
108 
109  if (inertIndex_ != -1)
110  {
111  Yi[inertIndex_] = scalar(1) - Yt;
112  Yi[inertIndex_].max(0);
113  }
114  else
115  {
116  forAll(Yi, i)
117  {
118  Yi[i] /= Yt;
119  Yi[i].max(0);
120  }
121  }
122 
124 }
125 
126 
127 template<class BasePhaseModel>
130 (
131  volScalarField& Yi
132 )
133 {
134  if
135  (
136  (inertIndex_ != -1)
137  && (
138  Yi.name()
139  == IOobject::groupName
140  (
141  this->thermo_->composition().species()[inertIndex_],
142  this->name()
143  )
144  )
145  )
146  {
147  return tmp<fvScalarMatrix>();
148  }
149 
150  const volScalarField& alpha = *this;
151  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi();
152  const volScalarField& rho = this->rho();
153 
154  return
155  (
156  fvm::ddt(alpha, rho, Yi)
157  + fvm::div(alphaRhoPhi, Yi, "div(" + alphaRhoPhi.name() + ",Yi)")
158  - fvm::Sp(this->continuityError(), Yi)
159 
161  (
163  *fvc::interpolate(this->turbulence().nut()*rho/Sc_),
164  Yi
165  )
166  ==
167  this->R(Yi)
168 
169  + fvc::ddt(residualAlpha_*rho, Yi)
170  - fvm::ddt(residualAlpha_*rho, Yi)
171  );
172 }
173 
174 
175 template<class BasePhaseModel>
178 {
179  return this->thermo_->composition().Y();
180 }
181 
182 
183 template<class BasePhaseModel>
186 {
187  return this->thermo_->composition().Y();
188 }
189 
190 
191 // ************************************************************************* //
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
fluid
multiphaseSystem & fluid
Definition: createFields.H:10
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:216
correctThermo
fluid correctThermo()
fvcDiv.H
Calculate the divergence of the given field.
inertSpecie
const word inertSpecie(thermo.lookup("inertSpecie"))
Foam::MULES::interpolate
tmp< surfaceScalarField > interpolate(const RhoType &rho)
Definition: IMULESTemplates.C:40
phaseSystem.H
Foam::fvc::div
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
fvMatrix.H
MultiComponentPhaseModel.H
Foam::MultiComponentPhaseModel::~MultiComponentPhaseModel
virtual ~MultiComponentPhaseModel()
Destructor.
fvmDiv.H
Calculate the matrix for the divergence of the given field and flux.
Foam::fvc::laplacian
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
Foam::MultiComponentPhaseModel::YiEqn
virtual tmp< fvScalarMatrix > YiEqn(volScalarField &Yi)
Return the species fraction equation.
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:59
Yt
volScalarField Yt(0.0 *Y[0])
fvmLaplacian.H
Calculate the matrix for the laplacian of the field.
Foam::MultiComponentPhaseModel::correctThermo
virtual void correctThermo()
Correct the thermodynamics.
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:41
Foam::PtrList
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:61
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
fvmSup.H
Calculate the matrix for implicit and explicit sources.
rho
rho
Definition: pEqn.H:3
Foam::MultiComponentPhaseModel::Y
virtual const PtrList< volScalarField > & Y() const
Constant access the species mass fractions.
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:52
timeName
word timeName
Definition: getTimeIndex.H:3
Foam::word::null
static const word null
An empty word.
Definition: word.H:77
fvcDdt.H
Calculate the first temporal derivative.
Foam::fvc::ddt
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
Foam::fvc::Sp
tmp< GeometricField< Type, fvPatchField, volMesh > > Sp(const volScalarField &sp, const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcSup.C:67
Foam::MultiComponentPhaseModel::MultiComponentPhaseModel
MultiComponentPhaseModel(const phaseSystem &fluid, const word &phaseName, const label index)
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Y
PtrList< volScalarField > & Y
Definition: createFields.H:36
fvmDdt.H
Calulate the matrix for the first temporal derivative.