InterfaceCompositionPhaseChangePhaseSystem.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  Copyright (C) 2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
30 #include "interfaceCompositionModel.H"
31 #include "massTransferModel.H"
32 
33 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
34 
35 template<class BasePhaseSystem>
38 (
39  const phasePairKey& key
40 ) const
41 {
42  tmp<volScalarField> tIDmdt = phaseSystem::dmdt(key);
43 
44  const phasePair unorderedPair
45  (
46  this->phases()[key.first()],
47  this->phases()[key.second()]
48  );
49 
50  forAllConstIter(phasePair, unorderedPair, iter)
51  {
52  const phaseModel& phase = iter();
53  const phaseModel& otherPhase = iter.otherPhase();
54  const phasePair pair(phase, otherPhase, true);
55 
56  if (interfaceCompositionModels_.found(pair))
57  {
58  const scalar iDmdtSign = Pair<word>::compare(pair, key);
59 
60  for
61  (
62  const word& member
63  : interfaceCompositionModels_[pair]->species()
64  )
65  {
66  tIDmdt.ref() +=
67  iDmdtSign
68  *(
69  *(*iDmdtSu_[pair])[member]
70  + *(*iDmdtSp_[pair])[member]*phase.Y(member)
71  );
72  }
73  }
74  }
75 
76  return tIDmdt;
77 }
78 
79 
80 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
81 
82 template<class BasePhaseSystem>
85 (
86  const fvMesh& mesh
87 )
88 :
89  BasePhaseSystem(mesh),
90  nInterfaceCorrectors_
91  (
92  this->template getOrDefault<label>("nInterfaceCorrectors", 1)
93  )
94 {
95  this->generatePairsAndSubModels
96  (
97  "interfaceComposition",
98  interfaceCompositionModels_
99  );
100 
101  this->generatePairsAndSubModels
102  (
103  "massTransfer",
104  massTransferModels_,
105  false
106  );
107 
108  // Check that models have been specified in the correct combinations
110  (
112  interfaceCompositionModels_,
113  interfaceCompositionModelIter
114  )
115  {
116  const phasePair& pair =
117  this->phasePairs_[interfaceCompositionModelIter.key()];
118  const phaseModel& phase = pair.phase1();
119  const phaseModel& otherPhase = pair.phase2();
120 
121  if (!pair.ordered())
122  {
124  << "An interfacial composition model is specified for the "
125  << "unordered " << pair << " pair. Composition models only "
126  << "apply to ordered pairs. An entry for a "
127  << phasePairKey("A", "B", true) << " pair means a model for "
128  << "the A side of the A-B interface; i.e., \"A in the presence "
129  << "of B\""
130  << exit(FatalError);
131  }
132 
133 
134  const phasePairKey key(phase.name(), otherPhase.name());
135 
136  if (!this->phasePairs_.found(key))
137  {
139  << "A mass transfer model the " << key << " pair is not "
140  << "specified. This is required by the corresponding interface "
141  << "composition model."
142  << exit(FatalError);
143  }
144 
145  const phasePair& uoPair = this->phasePairs_[key];
146 
147  if (!massTransferModels_[uoPair][uoPair.index(phase)])
148  {
150  << "A mass transfer model for the " << pair.phase1().name()
151  << " side of the " << uoPair << " pair is not "
152  << "specified. This is required by the corresponding interface "
153  << "composition model."
154  << exit(FatalError);
155  }
156  }
158  (
160  massTransferModels_,
161  massTransferModelIter
162  )
163  {
164  const phasePair& pair =
165  this->phasePairs_[massTransferModelIter.key()];
166 
167  if (!this->heatTransferModels_.found(pair))
168  {
170  << "A heat transfer model for " << pair << " pair is not "
171  << "specified. This is required by the corresponding species "
172  << "transfer model"
173  << exit(FatalError);
174  }
175  }
176 
177  // Generate mass transfer fields, initially assumed to be zero
179  (
181  interfaceCompositionModels_,
182  interfaceCompositionModelIter
183  )
184  {
185  const interfaceCompositionModel& compositionModel =
186  interfaceCompositionModelIter();
187 
188  const phasePair& pair =
189  this->phasePairs_[interfaceCompositionModelIter.key()];
190 
191  iDmdtSu_.set(pair, new HashPtrTable<volScalarField>());
192  iDmdtSp_.set(pair, new HashPtrTable<volScalarField>());
193 
194  for (const word& member : compositionModel.species())
195  {
196  iDmdtSu_[pair]->set
197  (
198  member,
199  new volScalarField
200  (
201  IOobject
202  (
203  IOobject::groupName("iDmdtSu", pair.name()),
204  this->mesh().time().timeName(),
205  this->mesh()
206  ),
207  this->mesh(),
209  )
210  );
211 
212  iDmdtSp_[pair]->set
213  (
214  member,
215  new volScalarField
216  (
217  IOobject
218  (
219  IOobject::groupName("iDmdtSp", pair.name()),
220  this->mesh().time().timeName(),
221  this->mesh()
222  ),
223  this->mesh(),
225  )
226  );
227  }
228  }
229 }
230 
231 
232 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
233 
234 template<class BasePhaseSystem>
237 {}
238 
239 
240 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
241 
242 template<class BasePhaseSystem>
245 (
246  const phasePairKey& key
247 ) const
248 {
249  return BasePhaseSystem::dmdt(key) + this->iDmdt(key);
250 }
251 
252 
253 template<class BasePhaseSystem>
256 {
257  PtrList<volScalarField> dmdts(BasePhaseSystem::dmdts());
258 
260  (
262  interfaceCompositionModels_,
263  interfaceCompositionModelIter
264  )
265  {
266  const interfaceCompositionModel& compositionModel =
267  interfaceCompositionModelIter();
268 
269  const phasePair& pair =
270  this->phasePairs_[interfaceCompositionModelIter.key()];
271  const phaseModel& phase = pair.phase1();
272  const phaseModel& otherPhase = pair.phase2();
273 
274  for (const word& member : compositionModel.species())
275  {
276  const volScalarField iDmdt
277  (
278  *(*iDmdtSu_[pair])[member]
279  + *(*iDmdtSp_[pair])[member]*phase.Y(member)
280  );
281 
282  this->addField(phase, "dmdt", iDmdt, dmdts);
283  this->addField(otherPhase, "dmdt", - iDmdt, dmdts);
284  }
285  }
286 
287  return dmdts;
288 }
289 
290 
291 template<class BasePhaseSystem>
295 {
298 
299  phaseSystem::massTransferTable& eqns = eqnsPtr();
300 
301  // Sum up the contribution from each interface composition model
303  (
305  interfaceCompositionModels_,
306  interfaceCompositionModelIter
307  )
308  {
309  const interfaceCompositionModel& compositionModel =
310  interfaceCompositionModelIter();
311 
312  const phasePair& pair =
313  this->phasePairs_[interfaceCompositionModelIter.key()];
314  const phaseModel& phase = pair.phase1();
315  const phaseModel& otherPhase = pair.phase2();
316  const phasePair& unorderedPair =
317  this->phasePairs_[phasePair(phase, otherPhase)];
318 
319  const volScalarField& Tf(*this->Tf_[unorderedPair]);
320 
321  const volScalarField K
322  (
323  massTransferModels_[unorderedPair][unorderedPair.index(phase)]->K()
324  );
325 
326  for
327  (
328  const word& member
329  : compositionModel.species()
330  )
331  {
332  const word name(IOobject::groupName(member, phase.name()));
333  const word otherName
334  (
335  IOobject::groupName(member, otherPhase.name())
336  );
337 
338  const volScalarField KD(K*compositionModel.D(member));
339 
340  const volScalarField Yf(compositionModel.Yf(member, Tf));
341 
342  *(*iDmdtSu_[pair])[member] = phase.rho()*KD*Yf;
343  *(*iDmdtSp_[pair])[member] = - phase.rho()*KD;
344 
345  const fvScalarMatrix eqn
346  (
347  *(*iDmdtSu_[pair])[member]
348  + fvm::Sp(*(*iDmdtSp_[pair])[member], phase.Y(member))
349  );
350 
351  const volScalarField iDmdt
352  (
353  *(*iDmdtSu_[pair])[member]
354  + *(*iDmdtSp_[pair])[member]*phase.Y(member)
355  );
356 
357  // Implicit transport through this phase
358  *eqns[name] += eqn;
359 
360  // Explicit transport out of the other phase
361  if (eqns.found(otherName))
362  {
363  *eqns[otherName] -= iDmdt;
364  }
365  }
366  }
367 
368  return eqnsPtr;
369 }
370 
371 
372 template<class BasePhaseSystem>
375 {
376  // This loop solves for the interface temperatures, Tf, and updates the
377  // interface composition models.
378  //
379  // The rate of heat transfer to the interface must equal the latent heat
380  // consumed at the interface, i.e.:
381  //
382  // H1*(T1 - Tf) + H2*(T2 - Tf) == mDotL
383  // == K*rho*(Yfi - Yi)*Li
384  //
385  // Yfi is likely to be a strong non-linear (typically exponential) function
386  // of Tf, so the solution for the temperature is newton-accelerated
387 
389  (
390  typename BasePhaseSystem::heatTransferModelTable,
391  this->heatTransferModels_,
392  heatTransferModelIter
393  )
394  {
395  const phasePair& pair =
396  this->phasePairs_[heatTransferModelIter.key()];
397 
398  const phasePairKey key12(pair.first(), pair.second(), true);
399  const phasePairKey key21(pair.second(), pair.first(), true);
400 
401  const volScalarField H1(heatTransferModelIter().first()->K());
402  const volScalarField H2(heatTransferModelIter().second()->K());
403  const dimensionedScalar HSmall("small", heatTransferModel::dimK, SMALL);
404 
405  volScalarField& Tf = *this->Tf_[pair];
406 
407  for (label i = 0; i < nInterfaceCorrectors_; ++ i)
408  {
409  volScalarField mDotL
410  (
411  IOobject
412  (
413  "mDotL",
414  this->mesh().time().timeName(),
415  this->mesh()
416  ),
417  this->mesh(),
419  );
420  volScalarField mDotLPrime
421  (
422  IOobject
423  (
424  "mDotLPrime",
425  this->mesh().time().timeName(),
426  this->mesh()
427  ),
428  this->mesh(),
429  dimensionedScalar(mDotL.dimensions()/dimTemperature)
430  );
431 
432  // Add latent heats from forward and backward models
433  if (this->interfaceCompositionModels_.found(key12))
434  {
435  this->interfaceCompositionModels_[key12]->addMDotL
436  (
437  massTransferModels_[pair].first()->K(),
438  Tf,
439  mDotL,
440  mDotLPrime
441  );
442  }
443  if (this->interfaceCompositionModels_.found(key21))
444  {
445  this->interfaceCompositionModels_[key21]->addMDotL
446  (
447  massTransferModels_[pair].second()->K(),
448  Tf,
449  mDotL,
450  mDotLPrime
451  );
452  }
453 
454  // Update the interface temperature by applying one step of newton's
455  // method to the interface relation
456  Tf -=
457  (
458  H1*(Tf - pair.phase1().thermo().T())
459  + H2*(Tf - pair.phase2().thermo().T())
460  + mDotL
461  )
462  /(
463  max(H1 + H2 + mDotLPrime, HSmall)
464  );
465 
467 
468  Info<< "Tf." << pair.name()
469  << ": min = " << min(Tf.primitiveField())
470  << ", mean = " << average(Tf.primitiveField())
471  << ", max = " << max(Tf.primitiveField())
472  << endl;
473 
474  // Update the interface compositions
475  if (this->interfaceCompositionModels_.found(key12))
476  {
477  this->interfaceCompositionModels_[key12]->update(Tf);
478  }
479  if (this->interfaceCompositionModels_.found(key21))
480  {
481  this->interfaceCompositionModels_[key21]->update(Tf);
482  }
483  }
484  }
485 }
486 
487 
488 template<class BasePhaseSystem>
490 {
491  if (BasePhaseSystem::read())
492  {
493  bool readOK = true;
494 
495  // Models ...
496 
497  return readOK;
498  }
499  else
500  {
501  return false;
502  }
503 }
504 
505 
506 // ************************************************************************* //
Foam::Pair::second
const T & second() const noexcept
Return second element, which is also the last element.
Definition: PairI.H:122
Foam::phaseModel
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:54
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::phasePair
Description for mass transfer between a pair of phases. The direction of the mass transfer is from th...
Definition: phasePair.H:53
Foam::InterfaceCompositionPhaseChangePhaseSystem::read
virtual bool read()
Read base phaseProperties dictionary.
Definition: InterfaceCompositionPhaseChangePhaseSystem.C:489
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::phase
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phase.H:54
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::dimEnergy
const dimensionSet dimEnergy
Foam::dimDensity
const dimensionSet dimDensity
Foam::interfaceCompositionModel::D
virtual tmp< volScalarField > D(const word &speciesName) const =0
Mass diffusivity.
Foam::GeometricField::primitiveField
const Internal::FieldType & primitiveField() const
Return a const-reference to the internal field.
Definition: GeometricFieldI.H:53
Sp
zeroField Sp
Definition: alphaSuSp.H:2
Foam::glTF::key
auto key(const Type &t) -> typename std::enable_if< std::is_enum< Type >::value, typename std::underlying_type< Type >::type >::type
Definition: foamGltfBase.H:108
Foam::phasePair::index
label index(const phaseModel &phase) const
Return the index of the given phase. Generates a FatalError if.
Definition: phasePairI.H:72
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::phaseModel::thermo
virtual const rhoThermo & thermo() const =0
Access const to phase thermo.
Foam::InterfaceCompositionPhaseChangePhaseSystem::correctInterfaceThermo
virtual void correctInterfaceThermo()
Correct the interface temperatures.
Definition: InterfaceCompositionPhaseChangePhaseSystem.C:374
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
forAllConstIter
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object.
Definition: stdFoam.H:344
Foam::interfaceCompositionModel
Generic base class for interface models. Mass transfer models are interface models between two thermo...
Definition: interfaceCompositionModel.H:59
K
CGAL::Exact_predicates_exact_constructions_kernel K
Definition: CGALTriangulation3DKernel.H:58
Foam::blockMeshTools::read
void read(Istream &, label &val, const dictionary &)
In-place read with dictionary lookup.
Definition: blockMeshTools.C:57
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
Foam::phasePairKey::ordered
bool ordered() const noexcept
Return the ordered flag.
Definition: phasePairKey.H:98
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:227
massTransferModel.H
massTransfer
phaseSystem::massTransferTable & massTransfer(massTransferPtr())
Foam::fac::average
tmp< GeometricField< Type, faPatchField, areaMesh > > average(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Area-weighted average a edgeField creating a areaField.
Definition: facAverage.C:47
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:59
Foam::phasePairKey
An ordered or unorder pair of phase names. Typically specified as follows.
Definition: phasePairKey.H:65
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
timeName
word timeName
Definition: getTimeIndex.H:3
Foam::FatalError
error FatalError
Foam::InterfaceCompositionPhaseChangePhaseSystem::dmdt
virtual tmp< volScalarField > dmdt(const phasePairKey &key) const
Return the mass transfer rate for a pair.
Definition: InterfaceCompositionPhaseChangePhaseSystem.C:245
Foam::interfaceCompositionModel::species
const hashedWordList & species() const
The transferring species names.
Definition: interfaceCompositionModel.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::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::phase::name
const word & name() const
Definition: phase.H:111
Foam::InterfaceCompositionPhaseChangePhaseSystem::InterfaceCompositionPhaseChangePhaseSystem
InterfaceCompositionPhaseChangePhaseSystem(const fvMesh &)
Construct from fvMesh.
Definition: InterfaceCompositionPhaseChangePhaseSystem.C:85
Foam::GeometricField::correctBoundaryConditions
void correctBoundaryConditions()
Correct boundary field.
Definition: GeometricField.C:940
Foam::HashTable< autoPtr< interfaceCompositionModel >, phasePairKey, phasePairKey::hash >
Foam::phaseModel::name
const word & name() const
Definition: phaseModel.H:140
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
Foam::InterfaceCompositionPhaseChangePhaseSystem::dmdts
virtual PtrList< volScalarField > dmdts() const
Return the mass transfer rates for each phase.
Definition: InterfaceCompositionPhaseChangePhaseSystem.C:255
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::Pair
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: Pair.H:54
Foam::basicThermo::T
virtual const volScalarField & T() const
Temperature [K].
Definition: basicThermo.C:614
Foam::HashPtrTable
A HashTable of pointers to objects of type <T>, with deallocation management of the pointers.
Definition: HashPtrTable.H:54
Foam::interfaceCompositionModel::Yf
virtual tmp< volScalarField > Yf(const word &speciesName, const volScalarField &Tf) const =0
Interface mass fraction.
Foam::phaseModel::otherPhase
const phaseModel & otherPhase() const
Return the other phase in this two-phase system.
Definition: phaseModel.C:218
Foam::phasePair::name
virtual word name() const
Pair name.
Definition: phasePair.C:69
Foam::InterfaceCompositionPhaseChangePhaseSystem::~InterfaceCompositionPhaseChangePhaseSystem
virtual ~InterfaceCompositionPhaseChangePhaseSystem()
Destructor.
Definition: InterfaceCompositionPhaseChangePhaseSystem.C:236
Foam::InterfaceCompositionPhaseChangePhaseSystem::massTransfer
virtual autoPtr< phaseSystem::massTransferTable > massTransfer() const
Return the mass transfer matrices.
Definition: InterfaceCompositionPhaseChangePhaseSystem.C:294
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:68
Foam::phasePair::phase2
const phaseModel & phase2() const
Definition: phasePairI.H:36
Foam::dimTemperature
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:54
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::InterfaceCompositionPhaseChangePhaseSystem::iDmdt
virtual tmp< volScalarField > iDmdt(const phasePairKey &key) const
Return the interfacial mass transfer rate for a pair for a pair.
Definition: InterfaceCompositionPhaseChangePhaseSystem.C:38
Foam::dimVolume
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:60
Foam::phasePair::phase1
const phaseModel & phase1() const
Definition: phasePairI.H:30
Foam::GeometricField< scalar, fvPatchField, volMesh >
InterfaceCompositionPhaseChangePhaseSystem.H
phases
multiphaseSystem::phaseModelList & phases
Definition: createFields.H:12
Foam::phase::rho
const dimensionedScalar & rho() const
Return const-access to phase1 density.
Definition: phase.H:140