Loading...
Searching...
No Matches
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 Copyright (C) 2020 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
30
31#include "phaseSystem.H"
32
33#include "fvmDdt.H"
34#include "fvmDiv.H"
35#include "fvmSup.H"
36#include "fvmLaplacian.H"
37#include "fvcDdt.H"
38#include "fvcDiv.H"
39#include "fvMatrix.H"
40
41// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
42
43template<class BasePhaseModel>
45(
46 const phaseSystem& fluid,
47 const word& phaseName,
48 const label index
49)
50:
51 BasePhaseModel(fluid, phaseName, index),
52 Sct_
53 (
54 "Sct",
55 dimless,
56 fluid.subDict(phaseName)
57 ),
59 (
60 "residualAlpha",
61 dimless,
62 fluid.mesh().solverDict("Yi")
63 ),
64 inertIndex_(-1)
65{
67 if
68 (
69 this->thermo_->readIfPresent("inertSpecie", inertSpecie)
70 && !inertSpecie.empty()
71 )
72 {
73 inertIndex_ = this->thermo_->composition().species().find(inertSpecie);
74 }
75
76 PtrList<volScalarField>& Y = this->thermo_->composition().Y();
77
78 forAll(Y, i)
79 {
80 if (i != inertIndex_ && this->thermo_->composition().active(i))
81 {
82 const label j = YActive_.size();
83 YActive_.resize(j + 1);
84 YActive_.set(j, &Y[i]);
85 }
86 }
87}
88
89
90// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
91
92template<class BasePhaseModel>
94{}
95
96
97// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
98
99template<class BasePhaseModel>
101{
103 (
105 (
106 IOobject::groupName("Yt", this->name()),
107 this->fluid().mesh().time().timeName(),
108 this->fluid().mesh()
109 ),
110 this->fluid().mesh(),
112 );
113
114 PtrList<volScalarField>& Yi = YRef();
115
116 forAll(Yi, i)
117 {
118 if (i != inertIndex_)
119 {
120 Yt += Yi[i];
121 }
122 }
123
124 if (inertIndex_ != -1)
125 {
126 Yi[inertIndex_] = scalar(1) - Yt;
127 Yi[inertIndex_].clamp_min(0);
128 }
129 else
130 {
131 forAll(Yi, i)
132 {
133 Yi[i] /= Yt;
134 Yi[i].clamp_min(0);
135 }
137
138 BasePhaseModel::correctThermo();
139}
140
141
142template<class BasePhaseModel>
145 return false;
146}
147
148
149template<class BasePhaseModel>
152{
153 const volScalarField& alpha = *this;
154 const surfaceScalarField alphaRhoPhi(this->alphaRhoPhi());
155 const tmp<volScalarField> trho(this->thermo().rho());
156 const volScalarField& rho = trho();
157
158 return
159 (
160 fvm::ddt(alpha, rho, Yi)
161 + fvm::div(alphaRhoPhi, Yi, "div(" + alphaRhoPhi.name() + ",Yi)")
162
164 (
166 *fvc::interpolate(this->muEff()/Sct_),
167 Yi
168 )
169 ==
170 alpha*this->R(Yi)
171
172 + fvc::ddt(residualAlpha_*rho, Yi)
173 - fvm::ddt(residualAlpha_*rho, Yi)
174 );
175}
176
177
178template<class BasePhaseModel>
182 return this->thermo_->composition().Y();
183}
184
185
186template<class BasePhaseModel>
190 return this->thermo_->composition().Y(name);
191}
192
193
194template<class BasePhaseModel>
198 return this->thermo_->composition().Y();
199}
200
201
202template<class BasePhaseModel>
206 return YActive_;
207}
208
209
210template<class BasePhaseModel>
213{
214 return YActive_;
215}
216
217
218// ************************************************************************* //
volScalarField Yt(0.0 *Y[0])
twoPhaseSystem & fluid
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const word & name() const noexcept
Return the object name.
Definition IOobjectI.H:205
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
virtual void correctThermo()
Correct the thermodynamics.
virtual const UPtrList< volScalarField > & YActive() const
Return the active species mass fractions.
virtual PtrList< volScalarField > & YRef()
Access the species mass fractions.
virtual tmp< fvScalarMatrix > YiEqn(volScalarField &Yi)
Return the species fraction equation.
virtual bool pure() const
Return whether the phase is pure (i.e., not multi-component).
dimensionedScalar residualAlpha_
Residual phase fraction.
MultiComponentPhaseModel(const multiphaseInterSystem &fluid, const word &phaseName)
label inertIndex_
Inert species index.
virtual ~MultiComponentPhaseModel()=default
Destructor.
UPtrList< volScalarField > YActive_
Pointer list to active species.
virtual UPtrList< volScalarField > & YActiveRef()
Access the active species mass fractions.
virtual const PtrList< volScalarField > & Y() const
Constant access the species mass fractions.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
A list of pointers to objects of type <T>, without allocation/deallocation management of the pointers...
Definition UPtrList.H:101
Class to represent a system of phases and model interfacial transfers between them.
Definition phaseSystem.H:72
A class for managing temporary objects.
Definition tmp.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
dynamicFvMesh & mesh
Calculate the first temporal derivative.
Calculate the divergence of the given field.
Calculate the matrix for the first temporal derivative.
Calculate the matrix for the divergence of the given field and flux.
Calculate the matrix for the laplacian of the field.
Calculate the finiteVolume matrix for implicit and explicit sources.
auto & name
word timeName
Definition getTimeIndex.H:3
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.
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition fvcDdt.C:40
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition fvmDiv.C:41
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition fvmDdt.C:41
const dimensionSet dimless
Dimensionless.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
volScalarField & alpha
tmp< volScalarField > trho
psiReactionThermo & thermo
const word inertSpecie(thermo.get< word >("inertSpecie"))
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299