externalTemperatureFvPatchScalarField.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) 2011-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 
28 #include "volFields.H"
31 
32 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
33 
35 (
37  const scalar d
38 ) const
39 {
40  if (!tf.valid())
41  {
42  tf = new scalarField(size(), d);
43  }
44  else
45  {
46  tf.ref() += d;
47  }
48 }
49 
50 
52 (
54  const tmp<scalarField>& tdf
55 ) const
56 {
57  if (!tdf.valid())
58  {
59  return;
60  }
61 
62  if (!tf.valid())
63  {
64  tf = tdf.ptr();
65  }
66  else
67  {
68  tf.ref() += tdf;
69  }
70 }
71 
72 
74 (
76  tmp<scalarField>& sumKappaTcByDelta,
77  tmp<scalarField>& sumKappaByDelta,
79  tmp<scalarField>& sumq
80 ) const
81 {
82  const thermophysicalTransportModel& ttm =
83  patch().boundaryMesh().mesh()
84  .lookupType<thermophysicalTransportModel>();
85 
86  kappa = ttm.kappaEff(patch().index());
87 
88  T = tmp<scalarField>(*this);
89 
90  plusEqOp(sumq, ttm.qCorr(patch().index()));
91 }
92 
93 
94 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
95 
98 (
99  const fvPatch& p,
101  const dictionary& dict
102 )
103 :
104  mixedFvPatchScalarField(p, iF, dict, false),
105  haveQ_(dict.found("Q")),
106  Q_
107  (
108  haveQ_
109  ? Function1<scalar>::New("Q", db().time().userUnits(), dimPower, dict)
110  : autoPtr<Function1<scalar>>()
111  ),
112  haveq_(dict.found("q")),
113  q_
114  (
115  haveq_
116  ? Function1<scalar>::New
117  (
118  "q",
119  db().time().userUnits(),
121  dict
122  )
123  : autoPtr<Function1<scalar>>()
124  ),
125  haveh_(dict.found("h")),
126  h_
127  (
128  haveh_
129  ? Function1<scalar>::New
130  (
131  "h",
132  db().time().userUnits(),
134  dict
135  ).ptr()
136  : nullptr
137  ),
138  haveEmissivity_(dict.found("emissivity")),
139  emissivity_
140  (
141  haveEmissivity_
142  ? dict.lookup<scalar>("emissivity", unitFraction)
143  : NaN
144  ),
145  haveLayers_(dict.found("thicknessLayers") || dict.found("kappaLayers")),
146  thicknessLayers_
147  (
148  haveLayers_
149  ? dict.lookup<scalarList>("thicknessLayers", dimLength)
150  : scalarList()
151  ),
152  kappaLayers_
153  (
154  haveLayers_
155  ? dict.lookup<scalarList>("kappaLayers", dimThermalConductivity)
156  : scalarList()
157  ),
158  Ta_
159  (
160  haveh_ || haveEmissivity_
161  ? Function1<scalar>::New
162  (
163  "Ta",
164  db().time().userUnits(),
166  dict
167  ).ptr()
168  : nullptr
169  ),
170  relax_(dict.lookupOrDefault<scalar>("relaxation", unitFraction, 1)),
171  qrName_(dict.lookupOrDefault<word>("qr", word::null)),
172  qrRelax_(dict.lookupOrDefault<scalar>("qrRelaxation", unitFraction, 1)),
173  qrPrevious_
174  (
175  qrName_ != word::null
176  ? dict.found("qrPrevious")
177  ? scalarField("qrPrevious", dimPower/dimArea, dict, p.size())
178  : scalarField(p.size(), 0)
179  : scalarField()
180  )
181 {
183  (
184  scalarField("value", iF.dimensions(), dict, p.size())
185  );
186 
187  if (haveEmissivity_ && (emissivity_ < 0 || emissivity_ > 1))
188  {
190  << "Emissivity must be in the range 0 to 1"
191  << exit(FatalIOError);
192  }
193 
194  if (thicknessLayers_.size() != kappaLayers_.size())
195  {
197  << "If either thicknessLayers or kappaLayers is specified, then "
198  << "both must be specified and be lists of the same length "
199  << exit(FatalIOError);
200  }
201 
202  if (haveEmissivity_ && haveLayers_)
203  {
205  << "Emissivity and thicknessLayers/kappaLayers are incompatible"
206  << exit(FatalIOError);
207  }
208 
209  if (dict.found("refValue"))
210  {
211  // Full restart
212  refValue() = scalarField("refValue", iF.dimensions(), dict, p.size());
213  refGrad() =
215  (
216  "refGradient",
217  iF.dimensions()/dimLength,
218  dict,
219  p.size()
220  );
221  valueFraction() =
222  scalarField("valueFraction", unitFraction, dict, p.size());
223  }
224  else
225  {
226  // Start from user entered data. Assume fixedValue.
227  refValue() = *this;
228  refGrad() = 0;
229  valueFraction() = 1;
230  }
231 }
232 
233 
236 (
238  const fvPatch& p,
240  const fieldMapper& mapper
241 )
242 :
243  mixedFvPatchScalarField(ptf, p, iF, mapper),
244  haveQ_(ptf.haveQ_),
245  Q_(ptf.Q_, false),
246  haveq_(ptf.haveq_),
247  q_(ptf.q_, false),
248  haveh_(ptf.haveh_),
249  h_(ptf.h_, false),
250  haveEmissivity_(ptf.haveEmissivity_),
251  emissivity_(ptf.emissivity_),
252  haveLayers_(ptf.haveLayers_),
253  thicknessLayers_(ptf.thicknessLayers_),
254  kappaLayers_(ptf.kappaLayers_),
255  Ta_(ptf.Ta_, false),
256  relax_(ptf.relax_),
257  qrName_(ptf.qrName_),
258  qrRelax_(ptf.qrRelax_),
259  qrPrevious_
260  (
261  qrName_ != word::null
262  ? mapper(ptf.qrPrevious_)()
263  : scalarField()
264  )
265 {}
266 
267 
270 (
273 )
274 :
275  mixedFvPatchScalarField(tppsf, iF),
276  haveQ_(tppsf.haveQ_),
277  Q_(tppsf.Q_, false),
278  haveq_(tppsf.haveq_),
279  q_(tppsf.q_, false),
280  haveh_(tppsf.haveh_),
281  h_(tppsf.h_, false),
282  haveEmissivity_(tppsf.haveEmissivity_),
283  emissivity_(tppsf.emissivity_),
284  haveLayers_(tppsf.haveLayers_),
285  thicknessLayers_(tppsf.thicknessLayers_),
286  kappaLayers_(tppsf.kappaLayers_),
287  Ta_(tppsf.Ta_, false),
288  relax_(tppsf.relax_),
289  qrName_(tppsf.qrName_),
290  qrRelax_(tppsf.qrRelax_),
291  qrPrevious_(tppsf.qrPrevious_)
292 {}
293 
294 
295 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
296 
298 (
299  const fvPatchScalarField& ptf,
300  const fieldMapper& mapper
301 )
302 {
303  mixedFvPatchScalarField::map(ptf, mapper);
304 
306  refCast<const externalTemperatureFvPatchScalarField>(ptf);
307 
308  if (qrName_ != word::null)
309  {
310  mapper(qrPrevious_, tiptf.qrPrevious_);
311  }
312 }
313 
314 
316 (
317  const fvPatchScalarField& ptf
318 )
319 {
320  mixedFvPatchScalarField::reset(ptf);
321 
323  refCast<const externalTemperatureFvPatchScalarField>(ptf);
324 
325  if (qrName_ != word::null)
326  {
327  qrPrevious_.reset(tiptf.qrPrevious_);
328  }
329 }
330 
331 
332 
334 {
335  if (updated())
336  {
337  return;
338  }
339 
340  const scalar t = db().time().value();
341 
342  // Get thermal conductivities, the patch temperature and any explicit
343  // correction heat flux
345  tmp<scalarField> sumKappaTcByDelta, sumKappaByDelta;
347  tmp<scalarField> sumq;
348  getKappa(kappa, sumKappaTcByDelta, sumKappaByDelta, T, sumq);
349 
350  // Add any user specified heat fluxes
351  if (haveQ_)
352  {
353  plusEqOp(sumq, Q_->value(t)/gSum(patch().magSf()));
354  }
355  if (haveq_)
356  {
357  plusEqOp(sumq, q_->value(t));
358  }
359 
360  // Add the (relaxed) radiative heat flux
361  if (qrName_ != word::null)
362  {
363  const fvPatchScalarField& qrCurrent =
364  patch().lookupPatchField<volScalarField, scalar>(qrName_);
365 
366  const scalarField qr(qrRelax_*qrCurrent + (1 - qrRelax_)*qrPrevious_);
367 
368  qrPrevious_ = qr;
369 
370  plusEqOp(sumq, qr);
371  }
372 
373  // Evaluate the ambient temperature
374  const scalar Ta = haveh_ || haveEmissivity_ ? Ta_->value(t) : NaN;
375 
376  // Evaluate the combined convective and radiative heat transfer coefficient
377  tmp<scalarField> hEff;
378  if (haveh_)
379  {
380  plusEqOp(hEff, h_->value(t));
381  }
382  if (haveEmissivity_)
383  {
384  plusEqOp
385  (
386  hEff,
387  emissivity_
389  *(sqr(Ta) + sqr(T()))
390  *(Ta + T())
391  );
392  }
393 
394  // Determine the (reciprocal of the) heat transfer coefficient for the
395  // layer resistances and combine with the convective and radiative heat
396  // transfer coefficients to create a complete effective coefficient
397  if (hEff.valid() && haveLayers_)
398  {
399  scalar oneByHLayers = 0;
400  if (haveLayers_)
401  {
402  forAll(thicknessLayers_, layeri)
403  {
404  oneByHLayers += thicknessLayers_[layeri]/kappaLayers_[layeri];
405  }
406  }
407 
408  hEff = 1/(1/hEff + oneByHLayers);
409  }
410 
411  // If we have a heat transfer coefficient then add it to the kappa sums
412  if (hEff.valid())
413  {
414  plusEqOp(sumKappaByDelta, hEff());
415  plusEqOp(sumKappaTcByDelta, hEff*Ta);
416  }
417 
418  // Set the mixed parameters
419  const scalarField kappaByDelta(kappa*patch().deltaCoeffs());
420  tmp<scalarField> kappaPlusSumKappaByDelta
421  (
422  sumKappaByDelta.valid()
423  ? kappaByDelta + sumKappaByDelta()
424  : tmp<scalarField>(kappaByDelta)
425  );
426 
427  // ... value fraction
428  if (sumKappaByDelta.valid())
429  {
430  valueFraction() = sumKappaByDelta()/kappaPlusSumKappaByDelta();
431  }
432  else
433  {
434  valueFraction() = Zero;
435  }
436 
437  // ... reference value
438  tmp<scalarField> trefValue;
439  if (sumKappaByDelta.valid())
440  {
441  plusEqOp
442  (
443  trefValue,
444  max(sumKappaTcByDelta, small*kappaByDelta*patchInternalField())
445  /max(sumKappaByDelta, small*kappaByDelta)
446  );
447  }
448  if (sumq.valid())
449  {
450  plusEqOp(trefValue, sumq()/kappaPlusSumKappaByDelta());
451  }
452  if (trefValue.valid())
453  {
454  refValue() = trefValue;
455  }
456  else
457  {
458  refValue() = Zero;
459  }
460 
461  // ... and reference gradient
462  if (sumq.valid())
463  {
464  refGrad() = sumq*patch().deltaCoeffs()/kappaPlusSumKappaByDelta();
465  }
466  else
467  {
468  refGrad() = Zero;
469  }
470 
471  // Modify the mixed parameters for under-relaxation
472  if (relax_ != 1)
473  {
474  const scalarField f(valueFraction());
475  valueFraction() = 1 - relax_*(1 - f);
476  refValue() = (f*relax_*refValue() + (1 - relax_)*T)/valueFraction();
477  }
478 
479  mixedFvPatchScalarField::updateCoeffs();
480 
481  if (debug)
482  {
483  const scalar Q = gSum(kappa*patch().magSf()*snGrad());
484 
485  Info<< patch().boundaryMesh().mesh().name() << ':'
486  << patch().name() << ':'
487  << this->internalField().name() << " :"
488  << " heat transfer rate:" << Q
489  << " walltemperature "
490  << " min:" << gMin(*this)
491  << " max:" << gMax(*this)
492  << " avg:" << gAverage(*this)
493  << endl;
494  }
495 }
496 
497 
499 (
500  Ostream& os
501 ) const
502 {
504 
505  if (haveQ_)
506  {
507  writeEntry(os, db().time().userUnits(), dimPower, Q_());
508  }
509 
510  if (haveq_)
511  {
512  writeEntry(os, db().time().userUnits(), dimPower/dimArea, q_());
513  }
514 
515  if (haveh_)
516  {
517  writeEntry
518  (
519  os,
520  db().time().userUnits(),
522  h_()
523  );
524  }
525 
526  if (haveEmissivity_)
527  {
528  writeEntry(os, "emissivity", emissivity_);
529  }
530 
531  if (haveLayers_)
532  {
533  writeEntry(os, "thicknessLayers", thicknessLayers_);
534  writeEntry(os, "kappaLayers", kappaLayers_);
535  }
536 
537  if (haveh_ || haveEmissivity_)
538  {
539  writeEntry(os, db().time().userUnits(), dimTemperature, Ta_());
540  }
541 
542  writeEntryIfDifferent(os, "relaxation", scalar(1), relax_);
543 
544  if (qrName_ != word::null)
545  {
546  writeEntry(os, "qr", qrName_);
547  writeEntry(os, "qrRelaxation", qrRelax_);
548  writeEntry(os, "qrPrevious", qrPrevious_);
549  }
550 
551  writeEntry(os, "refValue", refValue());
552  writeEntry(os, "refGradient", refGrad());
553  writeEntry(os, "valueFraction", valueFraction());
554  writeEntry(os, "value", *this);
555 }
556 
557 
558 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
559 
560 namespace Foam
561 {
563  (
566  );
567 
569  (
572  patchMapper,
573  externalWallHeatFluxTemperature,
574  "externalWallHeatFluxTemperature"
575  );
576 
578  (
581  dictionary,
582  externalWallHeatFluxTemperature,
583  "externalWallHeatFluxTemperature"
584  );
585 }
586 
587 
588 // ************************************************************************* //
bool found
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const dimensionSet & dimensions() const
Return dimensions.
Run-time selectable general function of one variable.
Definition: Function1.H:125
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
This boundary condition applies a heat flux condition to temperature on an external wall....
virtual void getKappa(scalarField &kappa, tmp< scalarField > &sumKappaTcByDelta, tmp< scalarField > &sumKappaByDelta, tmp< scalarField > &T, tmp< scalarField > &sumq) const
Get the patch kappa, kappa*Tc/delta, kappa/delta,.
void plusEqOp(tmp< scalarField > &tf, const scalar d) const
Plus-equals op for a tmp field. Will initialise the tmp if empty.
virtual void reset(const fvPatchScalarField &)
Reset the fvPatchField to the given fvPatchField.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void map(const fvPatchScalarField &, const fieldMapper &)
Map the given fvPatchField onto this fvPatchField.
externalTemperatureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
Abstract base class for field mapping.
Definition: fieldMapper.H:48
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:91
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:237
friend Ostream & operator(Ostream &, const fvPatchField< Type > &)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
Abstract base class for all fluid and solid thermophysical transport models.
virtual tmp< scalarField > qCorr(const label patchi) const =0
Return the patch heat flux correction [W/m^2].
virtual tmp< volScalarField > kappaEff() const =0
Effective thermal turbulent conductivity.
A class for managing temporary objects.
Definition: tmp.H:55
bool valid() const
Is this temporary object valid,.
Definition: tmpI.H:183
T * ptr() const
Return tmp pointer for reuse.
Definition: tmpI.H:221
A class for handling words, derived from string.
Definition: word.H:62
static const word null
An empty word.
Definition: word.H:77
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:346
const tensorField & tf
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m^2/K^4].
tmp< SurfaceField< Type > > snGrad(const VolField< Type > &vf, const word &name)
Definition: fvcSnGrad.C:45
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static const zero Zero
Definition: zero.H:97
Type gSum(const FieldField< Field, Type > &f)
void writeEntryIfDifferent(Ostream &os, const word &entryName, const EntryType &value1, const EntryType &value2)
Helper function to write the keyword and entry only if the.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
const dimensionSet dimPower
void T(LagrangianPatchField< Type > &f, const LagrangianPatchField< Type > &f1)
messageStream Info
const dimensionSet dimLength
const dimensionSet dimTemperature
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
addBackwardCompatibleToRunTimeSelectionTable(fvPatchScalarField, coupledTemperatureFvPatchScalarField, patchMapper, turbulentTemperatureCoupledBaffleMixed, "compressible::turbulentTemperatureCoupledBaffleMixed")
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:62
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
IOerror FatalIOError
void sqr(LagrangianPatchField< typename outerProduct< Type, Type >::type > &f, const LagrangianPatchField< Type > &f1)
Type gAverage(const FieldField< Field, Type > &f)
const dimensionSet dimArea
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Type gMin(const FieldField< Field, Type > &f)
const dimensionSet dimThermalConductivity
tmp< DimensionedField< TypeR, GeoMesh, Field > > New(const tmp< DimensionedField< TypeR, GeoMesh, Field >> &tdf1, const word &name, const dimensionSet &dimensions)
Type gMax(const FieldField< Field, Type > &f)
const unitConversion unitFraction
labelList f(nPoints)
dictionary dict
volScalarField & p