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 
34 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
35 
36 template<class BasePhaseSystem>
39 (
40  const phasePairKey& key
41 ) const
42 {
43  tmp<volScalarField> tIDmdt = phaseSystem::dmdt(key);
44 
45  const phasePair unorderedPair
46  (
47  this->phases()[key.first()],
48  this->phases()[key.second()]
49  );
50 
51  forAllConstIter(phasePair, unorderedPair, iter)
52  {
53  const phaseModel& phase = iter();
54  const phaseModel& otherPhase = iter.otherPhase();
55  const phasePair pair(phase, otherPhase, true);
56 
57  if (interfaceCompositionModels_.found(pair))
58  {
59  const scalar iDmdtSign = Pair<word>::compare(pair, key);
60 
62  (
64  interfaceCompositionModels_[pair]->species(),
65  memberIter
66  )
67  {
68  const word& member = *memberIter;
69 
70  const word name(IOobject::groupName(member, phase.name()));
71  const word otherName
72  (
73  IOobject::groupName(member, otherPhase.name())
74  );
75 
76  tIDmdt.ref() +=
77  iDmdtSign
78  *(
79  *(*iDmdtSu_[pair])[member]
80  + *(*iDmdtSp_[pair])[member]*phase.Y(member)
81  );
82  }
83  }
84  }
85 
86  return tIDmdt;
87 }
88 
89 
90 
91 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
92 
93 template<class BasePhaseSystem>
96 (
97  const fvMesh& mesh
98 )
99 :
100  BasePhaseSystem(mesh),
101  nInterfaceCorrectors_
102  (
103  this->template getOrDefault<label>("nInterfaceCorrectors", 1)
104  )
105 {
106  this->generatePairsAndSubModels
107  (
108  "interfaceComposition",
109  interfaceCompositionModels_
110  );
111 
112  this->generatePairsAndSubModels
113  (
114  "massTransfer",
115  massTransferModels_,
116  false
117  );
118 
119  // Check that models have been specified in the correct combinations
121  (
123  interfaceCompositionModels_,
124  interfaceCompositionModelIter
125  )
126  {
127  const phasePair& pair =
128  this->phasePairs_[interfaceCompositionModelIter.key()];
129  const phaseModel& phase = pair.phase1();
130  const phaseModel& otherPhase = pair.phase2();
131 
132  if (!pair.ordered())
133  {
135  << "An interfacial composition model is specified for the "
136  << "unordered " << pair << " pair. Composition models only "
137  << "apply to ordered pairs. An entry for a "
138  << phasePairKey("A", "B", true) << " pair means a model for "
139  << "the A side of the A-B interface; i.e., \"A in the presence "
140  << "of B\""
141  << exit(FatalError);
142  }
143 
144 
145  const phasePairKey key(phase.name(), otherPhase.name());
146 
147  if (!this->phasePairs_.found(key))
148  {
150  << "A mass transfer model the " << key << " pair is not "
151  << "specified. This is required by the corresponding interface "
152  << "composition model."
153  << exit(FatalError);
154  }
155 
156  const phasePair& uoPair = this->phasePairs_[key];
157 
158  if (!massTransferModels_[uoPair][uoPair.index(phase)].valid())
159  {
161  << "A mass transfer model for the " << pair.phase1().name()
162  << " side of the " << uoPair << " pair is not "
163  << "specified. This is required by the corresponding interface "
164  << "composition model."
165  << exit(FatalError);
166  }
167  }
169  (
171  massTransferModels_,
172  massTransferModelIter
173  )
174  {
175  const phasePair& pair =
176  this->phasePairs_[massTransferModelIter.key()];
177 
178  if (!this->heatTransferModels_.found(pair))
179  {
181  << "A heat transfer model for " << pair << " pair is not "
182  << "specified. This is required by the corresponding species "
183  << "transfer model"
184  << exit(FatalError);
185  }
186  }
187 
188  // Generate mass transfer fields, initially assumed to be zero
190  (
192  interfaceCompositionModels_,
193  interfaceCompositionModelIter
194  )
195  {
196  const interfaceCompositionModel& compositionModel =
197  interfaceCompositionModelIter();
198 
199  const phasePair& pair =
200  this->phasePairs_[interfaceCompositionModelIter.key()];
201 
202  iDmdtSu_.set(pair, new HashPtrTable<volScalarField>());
203  iDmdtSp_.set(pair, new HashPtrTable<volScalarField>());
204 
205  forAllConstIter(hashedWordList, compositionModel.species(), memberIter)
206  {
207  const word& member = *memberIter;
208 
209  iDmdtSu_[pair]->set
210  (
211  member,
212  new volScalarField
213  (
214  IOobject
215  (
216  IOobject::groupName("iDmdtSu", pair.name()),
217  this->mesh().time().timeName(),
218  this->mesh()
219  ),
220  this->mesh(),
222  )
223  );
224 
225  iDmdtSp_[pair]->set
226  (
227  member,
228  new volScalarField
229  (
230  IOobject
231  (
232  IOobject::groupName("iDmdtSp", pair.name()),
233  this->mesh().time().timeName(),
234  this->mesh()
235  ),
236  this->mesh(),
238  )
239  );
240  }
241  }
242 }
243 
244 
245 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
246 
247 template<class BasePhaseSystem>
250 {}
251 
252 
253 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
254 
255 template<class BasePhaseSystem>
258 (
259  const phasePairKey& key
260 ) const
261 {
262  return BasePhaseSystem::dmdt(key) + this->iDmdt(key);
263 }
264 
265 
266 template<class BasePhaseSystem>
269 {
270  PtrList<volScalarField> dmdts(BasePhaseSystem::dmdts());
271 
273  (
275  interfaceCompositionModels_,
276  interfaceCompositionModelIter
277  )
278  {
279  const interfaceCompositionModel& compositionModel =
280  interfaceCompositionModelIter();
281 
282  const phasePair& pair =
283  this->phasePairs_[interfaceCompositionModelIter.key()];
284  const phaseModel& phase = pair.phase1();
285  const phaseModel& otherPhase = pair.phase2();
286 
287  forAllConstIter(hashedWordList, compositionModel.species(), memberIter)
288  {
289  const word& member = *memberIter;
290 
291  const word name(IOobject::groupName(member, phase.name()));
292  const word otherName
293  (
294  IOobject::groupName(member, otherPhase.name())
295  );
296 
297  const volScalarField iDmdt
298  (
299  *(*iDmdtSu_[pair])[member]
300  + *(*iDmdtSp_[pair])[member]*phase.Y(member)
301  );
302 
303  this->addField(phase, "dmdt", iDmdt, dmdts);
304  this->addField(otherPhase, "dmdt", - iDmdt, dmdts);
305  }
306  }
307 
308  return dmdts;
309 }
310 
311 
312 template<class BasePhaseSystem>
316 {
319 
320  phaseSystem::massTransferTable& eqns = eqnsPtr();
321 
322  // Sum up the contribution from each interface composition model
324  (
326  interfaceCompositionModels_,
327  interfaceCompositionModelIter
328  )
329  {
330  const interfaceCompositionModel& compositionModel =
331  interfaceCompositionModelIter();
332 
333  const phasePair& pair =
334  this->phasePairs_[interfaceCompositionModelIter.key()];
335  const phaseModel& phase = pair.phase1();
336  const phaseModel& otherPhase = pair.phase2();
337  const phasePair& unorderedPair =
338  this->phasePairs_[phasePair(phase, otherPhase)];
339 
340  const volScalarField& Tf(*this->Tf_[unorderedPair]);
341 
342  const volScalarField K
343  (
344  massTransferModels_[unorderedPair][unorderedPair.index(phase)]->K()
345  );
346 
347  forAllConstIter(hashedWordList, compositionModel.species(), memberIter)
348  {
349  const word& member = *memberIter;
350 
351  const word name(IOobject::groupName(member, phase.name()));
352  const word otherName
353  (
354  IOobject::groupName(member, otherPhase.name())
355  );
356 
357  const volScalarField KD(K*compositionModel.D(member));
358 
359  const volScalarField Yf(compositionModel.Yf(member, Tf));
360 
361  *(*iDmdtSu_[pair])[member] = phase.rho()*KD*Yf;
362  *(*iDmdtSp_[pair])[member] = - phase.rho()*KD;
363 
364  const fvScalarMatrix eqn
365  (
366  *(*iDmdtSu_[pair])[member]
367  + fvm::Sp(*(*iDmdtSp_[pair])[member], phase.Y(member))
368  );
369 
370  const volScalarField iDmdt
371  (
372  *(*iDmdtSu_[pair])[member]
373  + *(*iDmdtSp_[pair])[member]*phase.Y(member)
374  );
375 
376  // Implicit transport through this phase
377  *eqns[name] += eqn;
378 
379  // Explicit transport out of the other phase
380  if (eqns.found(otherName))
381  {
382  *eqns[otherName] -= iDmdt;
383  }
384  }
385  }
386 
387  return eqnsPtr;
388 }
389 
390 
391 template<class BasePhaseSystem>
394 {
395  // This loop solves for the interface temperatures, Tf, and updates the
396  // interface composition models.
397  //
398  // The rate of heat transfer to the interface must equal the latent heat
399  // consumed at the interface, i.e.:
400  //
401  // H1*(T1 - Tf) + H2*(T2 - Tf) == mDotL
402  // == K*rho*(Yfi - Yi)*Li
403  //
404  // Yfi is likely to be a strong non-linear (typically exponential) function
405  // of Tf, so the solution for the temperature is newton-accelerated
406 
408  (
409  typename BasePhaseSystem::heatTransferModelTable,
410  this->heatTransferModels_,
411  heatTransferModelIter
412  )
413  {
414  const phasePair& pair =
415  this->phasePairs_[heatTransferModelIter.key()];
416 
417  const phasePairKey key12(pair.first(), pair.second(), true);
418  const phasePairKey key21(pair.second(), pair.first(), true);
419 
420  const volScalarField H1(heatTransferModelIter().first()->K());
421  const volScalarField H2(heatTransferModelIter().second()->K());
422  const dimensionedScalar HSmall("small", heatTransferModel::dimK, SMALL);
423 
424  volScalarField& Tf = *this->Tf_[pair];
425 
426  for (label i = 0; i < nInterfaceCorrectors_; ++ i)
427  {
428  volScalarField mDotL
429  (
430  IOobject
431  (
432  "mDotL",
433  this->mesh().time().timeName(),
434  this->mesh()
435  ),
436  this->mesh(),
438  );
439  volScalarField mDotLPrime
440  (
441  IOobject
442  (
443  "mDotLPrime",
444  this->mesh().time().timeName(),
445  this->mesh()
446  ),
447  this->mesh(),
448  dimensionedScalar(mDotL.dimensions()/dimTemperature)
449  );
450 
451  // Add latent heats from forward and backward models
452  if (this->interfaceCompositionModels_.found(key12))
453  {
454  this->interfaceCompositionModels_[key12]->addMDotL
455  (
456  massTransferModels_[pair].first()->K(),
457  Tf,
458  mDotL,
459  mDotLPrime
460  );
461  }
462  if (this->interfaceCompositionModels_.found(key21))
463  {
464  this->interfaceCompositionModels_[key21]->addMDotL
465  (
466  massTransferModels_[pair].second()->K(),
467  Tf,
468  mDotL,
469  mDotLPrime
470  );
471  }
472 
473  // Update the interface temperature by applying one step of newton's
474  // method to the interface relation
475  Tf -=
476  (
477  H1*(Tf - pair.phase1().thermo().T())
478  + H2*(Tf - pair.phase2().thermo().T())
479  + mDotL
480  )
481  /(
482  max(H1 + H2 + mDotLPrime, HSmall)
483  );
484 
486 
487  Info<< "Tf." << pair.name()
488  << ": min = " << min(Tf.primitiveField())
489  << ", mean = " << average(Tf.primitiveField())
490  << ", max = " << max(Tf.primitiveField())
491  << endl;
492 
493  // Update the interface compositions
494  if (this->interfaceCompositionModels_.found(key12))
495  {
496  this->interfaceCompositionModels_[key12]->update(Tf);
497  }
498  if (this->interfaceCompositionModels_.found(key21))
499  {
500  this->interfaceCompositionModels_[key21]->update(Tf);
501  }
502  }
503  }
504 }
505 
506 
507 template<class BasePhaseSystem>
509 {
510  if (BasePhaseSystem::read())
511  {
512  bool readOK = true;
513 
514  // Models ...
515 
516  return readOK;
517  }
518  else
519  {
520  return false;
521  }
522 }
523 
524 
525 // ************************************************************************* //
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:57
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::phasePair
Description for mass transfer between a pair of phases. The direction of the mass transfer is from th...
Definition: phasePair.H:51
Foam::InterfaceCompositionPhaseChangePhaseSystem::read
virtual bool read()
Read base phaseProperties dictionary.
Definition: InterfaceCompositionPhaseChangePhaseSystem.C:508
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
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:59
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::phaseModel::otherPhase
const phaseModel & otherPhase() const
Return the other phase in this two-phase system.
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:350
Foam::phaseModel::thermo
virtual const rhoThermo & thermo() const =0
Return the thermophysical model.
Foam::InterfaceCompositionPhaseChangePhaseSystem::correctInterfaceThermo
virtual void correctInterfaceThermo()
Correct the interface temperatures.
Definition: InterfaceCompositionPhaseChangePhaseSystem.C:393
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 composition models. These models describe the composition in phase 1...
Definition: interfaceCompositionModel.H:58
Foam::interfaceCompositionModel::species
const hashedWordList & species() const
Return the transferring species names.
Definition: interfaceCompositionModel.C:56
K
CGAL::Exact_predicates_exact_constructions_kernel K
Definition: CGALTriangulation3DKernel.H:58
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:54
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:258
massTransferModel.H
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::hashedWordList
A wordList with hashed named lookup, which can be faster in some situations than using the normal lis...
Definition: hashedWordList.H:54
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:43
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:62
Foam::phasePairKey
Definition: phasePairKey.H:59
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
Foam::blockMeshTools::read
void read(Istream &, label &, const dictionary &)
In-place read with dictionary lookup.
Definition: blockMeshTools.C:33
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:258
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::dimensioned< scalar >
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:84
Foam::phasePairKey::ordered
bool ordered() const
Return the ordered flag.
Definition: phasePairKey.C:60
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:96
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
Return the name of this phase.
Definition: phaseModel.C:95
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:268
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:372
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:524
Foam::HashPtrTable
A HashTable of pointers to objects of type <T>.
Definition: HashPtrTable.H:54
Foam::interfaceCompositionModel::Yf
virtual tmp< volScalarField > Yf(const word &speciesName, const volScalarField &Tf) const =0
Interface mass fraction.
phases
multiphaseSystem::phaseModelList & phases
Definition: createFields.H:12
Foam::phasePair::name
virtual word name() const
Pair name.
Definition: phasePair.C:92
Foam::InterfaceCompositionPhaseChangePhaseSystem::~InterfaceCompositionPhaseChangePhaseSystem
virtual ~InterfaceCompositionPhaseChangePhaseSystem()
Destructor.
Definition: InterfaceCompositionPhaseChangePhaseSystem.C:249
Foam::InterfaceCompositionPhaseChangePhaseSystem::massTransfer
virtual autoPtr< phaseSystem::massTransferTable > massTransfer() const
Return the mass transfer matrices.
Definition: InterfaceCompositionPhaseChangePhaseSystem.C:315
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:76
Foam::phasePair::phase2
const phaseModel & phase2() const
Return phase 2.
Definition: phasePairI.H:36
Foam::dimTemperature
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:55
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:39
Foam::dimVolume
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:61
massTransfer
phaseSystem::massTransferTable & massTransfer(massTransferPtr())
Foam::phasePair::phase1
const phaseModel & phase1() const
Return phase 1.
Definition: phasePairI.H:30
Foam::GeometricField< scalar, fvPatchField, volMesh >
InterfaceCompositionPhaseChangePhaseSystem.H
Foam::phase::rho
const dimensionedScalar & rho() const
Return const-access to phase1 density.
Definition: phase.H:140