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-2024 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 
33 
34 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
35 
37 (
39  scalarField& sumKappaTByDelta,
40  scalarField& sumKappaByDelta,
41  scalarField& Tref,
42  scalarField& Tw,
43  scalarField& sumq,
44  scalarField& qByKappa
45 ) const
46 {
47  const thermophysicalTransportModel& ttm =
48  patch().boundaryMesh().mesh()
49  .lookupType<thermophysicalTransportModel>();
50 
51  kappa = ttm.kappaEff(patch().index());
52 
53  tmp<scalarField> qCorr(ttm.qCorr(patch().index()));
54 
55  if (qCorr.valid())
56  {
57  sumq += qCorr;
58  }
59 
60  qByKappa = sumq/kappa;
61 
62  Tw = *this;
63 }
64 
65 
66 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
67 
70 (
71  const fvPatch& p,
73  const dictionary& dict
74 )
75 :
76  mixedFvPatchScalarField(p, iF, dict, false),
77  haveQ_(dict.found("Q")),
78  Q_
79  (
80  haveQ_
81  ? Function1<scalar>::New("Q", db().time().userUnits(), dimPower, dict)
82  : autoPtr<Function1<scalar>>()
83  ),
84  haveq_(dict.found("q")),
85  q_
86  (
87  haveq_
88  ? Function1<scalar>::New
89  (
90  "q",
91  db().time().userUnits(),
93  dict
94  )
95  : autoPtr<Function1<scalar>>()
96  ),
97  haveh_(dict.found("h")),
98  h_
99  (
100  haveh_
101  ? Function1<scalar>::New
102  (
103  "h",
104  db().time().userUnits(),
106  dict
107  ).ptr()
108  : nullptr
109  ),
110  Ta_
111  (
112  haveh_
113  ? Function1<scalar>::New
114  (
115  "Ta",
116  db().time().userUnits(),
118  dict
119  ).ptr()
120  : nullptr
121  ),
122  emissivity_(dict.lookupOrDefault<scalar>("emissivity", unitFraction, 0)),
123  thicknessLayers_
124  (
125  dict.lookupOrDefault<scalarList>
126  (
127  "thicknessLayers",
128  dimLength,
129  scalarList()
130  )
131  ),
132  kappaLayers_
133  (
134  dict.lookupOrDefault<scalarList>
135  (
136  "kappaLayers",
138  scalarList()
139  )
140  ),
141  relax_(dict.lookupOrDefault<scalar>("relaxation", unitFraction, 1)),
142  qrName_(dict.lookupOrDefault<word>("qr", word::null)),
143  qrRelax_(dict.lookupOrDefault<scalar>("qrRelaxation", unitFraction, 1)),
144  qrPrevious_
145  (
146  qrName_ != word::null
147  ? dict.found("qrPrevious")
148  ? scalarField("qrPrevious", dimPower/dimArea, dict, p.size())
149  : scalarField(0, p.size())
150  : scalarField()
151  )
152 {
154  (
155  scalarField("value", iF.dimensions(), dict, p.size())
156  );
157 
158  if (!haveQ_ && !haveq_ && !haveh_)
159  {
161  << "One or more of Q (heat power), q (heat flux), and h (heat "
162  << "transfer coefficient) must be specified"
163  << exit(FatalIOError);
164  }
165 
166  if (thicknessLayers_.size() != kappaLayers_.size())
167  {
169  << "If either thicknessLayers or kappaLayers is specified, then "
170  << "both must be specified and be lists of the same length "
171  << exit(FatalIOError);
172  }
173 
174  if (dict.found("refValue"))
175  {
176  // Full restart
177  refValue() = scalarField("refValue", iF.dimensions(), dict, p.size());
178  refGrad() =
180  (
181  "refGradient",
182  iF.dimensions()/dimLength,
183  dict,
184  p.size()
185  );
186  valueFraction() =
187  scalarField("valueFraction", unitFraction, dict, p.size());
188  }
189  else
190  {
191  // Start from user entered data. Assume fixedValue.
192  refValue() = *this;
193  refGrad() = 0;
194  valueFraction() = 1;
195  }
196 }
197 
198 
201 (
203  const fvPatch& p,
205  const fieldMapper& mapper
206 )
207 :
208  mixedFvPatchScalarField(ptf, p, iF, mapper),
209  haveQ_(ptf.haveQ_),
210  Q_(ptf.Q_, false),
211  haveq_(ptf.haveq_),
212  q_(ptf.q_, false),
213  haveh_(ptf.haveh_),
214  h_(ptf.h_, false),
215  Ta_(ptf.Ta_, false),
216  emissivity_(ptf.emissivity_),
217  thicknessLayers_(ptf.thicknessLayers_),
218  kappaLayers_(ptf.kappaLayers_),
219  relax_(ptf.relax_),
220  qrName_(ptf.qrName_),
221  qrRelax_(ptf.qrRelax_),
222  qrPrevious_
223  (
224  qrName_ != word::null
225  ? mapper(ptf.qrPrevious_)()
226  : scalarField()
227  )
228 {}
229 
230 
233 (
236 )
237 :
238  mixedFvPatchScalarField(tppsf, iF),
239  haveQ_(tppsf.haveQ_),
240  Q_(tppsf.Q_, false),
241  haveq_(tppsf.haveq_),
242  q_(tppsf.q_, false),
243  haveh_(tppsf.haveh_),
244  h_(tppsf.h_, false),
245  Ta_(tppsf.Ta_, false),
246  emissivity_(tppsf.emissivity_),
247  thicknessLayers_(tppsf.thicknessLayers_),
248  kappaLayers_(tppsf.kappaLayers_),
249  relax_(tppsf.relax_),
250  qrName_(tppsf.qrName_),
251  qrRelax_(tppsf.qrRelax_),
252  qrPrevious_(tppsf.qrPrevious_)
253 {}
254 
255 
256 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
257 
259 (
260  const fvPatchScalarField& ptf,
261  const fieldMapper& mapper
262 )
263 {
264  mixedFvPatchScalarField::map(ptf, mapper);
265 
267  refCast<const externalTemperatureFvPatchScalarField>(ptf);
268 
269  if (qrName_ != word::null)
270  {
271  mapper(qrPrevious_, tiptf.qrPrevious_);
272  }
273 }
274 
275 
277 (
278  const fvPatchScalarField& ptf
279 )
280 {
281  mixedFvPatchScalarField::reset(ptf);
282 
284  refCast<const externalTemperatureFvPatchScalarField>(ptf);
285 
286  if (qrName_ != word::null)
287  {
288  qrPrevious_.reset(tiptf.qrPrevious_);
289  }
290 }
291 
292 
294 {
295  if (updated())
296  {
297  return;
298  }
299 
300  // Get the radiative heat flux and relax
301  scalarField qr(size(), 0);
302  if (qrName_ != word::null)
303  {
304  qr =
305  qrRelax_
306  *patch().lookupPatchField<volScalarField, scalar>(qrName_)
307  + (1 - qrRelax_)*qrPrevious_;
308  qrPrevious_ = qr;
309  }
310 
311  // Compute the total non-convective heat flux
312  scalarField sumq(qr);
313  if (haveQ_)
314  {
315  sumq += Q_->value(db().time().value())/gSum(patch().magSf());
316  }
317  if (haveq_)
318  {
319  sumq += q_->value(db().time().value());
320  }
321 
322  scalarField kappa(size(), 0);
323  scalarField sumKappaTByDelta(size(), 0);
324  scalarField sumKappaByDelta(size(), 0);
325  scalarField Tref(*this);
326  scalarField Tw(*this);
327  scalarField qByKappa(size(), 0);
328  getKappa
329  (
330  kappa,
331  sumKappaTByDelta,
332  sumKappaByDelta,
333  Tref,
334  Tw,
335  sumq,
336  qByKappa
337  );
338 
339  // Add optional external convective heat transfer contribution
340  if (haveh_)
341  {
342  scalar totalSolidRes = 0;
343  if (thicknessLayers_.size())
344  {
345  forAll(thicknessLayers_, iLayer)
346  {
347  const scalar l = thicknessLayers_[iLayer];
348  if (kappaLayers_[iLayer] > 0)
349  {
350  totalSolidRes += l/kappaLayers_[iLayer];
351  }
352  }
353  }
354 
355  const scalar h = h_->value(this->db().time().value());
356  const scalar Ta = Ta_->value(this->db().time().value());
357 
358  const scalarField hp
359  (
360  1
361  /(
362  1
363  /(
364  (emissivity_ > 0)
365  ? (
366  h
367  + emissivity_*sigma.value()
368  *((pow3(Ta) + pow3(Tw)) + Ta*Tw*(Ta + Tw))
369  )()
370  : scalarField(size(), h)
371  ) + totalSolidRes
372  )
373  );
374 
375  sumKappaByDelta += hp;
376  sumKappaTByDelta += hp*Ta;
377 
378  refValue() = sumKappaTByDelta/sumKappaByDelta;
379  }
380  else
381  {
382  refValue() = Tref;
383  }
384 
385  valueFraction() =
386  sumKappaByDelta/(kappa*patch().deltaCoeffs() + sumKappaByDelta);
387 
388  refGrad() = qByKappa;
389 
390  // Modify mixed parameters for under-relaxation
391  if (relax_ != 1)
392  {
393  const scalarField f(valueFraction());
394  valueFraction() = 1 - relax_*(1 - f);
395  refValue() = (f*relax_*refValue() + (1 - relax_)*Tw)/valueFraction();
396  // refGrad() = No change
397  }
398 
399  mixedFvPatchScalarField::updateCoeffs();
400 
401  if (debug)
402  {
403  const scalar Q = gSum(kappa*patch().magSf()*snGrad());
404 
405  Info<< patch().boundaryMesh().mesh().name() << ':'
406  << patch().name() << ':'
407  << this->internalField().name() << " :"
408  << " heat transfer rate:" << Q
409  << " walltemperature "
410  << " min:" << gMin(*this)
411  << " max:" << gMax(*this)
412  << " avg:" << gAverage(*this)
413  << endl;
414  }
415 }
416 
417 
419 (
420  Ostream& os
421 ) const
422 {
424 
425  if (haveQ_)
426  {
427  writeEntry(os, db().time().userUnits(), dimPower, Q_());
428  }
429 
430  if (haveq_)
431  {
432  writeEntry(os, db().time().userUnits(), dimPower/dimArea, q_());
433  }
434 
435  if (haveh_)
436  {
437  writeEntry
438  (
439  os,
440  db().time().userUnits(),
442  h_()
443  );
444  writeEntry(os, db().time().userUnits(), dimTemperature, Ta_());
445  writeEntryIfDifferent(os, "emissivity", scalar(0), emissivity_);
447  (
448  os,
449  "thicknessLayers",
450  scalarList(),
451  thicknessLayers_
452  );
454  (
455  os,
456  "kappaLayers",
457  scalarList(),
458  kappaLayers_
459  );
460  }
461 
462  writeEntryIfDifferent(os, "relaxation", scalar(1), relax_);
463 
464  if (qrName_ != word::null)
465  {
466  writeEntry(os, "qr", qrName_);
467  writeEntry(os, "qrRelaxation", qrRelax_);
468  writeEntry(os, "qrPrevious", qrPrevious_);
469  }
470 
471  writeEntry(os, "refValue", refValue());
472  writeEntry(os, "refGradient", refGrad());
473  writeEntry(os, "valueFraction", valueFraction());
474  writeEntry(os, "value", *this);
475 }
476 
477 
478 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
479 
480 namespace Foam
481 {
483  (
486  );
487 
489  (
492  patchMapper,
493  externalWallHeatFluxTemperature,
494  "externalWallHeatFluxTemperature"
495  );
496 
498  (
501  dictionary,
502  externalWallHeatFluxTemperature,
503  "externalWallHeatFluxTemperature"
504  );
505 }
506 
507 
508 // ************************************************************************* //
bool found
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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 keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
const Type & value() const
Return const reference to value.
This boundary condition applies a heat flux condition to temperature on an external wall....
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.
virtual void getKappa(scalarField &kappa, scalarField &sumKappaTByDelta, scalarField &sumKappaByDelta, scalarField &Tref, scalarField &Tw, scalarField &sumq, scalarField &qByKappa) const
Get the patch kappa, kappa*Tc/delta, kappa/delta,.
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:88
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:229
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:167
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 dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m^2/K^4].
const dimensionedScalar h
Planck constant.
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
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:257
dimensionedScalar pow3(const dimensionedScalar &ds)
const dimensionSet dimPower
messageStream Info
const dimensionSet dimLength
const dimensionSet dimTemperature
List< scalar > scalarList
A List of scalars.
Definition: scalarList.H:50
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:64
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
IOerror FatalIOError
Type gAverage(const FieldField< Field, Type > &f)
const dimensionSet dimArea
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Type gMin(const FieldField< Field, Type > &f)
const dimensionSet dimThermalConductivity
Type gMax(const FieldField< Field, Type > &f)
const unitConversion unitFraction
labelList f(nPoints)
dictionary dict
volScalarField & p