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-2023 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 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35 
38 (
39  const fvPatch& p,
41  const dictionary& dict
42 )
43 :
44  mixedFvPatchScalarField(p, iF, dict, false),
45  haveQ_(dict.found("Q")),
46  Q_(haveQ_ ? dict.lookup<scalar>("Q") : NaN),
47  haveq_(dict.found("q")),
48  q_(haveq_ ? scalarField("q", dict, p.size()) : scalarField()),
49  haveh_(dict.found("h")),
50  h_(haveh_ ? scalarField("h", dict, p.size()) : scalarField()),
51  Ta_(haveh_ ? Function1<scalar>::New("Ta", dict).ptr() : nullptr),
52  emissivity_(dict.lookupOrDefault<scalar>("emissivity", 0)),
53  thicknessLayers_
54  (
55  dict.lookupOrDefault<scalarList>("thicknessLayers", scalarList())
56  ),
57  kappaLayers_
58  (
59  dict.lookupOrDefault<scalarList>("kappaLayers", scalarList())
60  ),
61  relaxation_(dict.lookupOrDefault<scalar>("relaxation", 1)),
62  qrName_(dict.lookupOrDefault<word>("qr", word::null)),
63  qrRelaxation_(dict.lookupOrDefault<scalar>("qrRelaxation", 1)),
64  qrPrevious_
65  (
66  qrName_ != word::null
67  ? dict.found("qrPrevious")
68  ? scalarField("qrPrevious", dict, p.size())
69  : scalarField(0, p.size())
70  : scalarField()
71  )
72 {
74 
75  if (!haveQ_ && !haveq_ && !haveh_)
76  {
78  << "One or more of Q (heat power), q (heat flux), and h (heat "
79  << "transfer coefficient) must be specified"
80  << exit(FatalIOError);
81  }
82 
83  if (thicknessLayers_.size() != kappaLayers_.size())
84  {
86  << "If either thicknessLayers or kappaLayers is specified, then "
87  << "both must be specified and be lists of the same length "
88  << exit(FatalIOError);
89  }
90 
91  if (dict.found("refValue"))
92  {
93  // Full restart
94  refValue() = scalarField("refValue", dict, p.size());
95  refGrad() = scalarField("refGradient", dict, p.size());
96  valueFraction() = scalarField("valueFraction", dict, p.size());
97  }
98  else
99  {
100  // Start from user entered data. Assume fixedValue.
101  refValue() = *this;
102  refGrad() = 0;
103  valueFraction() = 1;
104  }
105 }
106 
107 
110 (
112  const fvPatch& p,
114  const fvPatchFieldMapper& mapper
115 )
116 :
117  mixedFvPatchScalarField(ptf, p, iF, mapper),
118  haveQ_(ptf.haveQ_),
119  Q_(ptf.Q_),
120  haveq_(ptf.haveq_),
121  q_(haveq_ ? mapper(ptf.q_)() : scalarField()),
122  haveh_(ptf.haveh_),
123  h_(haveh_ ? mapper(ptf.h_)() : scalarField()),
124  Ta_(ptf.Ta_, false),
125  emissivity_(ptf.emissivity_),
126  thicknessLayers_(ptf.thicknessLayers_),
127  kappaLayers_(ptf.kappaLayers_),
128  relaxation_(ptf.relaxation_),
129  qrName_(ptf.qrName_),
130  qrRelaxation_(ptf.qrRelaxation_),
131  qrPrevious_
132  (
133  qrName_ != word::null
134  ? mapper(ptf.qrPrevious_)()
135  : scalarField()
136  )
137 {}
138 
139 
142 (
145 )
146 :
147  mixedFvPatchScalarField(tppsf, iF),
148  haveQ_(tppsf.haveQ_),
149  Q_(tppsf.Q_),
150  haveq_(tppsf.haveq_),
151  q_(tppsf.q_),
152  haveh_(tppsf.haveh_),
153  h_(tppsf.h_),
154  Ta_(tppsf.Ta_, false),
155  emissivity_(tppsf.emissivity_),
156  thicknessLayers_(tppsf.thicknessLayers_),
157  kappaLayers_(tppsf.kappaLayers_),
158  relaxation_(tppsf.relaxation_),
159  qrName_(tppsf.qrName_),
160  qrRelaxation_(tppsf.qrRelaxation_),
161  qrPrevious_(tppsf.qrPrevious_)
162 {}
163 
164 
165 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
166 
168 (
169  const fvPatchScalarField& ptf,
170  const fvPatchFieldMapper& mapper
171 )
172 {
173  mixedFvPatchScalarField::map(ptf, mapper);
174 
176  refCast<const externalTemperatureFvPatchScalarField>(ptf);
177 
178  if (haveq_)
179  {
180  mapper(q_, tiptf.q_);
181  }
182 
183  if (haveh_)
184  {
185  mapper(h_, tiptf.h_);
186  }
187 
188  if (qrName_ != word::null)
189  {
190  mapper(qrPrevious_, tiptf.qrPrevious_);
191  }
192 }
193 
194 
196 (
197  const fvPatchScalarField& ptf
198 )
199 {
200  mixedFvPatchScalarField::reset(ptf);
201 
203  refCast<const externalTemperatureFvPatchScalarField>(ptf);
204 
205  if (haveq_)
206  {
207  q_.reset(tiptf.q_);
208  }
209 
210  if (haveh_)
211  {
212  h_.reset(tiptf.h_);
213  }
214 
215  if (qrName_ != word::null)
216  {
217  qrPrevious_.reset(tiptf.qrPrevious_);
218  }
219 }
220 
221 
223 {
224  if (updated())
225  {
226  return;
227  }
228 
229  const scalarField& Tp(*this);
230 
231  // Store current valueFraction and refValue for relaxation
232  const scalarField valueFraction0(valueFraction());
233  const scalarField refValue0(refValue());
234 
235  // Get the radiative heat flux and relax
236  scalarField qr(Tp.size(), 0);
237  if (qrName_ != word::null)
238  {
239  qr =
240  qrRelaxation_
241  *patch().lookupPatchField<volScalarField, scalar>(qrName_)
242  + (1 - qrRelaxation_)*qrPrevious_;
243  qrPrevious_ = qr;
244  }
245 
246  // Compute the total non-convective heat flux
247  scalarField qTot(qr);
248  if (haveQ_)
249  {
250  qTot += Q_/gSum(patch().magSf());
251  }
252  if (haveq_)
253  {
254  qTot += q_;
255  }
256 
257  const thermophysicalTransportModel& ttm =
258  patch().boundaryMesh().mesh()
259  .lookupType<thermophysicalTransportModel>();
260 
261  const scalarField kappa(ttm.kappaEff(patch().index()));
262  tmp<scalarField> qCorr(ttm.qCorr(patch().index()));
263 
264  if (qCorr.valid())
265  {
266  qTot += qCorr;
267  }
268 
269  // Evaluate
270  if (!haveh_)
271  {
272  refGrad() = qTot/kappa;
273  refValue() = Tp;
274  valueFraction() = 0;
275  }
276  else
277  {
278  scalar totalSolidRes = 0;
279  if (thicknessLayers_.size())
280  {
281  forAll(thicknessLayers_, iLayer)
282  {
283  const scalar l = thicknessLayers_[iLayer];
284  if (kappaLayers_[iLayer] > 0)
285  {
286  totalSolidRes += l/kappaLayers_[iLayer];
287  }
288  }
289  }
290 
291  const scalar Ta = Ta_->value(this->db().time().userTimeValue());
292 
293  const scalarField hp
294  (
295  1
296  /(
297  1
298  /(
299  (emissivity_ > 0)
300  ? (
301  h_
302  + emissivity_*sigma.value()
303  *((pow3(Ta) + pow3(Tp)) + Ta*Tp*(Ta + Tp))
304  )()
305  : h_
306  ) + totalSolidRes
307  )
308  );
309 
310  const scalarField hpTa(hp*Ta);
311 
312  const scalarField kappaDeltaCoeffs
313  (
314  kappa*patch().deltaCoeffs()
315  );
316 
317  refGrad() = 0;
318  forAll(Tp, i)
319  {
320  if (qTot[i] < 0)
321  {
322  const scalar hpmqTot = hp[i] - qTot[i]/Tp[i];
323  refValue()[i] = hpTa[i]/hpmqTot;
324  valueFraction()[i] = hpmqTot/(hpmqTot + kappaDeltaCoeffs[i]);
325  }
326  else
327  {
328  refValue()[i] = (hpTa[i] + qTot[i])/hp[i];
329  valueFraction()[i] = hp[i]/(hp[i] + kappaDeltaCoeffs[i]);
330  }
331  }
332  }
333 
334  // Relax
335  valueFraction() =
336  relaxation_*valueFraction() + (1 - relaxation_)*valueFraction0;
337  refValue() =
338  relaxation_*refValue() + (1 - relaxation_)*refValue0;
339 
340  mixedFvPatchScalarField::updateCoeffs();
341 
342  if (debug)
343  {
344  const scalar Q = gSum(kappa*patch().magSf()*snGrad());
345 
346  Info<< patch().boundaryMesh().mesh().name() << ':'
347  << patch().name() << ':'
348  << this->internalField().name() << " :"
349  << " heat transfer rate:" << Q
350  << " walltemperature "
351  << " min:" << gMin(*this)
352  << " max:" << gMax(*this)
353  << " avg:" << gAverage(*this)
354  << endl;
355  }
356 }
357 
358 
360 (
361  Ostream& os
362 ) const
363 {
365 
366  if (haveQ_)
367  {
368  writeEntry(os, "Q", Q_);
369  }
370 
371  if (haveq_)
372  {
373  writeEntry(os, "q", q_);
374  }
375 
376  if (haveh_)
377  {
378  writeEntry(os, "h", h_);
379  writeEntry(os, Ta_());
380  writeEntryIfDifferent(os, "emissivity", scalar(0), emissivity_);
382  (
383  os,
384  "thicknessLayers",
385  scalarList(),
386  thicknessLayers_
387  );
389  (
390  os,
391  "kappaLayers",
392  scalarList(),
393  kappaLayers_
394  );
395  }
396 
397  writeEntryIfDifferent(os, "relaxation", scalar(1), relaxation_);
398 
399  if (qrName_ != word::null)
400  {
401  writeEntry(os, "qr", qrName_);
402  writeEntry(os, "qrRelaxation", qrRelaxation_);
403  writeEntry(os, "qrPrevious", qrPrevious_);
404  }
405 
406  writeEntry(os, "refValue", refValue());
407  writeEntry(os, "refGradient", refGrad());
408  writeEntry(os, "valueFraction", valueFraction());
409  writeEntry(os, "value", *this);
410 }
411 
412 
413 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
414 
415 namespace Foam
416 {
418  (
421  );
422 
424  (
427  patchMapper,
428  externalWallHeatFluxTemperature,
429  "externalWallHeatFluxTemperature"
430  );
431 
433  (
436  dictionary,
437  externalWallHeatFluxTemperature,
438  "externalWallHeatFluxTemperature"
439  );
440 }
441 
442 
443 // ************************************************************************* //
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...
Run-time selectable general function of one variable.
Definition: Function1.H:64
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
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
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 fvPatchFieldMapper &)
Map the given fvPatchField onto this fvPatchField.
externalTemperatureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
Foam::fvPatchFieldMapper.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:87
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:231
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:251
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:318
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
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
static Type NaN()
Return a primitive with all components set to NaN.
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:251
dimensionedScalar pow3(const dimensionedScalar &ds)
messageStream Info
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:61
IOerror FatalIOError
Type gAverage(const FieldField< Field, Type > &f)
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Type gMin(const FieldField< Field, Type > &f)
Type gMax(const FieldField< Field, Type > &f)
dictionary dict
volScalarField & p