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-2020 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 bool referencePhase,
46  const label index
47 )
48 :
49  BasePhaseModel(fluid, phaseName, referencePhase, index),
50  residualAlpha_
51  (
52  "residualAlpha",
53  dimless,
54  fluid.mesh().solverDict("Yi")
55  ),
56  inertIndex_(-1)
57 {
58  const word inertSpecie
59  (
60  this->thermo_->lookupOrDefault("inertSpecie", word::null)
61  );
62 
63  if (inertSpecie != word::null)
64  {
65  inertIndex_ = this->thermo_->composition().species()[inertSpecie];
66  }
67 
68  PtrList<volScalarField>& Y = this->thermo_->composition().Y();
69 
70  forAll(Y, i)
71  {
72  if (i != inertIndex_ && this->thermo_->composition().active(i))
73  {
74  const label j = YActive_.size();
75  YActive_.resize(j + 1);
76  YActive_.set(j, &Y[i]);
77  }
78  }
79 }
80 
81 
82 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
83 
84 template<class BasePhaseModel>
86 {}
87 
88 
89 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
90 
91 template<class BasePhaseModel>
93 {
95  (
96  IOobject
97  (
98  IOobject::groupName("Yt", this->name()),
99  this->fluid().mesh().time().timeName(),
100  this->fluid().mesh()
101  ),
102  this->fluid().mesh(),
104  );
105 
106  PtrList<volScalarField>& Yi = YRef();
107 
108  forAll(Yi, i)
109  {
110  if (i != inertIndex_)
111  {
112  Yi[i].max(0);
113  Yt += Yi[i];
114  }
115  }
116 
117  if (inertIndex_ != -1)
118  {
119  Yi[inertIndex_] = scalar(1) - Yt;
120  Yi[inertIndex_].max(0);
121  }
122  else
123  {
124  forAll(Yi, i)
125  {
126  Yi[i] /= Yt;
127  }
128  }
129 
131 }
132 
133 
134 template<class BasePhaseModel>
136 {
137  return false;
138 }
139 
140 
141 template<class BasePhaseModel>
144 {
145  const volScalarField& alpha = *this;
146  const surfaceScalarField alphaRhoPhi(this->alphaRhoPhi());
147  const volScalarField& rho = this->thermo().rho();
148 
149  return
150  (
151  fvm::ddt(alpha, rho, Yi)
152  + fvm::div(alphaRhoPhi, Yi, "div(" + alphaRhoPhi.name() + ",Yi)")
153  + this->divj(Yi)
154  ==
155  alpha*this->R(Yi)
156 
157  + fvc::ddt(residualAlpha_*rho, Yi)
158  - fvm::ddt(residualAlpha_*rho, Yi)
159  );
160 }
161 
162 
163 template<class BasePhaseModel>
166 {
167  return this->thermo_->composition().Y();
168 }
169 
170 
171 template<class BasePhaseModel>
174 {
175  return this->thermo_->composition().Y(name);
176 }
177 
178 
179 template<class BasePhaseModel>
182 {
183  return this->thermo_->composition().Y();
184 }
185 
186 
187 template<class BasePhaseModel>
190 {
191  return YActive_;
192 }
193 
194 
195 template<class BasePhaseModel>
198 {
199  return YActive_;
200 }
201 
202 
203 // ************************************************************************* //
virtual PtrList< volScalarField > & YRef()
Access the species mass fractions.
virtual void correctSpecies()
Correct the species 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
virtual tmp< fvScalarMatrix > YiEqn(volScalarField &Yi)
Return the species fraction equation.
Calculate the matrix for the laplacian of the field.
static const thermoModel & thermo(const TransportType &tt)
MultiComponentPhaseModel(const phaseSystem &fluid, const word &phaseName, const bool referencePhase, const label index)
fluid correctSpecies()
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:57
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)
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
phaseSystem & fluid
Definition: createFields.H:11
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.
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:46
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
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.
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
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])