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 -------------------------------------------------------------------------------
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 using namespace Foam::multiphaseInter;
37 
38 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
39 
40 template<class Thermo, class OtherThermo>
41 template<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 
57 template<class Thermo, class OtherThermo>
58 template<class ThermoType>
61 (
62  const word& speciesName,
63  const pureMixture<ThermoType>& globalThermo
64 ) const
65 {
66  return globalThermo.cellMixture(0);
67 }
68 
69 
70 template<class Thermo, class OtherThermo>
71 template<class ThermoType>
74 (
75  const word& speciesName,
77 ) const
78 {
79  const fvMesh& mesh = fromThermo_.p().mesh();
80 
81  auto tY = tmp<volScalarField>::New
82  (
83  IOobject
84  (
85  "tY",
86  mesh.time().timeName(),
87  mesh,
90  ),
91  mesh,
93  zeroGradientFvPatchScalarField::typeName
94  );
95 
96  auto& Ys = tY.ref();
97 
98  Ys = mixture.Y(speciesName);
99 
100  return tY;
101 }
102 
103 
104 template<class Thermo, class OtherThermo>
105 template<class ThermoType>
108 (
109  const word& speciesName,
111 ) const
112 {
113  const fvMesh& mesh = fromThermo_.p().mesh();
114 
116  (
117  IOobject
118  (
119  "tY",
120  mesh.time().timeName(),
121  mesh,
124  ),
125  mesh,
126  dimensionedScalar("one", dimless, scalar(1)),
127  zeroGradientFvPatchScalarField::typeName
128  );
129 }
130 
131 
132 template<class Thermo, class OtherThermo>
133 template<class ThermoType>
136 (
138 ) const
139 {
140  const fvMesh& mesh = fromThermo_.p().mesh();
141 
143  (
144  IOobject
145  (
146  "tM",
147  mesh.time().timeName(),
148  mesh,
151  ),
152  mesh,
154  (
155  "Mw",
157  1e-3*mixture.cellMixture(0).W()
158  ),
159  zeroGradientFvPatchScalarField::typeName
160  );
161 }
162 
163 
164 template<class Thermo, class OtherThermo>
165 template<class ThermoType>
168 (
170 ) const
171 {
172  return refCast<const basicSpecieMixture>(mixture).W();
173 }
174 
175 
176 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
177 
178 template<class Thermo, class OtherThermo>
180 (
181  const dictionary& dict,
182  const phasePair& pair
183 )
184 :
186  fromThermo_
187  (
188  pair.from().mesh().lookupObject<Thermo>
189  (
190  IOobject::groupName
191  (
193  pair.from().name()
194  )
195  )
196  ),
197  toThermo_
198  (
199  pair.to().mesh().lookupObject<OtherThermo>
200  (
201  IOobject::groupName
202  (
204  pair.to().name()
205  )
206  )
207  ),
208  Le_("Le", dimless, 1.0, dict)
209 {}
211 
212 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
213 
214 template<class Thermo, class OtherThermo>
217 (
218  const word& speciesName
219 ) const
220 {
221  const typename OtherThermo::thermoType& toThermoType =
222  getLocalThermo
223  (
224  speciesName,
225  toThermo_
226  );
227 
228  const volScalarField& p = toThermo_.p();
229 
230  const volScalarField& T = toThermo_.T();
231 
232  auto tmpD = tmp<volScalarField>::New
233  (
234  IOobject
235  (
236  IOobject::groupName("D", pair_.name()),
237  p.time().timeName(),
238  p.mesh()
239  ),
240  p.mesh(),
242  );
243 
244  auto& D = tmpD.ref();
245 
246  forAll(p, celli)
247  {
248  D[celli] =
249  toThermoType.alphah(p[celli], T[celli])
250  /toThermoType.rho(p[celli], T[celli]);
251  }
252 
253  D /= Le_;
254  D.correctBoundaryConditions();
255 
256  return tmpD;
257 }
258 
259 
260 template<class Thermo, class OtherThermo>
263 (
264  const word& speciesName
265 ) const
266 {
267  const typename Thermo::thermoType& fromThermoType =
268  getLocalThermo
269  (
270  speciesName,
271  fromThermo_
272  );
273 
274  const volScalarField& p(fromThermo_.p());
275 
276  const volScalarField& T(fromThermo_.T());
277 
278  auto tmpD = tmp<volScalarField>::New
279  (
280  IOobject
281  (
282  IOobject::groupName("D", pair_.name()),
283  p.time().timeName(),
284  p.mesh()
285  ),
286  p.mesh(),
288  );
289 
290  auto& D = tmpD.ref();
291 
292  forAll(p, celli)
293  {
294  D[celli] =
295  fromThermoType.alphah(p[celli], T[celli])
296  /fromThermoType.rho(p[celli], T[celli]);
297  }
298 
299  D /= Le_;
300  D.correctBoundaryConditions();
301 
302  return tmpD;
303 }
304 
305 
306 template<class Thermo, class OtherThermo>
309 (
310  const word& speciesName,
311  const volScalarField& Tf
312 ) const
313 {
314  const typename Thermo::thermoType& fromThermo =
315  getLocalThermo(speciesName, fromThermo_);
316  const typename OtherThermo::thermoType& toThermo =
317  getLocalThermo(speciesName, toThermo_);
318 
319  const volScalarField& p(fromThermo_.p());
320 
321  auto tmpL = tmp<volScalarField>::New
322  (
323  IOobject
324  (
325  IOobject::groupName("L", pair_.name()),
326  p.time().timeName(),
327  p.mesh()
328  ),
329  p.mesh(),
331  zeroGradientFvPatchScalarField::typeName
332  );
333 
334  auto& L = tmpL.ref();
335 
336  // from Thermo (from) to Thermo (to)
337  forAll(p, celli)
338  {
339  L[celli] = fromThermo.Hc() - toThermo.Hc();
340  }
341 
342  L.correctBoundaryConditions();
343 
344  return tmpL;
345 }
346 
347 
348 template<class Thermo, class OtherThermo>
351 (
352  const word& speciesName,
353  const volScalarField& Tf
354 ) const
355 {
357  return nullptr;
358 }
359 
360 
361 template<class Thermo, class OtherThermo>
364 (
365  const word& speciesName,
366  const volScalarField& Tf
367 ) const
368 {
370  return nullptr;
371 }
372 
373 
374 // ************************************************************************* //
dictionary dict
Abstract base-class for fluid and solid thermodynamic properties.
Definition: basicThermo.H:59
tmp< volScalarField > MwMixture(const pureMixture< ThermoType > &thermo) const
Return moleculas weight of the mixture for pureMixture [Kg/mol].
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:120
const vector L(dict.get< vector >("L"))
Info<< "Creating temperaturePhaseChangeTwoPhaseMixture\"<< endl;autoPtr< temperaturePhaseChangeTwoPhaseMixture > mixture
Definition: createFields.H:39
const word dictName("faMeshDefinition")
const speciesTable & species() const
Return the table of species.
Ignore writing from objectRegistry::writeObject()
const dimensionSet dimless
Dimensionless.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:413
const pureMixture< ThermoType >::thermoType & getLocalThermo(const word &speciesName, const pureMixture< ThermoType > &globalThermo) const
Get a reference to the local thermo for a pure mixture.
virtual tmp< volScalarField > Dfrom(const word &speciesName) const
Specie mass diffusivity for pure mixture.
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
const ThermoType & cellMixture(const label) const
Definition: pureMixture.H:92
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:52
tmp< volScalarField > getSpecieMassFraction(const word &speciesName, const pureMixture< ThermoType > &thermo) const
Return mass fraction for a pureMixture equal to one.
A class for handling words, derived from Foam::string.
Definition: word.H:63
Foam::multiComponentMixture.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:203
virtual tmp< volScalarField > L(const word &speciesName, const volScalarField &Tf) const
Latent heat (to - from)(thermo - otherThermo)
const Mesh & mesh() const noexcept
Return mesh.
static word timeName(const scalar t, const int precision=precision_)
Return time name of given scalar time formatted with the given precision.
Definition: Time.C:760
virtual tmp< volScalarField > dY(const word &speciesName, const volScalarField &Tf) const
Mass fraction difference between the interface and the field.
const dimensionSet dimMoles(0, 0, 0, 0, 1, 0, 0)
Definition: dimensionSets.H:53
virtual tmp< volScalarField > Dto(const word &speciesName) const
Specie mass diffusivity for specie in a multicomponent.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
const dimensionSet dimEnergy
ThermoType thermoType
The type of thermodynamics this mixture is instantiated for.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const Time & time() const
Return Time associated with the objectRegistry.
Definition: IOobject.C:447
ThermoType thermoType
The type of thermodynamics this mixture is instantiated for.
Definition: pureMixture.H:66
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:79
Nothing to be read.
virtual tmp< volScalarField > Yf(const word &speciesName, const volScalarField &Tf) const
Reference mass fraction for species based models.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
const dimensionedScalar & D
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:666
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:166
const ThermoType & getLocalThermo(const label speciei) const
Return thermo based on index.
Description for mass transfer between a pair of phases. The direction of the mass transfer is from th...
Definition: phasePair.H:49
InterfaceCompositionModel(const dictionary &dict, const phasePair &pair)
Construct from components.
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:157