All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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-2022 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 
27 #include "heatTransferModel.H"
29 #include "fvcVolumeIntegrate.H"
30 #include "fvmSup.H"
31 #include "rhoReactionThermo.H"
32 
33 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
34 
35 template<class BasePhaseSystem>
37 (
38  PtrList<volScalarField>& dmdts
39 ) const
40 {
41  forAllConstIter(phaseSystem::dmdtfTable, dmdtfs_, dmdtfIter)
42  {
43  const phaseInterface interface(*this, dmdtfIter.key());
44 
45  addField(interface.phase1(), "dmdt", *dmdtfIter(), dmdts);
46  addField(interface.phase2(), "dmdt", - *dmdtfIter(), dmdts);
47  }
48 
49  forAllConstIter(phaseSystem::dmdtfTable, nDmdtfs_, nDmdtfIter)
50  {
51  const phaseInterface interface(*this, nDmdtfIter.key());
52 
53  addField(interface.phase1(), "dmdt", *nDmdtfIter(), dmdts);
54  addField(interface.phase2(), "dmdt", - *nDmdtfIter(), dmdts);
55  }
56 }
57 
58 
59 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
60 
61 template<class BasePhaseSystem>
64 (
65  const fvMesh& mesh
66 )
67 :
68  BasePhaseSystem(mesh),
69  volatile_(this->template lookupOrDefault<word>("volatile", "none")),
70  dmdt0s_(this->phases().size()),
71  pressureImplicit_
72  (
73  this->template lookupOrDefault<Switch>("pressureImplicit", true)
74  )
75 {
76  this->generateInterfacialModels(saturationModels_);
77 
78  // Check that models have been specified in the correct combinations
80  (
81  saturationModelTable,
82  saturationModels_,
83  saturationModelIter
84  )
85  {
86  const phaseInterface& interface = saturationModelIter()->interface();
87  const phaseModel& phase1 = interface.phase1();
88  const phaseModel& phase2 = interface.phase2();
89 
90  this->template validateMassTransfer<saturationModel>(interface);
91 
92  if
93  (
94  !this->heatTransferModels_.found(interface)
95  || !this->heatTransferModels_[interface]->haveModelInThe(phase1)
96  || !this->heatTransferModels_[interface]->haveModelInThe(phase2)
97  )
98  {
100  << "A heat transfer model for both sides of the "
101  << interface.name() << " interface is not specified. This is "
102  << "required by the corresponding saturation model"
103  << exit(FatalError);
104  }
105  }
106 
107  // Generate interfacial mass transfer fields, initially assumed to be zero
109  (
110  saturationModelTable,
111  saturationModels_,
112  saturationModelIter
113  )
114  {
115  const phaseInterface& interface = saturationModelIter()->interface();
116 
117  dmdtfs_.insert
118  (
119  interface,
120  new volScalarField
121  (
122  IOobject
123  (
125  (
126  "thermalPhaseChange:dmdtf",
127  interface.name()
128  ),
129  this->mesh().time().timeName(),
130  this->mesh(),
133  ),
134  this->mesh(),
136  )
137  );
138 
139  d2mdtdpfs_.insert
140  (
141  interface,
142  new volScalarField
143  (
144  IOobject
145  (
147  (
148  "thermalPhaseChange:d2mdtdpf",
149  interface.name()
150  ),
151  this->mesh().time().timeName(),
152  this->mesh(),
155  ),
156  this->mesh(),
158  )
159  );
160 
161  Tfs_.insert
162  (
163  interface,
164  new volScalarField
165  (
166  IOobject
167  (
169  (
170  "thermalPhaseChange:Tf",
171  interface.name()
172  ),
173  this->mesh().time().timeName(),
174  this->mesh(),
177  ),
178  (
179  interface.phase1().thermo().T()
180  + interface.phase2().thermo().T()
181  )/2
182  )
183  );
184 
185  nDmdtfs_.insert
186  (
187  interface,
188  new volScalarField
189  (
190  IOobject
191  (
193  (
194  "thermalPhaseChange:nucleation:dmdtf",
195  interface.name()
196  ),
197  this->mesh().time().timeName(),
198  this->mesh(),
201  ),
202  this->mesh(),
204  )
205  );
206  }
207 }
208 
209 
210 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
211 
212 template<class BasePhaseSystem>
215 {}
216 
217 
218 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
219 
220 template<class BasePhaseSystem>
223 (
224  const phaseInterfaceKey& key
225 ) const
226 {
227  return saturationModels_[key];
228 }
229 
230 
231 template<class BasePhaseSystem>
234 (
235  const phaseInterfaceKey& key
236 ) const
237 {
238  tmp<volScalarField> tDmdtf = BasePhaseSystem::dmdtf(key);
239 
240  if (dmdtfs_.found(key))
241  {
242  tDmdtf.ref() += *dmdtfs_[key];
243  }
244 
245  if (nDmdtfs_.found(key))
246  {
247  tDmdtf.ref() += *nDmdtfs_[key];
248  }
249 
250  return tDmdtf;
251 }
252 
253 
254 template<class BasePhaseSystem>
257 {
258  PtrList<volScalarField> dmdts(BasePhaseSystem::dmdts());
259 
260  addDmdts(dmdts);
261 
262  return dmdts;
263 }
264 
265 
266 template<class BasePhaseSystem>
269 {
270  PtrList<volScalarField> d2mdtdps(BasePhaseSystem::d2mdtdps());
271 
272  forAllConstIter(phaseSystem::dmdtfTable, d2mdtdpfs_, d2mdtdpfIter)
273  {
274  const phaseInterface interface(*this, d2mdtdpfIter.key());
275 
276  addField(interface.phase1(), "d2mdtdp", *d2mdtdpfIter(), d2mdtdps);
277  addField(interface.phase2(), "d2mdtdp", - *d2mdtdpfIter(), d2mdtdps);
278  }
279 
280  return d2mdtdps;
281 }
282 
283 
284 template<class BasePhaseSystem>
288 {
289  autoPtr<phaseSystem::momentumTransferTable> eqnsPtr =
291 
292  phaseSystem::momentumTransferTable& eqns = eqnsPtr();
293 
294  this->addDmdtUfs(dmdtfs_, eqns);
295  this->addDmdtUfs(nDmdtfs_, eqns);
296 
297  return eqnsPtr;
298 }
299 
300 
301 template<class BasePhaseSystem>
305 {
306  autoPtr<phaseSystem::momentumTransferTable> eqnsPtr =
307  BasePhaseSystem::momentumTransferf();
308 
309  phaseSystem::momentumTransferTable& eqns = eqnsPtr();
310 
311  this->addDmdtUfs(dmdtfs_, eqns);
312  this->addDmdtUfs(nDmdtfs_, eqns);
313 
314  return eqnsPtr;
315 }
316 
317 
318 template<class BasePhaseSystem>
321 {
322  autoPtr<phaseSystem::heatTransferTable> eqnsPtr =
323  BasePhaseSystem::heatTransfer();
324 
325  phaseSystem::heatTransferTable& eqns = eqnsPtr();
326 
327  // Create temperatures at which to evaluate nucleation mass transfers
329  forAllConstIter(phaseSystem::dmdtfTable, nDmdtfs_, nDmdtfIter)
330  {
331  const phaseInterface interface(*this, nDmdtfIter.key());
332  const saturationModel& satModel = this->saturation(nDmdtfIter.key());
333 
334  Tns.insert
335  (
336  interface,
337  satModel.Tsat(interface.phase1().thermo().p()).ptr()
338  );
339  }
340 
341  // Mass transfer terms
342  if (volatile_ != "none")
343  {
344  {
345  phaseSystem::dmidtfTable dmidtfs;
346 
347  forAllConstIter(phaseSystem::dmdtfTable, dmdtfs_, dmdtfIter)
348  {
349  const phaseInterface interface(*this, dmdtfIter.key());
350 
351  dmidtfs.insert(interface, new HashPtrTable<volScalarField>());
352  dmidtfs[interface]->insert
353  (
354  volatile_,
355  new volScalarField(*dmdtfIter())
356  );
357  }
358 
359  this->addDmidtHefs
360  (
361  dmidtfs,
362  Tfs_,
363  latentHeatScheme::upwind,
364  latentHeatTransfer::heat,
365  eqns
366  );
367  }
368  {
369  phaseSystem::dmidtfTable nDmidtfs;
370 
371  forAllConstIter(phaseSystem::dmdtfTable, nDmdtfs_, nDmdtfIter)
372  {
373  const phaseInterface interface(*this, nDmdtfIter.key());
374 
375  nDmidtfs.insert(interface, new HashPtrTable<volScalarField>());
376  nDmidtfs[interface]->insert
377  (
378  volatile_,
379  new volScalarField(*nDmdtfIter())
380  );
381  }
382 
383  this->addDmidtHefs
384  (
385  nDmidtfs,
386  Tns,
387  0,
388  latentHeatScheme::upwind,
389  eqns
390  );
391  }
392  }
393  else
394  {
395  this->addDmdtHefs
396  (
397  dmdtfs_,
398  Tfs_,
399  latentHeatScheme::upwind,
400  latentHeatTransfer::heat,
401  eqns
402  );
403  this->addDmdtHefs
404  (
405  nDmdtfs_,
406  Tns,
407  0,
408  latentHeatScheme::upwind,
409  eqns
410  );
411  }
412 
413  // Lagging
414  {
415  PtrList<volScalarField> dmdts(this->phases().size());
416 
417  addDmdts(dmdts);
418 
420  {
421  const phaseModel& phase = phaseIter();
422 
423  if (dmdt0s_.set(phase.index()))
424  {
425  *eqns[phase.name()] +=
426  fvm::Sp
427  (
428  dmdt0s_[phase.index()] - dmdts[phase.index()],
429  phase.thermo().he()
430  );
431  }
432  }
433  }
434 
435  return eqnsPtr;
436 }
437 
438 
439 template<class BasePhaseSystem>
442 {
443  autoPtr<phaseSystem::specieTransferTable> eqnsPtr =
445 
446  phaseSystem::specieTransferTable& eqns = eqnsPtr();
447 
448  if (volatile_ != "none")
449  {
450  {
451  phaseSystem::dmidtfTable dmidtfs;
452 
453  forAllConstIter(phaseSystem::dmdtfTable, dmdtfs_, dmdtfIter)
454  {
455  const phaseInterface interface(*this, dmdtfIter.key());
456 
457  dmidtfs.insert(interface, new HashPtrTable<volScalarField>());
458  dmidtfs[interface]->insert
459  (
460  volatile_,
461  new volScalarField(*dmdtfIter())
462  );
463  }
464 
465  this->addDmidtYf(dmidtfs, eqns);
466  }
467 
468  {
469  phaseSystem::dmidtfTable nDmidtfs;
470 
471  forAllConstIter(phaseSystem::dmdtfTable, nDmdtfs_, nDmdtfIter)
472  {
473  const phaseInterface interface(*this, nDmdtfIter.key());
474 
475  nDmidtfs.insert(interface, new HashPtrTable<volScalarField>());
476  nDmidtfs[interface]->insert
477  (
478  volatile_,
479  new volScalarField(*nDmdtfIter())
480  );
481  }
482 
483  this->addDmidtYf(nDmidtfs, eqns);
484  }
485  }
486  else
487  {
488  this->addDmdtYfs(dmdtfs_, eqns);
489  this->addDmdtYfs(nDmdtfs_, eqns);
490  }
491 
492  return eqnsPtr;
493 }
494 
495 
496 template<class BasePhaseSystem>
499 {
500  dmdt0s_ = PtrList<volScalarField>(this->phases().size());
501 
502  addDmdts(dmdt0s_);
503 
504  BasePhaseSystem::correctContinuityError();
505 }
506 
507 
508 template<class BasePhaseSystem>
509 void
511 {
512  typedef compressible::alphatPhaseChangeWallFunctionFvPatchScalarField
513  alphatPhaseChangeWallFunction;
514 
516  (
517  saturationModelTable,
518  saturationModels_,
519  saturationModelIter
520  )
521  {
522  const phaseInterface& interface = saturationModelIter()->interface();
523  const phaseModel& phase1 = interface.phase1();
524  const phaseModel& phase2 = interface.phase2();
525  const rhoThermo& thermo1 = phase1.thermo();
526  const rhoThermo& thermo2 = phase2.thermo();
527  const volScalarField& T1(thermo1.T());
528  const volScalarField& T2(thermo2.T());
529 
530  const sidedBlendedHeatTransferModel& heatTransferModel =
531  this->heatTransferModels_[interface];
532 
533  // Interfacial mass transfer update
534  {
535  volScalarField& dmdtf(*this->dmdtfs_[interface]);
536  volScalarField& Tf(*this->Tfs_[interface]);
537 
538  const volScalarField Tsat(saturationModelIter()->Tsat(thermo1.p()));
539 
540  const volScalarField L
541  (
542  volatile_ != "none"
543  ? this->Li
544  (
545  interface,
546  volatile_,
547  dmdtf,
548  Tsat,
549  latentHeatScheme::upwind
550  )
551  : this->L
552  (
553  interface,
554  dmdtf,
555  Tsat,
556  latentHeatScheme::upwind
557  )
558  );
559 
560  volScalarField H1(heatTransferModel.modelInThe(phase1).K(0));
561  volScalarField H2(heatTransferModel.modelInThe(phase2).K(0));
562 
563  volScalarField dmdtfNew((H1*(Tsat - T1) + H2*(Tsat - T2))/L);
564 
565  if (volatile_ != "none")
566  {
567  dmdtfNew *=
568  neg0(dmdtfNew)*phase1.Y(volatile_)
569  + pos(dmdtfNew)*phase2.Y(volatile_);
570  }
571 
572  if (pressureImplicit_)
573  {
574  volScalarField& d2mdtdpf(*this->d2mdtdpfs_[interface]);
575 
576  const dimensionedScalar dp(rootSmall*thermo1.p().average());
577 
578  const volScalarField dTsatdp
579  (
580  (
581  saturationModelIter()->Tsat(thermo1.p() + dp/2)
582  - saturationModelIter()->Tsat(thermo1.p() - dp/2)
583  )/dp
584  );
585 
586  d2mdtdpf = (H1 + H2)*dTsatdp/L;
587 
588  if (volatile_ != "none")
589  {
590  d2mdtdpf *=
591  neg0(dmdtfNew)*phase1.Y(volatile_)
592  + pos(dmdtfNew)*phase2.Y(volatile_);
593  }
594  }
595 
596  H1 = heatTransferModel.modelInThe(phase1).K();
597  H2 = heatTransferModel.modelInThe(phase2).K();
598 
599  // Limit the H[12] to avoid /0
600  H1.max(small);
601  H2.max(small);
602 
603  Tf = (H1*T1 + H2*T2 + dmdtfNew*L)/(H1 + H2);
604 
605  Info<< Tf.name()
606  << ": min = " << min(Tf.primitiveField())
607  << ", mean = " << average(Tf.primitiveField())
608  << ", max = " << max(Tf.primitiveField())
609  << endl;
610 
611  const scalar dmdtfRelax =
612  this->mesh().solution().fieldRelaxationFactor(dmdtf.member());
613 
614  dmdtf = (1 - dmdtfRelax)*dmdtf + dmdtfRelax*dmdtfNew;
615 
616  Info<< dmdtf.name()
617  << ": min = " << min(dmdtf.primitiveField())
618  << ", mean = " << average(dmdtf.primitiveField())
619  << ", max = " << max(dmdtf.primitiveField())
620  << ", integral = " << fvc::domainIntegrate(dmdtf).value()
621  << endl;
622  }
623 
624  // Nucleation mass transfer update
625  {
626  volScalarField& nDmdtf(*this->nDmdtfs_[interface]);
627 
628  nDmdtf = Zero;
629 
630  bool wallBoilingActive = false;
631 
632  forAllConstIter(phaseInterface, interface, interfaceIter)
633  {
634  const phaseModel& phase = interfaceIter();
635 
636  const word alphatName =
637  IOobject::groupName("alphat", phase.name());
638 
639  if (phase.mesh().foundObject<volScalarField>(alphatName))
640  {
641  const volScalarField& alphat =
642  phase.mesh().lookupObject<volScalarField>(alphatName);
643 
644  forAll(alphat.boundaryField(), patchi)
645  {
646  const fvPatchScalarField& alphatp =
647  alphat.boundaryField()[patchi];
648 
649  if (isA<alphatPhaseChangeWallFunction>(alphatp))
650  {
651  const alphatPhaseChangeWallFunction& alphatw =
652  refCast<const alphatPhaseChangeWallFunction>
653  (
654  alphatp
655  );
656 
657  if (alphatw.activeInterface(interface))
658  {
659  wallBoilingActive = true;
660 
661  const scalarField& patchDmdtf =
662  alphatw.dmdtf(interface);
663 
664  const scalar sign =
665  interfaceIter.index() == 0 ? +1 : -1;
666 
667  forAll(patchDmdtf, facei)
668  {
669  const label celli =
670  alphatw.patch().faceCells()[facei];
671  nDmdtf[celli] -= sign*patchDmdtf[facei];
672  }
673  }
674  }
675  }
676  }
677  }
678 
679  if (wallBoilingActive)
680  {
681  Info<< nDmdtf.name()
682  << ": min = " << min(nDmdtf.primitiveField())
683  << ", mean = " << average(nDmdtf.primitiveField())
684  << ", max = " << max(nDmdtf.primitiveField())
685  << ", integral = " << fvc::domainIntegrate(nDmdtf).value()
686  << endl;
687  }
688  }
689  }
690 }
691 
692 
693 template<class BasePhaseSystem>
695 {
696  if (BasePhaseSystem::read())
697  {
698  bool readOK = true;
699 
700  // Models ...
701 
702  return readOK;
703  }
704  else
705  {
706  return false;
707  }
708 }
709 
710 
711 // ************************************************************************* //
dimensionedScalar sign(const dimensionedScalar &ds)
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
virtual PtrList< volScalarField > dmdts() const
Return the mass transfer rates for each phase.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
ThermalPhaseChangePhaseSystem(const fvMesh &)
Construct from fvMesh.
HashPtrTable< fvVectorMatrix > momentumTransferTable
Definition: phaseSystem.H:76
const dimensionSet dimPressure
virtual autoPtr< phaseSystem::specieTransferTable > specieTransfer() const
Return the specie transfer matrices.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
HashPtrTable< HashPtrTable< volScalarField >, phaseInterfaceKey, phaseInterfaceKey::hash > dmidtfTable
Definition: phaseSystem.H:102
virtual PtrList< volScalarField > dmdts() const
Return the mass transfer rates for each phase.
HashPtrTable< fvScalarMatrix > heatTransferTable
Definition: phaseSystem.H:78
dimensioned< Type > domainIntegrate(const GeometricField< Type, fvPatchField, volMesh > &vf)
HashPtrTable< fvScalarMatrix > specieTransferTable
Definition: phaseSystem.H:80
void read(Istream &, label &, const dictionary &)
In-place read with dictionary lookup.
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
fvMesh & mesh
Class to provide interfacial heat and mass transfer between a number of phases according the interfac...
bool insert(const Key &, const T &newElmt)
Insert a new hashedEntry.
Definition: HashTableI.H:80
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
dimensionedScalar pos(const dimensionedScalar &ds)
const dimensionSet dimTime
virtual PtrList< volScalarField > d2mdtdps() const
Return the mass transfer linearisation coeffs for each phase.
virtual void correctInterfaceThermo()
Correct the interface thermodynamics.
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:113
fvPatchField< scalar > fvPatchScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
static word groupName(Name name, const word &group)
const dimensionSet dimDensity
dimensionedScalar neg0(const dimensionedScalar &ds)
PtrListDictionary< phaseModel > phaseModelList
Definition: phaseSystem.H:82
virtual autoPtr< phaseSystem::heatTransferTable > heatTransfer() const
Return the heat transfer matrices.
virtual autoPtr< phaseSystem::momentumTransferTable > momentumTransferf()
Return the momentum transfer matrices for the face-based algorithm.
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
Volume integrate volField creating a volField.
static const zero Zero
Definition: zero.H:97
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &df)
phaseSystem::specieTransferTable & specieTransfer(specieTransferPtr())
const Mesh & mesh() const
Return mesh.
virtual bool read()
Read base phaseProperties dictionary.
HashPtrTable< volScalarField, phaseInterfaceKey, phaseInterfaceKey::hash > dmdtfTable
Definition: phaseSystem.H:93
Model to describe the dependence of saturation pressure on temperature, and vice versa.
phaseSystem::momentumTransferTable & momentumTransfer(momentumTransferPtr())
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
virtual void correctContinuityError()
Store phase dmdts at the during the continuity error update.
phaseSystem::phaseModelList & phases
Definition: createFields.H:12
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
virtual tmp< volScalarField > dmdtf(const phaseInterfaceKey &key) const
Return the mass transfer rate for an interface.
A class for managing temporary objects.
Definition: PtrList.H:53
virtual tmp< volScalarField > dmdtf(const phaseInterfaceKey &key) const
Return the mass transfer rate for an interface.
virtual autoPtr< phaseSystem::momentumTransferTable > momentumTransfer()
Return the momentum transfer matrices for the cell-based algorithm.
Calculate the matrix for implicit and explicit sources.
PtrList< volScalarField > d2mdtdps(fluid.d2mdtdps())
const saturationModel & saturation(const phaseInterfaceKey &key) const
Return the saturationModel.