externalWallHeatFluxTemperatureFvPatchScalarField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2015-2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
31 #include "fvPatchFieldMapper.H"
32 #include "volFields.H"
34 
36 
37 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38 
39 const Foam::Enum
40 <
42 >
44 ({
45  { operationMode::fixedPower, "power" },
46  { operationMode::fixedHeatFlux, "flux" },
47  { operationMode::fixedHeatTransferCoeff, "coefficient" },
48 });
49 
50 
51 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
52 
55 (
56  const fvPatch& p,
58 )
59 :
60  mixedFvPatchScalarField(p, iF),
62  (
63  patch(),
64  "undefined",
65  "undefined",
66  "undefined-K",
67  "undefined-alpha"
68  ),
69  mode_(fixedHeatFlux),
70  Q_(nullptr),
71  q_(nullptr),
72  h_(nullptr),
73  Ta_(nullptr),
74  relaxation_(1),
75  emissivity_(0),
76  qrRelaxation_(1),
77  qrName_("undefined-qr"),
78  thicknessLayers_(),
79  kappaLayers_()
80 {
81  refValue() = 0;
82  refGrad() = 0;
83  valueFraction() = 1;
84 }
85 
86 
89 (
90  const fvPatch& p,
92  const dictionary& dict
93 )
94 :
95  mixedFvPatchScalarField(p, iF),
97  mode_(operationModeNames.get("mode", dict)),
98  Q_(nullptr),
99  q_(nullptr),
100  h_(nullptr),
101  Ta_(nullptr),
102  relaxation_(dict.getOrDefault<scalar>("relaxation", 1)),
103  emissivity_(dict.getOrDefault<scalar>("emissivity", 0)),
104  qrRelaxation_(dict.getOrDefault<scalar>("qrRelaxation", 1)),
105  qrName_(dict.getOrDefault<word>("qr", "none")),
106  thicknessLayers_(),
107  kappaLayers_()
108 {
109  switch (mode_)
110  {
111  case fixedPower:
112  {
113  Q_ = Function1<scalar>::New("Q", dict, &db());
114  break;
115  }
116  case fixedHeatFlux:
117  {
119  break;
120  }
121  case fixedHeatTransferCoeff:
122  {
124  Ta_ = Function1<scalar>::New("Ta", dict, &db());
125 
126  if (dict.readIfPresent("thicknessLayers", thicknessLayers_))
127  {
128  dict.readEntry("kappaLayers", kappaLayers_);
129 
130  if (thicknessLayers_.size() != kappaLayers_.size())
131  {
133  << "\n number of layers for thicknessLayers and "
134  << "kappaLayers must be the same"
135  << "\n for patch " << p.name()
136  << " of field " << internalField().name()
137  << " in file " << internalField().objectPath()
138  << exit(FatalIOError);
139  }
140  }
141 
142  break;
143  }
144  }
145 
146  fvPatchScalarField::operator=(scalarField("value", dict, p.size()));
147 
148  if (qrName_ != "none")
149  {
150  if (dict.found("qrPrevious"))
151  {
152  qrPrevious_ = scalarField("qrPrevious", dict, p.size());
153  }
154  else
155  {
156  qrPrevious_.resize(p.size(), Zero);
157  }
158  }
159 
160  if (dict.found("refValue"))
161  {
162  // Full restart
163  refValue() = scalarField("refValue", dict, p.size());
164  refGrad() = scalarField("refGradient", dict, p.size());
165  valueFraction() = scalarField("valueFraction", dict, p.size());
166  }
167  else
168  {
169  // Start from user entered data. Assume fixedValue.
170  refValue() = *this;
171  refGrad() = 0;
172  valueFraction() = 1;
173  }
174 }
175 
176 
179 (
181  const fvPatch& p,
183  const fvPatchFieldMapper& mapper
184 )
185 :
186  mixedFvPatchScalarField(rhs, p, iF, mapper),
188  mode_(rhs.mode_),
189  Q_(rhs.Q_.clone()),
190  q_(rhs.q_.clone(patch().patch())),
191  h_(rhs.h_.clone(patch().patch())),
192  Ta_(rhs.Ta_.clone()),
193  relaxation_(rhs.relaxation_),
194  emissivity_(rhs.emissivity_),
195  qrPrevious_(),
196  qrRelaxation_(rhs.qrRelaxation_),
197  qrName_(rhs.qrName_),
198  thicknessLayers_(rhs.thicknessLayers_),
199  kappaLayers_(rhs.kappaLayers_)
200 {
201  if (qrName_ != "none")
202  {
203  qrPrevious_.resize(mapper.size());
204  qrPrevious_.map(rhs.qrPrevious_, mapper);
205  }
206 }
207 
208 
211 (
213 )
214 :
215  mixedFvPatchScalarField(rhs),
217  mode_(rhs.mode_),
218  Q_(rhs.Q_.clone()),
219  q_(rhs.q_.clone(patch().patch())),
220  h_(rhs.h_.clone(patch().patch())),
221  Ta_(rhs.Ta_.clone()),
222  relaxation_(rhs.relaxation_),
223  emissivity_(rhs.emissivity_),
224  qrPrevious_(rhs.qrPrevious_),
225  qrRelaxation_(rhs.qrRelaxation_),
226  qrName_(rhs.qrName_),
227  thicknessLayers_(rhs.thicknessLayers_),
228  kappaLayers_(rhs.kappaLayers_)
229 {}
230 
231 
234 (
237 )
238 :
239  mixedFvPatchScalarField(rhs, iF),
241  mode_(rhs.mode_),
242  Q_(rhs.Q_.clone()),
243  q_(rhs.q_.clone(patch().patch())),
244  h_(rhs.h_.clone(patch().patch())),
245  Ta_(rhs.Ta_.clone()),
246  relaxation_(rhs.relaxation_),
247  emissivity_(rhs.emissivity_),
248  qrPrevious_(rhs.qrPrevious_),
249  qrRelaxation_(rhs.qrRelaxation_),
250  qrName_(rhs.qrName_),
251  thicknessLayers_(rhs.thicknessLayers_),
252  kappaLayers_(rhs.kappaLayers_)
253 {}
254 
255 
256 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
257 
259 (
260  const fvPatchFieldMapper& mapper
261 )
262 {
263  mixedFvPatchScalarField::autoMap(mapper);
264  temperatureCoupledBase::autoMap(mapper);
265 
266  if (q_)
267  {
268  q_->autoMap(mapper);
269  }
270  if (h_)
271  {
272  h_->autoMap(mapper);
273  }
274 
275  if (qrName_ != "none")
276  {
277  qrPrevious_.autoMap(mapper);
278  }
279 }
280 
281 
283 (
284  const fvPatchScalarField& ptf,
285  const labelList& addr
286 )
287 {
288  mixedFvPatchScalarField::rmap(ptf, addr);
289 
290  const auto& rhs =
291  refCast<const externalWallHeatFluxTemperatureFvPatchScalarField>(ptf);
292 
293  temperatureCoupledBase::rmap(rhs, addr);
294 
295 
296  if (q_)
297  {
298  q_->rmap(rhs.q_(), addr);
299  }
300  if (h_)
301  {
302  h_->rmap(rhs.h_(), addr);
303  }
304 
305  if (qrName_ != "none")
306  {
307  qrPrevious_.rmap(rhs.qrPrevious_, addr);
308  }
309 }
310 
311 
313 {
314  if (updated())
315  {
316  return;
317  }
318 
319  const scalarField& Tp(*this);
320 
321  const scalarField valueFraction0(valueFraction());
322  const scalarField refValue0(refValue());
323 
324  scalarField qr(Tp.size(), Zero);
325  if (qrName_ != "none")
326  {
327  qr =
328  qrRelaxation_
329  *patch().lookupPatchField<volScalarField, scalar>(qrName_)
330  + (1 - qrRelaxation_)*qrPrevious_;
331 
332  qrPrevious_ = qr;
333  }
334 
335  switch (mode_)
336  {
337  case fixedPower:
338  {
339  const scalar heatPower =
340  Q_->value(this->db().time().timeOutputValue());
341 
342  refGrad() = (heatPower/gSum(patch().magSf()) + qr)/kappa(Tp);
343  refValue() = 0;
344  valueFraction() = 0;
345 
346  break;
347  }
348  case fixedHeatFlux:
349  {
350  tmp<scalarField> heatFlux =
351  q_->value(this->db().time().timeOutputValue());
352 
353  refGrad() = (heatFlux + qr)/kappa(Tp);
354  refValue() = 0;
355  valueFraction() = 0;
356 
357  break;
358  }
360  {
361  tmp<scalarField> thtcCoeff =
362  (
363  h_->value(this->db().time().timeOutputValue()) + VSMALL
364  );
365  const auto& htcCoeff = thtcCoeff();
366 
367  scalar totalSolidRes = 0;
368  if (thicknessLayers_.size())
369  {
370  forAll(thicknessLayers_, iLayer)
371  {
372  const scalar l = thicknessLayers_[iLayer];
373  if (kappaLayers_[iLayer] > 0)
374  {
375  totalSolidRes += l/kappaLayers_[iLayer];
376  }
377  }
378  }
379  scalarField hp(1/(1/htcCoeff + totalSolidRes));
380 
381  const scalar Ta =
382  Ta_->value(this->db().time().timeOutputValue());
383 
384  scalarField hpTa(hp*Ta);
385 
386  if (emissivity_ > 0)
387  {
388  // Evaluate the radiative flux to the environment
389  // from the surface temperature ...
390  if (totalSolidRes > 0)
391  {
392  // ... including the effect of the solid wall thermal
393  // resistance
394  scalarField TpLambda(htcCoeff/(htcCoeff + 1/totalSolidRes));
395  scalarField Ts(TpLambda*Tp + (1 - TpLambda)*Ta);
396  scalarField lambdaTa4(pow4((1 - TpLambda)*Ta));
397 
398  hp += emissivity_*sigma.value()*(pow4(Ts) - lambdaTa4)/Tp;
399  hpTa += emissivity_*sigma.value()*(lambdaTa4 + pow4(Ta));
400  }
401  else
402  {
403  // ... if there is no solid wall thermal resistance use
404  // the current wall temperature
405  hp += emissivity_*sigma.value()*pow3(Tp);
406  hpTa += emissivity_*sigma.value()*pow4(Ta);
407  }
408  }
409 
410  const scalarField kappaDeltaCoeffs
411  (
412  this->kappa(Tp)*patch().deltaCoeffs()
413  );
414 
415  refGrad() = 0;
416 
417  forAll(Tp, i)
418  {
419  if (qr[i] < 0)
420  {
421  const scalar hpmqr = hp[i] - qr[i]/Tp[i];
422 
423  refValue()[i] = hpTa[i]/hpmqr;
424  valueFraction()[i] = hpmqr/(hpmqr + kappaDeltaCoeffs[i]);
425  }
426  else
427  {
428  refValue()[i] = (hpTa[i] + qr[i])/hp[i];
429  valueFraction()[i] = hp[i]/(hp[i] + kappaDeltaCoeffs[i]);
430  }
431  }
432 
433  break;
434  }
435  }
436 
437  valueFraction() =
438  relaxation_*valueFraction() + (1 - relaxation_)*valueFraction0;
439  refValue() = relaxation_*refValue() + (1 - relaxation_)*refValue0;
440 
441  mixedFvPatchScalarField::updateCoeffs();
442 
443  DebugInfo
444  << patch().boundaryMesh().mesh().name() << ':' << patch().name() << ':'
445  << internalField().name() << " :"
446  << " heat transfer rate:" << gSum(kappa(Tp)*patch().magSf()*snGrad())
447  << " wall temperature "
448  << " min:" << gMin(*this)
449  << " max:" << gMax(*this)
450  << " avg:" << gAverage(*this) << nl;
451 }
452 
453 
455 (
456  Ostream& os
457 ) const
458 {
460 
461  os.writeEntry("mode", operationModeNames[mode_]);
463 
464  if (Q_)
465  {
466  Q_->writeData(os);
467  }
468  if (q_)
469  {
470  q_->writeData(os);
471  }
472  if (h_)
473  {
474  h_->writeData(os);
475  }
476  if (Ta_)
477  {
478  Ta_->writeData(os);
479  }
480 
481  switch (mode_)
482  {
483  case fixedHeatTransferCoeff:
484  {
485  if (relaxation_ < 1)
486  {
487  os.writeEntry("relaxation", relaxation_);
488  }
489 
490  if (emissivity_ > 0)
491  {
492  os.writeEntry("emissivity", emissivity_);
493  }
494 
495  if (thicknessLayers_.size())
496  {
497  thicknessLayers_.writeEntry("thicknessLayers", os);
498  kappaLayers_.writeEntry("kappaLayers", os);
499  }
500 
501  break;
502  }
503 
504  default:
505  break;
506  }
507 
508  os.writeEntry("qr", qrName_);
509 
510  if (qrName_ != "none")
511  {
512  os.writeEntry("qrRelaxation", qrRelaxation_);
513 
514  qrPrevious_.writeEntry("qrPrevious", os);
515  }
516 
517  refValue().writeEntry("refValue", os);
518  refGrad().writeEntry("refGradient", os);
519  valueFraction().writeEntry("valueFraction", os);
520  writeEntry("value", os);
521 }
522 
523 
524 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
525 
526 namespace Foam
527 {
529  (
532  );
533 }
534 
535 // ************************************************************************* //
Foam::fvPatchField< scalar >
volFields.H
Foam::fvPatchField< scalar >::write
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:384
Foam::fvc::snGrad
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:47
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::Enum
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: IOstreamOption.H:57
Foam::externalWallHeatFluxTemperatureFvPatchScalarField::fixedHeatFlux
Heat flux [W/m2].
Definition: externalWallHeatFluxTemperatureFvPatchScalarField.H:224
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::externalWallHeatFluxTemperatureFvPatchScalarField::autoMap
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Definition: externalWallHeatFluxTemperatureFvPatchScalarField.C:259
Foam::gAverage
Type gAverage(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:604
Foam::FieldMapper::size
virtual label size() const =0
Foam::externalWallHeatFluxTemperatureFvPatchScalarField::rmap
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Definition: externalWallHeatFluxTemperatureFvPatchScalarField.C:283
Foam::temperatureCoupledBase
Common functions used in temperature coupled boundaries.
Definition: temperatureCoupledBase.H:135
Foam::FatalIOError
IOerror FatalIOError
fvPatchFieldMapper.H
Foam::gSum
Type gSum(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:594
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::pow4
dimensionedScalar pow4(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:100
Foam::Field< scalar >
Foam::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:89
externalWallHeatFluxTemperatureFvPatchScalarField.H
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
Foam::externalWallHeatFluxTemperatureFvPatchScalarField::fixedHeatTransferCoeff
Heat transfer coefficient [W/m^2/K].
Definition: externalWallHeatFluxTemperatureFvPatchScalarField.H:225
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Foam::externalWallHeatFluxTemperatureFvPatchScalarField::operationMode
operationMode
Operation mode enumeration.
Definition: externalWallHeatFluxTemperatureFvPatchScalarField.H:221
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
os
OBJstream os(runTime.globalPath()/outputName)
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::externalWallHeatFluxTemperatureFvPatchScalarField::fixedPower
Heat power [W].
Definition: externalWallHeatFluxTemperatureFvPatchScalarField.H:223
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::externalWallHeatFluxTemperatureFvPatchScalarField::externalWallHeatFluxTemperatureFvPatchScalarField
externalWallHeatFluxTemperatureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: externalWallHeatFluxTemperatureFvPatchScalarField.C:55
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
Foam::externalWallHeatFluxTemperatureFvPatchScalarField::operationModeNames
static const Enum< operationMode > operationModeNames
Definition: externalWallHeatFluxTemperatureFvPatchScalarField.H:228
Foam::externalWallHeatFluxTemperatureFvPatchScalarField::write
void write(Ostream &) const
Write.
Definition: externalWallHeatFluxTemperatureFvPatchScalarField.C:455
DebugInfo
#define DebugInfo
Report an information message using Foam::Info.
Definition: messageStream.H:382
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
physicoChemicalConstants.H
Foam::List< label >
Foam::temperatureCoupledBase::kappa
virtual tmp< scalarField > kappa(const scalarField &Tp) const
Given patch temperature calculate corresponding K field.
Definition: temperatureCoupledBase.C:210
Foam::constant::physicoChemical::sigma
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m2/K4].
Foam::Ostream::writeEntry
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:236
FatalIOErrorInFunction
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:473
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
Foam::temperatureCoupledBase::write
void write(Ostream &os) const
Write.
Definition: temperatureCoupledBase.C:512
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
sigma
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
Foam::makePatchTypeField
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Foam::gMin
Type gMin(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:593
Foam::externalWallHeatFluxTemperatureFvPatchScalarField
This boundary condition applies a heat flux condition to temperature on an external wall in one of th...
Definition: externalWallHeatFluxTemperatureFvPatchScalarField.H:211
Foam::externalWallHeatFluxTemperatureFvPatchScalarField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: externalWallHeatFluxTemperatureFvPatchScalarField.C:312
Foam::gMax
Type gMax(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:592
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54