InterfaceCompositionPhaseChangePhaseSystem.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2015-2019 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
28 #include "massTransferModel.H"
29 
30 
31 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
32 
33 template<class BasePhaseSystem>
36 (
37  const phasePairKey& key
38 ) const
39 {
40  tmp<volScalarField> tIDmdt = phaseSystem::dmdt(key);
41 
42  const phasePair unorderedPair
43  (
44  this->phases()[key.first()],
45  this->phases()[key.second()]
46  );
47 
48  forAllConstIter(phasePair, unorderedPair, iter)
49  {
50  const phaseModel& phase = iter();
51  const phaseModel& otherPhase = iter.otherPhase();
52  const phasePair pair(phase, otherPhase, true);
53 
54  if (interfaceCompositionModels_.found(pair))
55  {
56  const scalar iDmdtSign = Pair<word>::compare(pair, key);
57 
59  (
60  hashedWordList,
61  interfaceCompositionModels_[pair]->species(),
62  memberIter
63  )
64  {
65  const word& member = *memberIter;
66 
67  const word name(IOobject::groupName(member, phase.name()));
68  const word otherName
69  (
70  IOobject::groupName(member, otherPhase.name())
71  );
72 
73  tIDmdt.ref() +=
74  iDmdtSign
75  *(
76  *(*iDmdtSu_[pair])[member]
77  + *(*iDmdtSp_[pair])[member]*phase.Y(member)
78  );
79  }
80  }
81  }
82 
83  return tIDmdt;
84 }
85 
86 
87 
88 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
89 
90 template<class BasePhaseSystem>
93 (
94  const fvMesh& mesh
95 )
96 :
97  BasePhaseSystem(mesh),
98  nInterfaceCorrectors_
99  (
100  this->template lookupOrDefault<label>("nInterfaceCorrectors", 1)
101  )
102 {
104  (
105  "interfaceComposition",
106  interfaceCompositionModels_
107  );
108 
110  (
111  "massTransfer",
112  massTransferModels_,
113  false
114  );
115 
116  // Check that models have been specified in the correct combinations
118  (
119  interfaceCompositionModelTable,
120  interfaceCompositionModels_,
121  interfaceCompositionModelIter
122  )
123  {
124  const phasePair& pair =
125  this->phasePairs_[interfaceCompositionModelIter.key()];
126  const phaseModel& phase = pair.phase1();
127  const phaseModel& otherPhase = pair.phase2();
128 
129  if (!pair.ordered())
130  {
132  << "An interfacial composition model is specified for the "
133  << "unordered " << pair << " pair. Composition models only "
134  << "apply to ordered pairs. An entry for a "
135  << phasePairKey("A", "B", true) << " pair means a model for "
136  << "the A side of the A-B interface; i.e., \"A in the presence "
137  << "of B\""
138  << exit(FatalError);
139  }
140 
141 
142  const phasePairKey key(phase.name(), otherPhase.name());
143 
144  if (!this->phasePairs_.found(key))
145  {
147  << "A mass transfer model the " << key << " pair is not "
148  << "specified. This is required by the corresponding interface "
149  << "composition model."
150  << exit(FatalError);
151  }
152 
153  const phasePair& uoPair = this->phasePairs_[key];
154 
155  if (!massTransferModels_[uoPair][uoPair.index(phase)].valid())
156  {
158  << "A mass transfer model for the " << pair.phase1().name()
159  << " side of the " << uoPair << " pair is not "
160  << "specified. This is required by the corresponding interface "
161  << "composition model."
162  << exit(FatalError);
163  }
164  }
166  (
167  massTransferModelTable,
168  massTransferModels_,
169  massTransferModelIter
170  )
171  {
172  const phasePair& pair =
173  this->phasePairs_[massTransferModelIter.key()];
174 
175  if (!this->heatTransferModels_.found(pair))
176  {
178  << "A heat transfer model for " << pair << " pair is not "
179  << "specified. This is required by the corresponding species "
180  << "transfer model"
181  << exit(FatalError);
182  }
183  }
184 
185  // Generate mass transfer fields, initially assumed to be zero
187  (
188  interfaceCompositionModelTable,
189  interfaceCompositionModels_,
190  interfaceCompositionModelIter
191  )
192  {
193  const interfaceCompositionModel& compositionModel =
194  interfaceCompositionModelIter();
195 
196  const phasePair& pair =
197  this->phasePairs_[interfaceCompositionModelIter.key()];
198 
199  iDmdtSu_.insert(pair, new HashPtrTable<volScalarField>());
200  iDmdtSp_.insert(pair, new HashPtrTable<volScalarField>());
201 
202  forAllConstIter(hashedWordList, compositionModel.species(), memberIter)
203  {
204  const word& member = *memberIter;
205 
206  iDmdtSu_[pair]->insert
207  (
208  member,
209  new volScalarField
210  (
211  IOobject
212  (
213  IOobject::groupName("iDmdtSu", pair.name()),
214  this->mesh().time().timeName(),
215  this->mesh()
216  ),
217  this->mesh(),
219  )
220  );
221 
222  iDmdtSp_[pair]->insert
223  (
224  member,
225  new volScalarField
226  (
227  IOobject
228  (
229  IOobject::groupName("iDmdtSp", pair.name()),
230  this->mesh().time().timeName(),
231  this->mesh()
232  ),
233  this->mesh(),
235  )
236  );
237  }
238  }
239 }
240 
241 
242 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
243 
244 template<class BasePhaseSystem>
247 {}
248 
249 
250 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
251 
252 template<class BasePhaseSystem>
255 (
256  const phasePairKey& key
257 ) const
258 {
259  return BasePhaseSystem::dmdt(key) + this->iDmdt(key);
260 }
261 
262 
263 template<class BasePhaseSystem>
266 {
267  PtrList<volScalarField> dmdts(BasePhaseSystem::dmdts());
268 
270  (
271  interfaceCompositionModelTable,
272  interfaceCompositionModels_,
273  interfaceCompositionModelIter
274  )
275  {
276  const interfaceCompositionModel& compositionModel =
277  interfaceCompositionModelIter();
278 
279  const phasePair& pair =
280  this->phasePairs_[interfaceCompositionModelIter.key()];
281  const phaseModel& phase = pair.phase1();
282  const phaseModel& otherPhase = pair.phase2();
283 
284  forAllConstIter(hashedWordList, compositionModel.species(), memberIter)
285  {
286  const word& member = *memberIter;
287 
288  const word name(IOobject::groupName(member, phase.name()));
289  const word otherName
290  (
291  IOobject::groupName(member, otherPhase.name())
292  );
293 
294  const volScalarField iDmdt
295  (
296  *(*iDmdtSu_[pair])[member]
297  + *(*iDmdtSp_[pair])[member]*phase.Y(member)
298  );
299 
300  this->addField(phase, "dmdt", iDmdt, dmdts);
301  this->addField(otherPhase, "dmdt", - iDmdt, dmdts);
302  }
303  }
304 
305  return dmdts;
306 }
307 
308 
309 template<class BasePhaseSystem>
312 massTransfer() const
313 {
314  autoPtr<phaseSystem::massTransferTable> eqnsPtr =
316 
317  phaseSystem::massTransferTable& eqns = eqnsPtr();
318 
319  // Sum up the contribution from each interface composition model
321  (
322  interfaceCompositionModelTable,
323  interfaceCompositionModels_,
324  interfaceCompositionModelIter
325  )
326  {
327  const interfaceCompositionModel& compositionModel =
328  interfaceCompositionModelIter();
329 
330  const phasePair& pair =
331  this->phasePairs_[interfaceCompositionModelIter.key()];
332  const phaseModel& phase = pair.phase1();
333  const phaseModel& otherPhase = pair.phase2();
334  const phasePair& unorderedPair =
335  this->phasePairs_[phasePair(phase, otherPhase)];
336 
337  const volScalarField& Tf(*this->Tf_[unorderedPair]);
338 
339  const volScalarField K
340  (
341  massTransferModels_[unorderedPair][unorderedPair.index(phase)]->K()
342  );
343 
344  forAllConstIter(hashedWordList, compositionModel.species(), memberIter)
345  {
346  const word& member = *memberIter;
347 
348  const word name(IOobject::groupName(member, phase.name()));
349  const word otherName
350  (
351  IOobject::groupName(member, otherPhase.name())
352  );
353 
354  const volScalarField KD(K*compositionModel.D(member));
355 
356  const volScalarField Yf(compositionModel.Yf(member, Tf));
357 
358  *(*iDmdtSu_[pair])[member] = phase.rho()*KD*Yf;
359  *(*iDmdtSp_[pair])[member] = - phase.rho()*KD;
360 
361  const fvScalarMatrix eqn
362  (
363  *(*iDmdtSu_[pair])[member]
364  + fvm::Sp(*(*iDmdtSp_[pair])[member], phase.Y(member))
365  );
366 
367  const volScalarField iDmdt
368  (
369  *(*iDmdtSu_[pair])[member]
370  + *(*iDmdtSp_[pair])[member]*phase.Y(member)
371  );
372 
373  // Implicit transport through this phase
374  *eqns[name] += eqn;
375 
376  // Explicit transport out of the other phase
377  if (eqns.found(otherName))
378  {
379  *eqns[otherName] -= iDmdt;
380  }
381  }
382  }
383 
384  return eqnsPtr;
385 }
386 
387 
388 template<class BasePhaseSystem>
391 {
392  // This loop solves for the interface temperatures, Tf, and updates the
393  // interface composition models.
394  //
395  // The rate of heat transfer to the interface must equal the latent heat
396  // consumed at the interface, i.e.:
397  //
398  // H1*(T1 - Tf) + H2*(T2 - Tf) == mDotL
399  // == K*rho*(Yfi - Yi)*Li
400  //
401  // Yfi is likely to be a strong non-linear (typically exponential) function
402  // of Tf, so the solution for the temperature is newton-accelerated
403 
405  (
406  typename BasePhaseSystem::heatTransferModelTable,
407  this->heatTransferModels_,
408  heatTransferModelIter
409  )
410  {
411  const phasePair& pair =
412  this->phasePairs_[heatTransferModelIter.key()];
413 
414  const phasePairKey key12(pair.first(), pair.second(), true);
415  const phasePairKey key21(pair.second(), pair.first(), true);
416 
417  const volScalarField H1(heatTransferModelIter().first()->K());
418  const volScalarField H2(heatTransferModelIter().second()->K());
419  const dimensionedScalar HSmall("small", heatTransferModel::dimK, small);
420 
421  volScalarField& Tf = *this->Tf_[pair];
422 
423  for (label i = 0; i < nInterfaceCorrectors_; ++ i)
424  {
425  volScalarField mDotL
426  (
427  IOobject
428  (
429  "mDotL",
430  this->mesh().time().timeName(),
431  this->mesh()
432  ),
433  this->mesh(),
435  );
436  volScalarField mDotLPrime
437  (
438  IOobject
439  (
440  "mDotLPrime",
441  this->mesh().time().timeName(),
442  this->mesh()
443  ),
444  this->mesh(),
445  dimensionedScalar(mDotL.dimensions()/dimTemperature, 0)
446  );
447 
448  // Add latent heats from forward and backward models
449  if (this->interfaceCompositionModels_.found(key12))
450  {
451  this->interfaceCompositionModels_[key12]->addMDotL
452  (
453  massTransferModels_[pair].first()->K(),
454  Tf,
455  mDotL,
456  mDotLPrime
457  );
458  }
459  if (this->interfaceCompositionModels_.found(key21))
460  {
461  this->interfaceCompositionModels_[key21]->addMDotL
462  (
463  massTransferModels_[pair].second()->K(),
464  Tf,
465  mDotL,
466  mDotLPrime
467  );
468  }
469 
470  // Update the interface temperature by applying one step of newton's
471  // method to the interface relation
472  Tf -=
473  (
474  H1*(Tf - pair.phase1().thermo().T())
475  + H2*(Tf - pair.phase2().thermo().T())
476  + mDotL
477  )
478  /(
479  max(H1 + H2 + mDotLPrime, HSmall)
480  );
481 
482  Tf.correctBoundaryConditions();
483 
484  Info<< "Tf." << pair.name()
485  << ": min = " << min(Tf.primitiveField())
486  << ", mean = " << average(Tf.primitiveField())
487  << ", max = " << max(Tf.primitiveField())
488  << endl;
489 
490  // Update the interface compositions
491  if (this->interfaceCompositionModels_.found(key12))
492  {
493  this->interfaceCompositionModels_[key12]->update(Tf);
494  }
495  if (this->interfaceCompositionModels_.found(key21))
496  {
497  this->interfaceCompositionModels_[key21]->update(Tf);
498  }
499  }
500  }
501 }
502 
503 
504 template<class BasePhaseSystem>
506 {
507  if (BasePhaseSystem::read())
508  {
509  bool readOK = true;
510 
511  // Models ...
512 
513  return readOK;
514  }
515  else
516  {
517  return false;
518  }
519 }
520 
521 
522 // ************************************************************************* //
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
virtual ~InterfaceCompositionPhaseChangePhaseSystem()
Destructor.
virtual PtrList< volScalarField > dmdts() const
Return the mass transfer rates for each phase.
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:59
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
HashPtrTable< fvScalarMatrix > massTransferTable
Definition: phaseSystem.H:79
error FatalError
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
multiphaseSystem::phaseModelList & phases
Definition: createFields.H:12
InterfaceCompositionPhaseChangePhaseSystem(const fvMesh &)
Construct from fvMesh.
virtual tmp< volScalarField > iDmdt(const phasePairKey &key) const
Return the interfacial mass transfer rate for a pair for a pair.
phasePairTable phasePairs_
Phase pairs.
Definition: phaseSystem.H:142
static int compare(const Pair< word > &a, const Pair< word > &b)
Compare Pairs.
Definition: Pair.H:142
IOobject(const word &name, const fileName &instance, const objectRegistry &registry, readOption r=NO_READ, writeOption w=NO_WRITE, bool registerObject=true)
Construct from name, instance, registry, io options.
Definition: IOobject.C:209
phaseSystem::massTransferTable & massTransfer(massTransferPtr())
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
virtual tmp< volScalarField > dmdt(const phasePairKey &key) const
Return the mass transfer rate for a pair.
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:626
void read(Istream &, label &, const dictionary &)
In-place read with dictionary lookup.
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:239
CGAL::Exact_predicates_exact_constructions_kernel K
virtual void correctInterfaceThermo()
Correct the interface temperatures.
static const dimensionSet dimK
Coefficient dimensions.
bool insert(const Key &, const T &newElmt)
Insert a new hashedEntry.
Definition: HashTableI.H:80
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:113
static word groupName(Name name, const word &group)
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:52
word timeName
Definition: getTimeIndex.H:3
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &df)
const word & name() const
Name function is needed to disambiguate those inherited.
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
virtual PtrList< volScalarField > dmdts() const
Return the mass transfer rates for each phase.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
const dimensionSet dimEnergy
const dimensionSet dimDensity
void generatePairsAndSubModels(const word &modelName, HashTable< autoPtr< modelType >, phasePairKey, phasePairKey::hash > &models)
Generate pairs and sub-model tables.
void addField(const phaseModel &phase, const word &fieldName, tmp< GeoField > field, PtrList< GeoField > &fieldList) const
Add the field to a phase-indexed list, with the given name,.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const Time & time() const
Return time.
Definition: IOobject.C:360
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:70
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
messageStream Info
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
A class for managing temporary objects.
Definition: PtrList.H:53
T * first()
Return the first entry.
Definition: UILList.H:109
virtual autoPtr< phaseSystem::massTransferTable > massTransfer() const
Return the mass transfer matrices.
const fvMesh & mesh() const
Return the mesh.
Definition: phaseSystemI.H:28
virtual bool read()
Read base phaseProperties dictionary.