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,
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  scalar(1),
127  dimless,
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,
152  ),
153  mesh,
155  (
156  "Mw",
158  1e-3*mixture.cellMixture(0).W()
159  ),
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 {}
212 
213 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
214 
215 template<class Thermo, class OtherThermo>
218 (
219  const word& speciesName
220 ) const
221 {
222  const typename OtherThermo::thermoType& toThermoType =
223  getLocalThermo
224  (
225  speciesName,
226  toThermo_
227  );
228 
229  const volScalarField& p = toThermo_.p();
230 
231  const volScalarField& T = toThermo_.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  toThermoType.alphah(p[celli], T[celli])
251  /toThermoType.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
267 {
268  const typename Thermo::thermoType& fromThermoType =
269  getLocalThermo
270  (
271  speciesName,
272  fromThermo_
273  );
274 
275  const volScalarField& p(fromThermo_.p());
276 
277  const volScalarField& T(fromThermo_.T());
278 
279  auto tmpD = tmp<volScalarField>::New
280  (
281  IOobject
282  (
283  IOobject::groupName("D", pair_.name()),
284  p.time().timeName(),
285  p.mesh()
286  ),
287  p.mesh(),
289  );
290 
291  auto& D = tmpD.ref();
292 
293  forAll(p, celli)
294  {
295  D[celli] =
296  fromThermoType.alphah(p[celli], T[celli])
297  /fromThermoType.rho(p[celli], T[celli]);
298  }
299 
300  D /= Le_;
301  D.correctBoundaryConditions();
302 
303  return tmpD;
304 }
305 
306 
307 template<class Thermo, class OtherThermo>
310 (
311  const word& speciesName,
312  const volScalarField& Tf
313 ) const
314 {
315  const typename Thermo::thermoType& fromThermo =
316  getLocalThermo(speciesName, fromThermo_);
317  const typename OtherThermo::thermoType& toThermo =
318  getLocalThermo(speciesName, toThermo_);
319 
320  const volScalarField& p(fromThermo_.p());
321 
322  auto tmpL = tmp<volScalarField>::New
323  (
324  IOobject
325  (
326  IOobject::groupName("L", pair_.name()),
327  p.time().timeName(),
328  p.mesh()
329  ),
330  p.mesh(),
333  );
334 
335  auto& L = tmpL.ref();
336 
337  // from Thermo (from) to Thermo (to)
338  forAll(p, celli)
339  {
340  L[celli] = fromThermo.Hc() - toThermo.Hc();
341  }
342 
343  L.correctBoundaryConditions();
344 
345  return tmpL;
346 }
347 
348 
349 template<class Thermo, class OtherThermo>
352 (
353  const word& speciesName,
354  const volScalarField& Tf
355 ) const
356 {
358  return nullptr;
359 }
360 
361 
362 template<class Thermo, class OtherThermo>
365 (
366  const word& speciesName,
367  const volScalarField& Tf
368 ) const
369 {
371  return nullptr;
372 }
373 
374 
375 // ************************************************************************* //
static const word & zeroGradientType() noexcept
The type name for zeroGradient patch fields.
Definition: fvPatchField.H:221
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:129
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:421
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 expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
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:206
virtual tmp< volScalarField > L(const word &speciesName, const volScalarField &Tf) const
Latent heat (to - from)(thermo - otherThermo)
const Time & time() const noexcept
Return Time associated with the objectRegistry.
Definition: IOobject.C:456
const Mesh & mesh() const noexcept
Return mesh.
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition: Time.C:714
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.
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:78
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:686
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
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:127