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 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
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"
32 #include "multiComponentMixture.H"
33 #include "rhoThermo.H"
35 
36 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
37 
38 template<class Thermo, class OtherThermo>
39 template<class ThermoType>
42 (
43  const word& speciesName,
44  const multiComponentMixture<ThermoType>& globalThermo
45 ) const
46 {
47  return
48  globalThermo.getLocalThermo
49  (
50  globalThermo.species()
51  [
52  speciesName
53  ]
54  );
55 }
56 
57 
58 template<class Thermo, class OtherThermo>
59 template<class ThermoType>
62 (
63  const word& speciesName,
64  const pureMixture<ThermoType>& globalThermo
65 ) const
66 {
67  return globalThermo.cellMixture(0);
68 }
69 
70 
71 template<class Thermo, class OtherThermo>
72 template<class ThermoType>
75 (
76  const word& speciesName,
78 ) const
79 {
80  const fvMesh& mesh = fromThermo_.p().mesh();
81 
82  auto tY = tmp<volScalarField>::New
83  (
84  IOobject
85  (
86  "tY",
87  mesh.time().timeName(),
88  mesh,
89  IOobject::NO_READ,
90  IOobject::NO_WRITE
91  ),
92  mesh,
94  zeroGradientFvPatchScalarField::typeName
95  );
96 
97  auto& Ys = tY.ref();
98 
99  Ys = mixture.Y(speciesName);
100 
101  return tY;
102 }
103 
104 
105 template<class Thermo, class OtherThermo>
106 template<class ThermoType>
109 (
110  const word& speciesName,
112 ) const
113 {
114  const fvMesh& mesh = fromThermo_.p().mesh();
115 
117  (
118  IOobject
119  (
120  "tY",
121  mesh.time().timeName(),
122  mesh,
123  IOobject::NO_READ,
124  IOobject::NO_WRITE
125  ),
126  mesh,
127  dimensionedScalar("one", dimless, scalar(1)),
128  zeroGradientFvPatchScalarField::typeName
129  );
130 }
131 
132 
133 template<class Thermo, class OtherThermo>
134 template<class ThermoType>
137 (
139 ) const
140 {
141  const fvMesh& mesh = fromThermo_.p().mesh();
142 
144  (
145  IOobject
146  (
147  "tM",
148  mesh.time().timeName(),
149  mesh,
150  IOobject::NO_READ,
151  IOobject::NO_WRITE
152  ),
153  mesh,
155  (
156  "Mw",
158  1e-3*mixture.cellMixture(0).W()
159  ),
160  zeroGradientFvPatchScalarField::typeName
161  );
162 }
163 
164 
165 template<class Thermo, class OtherThermo>
166 template<class ThermoType>
169 (
171 ) const
172 {
173  return refCast<const basicSpecieMixture>(mixture).W();
174 }
175 
176 
177 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
178 
179 template<class Thermo, class OtherThermo>
181 (
182  const dictionary& dict,
183  const phasePair& pair
184 )
185 :
187  fromThermo_
188  (
189  pair.from().mesh().lookupObject<Thermo>
190  (
191  IOobject::groupName
192  (
194  pair.from().name()
195  )
196  )
197  ),
198  toThermo_
199  (
200  pair.to().mesh().lookupObject<OtherThermo>
201  (
202  IOobject::groupName
203  (
205  pair.to().name()
206  )
207  )
208  ),
209  Le_("Le", dimless, 1.0, dict)
210 {}
211 
212 
213 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
214 
215 template<class Thermo, class OtherThermo>
218 (
219  const word& speciesName
220 ) const
221 {
222  const typename Thermo::thermoType& fromThermo =
223  getLocalThermo
224  (
225  speciesName,
226  fromThermo_
227  );
228 
229  const volScalarField& p(fromThermo_.p());
230 
231  const volScalarField& T(fromThermo_.T());
232 
233  auto tmpD = tmp<volScalarField>::New
234  (
235  IOobject
236  (
237  IOobject::groupName("D", pair_.name()),
238  p.time().timeName(),
239  p.mesh()
240  ),
241  p.mesh(),
243  );
244 
245  auto& D = tmpD.ref();
246 
247  forAll(p, cellI)
248  {
249  D[cellI] =
250  fromThermo.alphah(p[cellI], T[cellI])
251  /fromThermo.rho(p[cellI], T[cellI]);
252  }
253 
254  D /= Le_;
255  D.correctBoundaryConditions();
256 
257  return tmpD;
258 }
259 
260 
261 template<class Thermo, class OtherThermo>
264 (
265  const word& speciesName,
266  const volScalarField& Tf
267 ) const
268 {
269  const typename Thermo::thermoType& fromThermo =
270  getLocalThermo(speciesName, fromThermo_);
271  const typename OtherThermo::thermoType& toThermo =
272  getLocalThermo(speciesName, toThermo_);
273 
274  const volScalarField& p(fromThermo_.p());
275 
276  auto tmpL = tmp<volScalarField>::New
277  (
278  IOobject
279  (
280  IOobject::groupName("L", pair_.name()),
281  p.time().timeName(),
282  p.mesh()
283  ),
284  p.mesh(),
286  zeroGradientFvPatchScalarField::typeName
287  );
288 
289  auto& L = tmpL.ref();
290 
291  // from Thermo (from) to Thermo (to)
292  forAll(p, cellI)
293  {
294  L[cellI] = fromThermo.Hc() - toThermo.Hc();
295  }
296 
297  L.correctBoundaryConditions();
298 
299  return tmpL;
300 }
301 
302 
303 template<class Thermo, class OtherThermo>
306 (
307  const word& speciesName,
308  const volScalarField& Tf
309 ) const
310 {
312  return nullptr;
313 }
314 
315 
316 template<class Thermo, class OtherThermo>
319 (
320  const word& speciesName,
321  const volScalarField& Tf
322 ) const
323 {
325  return nullptr;
326 }
327 
328 
329 // ************************************************************************* //
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
L
const vector L(dict.get< vector >("L"))
Foam::phasePair
Description for mass transfer between a pair of phases. The direction of the mass transfer is from th...
Definition: phasePair.H:53
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::pureMixture
Foam::pureMixture.
Definition: InterfaceCompositionModel.H:50
Foam::InterfaceCompositionModel::L
virtual tmp< volScalarField > L(const word &speciesName, const volScalarField &Tf) const
Latent heat (to - from)(thermo - otherThermo)
Definition: InterfaceCompositionModel.C:264
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::dimEnergy
const dimensionSet dimEnergy
Foam::InterfaceCompositionModel::getLocalThermo
const pureMixture< ThermoType >::thermoType & getLocalThermo(const word &speciesName, const pureMixture< ThermoType > &globalThermo) const
Get a reference to the local thermo for a pure mixture.
dictName
const word dictName("faMeshDefinition")
Foam::phasePair::to
virtual const phaseModel & to() const
To phase.
Definition: phasePair.C:59
Foam::dimMoles
const dimensionSet dimMoles(0, 0, 0, 0, 1, 0, 0)
Definition: dimensionSets.H:55
Foam::pureMixture::cellMixture
const ThermoType & cellMixture(const label) const
Definition: pureMixture.H:88
Foam::pureMixture::thermoType
ThermoType thermoType
The type of thermodynamics this mixture is instantiated for.
Definition: pureMixture.H:66
Foam::interfaceCompositionModel
Generic base class for interface models. Mass transfer models are interface models between two thermo...
Definition: interfaceCompositionModel.H:59
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::InterfaceCompositionModel::Yf
virtual tmp< volScalarField > Yf(const word &speciesName, const volScalarField &Tf) const
Reference mass fraction for species based models.
Definition: InterfaceCompositionModel.C:319
rhoThermo.H
NotImplemented
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:517
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
Foam::InterfaceCompositionModel::InterfaceCompositionModel
InterfaceCompositionModel(const dictionary &dict, const phasePair &pair)
Construct from components.
Definition: InterfaceCompositionModel.C:181
Foam::dimArea
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:59
Foam::multiComponentMixture
Foam::multiComponentMixture.
Definition: InterfaceCompositionModel.H:51
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
Foam::InterfaceCompositionModel::MwMixture
tmp< volScalarField > MwMixture(const pureMixture< ThermoType > &thermo) const
Return moleculas weight of the mixture for pureMixture [Kg/mol].
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::dimensioned< scalar >
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam::dimMass
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:51
T
const volScalarField & T
Definition: createFieldRefs.H:2
Foam::InterfaceCompositionModel::getSpecieMassFraction
tmp< volScalarField > getSpecieMassFraction(const word &speciesName, const pureMixture< ThermoType > &thermo) const
Return mass fraction for a pureMixture equal to one.
InterfaceCompositionModel.H
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
Foam::phaseModel::name
const word & name() const
Definition: phaseModel.H:140
Foam::multiComponentMixture::getLocalThermo
const ThermoType & getLocalThermo(const label speciei) const
Return thermo based on index.
Definition: multiComponentMixture.H:159
Foam::phasePair::from
virtual const phaseModel & from() const
From phase.
Definition: phasePair.C:49
multiComponentMixture.H
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
D
const dimensionedScalar & D
Definition: solveBulkSurfactant.H:4
Foam::multiComponentMixture::thermoType
ThermoType thermoType
The type of thermodynamics this mixture is instantiated for.
Definition: multiComponentMixture.H:89
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::InterfaceCompositionModel::D
virtual tmp< volScalarField > D(const word &speciesName) const
Mass diffusivity of the local thermo.
Definition: InterfaceCompositionModel.C:218
pureMixture.H
zeroGradientFvPatchFields.H
Foam::InterfaceCompositionModel::dY
virtual tmp< volScalarField > dY(const word &speciesName, const volScalarField &Tf) const
Mass fraction difference between the interface and the field.
Definition: InterfaceCompositionModel.C:306
Foam::mixture
Definition: mixture.H:54
Foam::dimless
const dimensionSet dimless
Dimensionless.
Definition: dimensionSets.C:189