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-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 "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.xfer();
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  }
251 
252  return eqnsPtr;
253 }
254 
255 
256 template<class BasePhaseSystem>
259 {
260  autoPtr<phaseSystem::massTransferTable> eqnsPtr =
262 
263  phaseSystem::massTransferTable& eqns = eqnsPtr();
264 
266  (
268  this->phasePairs_,
269  phasePairIter
270  )
271  {
272  const phasePair& pair(phasePairIter());
273 
274  if (pair.ordered())
275  {
276  continue;
277  }
278 
279  const phaseModel& phase = pair.phase1();
280  const phaseModel& otherPhase = pair.phase2();
281 
282  const PtrList<volScalarField>& Yi = phase.Y();
283 
284  forAll(Yi, i)
285  {
286  if (Yi[i].member() != volatile_)
287  {
288  continue;
289  }
290 
291  const word name
292  (
293  IOobject::groupName(volatile_, phase.name())
294  );
295 
296  const word otherName
297  (
298  IOobject::groupName(volatile_, otherPhase.name())
299  );
300 
301  // Note that the phase YiEqn does not contain a continuity error
302  // term, so these additions represent the entire mass transfer
303 
304  const volScalarField dmdt(this->iDmdt(pair) + this->wDmdt(pair));
305 
306  *eqns[name] += dmdt;
307  *eqns[otherName] -= dmdt;
308  }
309  }
310 
311  return eqnsPtr;
312 }
313 
314 
315 template<class BasePhaseSystem>
316 void
318 {
319  typedef compressible::alphatPhaseChangeWallFunctionFvPatchScalarField
320  alphatPhaseChangeWallFunction;
321 
323  (
324  typename BasePhaseSystem::heatTransferModelTable,
325  this->heatTransferModels_,
326  heatTransferModelIter
327  )
328  {
329  const phasePair& pair
330  (
331  this->phasePairs_[heatTransferModelIter.key()]
332  );
333 
334  const phaseModel& phase1 = pair.phase1();
335  const phaseModel& phase2 = pair.phase2();
336 
337  const volScalarField& T1(phase1.thermo().T());
338  const volScalarField& T2(phase2.thermo().T());
339 
340  const volScalarField& he1(phase1.thermo().he());
341  const volScalarField& he2(phase2.thermo().he());
342 
343  volScalarField& iDmdt(*this->iDmdt_[pair]);
344  volScalarField& Tf(*this->Tf_[pair]);
345 
346  volScalarField hef1(phase1.thermo().he(phase1.thermo().p(), Tf));
347  volScalarField hef2(phase2.thermo().he(phase2.thermo().p(), Tf));
348 
350  (
351  (neg0(iDmdt)*hef2 + pos(iDmdt)*he2)
352  - (pos0(iDmdt)*hef1 + neg(iDmdt)*he1)
353  );
354 
355  volScalarField iDmdtNew(iDmdt);
356 
357  if (phaseChange_)
358  {
359  volScalarField H1(heatTransferModelIter().first()->K(0));
360  volScalarField H2(heatTransferModelIter().second()->K(0));
361 
362  Tf = saturationModel_->Tsat(phase1.thermo().p());
363 
364  iDmdtNew = (H1*(Tf - T1) + H2*(Tf - T2))/L;
365  }
366  else
367  {
368  iDmdtNew == dimensionedScalar("0", iDmdt.dimensions(), 0);
369  }
370 
371  volScalarField H1(heatTransferModelIter().first()->K());
372  volScalarField H2(heatTransferModelIter().second()->K());
373 
374  // Limit the H[12] to avoid /0
375  H1.max(small);
376  H2.max(small);
377 
378  Tf = (H1*T1 + H2*T2 + iDmdtNew*L)/(H1 + H2);
379 
380  Info<< "Tf." << pair.name()
381  << ": min = " << min(Tf.primitiveField())
382  << ", mean = " << average(Tf.primitiveField())
383  << ", max = " << max(Tf.primitiveField())
384  << endl;
385 
386  scalar iDmdtRelax(this->mesh().fieldRelaxationFactor("iDmdt"));
387  iDmdt = (1 - iDmdtRelax)*iDmdt + iDmdtRelax*iDmdtNew;
388 
389  if (phaseChange_)
390  {
391  Info<< "iDmdt." << pair.name()
392  << ": min = " << min(iDmdt.primitiveField())
393  << ", mean = " << average(iDmdt.primitiveField())
394  << ", max = " << max(iDmdt.primitiveField())
395  << ", integral = " << fvc::domainIntegrate(iDmdt).value()
396  << endl;
397  }
398 
399  volScalarField& wDmdt(*this->wDmdt_[pair]);
400  volScalarField& wMDotL(*this->wMDotL_[pair]);
401  wDmdt *= 0;
402  wMDotL *= 0;
403 
404  bool wallBoilingActive = false;
405 
406  forAllConstIter(phasePair, pair, iter)
407  {
408  const phaseModel& phase = iter();
409  const phaseModel& otherPhase = iter.otherPhase();
410 
411  if
412  (
413  phase.mesh().foundObject<volScalarField>
414  (
415  "alphat." + phase.name()
416  )
417  )
418  {
419  const volScalarField& alphat =
420  phase.mesh().lookupObject<volScalarField>
421  (
422  "alphat." + phase.name()
423  );
424 
425  const fvPatchList& patches = this->mesh().boundary();
426  forAll(patches, patchi)
427  {
428  const fvPatch& currPatch = patches[patchi];
429 
430  if
431  (
432  isA<alphatPhaseChangeWallFunction>
433  (
434  alphat.boundaryField()[patchi]
435  )
436  )
437  {
438  const alphatPhaseChangeWallFunction& PCpatch =
439  refCast<const alphatPhaseChangeWallFunction>
440  (
441  alphat.boundaryField()[patchi]
442  );
443 
444  phasePairKey key(phase.name(), otherPhase.name());
445 
446  if (PCpatch.activePhasePair(key))
447  {
448  wallBoilingActive = true;
449 
450  const scalarField& patchDmdt =
451  PCpatch.dmdt(key);
452  const scalarField& patchMDotL =
453  PCpatch.mDotL(key);
454 
455  const scalar sign
456  (
457  Pair<word>::compare(pair, key)
458  );
459 
460  forAll(patchDmdt, facei)
461  {
462  const label faceCelli =
463  currPatch.faceCells()[facei];
464  wDmdt[faceCelli] -= sign*patchDmdt[facei];
465  wMDotL[faceCelli] -= sign*patchMDotL[facei];
466  }
467  }
468  }
469  }
470  }
471  }
472 
473  if (wallBoilingActive)
474  {
475  Info<< "wDmdt." << pair.name()
476  << ": min = " << min(wDmdt.primitiveField())
477  << ", mean = " << average(wDmdt.primitiveField())
478  << ", max = " << max(wDmdt.primitiveField())
479  << ", integral = " << fvc::domainIntegrate(wDmdt).value()
480  << endl;
481  }
482  }
483 }
484 
485 
486 template<class BasePhaseSystem>
488 {
489  if (BasePhaseSystem::read())
490  {
491  bool readOK = true;
492 
493  // Models ...
494 
495  return readOK;
496  }
497  else
498  {
499  return false;
500  }
501 }
502 
503 
504 // ************************************************************************* //
dimensionedScalar sign(const dimensionedScalar &ds)
A simple container for copying or transferring objects of type <T>.
Definition: Xfer.H:85
virtual Xfer< PtrList< volScalarField > > dmdts() const
Return the mass transfer rates for each phase.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
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 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.
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.
virtual Xfer< PtrList< volScalarField > > dmdts() const
Return the mass transfer rates for each phase.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
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
A class for managing temporary objects.
Definition: PtrList.H:53
Calculate the matrix for implicit and explicit sources.
dimensionedScalar negPart(const dimensionedScalar &ds)