Loading...
Searching...
No Matches
InterfaceCompositionModel.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) 2017-2022 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
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#include "phaseModel.H"
30#include "phasePair.H"
31#include "pureMixture.H"
33#include "rhoThermo.H"
35
36using namespace Foam::multiphaseInter;
37
38// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
39
40template<class Thermo, class OtherThermo>
41template<class ThermoType>
44(
45 const word& speciesName,
46 const multiComponentMixture<ThermoType>& globalThermo
47) const
48{
49 return
50 globalThermo.getLocalThermo
51 (
52 globalThermo.species().find(speciesName)
53 );
54}
55
56
57template<class Thermo, class OtherThermo>
58template<class ThermoType>
61(
62 const word& speciesName,
63 const pureMixture<ThermoType>& globalThermo
64) const
65{
66 return globalThermo.cellMixture(0);
67}
68
69
70template<class Thermo, class OtherThermo>
71template<class ThermoType>
74(
75 const word& speciesName,
77) const
78{
79 const fvMesh& mesh = fromThermo_.p().mesh();
80
81 auto tY = volScalarField::New
82 (
83 "tY",
85 mesh,
88 );
89 auto& Ys = tY.ref();
90
91 Ys = mixture.Y(speciesName);
92
93 return tY;
94}
95
96
97template<class Thermo, class OtherThermo>
98template<class ThermoType>
101(
102 const word& speciesName,
104) const
105{
106 const fvMesh& mesh = fromThermo_.p().mesh();
107
109 (
110 "tY",
112 mesh,
113 scalar(1),
114 dimless,
116 );
117}
118
119
120template<class Thermo, class OtherThermo>
121template<class ThermoType>
124(
126) const
127{
128 const fvMesh& mesh = fromThermo_.p().mesh();
129
131 (
132 "tM",
134 mesh,
136 (
137 "Mw",
139 1e-3*mixture.cellMixture(0).W()
140 ),
142 );
143}
144
145
146template<class Thermo, class OtherThermo>
147template<class ThermoType>
150(
152) const
153{
155}
156
157
158// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
159
160template<class Thermo, class OtherThermo>
162(
163 const dictionary& dict,
164 const phasePair& pair
165)
166:
168 fromThermo_
169 (
170 pair.from().mesh().lookupObject<Thermo>
171 (
172 IOobject::groupName
173 (
175 pair.from().name()
176 )
177 )
178 ),
179 toThermo_
180 (
181 pair.to().mesh().lookupObject<OtherThermo>
182 (
183 IOobject::groupName
184 (
186 pair.to().name()
187 )
188 )
189 ),
190 Le_("Le", dimless, 1.0, dict)
192
193
194// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
195
196template<class Thermo, class OtherThermo>
199(
200 const word& speciesName
201) const
202{
203 const typename OtherThermo::thermoType& toThermoType =
204 getLocalThermo
205 (
206 speciesName,
207 toThermo_
208 );
209
210 const volScalarField& p = toThermo_.p();
211
212 const volScalarField& T = toThermo_.T();
213
214 auto tmpD = volScalarField::New
215 (
216 IOobject::groupName("D", pair_.name()),
218 p.mesh(),
220 );
221 auto& D = tmpD.ref();
222
223 forAll(p, celli)
224 {
225 D[celli] =
226 toThermoType.alphah(p[celli], T[celli])
227 /toThermoType.rho(p[celli], T[celli]);
228 }
229
230 D /= Le_;
231 D.correctBoundaryConditions();
233 return tmpD;
234}
235
236
237template<class Thermo, class OtherThermo>
240(
241 const word& speciesName
242) const
243{
244 const typename Thermo::thermoType& fromThermoType =
245 getLocalThermo
246 (
247 speciesName,
248 fromThermo_
249 );
250
251 const volScalarField& p(fromThermo_.p());
252
253 const volScalarField& T(fromThermo_.T());
254
255 auto tmpD = volScalarField::New
256 (
257 IOobject::groupName("D", pair_.name()),
259 p.mesh(),
261 );
262 auto& D = tmpD.ref();
263
264 forAll(p, celli)
265 {
266 D[celli] =
267 fromThermoType.alphah(p[celli], T[celli])
268 /fromThermoType.rho(p[celli], T[celli]);
269 }
270
271 D /= Le_;
272 D.correctBoundaryConditions();
274 return tmpD;
275}
276
277
278template<class Thermo, class OtherThermo>
281(
282 const word& speciesName,
283 const volScalarField& Tf
284) const
285{
286 const typename Thermo::thermoType& fromThermo =
287 getLocalThermo(speciesName, fromThermo_);
288 const typename OtherThermo::thermoType& toThermo =
289 getLocalThermo(speciesName, toThermo_);
290
291 const volScalarField& p(fromThermo_.p());
292
293 auto tmpL = volScalarField::New
294 (
295 IOobject::groupName("L", pair_.name()),
297 p.mesh(),
300 );
301 auto& L = tmpL.ref();
302
303 // from Thermo (from) to Thermo (to)
304 forAll(p, celli)
305 {
306 L[celli] = fromThermo.Hc() - toThermo.Hc();
307 }
308
309 L.correctBoundaryConditions();
311 return tmpL;
312}
313
314
315template<class Thermo, class OtherThermo>
318(
319 const word& speciesName,
320 const volScalarField& Tf
321) const
322{
324 return nullptr;
325}
326
327
328template<class Thermo, class OtherThermo>
331(
332 const word& speciesName,
333 const volScalarField& Tf
334) const
335{
337 return nullptr;
338}
339
340
341// ************************************************************************* //
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
@ NO_REGISTER
Do not request registration (bool: false).
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
virtual tmp< volScalarField > Dto(const word &speciesName) const
Specie mass diffusivity for specie in a multicomponent.
const Thermo & fromThermo_
Thermo (from).
tmp< volScalarField > MwMixture(const pureMixture< ThermoType > &thermo) const
Return moleculas weight of the mixture for pureMixture [Kg/mol].
virtual tmp< volScalarField > D(const word &speciesName) const
Mass diffusivity.
virtual tmp< volScalarField > dY(const word &speciesName, const volScalarField &Tf) const
Mass fraction difference between the interface and the field.
tmp< volScalarField > getSpecieMassFraction(const word &speciesName, const pureMixture< ThermoType > &thermo) const
Return mass fraction for a pureMixture equal to one.
virtual tmp< volScalarField > Dfrom(const word &speciesName) const
Specie mass diffusivity for pure mixture.
InterfaceCompositionModel(const dictionary &dict, const phasePair &pair)
Construct from components.
virtual tmp< volScalarField > L(const word &speciesName, const volScalarField &Tf) const
Latent heat (to - from)(thermo - otherThermo).
const dimensionedScalar Le_
Lewis number.
const OtherThermo & toThermo_
Other Thermo (to).
virtual tmp< volScalarField > Yf(const word &speciesName, const volScalarField &Tf) const
Reference mass fraction for species based models.
const pureMixture< ThermoType >::thermoType & getLocalThermo(const word &speciesName, const pureMixture< ThermoType > &globalThermo) const
Get a reference to the local thermo for a pure mixture.
const speciesTable & species() const
Return the table of species.
Abstract base-class for fluid and solid thermodynamic properties.
Definition basicThermo.H:62
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
static const word & zeroGradientType() noexcept
The type name for zeroGradient patch fields.
label find(const word &val) const
Find index of the value (searches the hash).
Foam::multiComponentMixture.
ThermoType thermoType
The type of thermodynamics this mixture is instantiated for.
const ThermoType & getLocalThermo(const label speciei) const
Return thermo based on index.
interfaceCompositionModel(const dictionary &dict, const phasePair &pair)
Construct from a dictionary and a phase pair.
Description for mass transfer between a pair of phases. The direction of the mass transfer is from th...
Definition phasePair.H:52
Foam::pureMixture.
Definition pureMixture.H:50
ThermoType thermoType
The type of thermodynamics this mixture is instantiated for.
Definition pureMixture.H:66
const ThermoType & cellMixture(const label) const
Definition pureMixture.H:92
A class for managing temporary objects.
Definition tmp.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
const volScalarField & T
dynamicFvMesh & mesh
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition error.H:688
const word dictName("faMeshDefinition")
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
const dimensionSet dimless
Dimensionless.
const dimensionSet dimEnergy
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const dimensionSet dimMoles(0, 0, 0, 0, 1, 0, 0)
const dimensionSet dimArea(sqr(dimLength))
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
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.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
dictionary dict
const dimensionedScalar & D
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
const vector L(dict.get< vector >("L"))