externalWallHeatFluxTemperatureFvPatchScalarField.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-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 "volFields.H"
30 
32 
33 
34 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35 
38 (
39  const fvPatch& p,
41 )
42 :
43  mixedFvPatchScalarField(p, iF),
44  temperatureCoupledBase(patch()),
45  haveQ_(false),
46  Q_(NaN),
47  haveq_(false),
48  q_(),
49  haveh_(false),
50  h_(),
51  Ta_(),
52  emissivity_(0),
53  thicknessLayers_(),
54  kappaLayers_(),
55  relaxation_(1),
56  qrName_(word::null),
57  qrRelaxation_(1),
58  qrPrevious_()
59 {
60  refValue() = 0;
61  refGrad() = 0;
62  valueFraction() = 1;
63 }
64 
65 
68 (
69  const fvPatch& p,
71  const dictionary& dict
72 )
73 :
74  mixedFvPatchScalarField(p, iF),
75  temperatureCoupledBase(patch(), dict),
76  haveQ_(dict.found("Q")),
77  Q_(haveQ_ ? dict.lookup<scalar>("Q") : NaN),
78  haveq_(dict.found("q")),
79  q_(haveq_ ? scalarField("q", dict, p.size()) : scalarField()),
80  haveh_(dict.found("h")),
81  h_(haveh_ ? scalarField("h", dict, p.size()) : scalarField()),
82  Ta_(haveh_ ? Function1<scalar>::New("Ta", dict).ptr() : nullptr),
83  emissivity_(dict.lookupOrDefault<scalar>("emissivity", 0)),
84  thicknessLayers_
85  (
86  dict.lookupOrDefault<scalarList>("thicknessLayers", scalarList())
87  ),
88  kappaLayers_
89  (
90  dict.lookupOrDefault<scalarList>("kappaLayers", scalarList())
91  ),
92  relaxation_(dict.lookupOrDefault<scalar>("relaxation", 1)),
93  qrName_(dict.lookupOrDefault<word>("qr", word::null)),
94  qrRelaxation_(dict.lookupOrDefault<scalar>("qrRelaxation", 1)),
95  qrPrevious_
96  (
97  qrName_ != word::null
98  ? dict.found("qrPrevious")
99  ? scalarField("qrPrevious", dict, p.size())
100  : scalarField(0, p.size())
101  : scalarField()
102  )
103 {
104  fvPatchScalarField::operator=(scalarField("value", dict, p.size()));
105 
106  if (!haveQ_ && !haveq_ && !haveh_)
107  {
109  << "One or more of Q (heat power), q (heat flux), and h (heat "
110  << "transfer coefficient) must be specified"
111  << exit(FatalIOError);
112  }
113 
114  if (thicknessLayers_.size() != kappaLayers_.size())
115  {
117  << "If either thicknessLayers or kappaLayers is specified, then "
118  << "both must be specified and be lists of the same length "
119  << exit(FatalIOError);
120  }
121 
122  if (dict.found("refValue"))
123  {
124  // Full restart
125  refValue() = scalarField("refValue", dict, p.size());
126  refGrad() = scalarField("refGradient", dict, p.size());
127  valueFraction() = scalarField("valueFraction", dict, p.size());
128  }
129  else
130  {
131  // Start from user entered data. Assume fixedValue.
132  refValue() = *this;
133  refGrad() = 0;
134  valueFraction() = 1;
135  }
136 }
137 
138 
141 (
143  const fvPatch& p,
145  const fvPatchFieldMapper& mapper
146 )
147 :
148  mixedFvPatchScalarField(ptf, p, iF, mapper),
149  temperatureCoupledBase(patch(), ptf),
150  haveQ_(ptf.haveQ_),
151  Q_(ptf.Q_),
152  haveq_(ptf.haveq_),
153  q_(haveq_ ? mapper(ptf.q_)() : scalarField()),
154  haveh_(ptf.haveh_),
155  h_(haveh_ ? mapper(ptf.h_)() : scalarField()),
156  Ta_(ptf.Ta_, false),
157  emissivity_(ptf.emissivity_),
158  thicknessLayers_(ptf.thicknessLayers_),
159  kappaLayers_(ptf.kappaLayers_),
160  relaxation_(ptf.relaxation_),
161  qrName_(ptf.qrName_),
162  qrRelaxation_(ptf.qrRelaxation_),
163  qrPrevious_
164  (
165  qrName_ != word::null
166  ? mapper(ptf.qrPrevious_)()
167  : scalarField()
168  )
169 {}
170 
171 
174 (
177 )
178 :
179  mixedFvPatchScalarField(tppsf, iF),
180  temperatureCoupledBase(patch(), tppsf),
181  haveQ_(tppsf.haveQ_),
182  Q_(tppsf.Q_),
183  haveq_(tppsf.haveq_),
184  q_(tppsf.q_),
185  haveh_(tppsf.haveh_),
186  h_(tppsf.h_),
187  Ta_(tppsf.Ta_, false),
188  emissivity_(tppsf.emissivity_),
189  thicknessLayers_(tppsf.thicknessLayers_),
190  kappaLayers_(tppsf.kappaLayers_),
191  relaxation_(tppsf.relaxation_),
192  qrName_(tppsf.qrName_),
193  qrRelaxation_(tppsf.qrRelaxation_),
194  qrPrevious_(tppsf.qrPrevious_)
195 {}
196 
197 
198 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
199 
201 (
202  const fvPatchFieldMapper& m
203 )
204 {
205  mixedFvPatchScalarField::autoMap(m);
206 
207  if (haveq_)
208  {
209  m(q_, q_);
210  }
211 
212  if (haveh_)
213  {
214  m(h_, h_);
215  }
216 
217  if (qrName_ != word::null)
218  {
219  m(qrPrevious_, qrPrevious_);
220  }
221 }
222 
223 
225 (
226  const fvPatchScalarField& ptf,
227  const labelList& addr
228 )
229 {
230  mixedFvPatchScalarField::rmap(ptf, addr);
231 
233  refCast<const externalWallHeatFluxTemperatureFvPatchScalarField>(ptf);
234 
235  if (haveq_)
236  {
237  q_.rmap(tiptf.q_, addr);
238  }
239 
240  if (haveh_)
241  {
242  h_.rmap(tiptf.h_, addr);
243  }
244 
245  if (qrName_ != word::null)
246  {
247  qrPrevious_.rmap(tiptf.qrPrevious_, addr);
248  }
249 }
250 
251 
253 (
254  const fvPatchScalarField& ptf
255 )
256 {
257  mixedFvPatchScalarField::reset(ptf);
258 
260  refCast<const externalWallHeatFluxTemperatureFvPatchScalarField>(ptf);
261 
262  if (haveq_)
263  {
264  q_.reset(tiptf.q_);
265  }
266 
267  if (haveh_)
268  {
269  h_.reset(tiptf.h_);
270  }
271 
272  if (qrName_ != word::null)
273  {
274  qrPrevious_.reset(tiptf.qrPrevious_);
275  }
276 }
277 
278 
280 {
281  if (updated())
282  {
283  return;
284  }
285 
286  const scalarField& Tp(*this);
287 
288  // Store current valueFraction and refValue for relaxation
289  const scalarField valueFraction0(valueFraction());
290  const scalarField refValue0(refValue());
291 
292  // Get the radiative heat flux and relax
293  scalarField qr(Tp.size(), 0);
294  if (qrName_ != word::null)
295  {
296  qr =
297  qrRelaxation_
298  *patch().lookupPatchField<volScalarField, scalar>(qrName_)
299  + (1 - qrRelaxation_)*qrPrevious_;
300  qrPrevious_ = qr;
301  }
302 
303  // Compute the total non-convective heat flux
304  scalarField qTot(qr);
305  if (haveQ_)
306  {
307  qTot += Q_/gSum(patch().magSf());
308  }
309  if (haveq_)
310  {
311  qTot += q_;
312  }
313 
314  // Evaluate
315  if (!haveh_)
316  {
317  refGrad() = qTot/kappa(*this);
318  refValue() = Tp;
319  valueFraction() = 0;
320  }
321  else
322  {
323  scalar totalSolidRes = 0;
324  if (thicknessLayers_.size())
325  {
326  forAll(thicknessLayers_, iLayer)
327  {
328  const scalar l = thicknessLayers_[iLayer];
329  if (kappaLayers_[iLayer] > 0)
330  {
331  totalSolidRes += l/kappaLayers_[iLayer];
332  }
333  }
334  }
335 
336  const scalar Ta = Ta_->value(this->db().time().userTimeValue());
337 
338  const scalarField hp
339  (
340  1
341  /(
342  1
343  /(
344  (emissivity_ > 0)
345  ? (
346  h_
347  + emissivity_*sigma.value()
348  *((pow3(Ta) + pow3(Tp)) + Ta*Tp*(Ta + Tp))
349  )()
350  : h_
351  ) + totalSolidRes
352  )
353  );
354 
355  const scalarField hpTa(hp*Ta);
356 
357  const scalarField kappaDeltaCoeffs
358  (
359  this->kappa(*this)*patch().deltaCoeffs()
360  );
361 
362  refGrad() = 0;
363  forAll(Tp, i)
364  {
365  if (qTot[i] < 0)
366  {
367  const scalar hpmqTot = hp[i] - qTot[i]/Tp[i];
368  refValue()[i] = hpTa[i]/hpmqTot;
369  valueFraction()[i] = hpmqTot/(hpmqTot + kappaDeltaCoeffs[i]);
370  }
371  else
372  {
373  refValue()[i] = (hpTa[i] + qTot[i])/hp[i];
374  valueFraction()[i] = hp[i]/(hp[i] + kappaDeltaCoeffs[i]);
375  }
376  }
377  }
378 
379  // Relax
380  valueFraction() =
381  relaxation_*valueFraction() + (1 - relaxation_)*valueFraction0;
382  refValue() =
383  relaxation_*refValue() + (1 - relaxation_)*refValue0;
384 
385  mixedFvPatchScalarField::updateCoeffs();
386 
387  if (debug)
388  {
389  const scalar Q = gSum(kappa(*this)*patch().magSf()*snGrad());
390 
391  Info<< patch().boundaryMesh().mesh().name() << ':'
392  << patch().name() << ':'
393  << this->internalField().name() << " :"
394  << " heat transfer rate:" << Q
395  << " walltemperature "
396  << " min:" << gMin(*this)
397  << " max:" << gMax(*this)
398  << " avg:" << gAverage(*this)
399  << endl;
400  }
401 }
402 
403 
405 (
406  Ostream& os
407 ) const
408 {
410 
412 
413  if (haveQ_)
414  {
415  writeEntry(os, "Q", Q_);
416  }
417 
418  if (haveq_)
419  {
420  writeEntry(os, "q", q_);
421  }
422 
423  if (haveh_)
424  {
425  writeEntry(os, "h", h_);
426  writeEntry(os, Ta_());
427  writeEntryIfDifferent(os, "emissivity", scalar(0), emissivity_);
429  (
430  os,
431  "thicknessLayers",
432  scalarList(),
433  thicknessLayers_
434  );
436  (
437  os,
438  "kappaLayers",
439  scalarList(),
440  kappaLayers_
441  );
442  }
443 
444  writeEntryIfDifferent(os, "relaxation", scalar(1), relaxation_);
445 
446  if (qrName_ != word::null)
447  {
448  writeEntry(os, "qr", qrName_);
449  writeEntry(os, "qrRelaxation", qrRelaxation_);
450  writeEntry(os, "qrPrevious", qrPrevious_);
451  }
452 
453  writeEntry(os, "refValue", refValue());
454  writeEntry(os, "refGradient", refGrad());
455  writeEntry(os, "valueFraction", valueFraction());
456  writeEntry(os, "value", *this);
457 }
458 
459 
460 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
461 
462 namespace Foam
463 {
465  (
468  );
469 }
470 
471 // ************************************************************************* //
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Run-time selectable general function of one variable.
Definition: Function1.H:52
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:663
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
Type gMin(const FieldField< Field, Type > &f)
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
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:251
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:243
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m^2/K^4].
Type gSum(const FieldField< Field, Type > &f)
A class for handling words, derived from string.
Definition: word.H:59
Foam::fvPatchFieldMapper.
List< scalar > scalarList
A List of scalars.
Definition: scalarList.H:50
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
const Type & value() const
Return const reference to value.
static const word null
An empty word.
Definition: word.H:77
virtual label size() const
Return size.
Definition: fvPatch.H:157
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
volScalarField scalarField(fieldObject, mesh)
Type gMax(const FieldField< Field, Type > &f)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
dimensionedScalar pow3(const dimensionedScalar &ds)
Common functions used in temperature coupled boundaries.
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:318
tmp< scalarField > kappa(const fvPatchScalarField &Tp) const
Given patch temperature calculate corresponding K field.
Type gAverage(const FieldField< Field, Type > &f)
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
messageStream Info
externalWallHeatFluxTemperatureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
This boundary condition applies a heat flux condition to temperature on an external wall...
void write(Ostream &) const
Write.
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:45
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:864
virtual void reset(const fvPatchScalarField &)
Reset the fvPatchField to the given fvPatchField.
IOerror FatalIOError