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) 2015-2018 OpenFOAM Foundation
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"
34
35// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
36
37template<class Thermo, class OtherThermo>
38template<class ThermoType>
41(
42 const word& speciesName,
43 const multiComponentMixture<ThermoType>& globalThermo
44) const
45{
46 return
47 globalThermo.getLocalThermo
48 (
49 globalThermo.species().find(speciesName)
50 );
51}
52
53
54template<class Thermo, class OtherThermo>
55template<class ThermoType>
58(
59 const word& speciesName,
60 const pureMixture<ThermoType>& globalThermo
61) const
62{
63 return globalThermo.cellMixture(0);
64}
65
66
67// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
68
69template<class Thermo, class OtherThermo>
71(
72 const dictionary& dict,
73 const phasePair& pair
74)
75:
76 interfaceCompositionModel(dict, pair),
77 thermo_
78 (
79 pair.phase1().mesh().lookupObject<Thermo>
80 (
81 IOobject::groupName(basicThermo::dictName, pair.phase1().name())
82 )
83 ),
84 otherThermo_
85 (
86 pair.phase2().mesh().lookupObject<OtherThermo>
87 (
88 IOobject::groupName(basicThermo::dictName, pair.phase2().name())
89 )
90 ),
91 Le_("Le", dimless, dict)
92{}
93
94
95// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
96
97template<class Thermo, class OtherThermo>
98Foam::tmp<Foam::volScalarField>
100(
101 const word& speciesName,
102 const volScalarField& Tf
103) const
104{
105 return
106 Yf(speciesName, Tf)
107 - thermo_.composition().Y()
108 [
109 thermo_.composition().species().find(speciesName)
110 ];
111}
112
113
114template<class Thermo, class OtherThermo>
117(
118 const word& speciesName
119) const
120{
121 const typename Thermo::thermoType& localThermo =
122 getLocalThermo
123 (
124 speciesName,
125 thermo_
126 );
127
128 const volScalarField& p(thermo_.p());
129
130 const volScalarField& T(thermo_.T());
131
133 (
135 (
136 IOobject::groupName("D", pair_.name()),
137 p.mesh(),
139 )
140 );
141
142 volScalarField& D = tmpD.ref();
143
144 forAll(p, celli)
145 {
146 D[celli] =
147 localThermo.alphah(p[celli], T[celli])
148 /localThermo.rho(p[celli], T[celli]);
149 }
150
151 D /= Le_;
152
153 return tmpD;
154}
155
156
157template<class Thermo, class OtherThermo>
158Foam::tmp<Foam::volScalarField>
160(
161 const word& speciesName,
162 const volScalarField& Tf
163) const
164{
165 const typename Thermo::thermoType& localThermo =
166 getLocalThermo
167 (
168 speciesName,
169 thermo_
170 );
171 const typename OtherThermo::thermoType& otherLocalThermo =
172 getLocalThermo
173 (
174 speciesName,
175 otherThermo_
176 );
177
178 const volScalarField& p(thermo_.p());
179 const volScalarField& otherP(otherThermo_.p());
180
182 (
184 (
185 IOobject::groupName("L", pair_.name()),
186 p.mesh(),
188 )
189 );
190
191 volScalarField& L = tmpL.ref();
192
193 forAll(p, celli)
194 {
195 L[celli] =
196 localThermo.Ha(p[celli], Tf[celli])
197 - otherLocalThermo.Ha(otherP[celli], Tf[celli]);
199
200 return tmpL;
201}
202
203
204template<class Thermo, class OtherThermo>
206(
207 const volScalarField& K,
208 const volScalarField& Tf,
209 volScalarField& mDotL,
210 volScalarField& mDotLPrime
211) const
212{
213 for (const word& speciesName : this->speciesNames_)
214 {
215 volScalarField rhoKDL
216 (
217 thermo_.rhoThermo::rho()
218 *K
219 *D(speciesName)
220 *L(speciesName, Tf)
221 );
222
223 mDotL += rhoKDL*dY(speciesName, Tf);
224 mDotLPrime += rhoKDL*YfPrime(speciesName, Tf);
225 }
226}
227
228
229// ************************************************************************* //
CGAL::Exact_predicates_exact_constructions_kernel K
phaseModel & phase1
phaseModel & phase2
const Mesh & mesh() const noexcept
Return const reference to mesh.
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())
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
virtual void addMDotL(const volScalarField &K, const volScalarField &Tf, volScalarField &mDotL, volScalarField &mDotLPrime) const
Add latent heat flow rate to total.
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.
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 pureMixture< ThermoType >::thermoType & getLocalThermo(const word &speciesName, const pureMixture< ThermoType > &globalThermo) const
Get a reference to the local thermo for a pure mixture.
ThermoType thermoType
The type of thermodynamics this mixture is instantiated for.
ThermoType thermoType
The type of thermodynamics this mixture is instantiated for.
Definition pureMixture.H:66
A class for managing temporary objects.
Definition tmp.H:75
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition tmpI.H:235
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
const volScalarField & T
dynamicFvMesh & mesh
const word dictName("faMeshDefinition")
auto & name
const dimensionSet dimless
Dimensionless.
const dimensionSet dimEnergy
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const dimensionSet dimArea(sqr(dimLength))
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
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
const vector L(dict.get< vector >("L"))