heatTransferSystem.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-2025 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 
26 #include "heatTransferSystem.H"
27 
28 #include "fvmSup.H"
29 
30 #include "heatTransferModel.H"
32 
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
40 }
41 
42 
44 (
45  "heatTransfer"
46 );
47 
48 
49 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
50 
51 Foam::IOobject Foam::heatTransferSystem::io(const phaseSystem& fluid) const
52 {
53  typeIOobject<IOdictionary> result
54  (
56  fluid.mesh().time().constant(),
57  fluid.mesh(),
60  );
61 
62  // Check for and warn about heat transfer models being in the old
63  // location in constant/phaseProperties
64  if (!result.headerOk())
65  {
66  result.readOpt() = IOobject::NO_READ;
67 
68  if
69  (
70  fluid.found(modelName<heatTransferModel>())
71  || !fluid.thermalPhases().empty()
72  )
73  {
75  << "Specifying a heat transfer model entry - "
76  << modelName<heatTransferModel>() << " - in "
77  << fluid.relativeObjectPath() << " is deprecated. The contents "
78  << "of this entry should now be specified in "
79  << result.relativeObjectPath() << "." << endl;
80  }
81  }
82  else
83  {
84  if (fluid.found(modelName<heatTransferModel>()))
85  {
87  << "Heat transfer model entry - "
88  << modelName<heatTransferModel>() << " - in "
89  << fluid.relativeObjectPath() << " is no longer used. The "
90  << "contents of this entry are now read from "
91  << result.relativeObjectPath() << "." << endl;
92  }
93  }
94 
95  return result;
96 }
97 
98 
99 const Foam::dictionary& Foam::heatTransferSystem::modelsDict() const
100 {
101  const word key = modelName<heatTransferModel>();
102 
103  return
104  readOpt() == IOobject::NO_READ && fluid_.found(key)
105  ? fluid_.subDict(key)
106  : *this;
107 }
108 
109 
110 template<class ... Args>
111 Foam::Pair<Foam::tmp<Foam::volScalarField>> Foam::heatTransferSystem::Hs
112 (
113  const phaseModel& phase1,
114  const phaseModel& phase2,
115  Args ... args
116 ) const
117 {
118  auto error = [&](const word& why)
119  {
121  << why << " two-resistance heat transfer models found that "
122  << "provide a heat transfer coefficient between phases "
123  << phase1.name() << " and " << phase2.name() << exit(FatalError);
124  };
125 
127 
128  const phaseInterface interface(phase1, phase2);
129 
130  auto iter = sidedModels_.find(interface);
131 
132  if (iter != sidedModels_.end())
133  {
134  HsPtr.set
135  (
137  (
138  iter()->KinThe(phase1, args ...),
139  iter()->KinThe(phase2, args ...)
140  )
141  );
142  }
143 
144  const Foam::fvModels& fvModels = Foam::fvModels::New(fluid_.mesh());
145 
146  forAll(fvModels, i)
147  {
148  if (!isA<fv::twoResistanceHeatTransfer>(fvModels[i])) continue;
149 
150  const fv::twoResistanceHeatTransfer& heatTransferFvModel =
151  refCast<const fv::twoResistanceHeatTransfer>(fvModels[i]);
152 
154  heatTransferFvModel.Ks(phase1, phase2, args ...);
155 
156  if (!Hs.first().valid() || !Hs.second().valid()) continue;
157 
158  if (HsPtr.valid()) error("Multiple");
159 
160  HsPtr.set(new Pair<tmp<volScalarField>>(Hs.first(), Hs.second()));
161  }
162 
163  if (!HsPtr.valid()) error("No");
164 
165  return HsPtr();
166 }
167 
168 
169 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
170 
172 (
173  const phaseSystem& fluid
174 )
175 :
176  IOdictionary(io(fluid)),
177  fluid_(fluid),
178  models_(),
179  sidedModels_()
180 {
181  // If we have no entries and there are no thermal phases then this case
182  // does not need any heat transfer modelling and we can just quit without
183  // trying to construct anything
184  if
185  (
187  && !fluid_.found(modelName<heatTransferModel>())
188  && fluid.thermalPhases().empty()
189  ) return;
190 
191  modelsTable models
192  (
193  generateBlendedInterfacialModels<blendedHeatTransferModel>
194  (
195  fluid,
196  modelsDict(),
197  wordHashSet(),
198  true
199  )
200  );
201 
202  models_.transfer(models);
203 
204  sidedModelsTable sidedModels
205  (
206  generateBlendedInterfacialModels<blendedSidedHeatTransferModel>
207  (
208  fluid,
209  modelsDict(),
210  wordHashSet(),
211  true
212  )
213  );
214 
215  sidedModels_.transfer(sidedModels);
216 
217  forAllConstIter(modelsTable, models_, modelIter)
218  {
219  if (sidedModels_.found(modelIter.key()))
220  {
221  const phaseInterface interface(fluid_, modelIter.key());
222 
223  FatalIOErrorInFunction(modelsDict())
224  << "One-resistance and two-resistance heat transfer models "
225  << "both specified between phases "
226  << interface.phase1().name() << " and "
227  << interface.phase2().name() << exit(FatalIOError);
228  }
229  }
230 }
231 
232 
233 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
234 
236 {}
237 
238 
239 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
240 
241 Foam::Pair<Foam::tmp<Foam::volScalarField>> Foam::heatTransferSystem::Hs
242 (
243  const phaseModel& phase1,
244  const phaseModel& phase2
245 ) const
246 {
247  return Hs<>(phase1, phase2);
248 }
249 
250 
251 Foam::Pair<Foam::tmp<Foam::volScalarField>> Foam::heatTransferSystem::Hs
252 (
253  const phaseModel& phase1,
254  const phaseModel& phase2,
255  const scalar residualAlpha
256 ) const
257 {
258  return Hs<scalar>(phase1, phase2, residualAlpha);
259 }
260 
261 
264 {
266  (
268  );
269  HashPtrTable<fvScalarMatrix>& eqns = eqnsPtr();
270 
271  forAll(fluid_.phases(), phasei)
272  {
273  const phaseModel& phase = fluid_.phases()[phasei];
274 
275  eqns.insert
276  (
277  phase.name(),
278  new fvScalarMatrix(phase.thermo().he(), dimEnergy/dimTime)
279  );
280  }
281 
282  forAllConstIter(modelsTable, models_, modelIter)
283  {
284  const phaseInterface interface(fluid_, modelIter.key());
285 
286  const volScalarField H(modelIter()->K());
287 
288  forAllConstIter(phaseInterface, interface, iter)
289  {
290  const phaseModel& phase = iter();
291  const phaseModel& otherPhase = iter.otherPhase();
292 
293  const volScalarField& he = phase.thermo().he();
294  const volScalarField Cpv(phase.thermo().Cpv());
295 
296  const volScalarField Hstabilised
297  (
298  iter.otherPhase()
299  /max(iter.otherPhase(), iter.otherPhase().residualAlpha())
300  *H
301  );
302 
303  *eqns[phase.name()] +=
304  Hstabilised*(otherPhase.thermo().T() - phase.thermo().T())
305  + Hstabilised/Cpv*he - fvm::Sp(Hstabilised/Cpv, he);
306  }
307  }
308 
309  forAllConstIter(sidedModelsTable, sidedModels_, sidedModelIter)
310  {
311  const phaseInterface interface(fluid_, sidedModelIter.key());
312 
314  (
315  sidedModelIter()->KinThe(interface.phase1()),
316  sidedModelIter()->KinThe(interface.phase2())
317  );
318 
319  const volScalarField HEff
320  (
321  Hs.first()*Hs.second()/(Hs.first() + Hs.second())
322  );
323 
324  forAllConstIter(phaseInterface, interface, iter)
325  {
326  const phaseModel& phase = iter();
327  const phaseModel& otherPhase = iter.otherPhase();
328 
329  const volScalarField& he = phase.thermo().he();
330  const volScalarField Cpv(phase.thermo().Cpv());
331 
332  const volScalarField& H = Hs[iter.index()];
333 
334  *eqns[phase.name()] +=
335  HEff*(otherPhase.thermo().T() - phase.thermo().T())
336  + H/Cpv*he - fvm::Sp(H/Cpv, he);
337  }
338  }
339 
340  return eqnsPtr;
341 }
342 
343 
345 {
346  if (regIOobject::read())
347  {
348  bool readOK = true;
349 
350  // models ...
351 
352  return readOK;
353  }
354  else
355  {
356  return false;
357  }
358 }
359 
360 
361 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:476
Generic GeometricField class.
void transfer(HashTable< T, Key, Hash > &)
Transfer the contents of the argument table into this table.
Definition: HashTable.C:587
bool insert(const Key &, const T &newElmt)
Insert a new hashedEntry.
Definition: HashTableI.H:80
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:138
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:57
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
@ MUST_READ_IF_MODIFIED
Definition: IOobject.H:119
readOption readOpt() const
Definition: IOobject.H:357
An ordered pair of two objects of type <Type> with first() and second() elements.
Definition: Pair.H:66
const Type & second() const
Return second.
Definition: PairI.H:121
const Type & first() const
Return first.
Definition: PairI.H:107
bool empty() const
Return true if the UPtrList is empty (ie, size() is zero)
Definition: UPtrListI.H:36
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
bool valid() const
Return true if the autoPtr valid (ie, the pointer is set)
Definition: autoPtrI.H:83
void set(T *)
Set pointer to that given.
Definition: autoPtrI.H:99
virtual const volScalarField & T() const =0
Temperature [K].
virtual const volScalarField & Cpv() const =0
Heat capacity at constant pressure/volume [J/kg/K].
virtual const volScalarField & he() const =0
Enthalpy/Internal energy [J/kg].
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:539
static autoPtr< dictionary > New(Istream &)
Construct top-level dictionary on freestore from Istream.
Definition: dictionaryIO.C:103
Class to handle errors and exceptions in a simple, consistent stream-based manner.
Definition: error.H:71
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
Finite volume models.
Definition: fvModels.H:65
Model for heat transfer between two phases. Two heat transfer coefficients are used to calculate the ...
Class which provides interfacial momentum transfer between a number of phases. Drag,...
static const word propertiesName
Default name of the phase properties dictionary.
heatTransferSystem(const phaseSystem &)
Construct from a phase system.
virtual ~heatTransferSystem()
Destructor.
autoPtr< HashPtrTable< fvScalarMatrix > > heatTransfer() const
Return the heat transfer matrices.
virtual bool read()
Read base phaseProperties dictionary.
Class to represent an interface between phases. Derivations can further specify the configuration of ...
const phaseModel & phase1() const
Return phase 1.
const phaseModel & phase2() const
Return phase 2.
virtual const rhoThermo & thermo() const =0
Return the thermophysical model.
const word & name() const
Return the name of this phase.
Definition: phaseModel.C:145
Class to represent a system of phases.
Definition: phaseSystem.H:74
const phaseModelPartialList & thermalPhases() const
Return the models for phases that have variable temperature.
Definition: phaseSystemI.H:146
virtual bool read()
Read object.
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:346
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
Calculate the matrix for implicit and explicit sources.
K
Definition: pEqn.H:75
#define WarningInFunction
Report a warning using Foam::Warning.
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const VolField< Type > &)
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
const dimensionSet dimEnergy
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
const dimensionSet dimTime
defineTypeNameAndDebug(combustionModel, 0)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
IOerror FatalIOError
HashSet wordHashSet
A HashSet with word keys.
Definition: HashSet.H:210
error FatalError
thermo he()
Foam::argList args(argc, argv)