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-2018 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 
127  if (!pair.ordered())
128  {
130  << "An interfacial composition model is specified for the "
131  << "unordered " << pair << " pair. Composition models only "
132  << "apply to ordered pairs. A entry for an "
133  << phasePairKey("A", "B", true) << " pair means a model for "
134  << "the A side of the A-B interface; i.e., \"A in the presence "
135  << "of B\""
136  << exit(FatalError);
137  }
138 
139  const phasePairKey key(pair.phase1().name(), pair.phase2().name());
140 
141  if (!massTransferModels_[key][pair.index(pair.phase1())].valid())
142  {
144  << "A mass transfer model for the " << pair.phase1().name()
145  << " side of the " << key << " pair is not specified. This is "
146  << "required by the corresponding interface composition model."
147  << exit(FatalError);
148  }
149  }
151  (
152  massTransferModelTable,
153  massTransferModels_,
154  massTransferModelIter
155  )
156  {
157  const phasePair& pair =
158  this->phasePairs_[massTransferModelIter.key()];
159 
160  if (!this->heatTransferModels_.found(pair))
161  {
163  << "A heat transfer model for " << pair << " pair is not "
164  << "specified. This is required by the corresponding species "
165  << "transfer model"
166  << exit(FatalError);
167  }
168  }
169 
170  // Generate mass transfer fields, initially assumed to be zero
172  (
173  interfaceCompositionModelTable,
174  interfaceCompositionModels_,
175  interfaceCompositionModelIter
176  )
177  {
178  const interfaceCompositionModel& compositionModel =
179  interfaceCompositionModelIter();
180 
181  const phasePair& pair =
182  this->phasePairs_[interfaceCompositionModelIter.key()];
183  const phasePair unorderedPair(pair.phase1(), pair.phase2());
184 
185  iDmdtSu_.insert(pair, new HashPtrTable<volScalarField>());
186  iDmdtSp_.insert(pair, new HashPtrTable<volScalarField>());
187 
188  forAllConstIter(hashedWordList, compositionModel.species(), memberIter)
189  {
190  const word& member = *memberIter;
191 
192  iDmdtSu_[pair]->insert
193  (
194  member,
195  new volScalarField
196  (
197  IOobject
198  (
199  IOobject::groupName("iDmdtSu", pair.name()),
200  this->mesh().time().timeName(),
201  this->mesh()
202  ),
203  this->mesh(),
205  )
206  );
207 
208  iDmdtSp_[pair]->insert
209  (
210  member,
211  new volScalarField
212  (
213  IOobject
214  (
215  IOobject::groupName("iDmdtSp", pair.name()),
216  this->mesh().time().timeName(),
217  this->mesh()
218  ),
219  this->mesh(),
221  )
222  );
223  }
224  }
225 }
226 
227 
228 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
229 
230 template<class BasePhaseSystem>
233 {}
234 
235 
236 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
237 
238 template<class BasePhaseSystem>
241 (
242  const phasePairKey& key
243 ) const
244 {
245  return BasePhaseSystem::dmdt(key) + this->iDmdt(key);
246 }
247 
248 
249 template<class BasePhaseSystem>
252 {
253  PtrList<volScalarField> dmdts(BasePhaseSystem::dmdts());
254 
256  (
257  interfaceCompositionModelTable,
258  interfaceCompositionModels_,
259  interfaceCompositionModelIter
260  )
261  {
262  const interfaceCompositionModel& compositionModel =
263  interfaceCompositionModelIter();
264 
265  const phasePair& pair =
266  this->phasePairs_[interfaceCompositionModelIter.key()];
267  const phaseModel& phase = pair.phase1();
268  const phaseModel& otherPhase = pair.phase2();
269 
270  forAllConstIter(hashedWordList, compositionModel.species(), memberIter)
271  {
272  const word& member = *memberIter;
273 
274  const word name(IOobject::groupName(member, phase.name()));
275  const word otherName
276  (
277  IOobject::groupName(member, otherPhase.name())
278  );
279 
280  const volScalarField iDmdt
281  (
282  *(*iDmdtSu_[pair])[member]
283  + *(*iDmdtSp_[pair])[member]*phase.Y(member)
284  );
285 
286  this->addField(phase, "dmdt", iDmdt, dmdts);
287  this->addField(otherPhase, "dmdt", - iDmdt, dmdts);
288  }
289  }
290 
291  return dmdts.xfer();
292 }
293 
294 
295 template<class BasePhaseSystem>
298 massTransfer() const
299 {
300  autoPtr<phaseSystem::massTransferTable> eqnsPtr =
302 
303  phaseSystem::massTransferTable& eqns = eqnsPtr();
304 
305  // Sum up the contribution from each interface composition model
307  (
308  interfaceCompositionModelTable,
309  interfaceCompositionModels_,
310  interfaceCompositionModelIter
311  )
312  {
313  const interfaceCompositionModel& compositionModel =
314  interfaceCompositionModelIter();
315 
316  const phasePair& pair =
317  this->phasePairs_[interfaceCompositionModelIter.key()];
318  const phaseModel& phase = pair.phase1();
319  const phaseModel& otherPhase = pair.phase2();
320  const phasePair unorderedPair(phase, otherPhase);
321 
322  const volScalarField& Tf(*this->Tf_[unorderedPair]);
323 
324  const volScalarField K
325  (
326  massTransferModels_[unorderedPair][pair.index(phase)]->K()
327  );
328 
329  forAllConstIter(hashedWordList, compositionModel.species(), memberIter)
330  {
331  const word& member = *memberIter;
332 
333  const word name(IOobject::groupName(member, phase.name()));
334  const word otherName
335  (
336  IOobject::groupName(member, otherPhase.name())
337  );
338 
339  const volScalarField KD(K*compositionModel.D(member));
340 
341  const volScalarField Yf(compositionModel.Yf(member, Tf));
342 
343  *(*iDmdtSu_[pair])[member] = phase.rho()*KD*Yf;
344  *(*iDmdtSp_[pair])[member] = - phase.rho()*KD;
345 
346  const fvScalarMatrix eqn
347  (
348  *(*iDmdtSu_[pair])[member]
349  + fvm::Sp(*(*iDmdtSp_[pair])[member], phase.Y(member))
350  );
351 
352  const volScalarField iDmdt
353  (
354  *(*iDmdtSu_[pair])[member]
355  + *(*iDmdtSp_[pair])[member]*phase.Y(member)
356  );
357 
358  // Implicit transport through this phase
359  *eqns[name] += eqn;
360 
361  // Explicit transport out of the other phase
362  if (eqns.found(otherName))
363  {
364  *eqns[otherName] -= iDmdt;
365  }
366  }
367  }
368 
369  return eqnsPtr;
370 }
371 
372 
373 template<class BasePhaseSystem>
376 {
377  // This loop solves for the interface temperatures, Tf, and updates the
378  // interface composition models.
379  //
380  // The rate of heat transfer to the interface must equal the latent heat
381  // consumed at the interface, i.e.:
382  //
383  // H1*(T1 - Tf) + H2*(T2 - Tf) == mDotL
384  // == K*rho*(Yfi - Yi)*Li
385  //
386  // Yfi is likely to be a strong non-linear (typically exponential) function
387  // of Tf, so the solution for the temperature is newton-accelerated
388 
390  (
391  typename BasePhaseSystem::heatTransferModelTable,
392  this->heatTransferModels_,
393  heatTransferModelIter
394  )
395  {
396  const phasePair& pair =
397  this->phasePairs_[heatTransferModelIter.key()];
398 
399  const phasePairKey key12(pair.first(), pair.second(), true);
400  const phasePairKey key21(pair.second(), pair.first(), true);
401 
402  const volScalarField H1(heatTransferModelIter().first()->K());
403  const volScalarField H2(heatTransferModelIter().second()->K());
404  const dimensionedScalar HSmall("small", heatTransferModel::dimK, small);
405 
406  volScalarField& Tf = *this->Tf_[pair];
407 
408  for (label i = 0; i < nInterfaceCorrectors_; ++ i)
409  {
410  volScalarField mDotL
411  (
412  IOobject
413  (
414  "mDotL",
415  this->mesh().time().timeName(),
416  this->mesh()
417  ),
418  this->mesh(),
420  );
421  volScalarField mDotLPrime
422  (
423  IOobject
424  (
425  "mDotLPrime",
426  this->mesh().time().timeName(),
427  this->mesh()
428  ),
429  this->mesh(),
430  dimensionedScalar("zero", mDotL.dimensions()/dimTemperature, 0)
431  );
432 
433  // Add latent heats from forward and backward models
434  if (this->interfaceCompositionModels_.found(key12))
435  {
436  this->interfaceCompositionModels_[key12]->addMDotL
437  (
438  massTransferModels_[pair].first()->K(),
439  Tf,
440  mDotL,
441  mDotLPrime
442  );
443  }
444  if (this->interfaceCompositionModels_.found(key21))
445  {
446  this->interfaceCompositionModels_[key21]->addMDotL
447  (
448  massTransferModels_[pair].second()->K(),
449  Tf,
450  mDotL,
451  mDotLPrime
452  );
453  }
454 
455  // Update the interface temperature by applying one step of newton's
456  // method to the interface relation
457  Tf -=
458  (
459  H1*(Tf - pair.phase1().thermo().T())
460  + H2*(Tf - pair.phase2().thermo().T())
461  + mDotL
462  )
463  /(
464  max(H1 + H2 + mDotLPrime, HSmall)
465  );
466 
467  Tf.correctBoundaryConditions();
468 
469  Info<< "Tf." << pair.name()
470  << ": min = " << min(Tf.primitiveField())
471  << ", mean = " << average(Tf.primitiveField())
472  << ", max = " << max(Tf.primitiveField())
473  << endl;
474 
475  // Update the interface compositions
476  if (this->interfaceCompositionModels_.found(key12))
477  {
478  this->interfaceCompositionModels_[key12]->update(Tf);
479  }
480  if (this->interfaceCompositionModels_.found(key21))
481  {
482  this->interfaceCompositionModels_[key21]->update(Tf);
483  }
484  }
485  }
486 }
487 
488 
489 template<class BasePhaseSystem>
491 {
492  if (BasePhaseSystem::read())
493  {
494  bool readOK = true;
495 
496  // Models ...
497 
498  return readOK;
499  }
500  else
501  {
502  return false;
503  }
504 }
505 
506 
507 // ************************************************************************* //
A simple container for copying or transferring objects of type <T>.
Definition: Xfer.H:85
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
virtual ~InterfaceCompositionPhaseChangePhaseSystem()
Destructor.
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:243
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
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
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:367
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
virtual Xfer< PtrList< volScalarField > > dmdts() const
Return the mass transfer rates for each phase.
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:106
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.
virtual Xfer< PtrList< volScalarField > > dmdts() const
Return the mass transfer rates for each phase.