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-2019 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_(0),
71  q_(),
72  h_(),
73  Ta_(),
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_(0),
99  q_(),
100  h_(),
101  Ta_(),
102  relaxation_(dict.lookupOrDefault<scalar>("relaxation", 1)),
103  emissivity_(dict.lookupOrDefault<scalar>("emissivity", 0)),
104  qrRelaxation_(dict.lookupOrDefault<scalar>("qrRelaxation", 1)),
105  qrName_(dict.lookupOrDefault<word>("qr", "none")),
106  thicknessLayers_(),
107  kappaLayers_()
108 {
109  switch (mode_)
110  {
111  case fixedPower:
112  {
113  dict.readEntry("Q", Q_);
114 
115  break;
116  }
117  case fixedHeatFlux:
118  {
119  q_ = scalarField("q", dict, p.size());
120 
121  break;
122  }
123  case fixedHeatTransferCoeff:
124  {
125  h_ = scalarField("h", dict, p.size());
126  Ta_ = Function1<scalar>::New("Ta", dict);
127 
128  if (dict.readIfPresent("thicknessLayers", thicknessLayers_))
129  {
130  dict.readEntry("kappaLayers", kappaLayers_);
131 
132  if (thicknessLayers_.size() != kappaLayers_.size())
133  {
135  << "\n number of layers for thicknessLayers and "
136  << "kappaLayers must be the same"
137  << "\n for patch " << p.name()
138  << " of field " << internalField().name()
139  << " in file " << internalField().objectPath()
140  << exit(FatalIOError);
141  }
142  }
143 
144  break;
145  }
146  }
147 
148  fvPatchScalarField::operator=(scalarField("value", dict, p.size()));
149 
150  if (qrName_ != "none")
151  {
152  if (dict.found("qrPrevious"))
153  {
154  qrPrevious_ = scalarField("qrPrevious", dict, p.size());
155  }
156  else
157  {
158  qrPrevious_.setSize(p.size(), 0);
159  }
160  }
161 
162  if (dict.found("refValue"))
163  {
164  // Full restart
165  refValue() = scalarField("refValue", dict, p.size());
166  refGrad() = scalarField("refGradient", dict, p.size());
167  valueFraction() = scalarField("valueFraction", dict, p.size());
168  }
169  else
170  {
171  // Start from user entered data. Assume fixedValue.
172  refValue() = *this;
173  refGrad() = 0;
174  valueFraction() = 1;
175  }
176 }
177 
178 
181 (
183  const fvPatch& p,
185  const fvPatchFieldMapper& mapper
186 )
187 :
188  mixedFvPatchScalarField(ptf, p, iF, mapper),
190  mode_(ptf.mode_),
191  Q_(ptf.Q_),
192  q_(),
193  h_(),
194  Ta_(ptf.Ta_.clone()),
195  relaxation_(ptf.relaxation_),
196  emissivity_(ptf.emissivity_),
197  qrPrevious_(),
198  qrRelaxation_(ptf.qrRelaxation_),
199  qrName_(ptf.qrName_),
200  thicknessLayers_(ptf.thicknessLayers_),
201  kappaLayers_(ptf.kappaLayers_)
202 {
203  switch (mode_)
204  {
205  case fixedPower:
206  {
207  break;
208  }
209  case fixedHeatFlux:
210  {
211  q_.setSize(mapper.size());
212  q_.map(ptf.q_, mapper);
213 
214  break;
215  }
216  case fixedHeatTransferCoeff:
217  {
218  h_.setSize(mapper.size());
219  h_.map(ptf.h_, mapper);
220 
221  break;
222  }
223  }
224 
225  if (qrName_ != "none")
226  {
227  qrPrevious_.setSize(mapper.size());
228  qrPrevious_.map(ptf.qrPrevious_, mapper);
229  }
230 }
231 
232 
235 (
237 )
238 :
239  mixedFvPatchScalarField(ewhftpsf),
240  temperatureCoupledBase(ewhftpsf),
241  mode_(ewhftpsf.mode_),
242  Q_(ewhftpsf.Q_),
243  q_(ewhftpsf.q_),
244  h_(ewhftpsf.h_),
245  Ta_(ewhftpsf.Ta_.clone()),
246  relaxation_(ewhftpsf.relaxation_),
247  emissivity_(ewhftpsf.emissivity_),
248  qrPrevious_(ewhftpsf.qrPrevious_),
249  qrRelaxation_(ewhftpsf.qrRelaxation_),
250  qrName_(ewhftpsf.qrName_),
251  thicknessLayers_(ewhftpsf.thicknessLayers_),
252  kappaLayers_(ewhftpsf.kappaLayers_)
253 {}
254 
255 
258 (
261 )
262 :
263  mixedFvPatchScalarField(ewhftpsf, iF),
264  temperatureCoupledBase(patch(), ewhftpsf),
265  mode_(ewhftpsf.mode_),
266  Q_(ewhftpsf.Q_),
267  q_(ewhftpsf.q_),
268  h_(ewhftpsf.h_),
269  Ta_(ewhftpsf.Ta_.clone()),
270  relaxation_(ewhftpsf.relaxation_),
271  emissivity_(ewhftpsf.emissivity_),
272  qrPrevious_(ewhftpsf.qrPrevious_),
273  qrRelaxation_(ewhftpsf.qrRelaxation_),
274  qrName_(ewhftpsf.qrName_),
275  thicknessLayers_(ewhftpsf.thicknessLayers_),
276  kappaLayers_(ewhftpsf.kappaLayers_)
277 {}
278 
279 
280 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
281 
283 (
284  const fvPatchFieldMapper& m
285 )
286 {
287  mixedFvPatchScalarField::autoMap(m);
288 
289  switch (mode_)
290  {
291  case fixedPower:
292  {
293  break;
294  }
295  case fixedHeatFlux:
296  {
297  q_.autoMap(m);
298 
299  break;
300  }
301  case fixedHeatTransferCoeff:
302  {
303  h_.autoMap(m);
304 
305  break;
306  }
307  }
308 
309  if (qrName_ != "none")
310  {
311  qrPrevious_.autoMap(m);
312  }
313 }
314 
315 
317 (
318  const fvPatchScalarField& ptf,
319  const labelList& addr
320 )
321 {
322  mixedFvPatchScalarField::rmap(ptf, addr);
323 
325  refCast<const externalWallHeatFluxTemperatureFvPatchScalarField>(ptf);
326 
327  switch (mode_)
328  {
329  case fixedPower:
330  {
331  break;
332  }
333  case fixedHeatFlux:
334  {
335  q_.rmap(ewhftpsf.q_, addr);
336 
337  break;
338  }
339  case fixedHeatTransferCoeff:
340  {
341  h_.rmap(ewhftpsf.h_, addr);
342 
343  break;
344  }
345  }
346 
347  if (qrName_ != "none")
348  {
349  qrPrevious_.rmap(ewhftpsf.qrPrevious_, addr);
350  }
351 }
352 
353 
355 {
356  if (updated())
357  {
358  return;
359  }
360 
361  const scalarField& Tp(*this);
362 
363  scalarField qr(Tp.size(), Zero);
364  if (qrName_ != "none")
365  {
366  qr =
367  qrRelaxation_
368  *patch().lookupPatchField<volScalarField, scalar>(qrName_)
369  + (1 - qrRelaxation_)*qrPrevious_;
370 
371  qrPrevious_ = qr;
372  }
373 
374  switch (mode_)
375  {
376  case fixedPower:
377  {
378  refGrad() = (Q_/gSum(patch().magSf()) + qr)/kappa(Tp);
379  refValue() = 0;
380  valueFraction() = 0;
381 
382  break;
383  }
384  case fixedHeatFlux:
385  {
386  refGrad() = (q_ + qr)/kappa(Tp);
387  refValue() = 0;
388  valueFraction() = 0;
389 
390  break;
391  }
393  {
394  scalar totalSolidRes = 0;
395  if (thicknessLayers_.size())
396  {
397  forAll(thicknessLayers_, iLayer)
398  {
399  const scalar l = thicknessLayers_[iLayer];
400  if (kappaLayers_[iLayer] > 0)
401  {
402  totalSolidRes += l/kappaLayers_[iLayer];
403  }
404  }
405  }
406  scalarField hp(1/(1/h_ + totalSolidRes));
407 
408  const scalar Ta = Ta_->value(this->db().time().timeOutputValue());
409  scalarField hpTa(hp*Ta);
410 
411  if (emissivity_ > 0)
412  {
413  // Evaluate the radiative flux to the environment
414  // from the surface temperature ...
415  if (totalSolidRes > 0)
416  {
417  // ... including the effect of the solid wall thermal
418  // resistance
419  scalarField TpLambda(h_/(h_ + 1/totalSolidRes));
420  scalarField Ts(TpLambda*Tp + (1 - TpLambda)*Ta);
421  scalarField lambdaTa4(pow4((1 - TpLambda)*Ta));
422 
423  hp += emissivity_*sigma.value()*(pow4(Ts) - lambdaTa4)/Tp;
424  hpTa += emissivity_*sigma.value()*(lambdaTa4 + pow4(Ta));
425  }
426  else
427  {
428  // ... if there is no solid wall thermal resistance use
429  // the current wall temperature
430  hp += emissivity_*sigma.value()*pow3(Tp);
431  hpTa += emissivity_*sigma.value()*pow4(Ta);
432  }
433  }
434 
435  const scalarField kappaDeltaCoeffs
436  (
437  this->kappa(Tp)*patch().deltaCoeffs()
438  );
439 
440  refGrad() = 0;
441 
442  forAll(Tp, i)
443  {
444  if (qr[i] < 0)
445  {
446  const scalar hpmqr = hp[i] - qr[i]/Tp[i];
447 
448  refValue()[i] = hpTa[i]/hpmqr;
449  valueFraction()[i] = hpmqr/(hpmqr + kappaDeltaCoeffs[i]);
450  }
451  else
452  {
453  refValue()[i] = (hpTa[i] + qr[i])/hp[i];
454  valueFraction()[i] = hp[i]/(hp[i] + kappaDeltaCoeffs[i]);
455  }
456  }
457 
458  break;
459  }
460  }
461 
462  valueFraction() = relaxation_*valueFraction() + (1 - relaxation_);
463  refValue() = relaxation_*refValue() + (1 - relaxation_)*Tp;
464 
465  mixedFvPatchScalarField::updateCoeffs();
466 
467  DebugInfo
468  << patch().boundaryMesh().mesh().name() << ':' << patch().name() << ':'
469  << internalField().name() << " :"
470  << " heat transfer rate:" << gSum(kappa(Tp)*patch().magSf()*snGrad())
471  << " wall temperature "
472  << " min:" << gMin(*this)
473  << " max:" << gMax(*this)
474  << " avg:" << gAverage(*this) << nl;
475 }
476 
477 
479 (
480  Ostream& os
481 ) const
482 {
484 
485  os.writeEntry("mode", operationModeNames[mode_]);
487 
488  switch (mode_)
489  {
490  case fixedPower:
491  {
492  os.writeEntry("Q", Q_);
493  break;
494  }
495  case fixedHeatFlux:
496  {
497  q_.writeEntry("q", os);
498 
499  break;
500  }
501  case fixedHeatTransferCoeff:
502  {
503  h_.writeEntry("h", os);
504  Ta_->writeData(os);
505 
506  if (relaxation_ < 1)
507  {
508  os.writeEntry("relaxation", relaxation_);
509  }
510 
511  if (emissivity_ > 0)
512  {
513  os.writeEntry("emissivity", emissivity_);
514  }
515 
516  if (thicknessLayers_.size())
517  {
518  thicknessLayers_.writeEntry("thicknessLayers", os);
519  kappaLayers_.writeEntry("kappaLayers", os);
520  }
521 
522  break;
523  }
524  }
525 
526  os.writeEntry("qr", qrName_);
527 
528  if (qrName_ != "none")
529  {
530  os.writeEntry("qrRelaxation", qrRelaxation_);
531 
532  qrPrevious_.writeEntry("qrPrevious", os);
533  }
534 
535  refValue().writeEntry("refValue", os);
536  refGrad().writeEntry("refGradient", os);
537  valueFraction().writeEntry("valueFraction", os);
538  writeEntry("value", os);
539 }
540 
541 
542 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
543 
544 namespace Foam
545 {
547  (
550  );
551 }
552 
553 // ************************************************************************* //
Foam::fvPatchField< scalar >
volFields.H
Foam::fvPatchField< scalar >::write
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:364
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:51
Foam::externalWallHeatFluxTemperatureFvPatchScalarField::fixedHeatFlux
Fixed heat flux [W/m2].
Definition: externalWallHeatFluxTemperatureFvPatchScalarField.H:217
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::externalWallHeatFluxTemperatureFvPatchScalarField::autoMap
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Definition: externalWallHeatFluxTemperatureFvPatchScalarField.C:283
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:317
Foam::temperatureCoupledBase
Common functions used in temperature coupled boundaries.
Definition: temperatureCoupledBase.H:110
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:290
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:63
Foam::externalWallHeatFluxTemperatureFvPatchScalarField::fixedHeatTransferCoeff
Fixed heat transfer coefficient.
Definition: externalWallHeatFluxTemperatureFvPatchScalarField.H:218
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Foam::externalWallHeatFluxTemperatureFvPatchScalarField::operationMode
operationMode
Operation mode enumeration.
Definition: externalWallHeatFluxTemperatureFvPatchScalarField.H:214
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:121
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::externalWallHeatFluxTemperatureFvPatchScalarField::fixedPower
Fixed heat power [W].
Definition: externalWallHeatFluxTemperatureFvPatchScalarField.H:216
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:221
Foam::externalWallHeatFluxTemperatureFvPatchScalarField::write
void write(Ostream &) const
Write.
Definition: externalWallHeatFluxTemperatureFvPatchScalarField.C:479
DebugInfo
#define DebugInfo
Report an information message using Foam::Info.
Definition: messageStream.H:350
Foam::nl
constexpr char nl
Definition: Ostream.H:372
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:138
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:219
FatalIOErrorInFunction
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:375
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
Foam::temperatureCoupledBase::write
void write(Ostream &os) const
Write.
Definition: temperatureCoupledBase.C:426
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:204
Foam::externalWallHeatFluxTemperatureFvPatchScalarField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: externalWallHeatFluxTemperatureFvPatchScalarField.C:354
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