ThermalPhaseChangePhaseSystem.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 "fvcVolumeIntegrate.H"
29 #include "fvmSup.H"
30 
31 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
32 
33 template<class BasePhaseSystem>
36 (
37  const phasePairKey& key
38 ) const
39 {
40  if (!iDmdt_.found(key))
41  {
42  return phaseSystem::dmdt(key);
43  }
44 
45  const scalar dmdtSign(Pair<word>::compare(iDmdt_.find(key).key(), key));
46 
47  return dmdtSign**iDmdt_[key];
48 }
49 
50 
51 template<class BasePhaseSystem>
54 (
55  const phasePairKey& key
56 ) const
57 {
58  if (!wDmdt_.found(key))
59  {
60  return phaseSystem::dmdt(key);
61  }
62 
63  const scalar dmdtSign(Pair<word>::compare(wDmdt_.find(key).key(), key));
64 
65  return dmdtSign**wDmdt_[key];
66 }
67 
68 
69 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
70 
71 template<class BasePhaseSystem>
74 (
75  const fvMesh& mesh
76 )
77 :
78  BasePhaseSystem(mesh),
79  volatile_(this->template lookupOrDefault<word>("volatile", "none")),
80  saturationModel_
81  (
82  saturationModel::New(this->subDict("saturationModel"), mesh)
83  ),
84  phaseChange_(this->lookup("phaseChange"))
85 {
86 
88  (
90  this->phasePairs_,
91  phasePairIter
92  )
93  {
94  const phasePair& pair(phasePairIter());
95 
96  if (pair.ordered())
97  {
98  continue;
99  }
100 
101  // Initially assume no mass transfer
102  iDmdt_.insert
103  (
104  pair,
105  new volScalarField
106  (
107  IOobject
108  (
109  IOobject::groupName("iDmdt", pair.name()),
110  this->mesh().time().timeName(),
111  this->mesh(),
114  ),
115  this->mesh(),
117  )
118  );
119 
120  // Initially assume no mass transfer
121  wDmdt_.insert
122  (
123  pair,
124  new volScalarField
125  (
126  IOobject
127  (
128  IOobject::groupName("wDmdt", pair.name()),
129  this->mesh().time().timeName(),
130  this->mesh(),
133  ),
134  this->mesh(),
136  )
137  );
138 
139  // Initially assume no mass transfer
140  wMDotL_.insert
141  (
142  pair,
143  new volScalarField
144  (
145  IOobject
146  (
147  IOobject::groupName("wMDotL", pair.name()),
148  this->mesh().time().timeName(),
149  this->mesh(),
152  ),
153  this->mesh(),
155  )
156  );
157  }
158 }
159 
160 
161 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
162 
163 template<class BasePhaseSystem>
166 {}
167 
168 
169 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
170 
171 template<class BasePhaseSystem>
174 {
175  return saturationModel_();
176 }
177 
178 
179 template<class BasePhaseSystem>
182 (
183  const phasePairKey& key
184 ) const
185 {
186  return BasePhaseSystem::dmdt(key) + this->iDmdt(key) + this->wDmdt(key);
187 }
188 
189 
190 template<class BasePhaseSystem>
193 {
194  PtrList<volScalarField> dmdts(BasePhaseSystem::dmdts());
195 
196  forAllConstIter(iDmdtTable, iDmdt_, iDmdtIter)
197  {
198  const phasePair& pair = this->phasePairs_[iDmdtIter.key()];
199  const volScalarField& iDmdt = *iDmdtIter();
200 
201  this->addField(pair.phase1(), "dmdt", iDmdt, dmdts);
202  this->addField(pair.phase2(), "dmdt", - iDmdt, dmdts);
203  }
204 
205  forAllConstIter(wDmdtTable, wDmdt_, wDmdtIter)
206  {
207  const phasePair& pair = this->phasePairs_[wDmdtIter.key()];
208  const volScalarField& wDmdt = *wDmdtIter();
209 
210  this->addField(pair.phase1(), "dmdt", wDmdt, dmdts);
211  this->addField(pair.phase2(), "dmdt", - wDmdt, dmdts);
212  }
213 
214  return dmdts;
215 }
216 
217 
218 template<class BasePhaseSystem>
221 {
222  autoPtr<phaseSystem::heatTransferTable> eqnsPtr =
223  BasePhaseSystem::heatTransfer();
224 
225  phaseSystem::heatTransferTable& eqns = eqnsPtr();
226 
227  // Add boundary term
229  (
231  this->phasePairs_,
232  phasePairIter
233  )
234  {
235  if (this->wMDotL_.found(phasePairIter.key()))
236  {
237  const phasePair& pair(phasePairIter());
238 
239  if (pair.ordered())
240  {
241  continue;
242  }
243 
244  const phaseModel& phase1 = pair.phase1();
245  const phaseModel& phase2 = pair.phase2();
246 
247  *eqns[phase1.name()] += negPart(*this->wMDotL_[pair]);
248  *eqns[phase2.name()] -= posPart(*this->wMDotL_[pair]);
249 
250  if
251  (
252  phase1.thermo().he().member() == "e"
253  || phase2.thermo().he().member() == "e"
254  )
255  {
256  const volScalarField dmdt
257  (
258  this->iDmdt(pair) + this->wDmdt(pair)
259  );
260 
261  if (phase1.thermo().he().member() == "e")
262  {
263  *eqns[phase1.name()] +=
264  phase1.thermo().p()*dmdt/phase1.thermo().rho();
265  }
266 
267  if (phase2.thermo().he().member() == "e")
268  {
269  *eqns[phase2.name()] -=
270  phase2.thermo().p()*dmdt/phase2.thermo().rho();
271  }
272  }
273  }
274  }
275 
276  return eqnsPtr;
277 }
278 
279 
280 template<class BasePhaseSystem>
283 {
284  autoPtr<phaseSystem::massTransferTable> eqnsPtr =
286 
287  phaseSystem::massTransferTable& eqns = eqnsPtr();
288 
290  (
292  this->phasePairs_,
293  phasePairIter
294  )
295  {
296  const phasePair& pair(phasePairIter());
297 
298  if (pair.ordered())
299  {
300  continue;
301  }
302 
303  const phaseModel& phase = pair.phase1();
304  const phaseModel& otherPhase = pair.phase2();
305 
306  const PtrList<volScalarField>& Yi = phase.Y();
307 
308  forAll(Yi, i)
309  {
310  if (Yi[i].member() != volatile_)
311  {
312  continue;
313  }
314 
315  const word name
316  (
317  IOobject::groupName(volatile_, phase.name())
318  );
319 
320  const word otherName
321  (
322  IOobject::groupName(volatile_, otherPhase.name())
323  );
324 
325  // Note that the phase YiEqn does not contain a continuity error
326  // term, so these additions represent the entire mass transfer
327 
328  const volScalarField dmdt(this->iDmdt(pair) + this->wDmdt(pair));
329 
330  *eqns[name] += dmdt;
331  *eqns[otherName] -= dmdt;
332  }
333  }
334 
335  return eqnsPtr;
336 }
337 
338 
339 template<class BasePhaseSystem>
340 void
342 {
343  typedef compressible::alphatPhaseChangeWallFunctionFvPatchScalarField
344  alphatPhaseChangeWallFunction;
345 
347  (
348  typename BasePhaseSystem::heatTransferModelTable,
349  this->heatTransferModels_,
350  heatTransferModelIter
351  )
352  {
353  const phasePair& pair
354  (
355  this->phasePairs_[heatTransferModelIter.key()]
356  );
357 
358  const phaseModel& phase1 = pair.phase1();
359  const phaseModel& phase2 = pair.phase2();
360 
361  const volScalarField& T1(phase1.thermo().T());
362  const volScalarField& T2(phase2.thermo().T());
363 
364  const volScalarField& he1(phase1.thermo().he());
365  const volScalarField& he2(phase2.thermo().he());
366 
367  const volScalarField& p(phase1.thermo().p());
368 
369  volScalarField& iDmdt(*this->iDmdt_[pair]);
370  volScalarField& Tf(*this->Tf_[pair]);
371 
372  const volScalarField Tsat(saturationModel_->Tsat(phase1.thermo().p()));
373 
374  volScalarField hf1
375  (
376  he1.member() == "e"
377  ? phase1.thermo().he(p, Tsat) + p/phase1.rho()
378  : phase1.thermo().he(p, Tsat)
379  );
380  volScalarField hf2
381  (
382  he2.member() == "e"
383  ? phase2.thermo().he(p, Tsat) + p/phase2.rho()
384  : phase2.thermo().he(p, Tsat)
385  );
386 
387  volScalarField h1
388  (
389  he1.member() == "e"
390  ? he1 + p/phase1.rho()
391  : tmp<volScalarField>(he1)
392  );
393 
394  volScalarField h2
395  (
396  he2.member() == "e"
397  ? he2 + p/phase2.rho()
398  : tmp<volScalarField>(he2)
399  );
400 
402  (
403  (neg0(iDmdt)*hf2 + pos(iDmdt)*h2)
404  - (pos0(iDmdt)*hf1 + neg(iDmdt)*h1)
405  );
406 
407  volScalarField iDmdtNew(iDmdt);
408 
409  if (phaseChange_)
410  {
411  volScalarField H1(heatTransferModelIter().first()->K(0));
412  volScalarField H2(heatTransferModelIter().second()->K(0));
413 
414  iDmdtNew = (H1*(Tsat - T1) + H2*(Tsat - T2))/L;
415  }
416  else
417  {
418  iDmdtNew == dimensionedScalar(iDmdt.dimensions(), 0);
419  }
420 
421  volScalarField H1(heatTransferModelIter().first()->K());
422  volScalarField H2(heatTransferModelIter().second()->K());
423 
424  // Limit the H[12] to avoid /0
425  H1.max(small);
426  H2.max(small);
427 
428  Tf = (H1*T1 + H2*T2 + iDmdtNew*L)/(H1 + H2);
429 
430  Info<< "Tf." << pair.name()
431  << ": min = " << min(Tf.primitiveField())
432  << ", mean = " << average(Tf.primitiveField())
433  << ", max = " << max(Tf.primitiveField())
434  << endl;
435 
436  scalar iDmdtRelax(this->mesh().fieldRelaxationFactor("iDmdt"));
437  iDmdt = (1 - iDmdtRelax)*iDmdt + iDmdtRelax*iDmdtNew;
438 
439  if (phaseChange_)
440  {
441  Info<< "iDmdt." << pair.name()
442  << ": min = " << min(iDmdt.primitiveField())
443  << ", mean = " << average(iDmdt.primitiveField())
444  << ", max = " << max(iDmdt.primitiveField())
445  << ", integral = " << fvc::domainIntegrate(iDmdt).value()
446  << endl;
447  }
448 
449  volScalarField& wDmdt(*this->wDmdt_[pair]);
450  volScalarField& wMDotL(*this->wMDotL_[pair]);
451  wDmdt = Zero;
452  wMDotL = Zero;
453 
454  bool wallBoilingActive = false;
455 
456  forAllConstIter(phasePair, pair, iter)
457  {
458  const phaseModel& phase = iter();
459  const phaseModel& otherPhase = iter.otherPhase();
460 
461  if
462  (
463  phase.mesh().foundObject<volScalarField>
464  (
465  "alphat." + phase.name()
466  )
467  )
468  {
469  const volScalarField& alphat =
470  phase.mesh().lookupObject<volScalarField>
471  (
472  "alphat." + phase.name()
473  );
474 
475  const fvPatchList& patches = this->mesh().boundary();
476  forAll(patches, patchi)
477  {
478  const fvPatch& currPatch = patches[patchi];
479 
480  if
481  (
482  isA<alphatPhaseChangeWallFunction>
483  (
484  alphat.boundaryField()[patchi]
485  )
486  )
487  {
488  const alphatPhaseChangeWallFunction& PCpatch =
489  refCast<const alphatPhaseChangeWallFunction>
490  (
491  alphat.boundaryField()[patchi]
492  );
493 
494  phasePairKey key(phase.name(), otherPhase.name());
495 
496  if (PCpatch.activePhasePair(key))
497  {
498  wallBoilingActive = true;
499 
500  const scalarField& patchDmdt =
501  PCpatch.dmdt(key);
502  const scalarField& patchMDotL =
503  PCpatch.mDotL(key);
504 
505  const scalar sign
506  (
507  Pair<word>::compare(pair, key)
508  );
509 
510  forAll(patchDmdt, facei)
511  {
512  const label faceCelli =
513  currPatch.faceCells()[facei];
514  wDmdt[faceCelli] -= sign*patchDmdt[facei];
515  wMDotL[faceCelli] -= sign*patchMDotL[facei];
516  }
517  }
518  }
519  }
520  }
521  }
522 
523  if (wallBoilingActive)
524  {
525  Info<< "wDmdt." << pair.name()
526  << ": min = " << min(wDmdt.primitiveField())
527  << ", mean = " << average(wDmdt.primitiveField())
528  << ", max = " << max(wDmdt.primitiveField())
529  << ", integral = " << fvc::domainIntegrate(wDmdt).value()
530  << endl;
531  }
532  }
533 }
534 
535 
536 template<class BasePhaseSystem>
538 {
539  if (BasePhaseSystem::read())
540  {
541  bool readOK = true;
542 
543  // Models ...
544 
545  return readOK;
546  }
547  else
548  {
549  return false;
550  }
551 }
552 
553 
554 // ************************************************************************* //
dimensionedScalar sign(const dimensionedScalar &ds)
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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
HashPtrTable< fvScalarMatrix > massTransferTable
Definition: phaseSystem.H:79
HashTable< autoPtr< phasePair >, phasePairKey, phasePairKey::hash > phasePairTable
Definition: phaseSystem.H:87
virtual PtrList< volScalarField > dmdts() const
Return the mass transfer rates for each phase.
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
ThermalPhaseChangePhaseSystem(const fvMesh &)
Construct from fvMesh.
static int compare(const Pair< word > &a, const Pair< word > &b)
Compare Pairs.
Definition: Pair.H:142
phaseSystem::massTransferTable & massTransfer(massTransferPtr())
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
virtual PtrList< volScalarField > dmdts() const
Return the mass transfer rates for each phase.
virtual tmp< volScalarField > dmdt(const phasePairKey &key) const
Return the mass transfer rate for a pair.
virtual autoPtr< phaseSystem::massTransferTable > massTransfer() const
Return the mass transfer matrices.
HashPtrTable< fvScalarMatrix > heatTransferTable
Definition: phaseSystem.H:77
dimensioned< Type > domainIntegrate(const GeometricField< Type, fvPatchField, volMesh > &vf)
void read(Istream &, label &, const dictionary &)
In-place read with dictionary lookup.
CGAL::Exact_predicates_exact_constructions_kernel K
dimensionedScalar posPart(const dimensionedScalar &ds)
dimensionedScalar neg(const dimensionedScalar &ds)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
dimensionedScalar pos(const dimensionedScalar &ds)
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
stressControl lookup("compactNormalStress") >> compactNormalStress
virtual void correctInterfaceThermo()
Correct the interface thermodynamics.
dynamicFvMesh & mesh
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
static word groupName(Name name, const word &group)
dimensionedScalar neg0(const dimensionedScalar &ds)
virtual autoPtr< phaseSystem::heatTransferTable > heatTransfer() const
Return the heat transfer matrices.
Volume integrate volField creating a volField.
static const zero Zero
Definition: zero.H:97
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &df)
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
dimensionedScalar pos0(const dimensionedScalar &ds)
volScalarField & he2
Definition: EEqns.H:3
const Mesh & mesh() const
Return mesh.
virtual bool read()
Read base phaseProperties dictionary.
tmp< volScalarField > iDmdt(const phasePairKey &key) const
Return the interfacial mass transfer rate for a pair.
tmp< volScalarField > wDmdt(const phasePairKey &key) const
Return the boundary mass transfer rate for a pair.
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
const saturationModel & saturation() const
Return the saturationModel.
PtrList< fvPatch > fvPatchList
container classes for fvPatch
Definition: fvPatchList.H:45
label patchi
virtual ~ThermalPhaseChangePhaseSystem()
Destructor.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
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
virtual tmp< volScalarField > dmdt(const phasePairKey &key) const
Return the mass transfer rate for a pair.
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
volScalarField & p
A class for managing temporary objects.
Definition: PtrList.H:53
Calculate the matrix for implicit and explicit sources.
dimensionedScalar negPart(const dimensionedScalar &ds)