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-2017 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  (
139  Yi.name()
141  (
142  this->thermo_->composition().species()[inertIndex_],
143  this->name()
144  )
145  )
146  || (
147  !this->thermo_->composition().active
148  (
149  this->thermo_->composition().species()[Yi.member()]
150  )
151  )
152  )
153  )
154  {
155  return tmp<fvScalarMatrix>();
156  }
157 
158  const volScalarField& alpha = *this;
159  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi();
160  const volScalarField& rho = this->rho();
161 
162  return
163  (
164  fvm::ddt(alpha, rho, Yi)
165  + fvm::div(alphaRhoPhi, Yi, "div(" + alphaRhoPhi.name() + ",Yi)")
166  - fvm::Sp(this->continuityError(), Yi)
167 
169  (
170  fvc::interpolate(alpha)
171  *fvc::interpolate(this->turbulence().nut()*rho/Sc_),
172  Yi
173  )
174  ==
175  alpha*this->R(Yi)
176 
177  + fvc::ddt(residualAlpha_*rho, Yi)
178  - fvm::ddt(residualAlpha_*rho, Yi)
179  );
180 }
181 
182 
183 template<class BasePhaseModel>
186 {
187  return this->thermo_->composition().Y();
188 }
189 
190 
191 template<class BasePhaseModel>
194 {
195  return this->thermo_->composition().Y();
196 }
197 
198 
199 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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.
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
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
Constant access the species mass fractions.
static word groupName(Name name, const word &group)
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
Calulate the matrix for the first temporal derivative.
static const word null
An empty word.
Definition: word.H:77
word timeName
Definition: getTimeIndex.H:3
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
Info<< "Reading field U\"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);volScalarField rho(IOobject("rho", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), thermo.rho());volVectorField rhoU(IOobject("rhoU", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *U);volScalarField rhoE(IOobject("rhoE", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *(e+0.5 *magSqr(U)));surfaceScalarField pos(IOobject("pos", runTime.timeName(), mesh), mesh, dimensionedScalar("pos", dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar("neg", dimless, -1.0));surfaceScalarField phi("phi", fvc::flux(rhoU));Info<< "Creating turbulence model\"<< endl;autoPtr< compressible::turbulenceModel > turbulence(compressible::turbulenceModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
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:63
void max(const dimensioned< Type > &)
PtrList< volScalarField > & Y
A class for managing temporary objects.
Definition: PtrList.H:53
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
const word inertSpecie(thermo.lookup("inertSpecie"))
Calculate the matrix for implicit and explicit sources.
scalar nut
mixture correctThermo()
volScalarField Yt(0.0 *Y[0])