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 -------------------------------------------------------------------------------
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 "interfaceCompositionModel.H"
30 #include "massTransferModel.H"
31 
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 
61  (
63  interfaceCompositionModels_[pair]->species(),
64  memberIter
65  )
66  {
67  const word& member = *memberIter;
68 
69  const word name(IOobject::groupName(member, phase.name()));
70  const word otherName
71  (
72  IOobject::groupName(member, otherPhase.name())
73  );
74 
75  tIDmdt.ref() +=
76  iDmdtSign
77  *(
78  *(*iDmdtSu_[pair])[member]
79  + *(*iDmdtSp_[pair])[member]*phase.Y(member)
80  );
81  }
82  }
83  }
84 
85  return tIDmdt;
86 }
87 
88 
89 
90 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
91 
92 template<class BasePhaseSystem>
95 (
96  const fvMesh& mesh
97 )
98 :
99  BasePhaseSystem(mesh),
100  nInterfaceCorrectors_
101  (
102  this->template lookupOrDefault<label>("nInterfaceCorrectors", 1)
103  )
104 {
105  this->generatePairsAndSubModels
106  (
107  "interfaceComposition",
108  interfaceCompositionModels_
109  );
110 
111  this->generatePairsAndSubModels
112  (
113  "massTransfer",
114  massTransferModels_,
115  false
116  );
117 
118  // Check that models have been specified in the correct combinations
120  (
122  interfaceCompositionModels_,
123  interfaceCompositionModelIter
124  )
125  {
126  const phasePair& pair =
127  this->phasePairs_[interfaceCompositionModelIter.key()];
128  const phaseModel& phase = pair.phase1();
129  const phaseModel& otherPhase = pair.phase2();
130 
131  if (!pair.ordered())
132  {
134  << "An interfacial composition model is specified for the "
135  << "unordered " << pair << " pair. Composition models only "
136  << "apply to ordered pairs. An entry for a "
137  << phasePairKey("A", "B", true) << " pair means a model for "
138  << "the A side of the A-B interface; i.e., \"A in the presence "
139  << "of B\""
140  << exit(FatalError);
141  }
142 
143 
144  const phasePairKey key(phase.name(), otherPhase.name());
145 
146  if (!this->phasePairs_.found(key))
147  {
149  << "A mass transfer model the " << key << " pair is not "
150  << "specified. This is required by the corresponding interface "
151  << "composition model."
152  << exit(FatalError);
153  }
154 
155  const phasePair& uoPair = this->phasePairs_[key];
156 
157  if (!massTransferModels_[uoPair][uoPair.index(phase)].valid())
158  {
160  << "A mass transfer model for the " << pair.phase1().name()
161  << " side of the " << uoPair << " pair is not "
162  << "specified. This is required by the corresponding interface "
163  << "composition model."
164  << exit(FatalError);
165  }
166  }
168  (
170  massTransferModels_,
171  massTransferModelIter
172  )
173  {
174  const phasePair& pair =
175  this->phasePairs_[massTransferModelIter.key()];
176 
177  if (!this->heatTransferModels_.found(pair))
178  {
180  << "A heat transfer model for " << pair << " pair is not "
181  << "specified. This is required by the corresponding species "
182  << "transfer model"
183  << exit(FatalError);
184  }
185  }
186 
187  // Generate mass transfer fields, initially assumed to be zero
189  (
191  interfaceCompositionModels_,
192  interfaceCompositionModelIter
193  )
194  {
195  const interfaceCompositionModel& compositionModel =
196  interfaceCompositionModelIter();
197 
198  const phasePair& pair =
199  this->phasePairs_[interfaceCompositionModelIter.key()];
200 
201  iDmdtSu_.set(pair, new HashPtrTable<volScalarField>());
202  iDmdtSp_.set(pair, new HashPtrTable<volScalarField>());
203 
204  forAllConstIter(hashedWordList, compositionModel.species(), memberIter)
205  {
206  const word& member = *memberIter;
207 
208  iDmdtSu_[pair]->set
209  (
210  member,
211  new volScalarField
212  (
213  IOobject
214  (
215  IOobject::groupName("iDmdtSu", pair.name()),
216  this->mesh().time().timeName(),
217  this->mesh()
218  ),
219  this->mesh(),
221  )
222  );
223 
224  iDmdtSp_[pair]->set
225  (
226  member,
227  new volScalarField
228  (
229  IOobject
230  (
231  IOobject::groupName("iDmdtSp", pair.name()),
232  this->mesh().time().timeName(),
233  this->mesh()
234  ),
235  this->mesh(),
237  )
238  );
239  }
240  }
241 }
242 
243 
244 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
245 
246 template<class BasePhaseSystem>
249 {}
250 
251 
252 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
253 
254 template<class BasePhaseSystem>
257 (
258  const phasePairKey& key
259 ) const
260 {
261  return BasePhaseSystem::dmdt(key) + this->iDmdt(key);
262 }
263 
264 
265 template<class BasePhaseSystem>
268 {
269  PtrList<volScalarField> dmdts(BasePhaseSystem::dmdts());
270 
272  (
274  interfaceCompositionModels_,
275  interfaceCompositionModelIter
276  )
277  {
278  const interfaceCompositionModel& compositionModel =
279  interfaceCompositionModelIter();
280 
281  const phasePair& pair =
282  this->phasePairs_[interfaceCompositionModelIter.key()];
283  const phaseModel& phase = pair.phase1();
284  const phaseModel& otherPhase = pair.phase2();
285 
286  forAllConstIter(hashedWordList, compositionModel.species(), memberIter)
287  {
288  const word& member = *memberIter;
289 
290  const word name(IOobject::groupName(member, phase.name()));
291  const word otherName
292  (
293  IOobject::groupName(member, otherPhase.name())
294  );
295 
296  const volScalarField iDmdt
297  (
298  *(*iDmdtSu_[pair])[member]
299  + *(*iDmdtSp_[pair])[member]*phase.Y(member)
300  );
301 
302  this->addField(phase, "dmdt", iDmdt, dmdts);
303  this->addField(otherPhase, "dmdt", - iDmdt, dmdts);
304  }
305  }
306 
307  return dmdts;
308 }
309 
310 
311 template<class BasePhaseSystem>
315 {
318 
319  phaseSystem::massTransferTable& eqns = eqnsPtr();
320 
321  // Sum up the contribution from each interface composition model
323  (
325  interfaceCompositionModels_,
326  interfaceCompositionModelIter
327  )
328  {
329  const interfaceCompositionModel& compositionModel =
330  interfaceCompositionModelIter();
331 
332  const phasePair& pair =
333  this->phasePairs_[interfaceCompositionModelIter.key()];
334  const phaseModel& phase = pair.phase1();
335  const phaseModel& otherPhase = pair.phase2();
336  const phasePair& unorderedPair =
337  this->phasePairs_[phasePair(phase, otherPhase)];
338 
339  const volScalarField& Tf(*this->Tf_[unorderedPair]);
340 
341  const volScalarField K
342  (
343  massTransferModels_[unorderedPair][unorderedPair.index(phase)]->K()
344  );
345 
346  forAllConstIter(hashedWordList, compositionModel.species(), memberIter)
347  {
348  const word& member = *memberIter;
349 
350  const word name(IOobject::groupName(member, phase.name()));
351  const word otherName
352  (
353  IOobject::groupName(member, otherPhase.name())
354  );
355 
356  const volScalarField KD(K*compositionModel.D(member));
357 
358  const volScalarField Yf(compositionModel.Yf(member, Tf));
359 
360  *(*iDmdtSu_[pair])[member] = phase.rho()*KD*Yf;
361  *(*iDmdtSp_[pair])[member] = - phase.rho()*KD;
362 
363  const fvScalarMatrix eqn
364  (
365  *(*iDmdtSu_[pair])[member]
366  + fvm::Sp(*(*iDmdtSp_[pair])[member], phase.Y(member))
367  );
368 
369  const volScalarField iDmdt
370  (
371  *(*iDmdtSu_[pair])[member]
372  + *(*iDmdtSp_[pair])[member]*phase.Y(member)
373  );
374 
375  // Implicit transport through this phase
376  *eqns[name] += eqn;
377 
378  // Explicit transport out of the other phase
379  if (eqns.found(otherName))
380  {
381  *eqns[otherName] -= iDmdt;
382  }
383  }
384  }
385 
386  return eqnsPtr;
387 }
388 
389 
390 template<class BasePhaseSystem>
393 {
394  // This loop solves for the interface temperatures, Tf, and updates the
395  // interface composition models.
396  //
397  // The rate of heat transfer to the interface must equal the latent heat
398  // consumed at the interface, i.e.:
399  //
400  // H1*(T1 - Tf) + H2*(T2 - Tf) == mDotL
401  // == K*rho*(Yfi - Yi)*Li
402  //
403  // Yfi is likely to be a strong non-linear (typically exponential) function
404  // of Tf, so the solution for the temperature is newton-accelerated
405 
407  (
408  typename BasePhaseSystem::heatTransferModelTable,
409  this->heatTransferModels_,
410  heatTransferModelIter
411  )
412  {
413  const phasePair& pair =
414  this->phasePairs_[heatTransferModelIter.key()];
415 
416  const phasePairKey key12(pair.first(), pair.second(), true);
417  const phasePairKey key21(pair.second(), pair.first(), true);
418 
419  const volScalarField H1(heatTransferModelIter().first()->K());
420  const volScalarField H2(heatTransferModelIter().second()->K());
421  const dimensionedScalar HSmall("small", heatTransferModel::dimK, SMALL);
422 
423  volScalarField& Tf = *this->Tf_[pair];
424 
425  for (label i = 0; i < nInterfaceCorrectors_; ++ i)
426  {
427  volScalarField mDotL
428  (
429  IOobject
430  (
431  "mDotL",
432  this->mesh().time().timeName(),
433  this->mesh()
434  ),
435  this->mesh(),
437  );
438  volScalarField mDotLPrime
439  (
440  IOobject
441  (
442  "mDotLPrime",
443  this->mesh().time().timeName(),
444  this->mesh()
445  ),
446  this->mesh(),
447  dimensionedScalar(mDotL.dimensions()/dimTemperature)
448  );
449 
450  // Add latent heats from forward and backward models
451  if (this->interfaceCompositionModels_.found(key12))
452  {
453  this->interfaceCompositionModels_[key12]->addMDotL
454  (
455  massTransferModels_[pair].first()->K(),
456  Tf,
457  mDotL,
458  mDotLPrime
459  );
460  }
461  if (this->interfaceCompositionModels_.found(key21))
462  {
463  this->interfaceCompositionModels_[key21]->addMDotL
464  (
465  massTransferModels_[pair].second()->K(),
466  Tf,
467  mDotL,
468  mDotLPrime
469  );
470  }
471 
472  // Update the interface temperature by applying one step of newton's
473  // method to the interface relation
474  Tf -=
475  (
476  H1*(Tf - pair.phase1().thermo().T())
477  + H2*(Tf - pair.phase2().thermo().T())
478  + mDotL
479  )
480  /(
481  max(H1 + H2 + mDotLPrime, HSmall)
482  );
483 
485 
486  Info<< "Tf." << pair.name()
487  << ": min = " << min(Tf.primitiveField())
488  << ", mean = " << average(Tf.primitiveField())
489  << ", max = " << max(Tf.primitiveField())
490  << endl;
491 
492  // Update the interface compositions
493  if (this->interfaceCompositionModels_.found(key12))
494  {
495  this->interfaceCompositionModels_[key12]->update(Tf);
496  }
497  if (this->interfaceCompositionModels_.found(key21))
498  {
499  this->interfaceCompositionModels_[key21]->update(Tf);
500  }
501  }
502  }
503 }
504 
505 
506 template<class BasePhaseSystem>
508 {
509  if (BasePhaseSystem::read())
510  {
511  bool readOK = true;
512 
513  // Models ...
514 
515  return readOK;
516  }
517  else
518  {
519  return false;
520  }
521 }
522 
523 
524 // ************************************************************************* //
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:507
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:337
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:392
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:338
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
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:62
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:65
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:257
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:95
Foam::GeometricField::correctBoundaryConditions
void correctBoundaryConditions()
Correct boundary field.
Definition: GeometricField.C:909
Foam::HashTable< autoPtr< interfaceCompositionModel >, phasePairKey, phasePairKey::hash >
Foam::phaseModel::name
const word & name() const
Return the name of this phase.
Definition: phaseModel.C:94
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:267
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
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:248
Foam::InterfaceCompositionPhaseChangePhaseSystem::massTransfer
virtual autoPtr< phaseSystem::massTransferTable > massTransfer() const
Return the mass transfer matrices.
Definition: InterfaceCompositionPhaseChangePhaseSystem.C:314
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:38
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