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-2020 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 
35 namespace Foam
36 {
37  template<>
38  const char*
39  NamedEnum
40  <
42  3
43  >::names[] =
44  {
45  "power",
46  "flux",
47  "coefficient"
48  };
49 }
50 
51 const Foam::NamedEnum
52 <
54  3
56 
57 
58 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
59 
62 (
63  const fvPatch& p,
65 )
66 :
67  mixedFvPatchScalarField(p, iF),
68  temperatureCoupledBase(patch()),
69  mode_(fixedHeatFlux),
70  Q_(0),
71  Ta_(),
72  relaxation_(1),
73  emissivity_(0),
74  qrRelaxation_(1),
75  qrName_("undefined-qr"),
76  thicknessLayers_(),
77  kappaLayers_()
78 {
79  refValue() = 0;
80  refGrad() = 0;
81  valueFraction() = 1;
82 }
83 
84 
87 (
88  const fvPatch& p,
90  const dictionary& dict
91 )
92 :
93  mixedFvPatchScalarField(p, iF),
94  temperatureCoupledBase(patch(), dict),
95  mode_(operationModeNames.read(dict.lookup("mode"))),
96  Q_(0),
97  Ta_(),
98  relaxation_(dict.lookupOrDefault<scalar>("relaxation", 1)),
99  emissivity_(dict.lookupOrDefault<scalar>("emissivity", 0)),
100  qrRelaxation_(dict.lookupOrDefault<scalar>("qrRelaxation", 1)),
101  qrName_(dict.lookupOrDefault<word>("qr", "none")),
102  thicknessLayers_(),
103  kappaLayers_()
104 {
105  switch (mode_)
106  {
107  case fixedPower:
108  {
109  dict.lookup("Q") >> Q_;
110 
111  break;
112  }
113  case fixedHeatFlux:
114  {
115  q_ = scalarField("q", dict, p.size());
116 
117  break;
118  }
119  case fixedHeatTransferCoeff:
120  {
121  h_ = scalarField("h", dict, p.size());
122  Ta_ = Function1<scalar>::New("Ta", dict);
123 
124  if (dict.found("thicknessLayers"))
125  {
126  dict.lookup("thicknessLayers") >> thicknessLayers_;
127  dict.lookup("kappaLayers") >> kappaLayers_;
128  }
129 
130  break;
131  }
132  }
133 
134  fvPatchScalarField::operator=(scalarField("value", dict, p.size()));
135 
136  if (qrName_ != "none")
137  {
138  if (dict.found("qrPrevious"))
139  {
140  qrPrevious_ = scalarField("qrPrevious", dict, p.size());
141  }
142  else
143  {
144  qrPrevious_.setSize(p.size(), 0);
145  }
146  }
147 
148  if (dict.found("refValue"))
149  {
150  // Full restart
151  refValue() = scalarField("refValue", dict, p.size());
152  refGrad() = scalarField("refGradient", dict, p.size());
153  valueFraction() = scalarField("valueFraction", dict, p.size());
154  }
155  else
156  {
157  // Start from user entered data. Assume fixedValue.
158  refValue() = *this;
159  refGrad() = 0;
160  valueFraction() = 1;
161  }
162 }
163 
164 
167 (
169  const fvPatch& p,
171  const fvPatchFieldMapper& mapper
172 )
173 :
174  mixedFvPatchScalarField(ptf, p, iF, mapper),
175  temperatureCoupledBase(patch(), ptf),
176  mode_(ptf.mode_),
177  Q_(ptf.Q_),
178  Ta_(ptf.Ta_, false),
179  relaxation_(ptf.relaxation_),
180  emissivity_(ptf.emissivity_),
181  qrRelaxation_(ptf.qrRelaxation_),
182  qrName_(ptf.qrName_),
183  thicknessLayers_(ptf.thicknessLayers_),
184  kappaLayers_(ptf.kappaLayers_)
185 {
186  switch (mode_)
187  {
188  case fixedPower:
189  {
190  break;
191  }
192  case fixedHeatFlux:
193  {
194  mapper(q_, ptf.q_);
195  break;
196  }
197  case fixedHeatTransferCoeff:
198  {
199  mapper(h_, ptf.h_);
200  break;
201  }
202  }
203 
204  if (qrName_ != "none")
205  {
206  mapper(qrPrevious_, ptf.qrPrevious_);
207  }
208 }
209 
210 
213 (
215 )
216 :
217  mixedFvPatchScalarField(tppsf),
218  temperatureCoupledBase(tppsf),
219  mode_(tppsf.mode_),
220  Q_(tppsf.Q_),
221  q_(tppsf.q_),
222  h_(tppsf.h_),
223  Ta_(tppsf.Ta_, false),
224  relaxation_(tppsf.relaxation_),
225  emissivity_(tppsf.emissivity_),
226  qrPrevious_(tppsf.qrPrevious_),
227  qrRelaxation_(tppsf.qrRelaxation_),
228  qrName_(tppsf.qrName_),
229  thicknessLayers_(tppsf.thicknessLayers_),
230  kappaLayers_(tppsf.kappaLayers_)
231 {}
232 
233 
236 (
239 )
240 :
241  mixedFvPatchScalarField(tppsf, iF),
242  temperatureCoupledBase(patch(), tppsf),
243  mode_(tppsf.mode_),
244  Q_(tppsf.Q_),
245  q_(tppsf.q_),
246  h_(tppsf.h_),
247  Ta_(tppsf.Ta_, false),
248  relaxation_(tppsf.relaxation_),
249  emissivity_(tppsf.emissivity_),
250  qrPrevious_(tppsf.qrPrevious_),
251  qrRelaxation_(tppsf.qrRelaxation_),
252  qrName_(tppsf.qrName_),
253  thicknessLayers_(tppsf.thicknessLayers_),
254  kappaLayers_(tppsf.kappaLayers_)
255 {}
256 
257 
258 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
259 
261 (
262  const fvPatchFieldMapper& m
263 )
264 {
265  mixedFvPatchScalarField::autoMap(m);
266 
267  switch (mode_)
268  {
269  case fixedPower:
270  {
271  break;
272  }
273  case fixedHeatFlux:
274  {
275  m(q_, q_);
276 
277  break;
278  }
279  case fixedHeatTransferCoeff:
280  {
281  m(h_, h_);
282 
283  break;
284  }
285  }
286 
287  if (qrName_ != "none")
288  {
289  m(qrPrevious_, qrPrevious_);
290  }
291 }
292 
293 
295 (
296  const fvPatchScalarField& ptf,
297  const labelList& addr
298 )
299 {
300  mixedFvPatchScalarField::rmap(ptf, addr);
301 
303  refCast<const externalWallHeatFluxTemperatureFvPatchScalarField>(ptf);
304 
305  switch (mode_)
306  {
307  case fixedPower:
308  {
309  break;
310  }
311  case fixedHeatFlux:
312  {
313  q_.rmap(tiptf.q_, addr);
314 
315  break;
316  }
317  case fixedHeatTransferCoeff:
318  {
319  h_.rmap(tiptf.h_, addr);
320 
321  break;
322  }
323  }
324 
325  if (qrName_ != "none")
326  {
327  qrPrevious_.rmap(tiptf.qrPrevious_, addr);
328  }
329 }
330 
331 
333 {
334  if (updated())
335  {
336  return;
337  }
338 
339  const scalarField& Tp(*this);
340 
341  // Store current valueFraction and refValue for relaxation
342  const scalarField valueFraction0(valueFraction());
343  const scalarField refValue0(refValue());
344 
345  scalarField qr(Tp.size(), 0);
346  if (qrName_ != "none")
347  {
348  qr =
349  qrRelaxation_
350  *patch().lookupPatchField<volScalarField, scalar>(qrName_)
351  + (1 - qrRelaxation_)*qrPrevious_;
352 
353  qrPrevious_ = qr;
354  }
355 
356  switch (mode_)
357  {
358  case fixedPower:
359  {
360  refGrad() = (Q_/gSum(patch().magSf()) + qr)/kappa(*this);
361  refValue() = Tp;
362  valueFraction() = 0;
363 
364  break;
365  }
366  case fixedHeatFlux:
367  {
368  refGrad() = (q_ + qr)/kappa(*this);
369  refValue() = Tp;
370  valueFraction() = 0;
371 
372  break;
373  }
374  case fixedHeatTransferCoeff:
375  {
376  scalar totalSolidRes = 0;
377  if (thicknessLayers_.size())
378  {
379  forAll(thicknessLayers_, iLayer)
380  {
381  const scalar l = thicknessLayers_[iLayer];
382  if (kappaLayers_[iLayer] > 0)
383  {
384  totalSolidRes += l/kappaLayers_[iLayer];
385  }
386  }
387  }
388  scalarField hp(1/(1/h_ + totalSolidRes));
389 
390  const scalar Ta = Ta_->value(this->db().time().timeOutputValue());
391  scalarField hpTa(hp*Ta);
392 
393  if (emissivity_ > 0)
394  {
395  // Evaluate the radiative flux to the environment
396  // from the surface temperature ...
397  if (totalSolidRes > 0)
398  {
399  // ... including the effect of the solid wall thermal
400  // resistance
401  scalarField TpLambda(h_/(h_ + 1/totalSolidRes));
402  scalarField Ts(TpLambda*Tp + (1 - TpLambda)*Ta);
403  scalarField lambdaTa4(pow4((1 - TpLambda)*Ta));
404 
405  hp += emissivity_*sigma.value()*(pow4(Ts) - lambdaTa4)/Tp;
406  hpTa += emissivity_*sigma.value()*(lambdaTa4 + pow4(Ta));
407  }
408  else
409  {
410  // ... if there is no solid wall thermal resistance use
411  // the current wall temperature
412  hp += emissivity_*sigma.value()*pow3(Tp);
413  hpTa += emissivity_*sigma.value()*pow4(Ta);
414  }
415  }
416 
417  const scalarField kappaDeltaCoeffs
418  (
419  this->kappa(*this)*patch().deltaCoeffs()
420  );
421 
422  refGrad() = 0;
423 
424  forAll(Tp, i)
425  {
426  if (qr[i] < 0)
427  {
428  const scalar hpmqr = hp[i] - qr[i]/Tp[i];
429 
430  refValue()[i] = hpTa[i]/hpmqr;
431  valueFraction()[i] = hpmqr/(hpmqr + kappaDeltaCoeffs[i]);
432  }
433  else
434  {
435  refValue()[i] = (hpTa[i] + qr[i])/hp[i];
436  valueFraction()[i] = hp[i]/(hp[i] + kappaDeltaCoeffs[i]);
437  }
438  }
439 
440  break;
441  }
442  }
443 
444  valueFraction() =
445  relaxation_*valueFraction()
446  + (1 - relaxation_)*valueFraction0;
447 
448  refValue() = relaxation_*refValue() + (1 - relaxation_)*refValue0;
449 
450  mixedFvPatchScalarField::updateCoeffs();
451 
452  if (debug)
453  {
454  const scalar Q = gSum(kappa(*this)*patch().magSf()*snGrad());
455 
456  Info<< patch().boundaryMesh().mesh().name() << ':'
457  << patch().name() << ':'
458  << this->internalField().name() << " :"
459  << " heat transfer rate:" << Q
460  << " walltemperature "
461  << " min:" << gMin(*this)
462  << " max:" << gMax(*this)
463  << " avg:" << gAverage(*this)
464  << endl;
465  }
466 }
467 
468 
470 (
471  Ostream& os
472 ) const
473 {
475 
476  writeEntry(os, "mode", operationModeNames[mode_]);
478 
479  switch (mode_)
480  {
481  case fixedPower:
482  {
483  writeEntry(os, "Q", Q_);
484 
485  break;
486  }
487  case fixedHeatFlux:
488  {
489  writeEntry(os, "q", q_);
490 
491  break;
492  }
493  case fixedHeatTransferCoeff:
494  {
495  writeEntry(os, "h", h_);
496  writeEntry(os, Ta_());
497 
498  if (relaxation_ < 1)
499  {
500  writeEntry(os, "relaxation", relaxation_);
501  }
502 
503  if (emissivity_ > 0)
504  {
505  writeEntry(os, "emissivity", emissivity_);
506  }
507 
508  if (thicknessLayers_.size())
509  {
510  writeEntry(os, "thicknessLayers", thicknessLayers_);
511  writeEntry(os, "kappaLayers", kappaLayers_);
512  }
513 
514  break;
515  }
516  }
517 
518  writeEntry(os, "qr", qrName_);
519 
520  if (qrName_ != "none")
521  {
522  writeEntry(os, "qrRelaxation", qrRelaxation_);
523 
524  writeEntry(os, "qrPrevious", qrPrevious_);
525  }
526 
527  writeEntry(os, "refValue", refValue());
528  writeEntry(os, "refGradient", refGrad());
529  writeEntry(os, "valueFraction", valueFraction());
530  writeEntry(os, "value", *this);
531 }
532 
533 
534 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
535 
536 namespace Foam
537 {
539  (
542  );
543 }
544 
545 // ************************************************************************* //
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:667
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
const dimensionedScalar & kappa
Coulomb constant: default SI units: [N.m2/C2].
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
Type gMin(const FieldField< Field, Type > &f)
const dimensionedScalar & sigma
Stefan-Boltzmann constant: default SI units: [W/m^2/K^4].
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
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:61
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:280
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:51
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Type gSum(const FieldField< Field, Type > &f)
A class for handling words, derived from string.
Definition: word.H:59
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Foam::fvPatchFieldMapper.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
const Type & value() const
Return const reference to value.
virtual label size() const
Return size.
Definition: fvPatch.H:155
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
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)
void setSize(const label)
Reset size of List.
Definition: List.C:281
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,.
Type gAverage(const FieldField< Field, Type > &f)
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:295
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
dimensionedScalar pow4(const dimensionedScalar &ds)
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 in one of th...
void write(Ostream &) const
Write.
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:45
static autoPtr< Function1< Type > > New(const word &entryName, const dictionary &dict)
Selector.
Definition: Function1New.C:32
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:812