Loading...
Searching...
No Matches
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-2022 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
39const Foam::Enum
40<
42>
44({
45 { operationMode::fixedPower, "power" },
47 { operationMode::fixedHeatTransferCoeff, "coefficient" },
48});
49
50
51// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
52
55(
56 const fvPatch& p,
57 const DimensionedField<scalar, volMesh>& iF
58)
59:
60 mixedFvPatchScalarField(p, iF),
61 temperatureCoupledBase(patch()), // default method (fluidThermo)
62 mode_(fixedHeatFlux),
63 Q_(nullptr),
64 q_(nullptr),
65 h_(nullptr),
66 Ta_(nullptr),
67 relaxation_(1),
68 emissivity_(0),
69 qrRelaxation_(1),
70 qrName_("undefined-qr"),
71 thicknessLayers_(),
72 kappaLayers_()
74 refValue() = 0;
75 refGrad() = 0;
76 valueFraction() = 1;
77}
78
79
82(
83 const fvPatch& p,
85 const dictionary& dict
86)
87:
88 mixedFvPatchScalarField(p, iF),
90 mode_(operationModeNames.get("mode", dict)),
91 Q_(nullptr),
92 q_(nullptr),
93 h_(nullptr),
94 Ta_(nullptr),
95 relaxation_(dict.getOrDefault<scalar>("relaxation", 1)),
96 emissivity_(dict.getOrDefault<scalar>("emissivity", 0)),
97 qrRelaxation_(dict.getOrDefault<scalar>("qrRelaxation", 1)),
98 qrName_(dict.getOrDefault<word>("qr", "none")),
99 thicknessLayers_(),
100 kappaLayers_()
101{
102 switch (mode_)
103 {
104 case fixedPower:
105 {
106 Q_ = Function1<scalar>::New("Q", dict, &db());
107 break;
108 }
109 case fixedHeatFlux:
110 {
111 q_ = PatchFunction1<scalar>::New(patch().patch(), "q", dict);
112 break;
113 }
115 {
116 h_ = PatchFunction1<scalar>::New(patch().patch(), "h", dict);
117 Ta_ = Function1<scalar>::New("Ta", dict, &db());
118
119 if (dict.readIfPresent("thicknessLayers", thicknessLayers_))
120 {
121 dict.readEntry("kappaLayers", kappaLayers_);
122
123 if (thicknessLayers_.size() != kappaLayers_.size())
124 {
126 << "\n number of layers for thicknessLayers and "
127 << "kappaLayers must be the same"
128 << "\n for patch " << p.name()
129 << " of field " << internalField().name()
130 << " in file " << internalField().objectPath()
131 << exit(FatalIOError);
132 }
133 }
134
135 break;
136 }
137 }
138
139 if (qrName_ != "none")
140 {
141 if (dict.found("qrPrevious"))
142 {
143 qrPrevious_ = scalarField("qrPrevious", dict, p.size());
144 }
145 else
146 {
147 qrPrevious_.resize(p.size(), Zero);
148 }
149 }
150
151 this->readValueEntry(dict, IOobjectOption::MUST_READ);
152
153 if (this->readMixedEntries(dict))
154 {
155 // Full restart
156 }
157 else
158 {
159 // Start from user entered data. Assume fixedValue.
160 refValue() = *this;
161 refGrad() = 0;
162 valueFraction() = 1;
163 }
164}
165
166
169(
170 const externalWallHeatFluxTemperatureFvPatchScalarField& rhs,
171 const fvPatch& p,
172 const DimensionedField<scalar, volMesh>& iF,
173 const fvPatchFieldMapper& mapper
174)
175:
176 mixedFvPatchScalarField(rhs, p, iF, mapper),
177 temperatureCoupledBase(patch(), rhs),
178 mode_(rhs.mode_),
179 Q_(rhs.Q_.clone()),
180 q_(rhs.q_.clone(patch().patch())),
181 h_(rhs.h_.clone(patch().patch())),
182 Ta_(rhs.Ta_.clone()),
183 relaxation_(rhs.relaxation_),
184 emissivity_(rhs.emissivity_),
185 qrPrevious_(),
186 qrRelaxation_(rhs.qrRelaxation_),
187 qrName_(rhs.qrName_),
188 thicknessLayers_(rhs.thicknessLayers_),
189 kappaLayers_(rhs.kappaLayers_)
190{
191 if (qrName_ != "none")
193 qrPrevious_.resize(mapper.size());
194 qrPrevious_.map(rhs.qrPrevious_, mapper);
195 }
196}
197
198
201(
202 const externalWallHeatFluxTemperatureFvPatchScalarField& rhs
203)
204:
205 mixedFvPatchScalarField(rhs),
206 temperatureCoupledBase(rhs),
207 mode_(rhs.mode_),
208 Q_(rhs.Q_.clone()),
209 q_(rhs.q_.clone(patch().patch())),
210 h_(rhs.h_.clone(patch().patch())),
211 Ta_(rhs.Ta_.clone()),
212 relaxation_(rhs.relaxation_),
213 emissivity_(rhs.emissivity_),
214 qrPrevious_(rhs.qrPrevious_),
215 qrRelaxation_(rhs.qrRelaxation_),
216 qrName_(rhs.qrName_),
217 thicknessLayers_(rhs.thicknessLayers_),
218 kappaLayers_(rhs.kappaLayers_)
219{}
220
221
224(
227)
228:
229 mixedFvPatchScalarField(rhs, iF),
230 temperatureCoupledBase(patch(), rhs),
231 mode_(rhs.mode_),
232 Q_(rhs.Q_.clone()),
233 q_(rhs.q_.clone(patch().patch())),
234 h_(rhs.h_.clone(patch().patch())),
235 Ta_(rhs.Ta_.clone()),
236 relaxation_(rhs.relaxation_),
237 emissivity_(rhs.emissivity_),
238 qrPrevious_(rhs.qrPrevious_),
239 qrRelaxation_(rhs.qrRelaxation_),
240 qrName_(rhs.qrName_),
241 thicknessLayers_(rhs.thicknessLayers_),
242 kappaLayers_(rhs.kappaLayers_)
243{}
244
245
246// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
247
249(
250 const fvPatchFieldMapper& mapper
251)
252{
253 mixedFvPatchScalarField::autoMap(mapper);
255
256 if (q_)
257 {
258 q_->autoMap(mapper);
259 }
260 if (h_)
261 {
262 h_->autoMap(mapper);
263 }
264
265 if (qrName_ != "none")
266 {
267 qrPrevious_.autoMap(mapper);
268 }
269}
270
271
273(
274 const fvPatchScalarField& ptf,
275 const labelList& addr
276)
277{
278 mixedFvPatchScalarField::rmap(ptf, addr);
279
280 const auto& rhs =
282
284
285
286 if (q_)
287 {
288 q_->rmap(rhs.q_(), addr);
289 }
290 if (h_)
291 {
292 h_->rmap(rhs.h_(), addr);
293 }
294
295 if (qrName_ != "none")
296 {
297 qrPrevious_.rmap(rhs.qrPrevious_, addr);
298 }
299}
300
301
303{
304 if (updated())
305 {
306 return;
307 }
308
309 const scalarField& Tp(*this);
310
311 const scalarField valueFraction0(valueFraction());
312 const scalarField refValue0(refValue());
313
314 scalarField qr(Tp.size(), Zero);
315 if (qrName_ != "none")
316 {
317 qr = lerp
318 (
319 qrPrevious_,
320 patch().lookupPatchField<volScalarField>(qrName_),
321 qrRelaxation_
322 );
323 qrPrevious_ = qr;
324 }
325
326 switch (mode_)
327 {
328 case fixedPower:
329 {
330 const scalar heatPower =
331 Q_->value(this->db().time().timeOutputValue());
332
333 refGrad() = (heatPower/gSum(patch().magSf()) + qr)/kappa(Tp);
334 refValue() = 293.15; // prevents FPE, no impact on condition
335 valueFraction() = 0;
336
337 break;
338 }
339 case fixedHeatFlux:
340 {
341 tmp<scalarField> heatFlux =
342 q_->value(this->db().time().timeOutputValue());
343
344 refGrad() = (heatFlux + qr)/kappa(Tp);
345 refValue() = 293.15; // prevents FPE, no impact on condition
346 valueFraction() = 0;
347
348 break;
349 }
350 case fixedHeatTransferCoeff:
351 {
352 tmp<scalarField> thtcCoeff =
353 (
354 h_->value(this->db().time().timeOutputValue()) + VSMALL
355 );
356 const auto& htcCoeff = thtcCoeff();
357
358 scalar totalSolidRes = 0;
359 if (thicknessLayers_.size())
360 {
361 forAll(thicknessLayers_, iLayer)
362 {
363 const scalar l = thicknessLayers_[iLayer];
364 if (kappaLayers_[iLayer] > 0)
365 {
366 totalSolidRes += l/kappaLayers_[iLayer];
367 }
368 }
369 }
370
371 const scalar Ta =
372 Ta_->value(this->db().time().timeOutputValue());
373
374 scalarField hrad(Tp.size(), Zero);
375
376 if (emissivity_ > 0)
377 {
378 const scalar eSig(emissivity_*sigma.value());
379 // Evaluate the radiative flux to the environment
380 // from the surface temperature ...
381 if (totalSolidRes > 0)
382 {
383 // ... including the effect of the solid wall thermal
384 // resistance
385 scalarField TpLambda(htcCoeff/(htcCoeff + 1/totalSolidRes));
386 scalarField Ts(TpLambda*Ta + (1 - TpLambda)*Tp);
387 hrad = eSig*((pow3(Ta) + pow3(Ts)) + Ta*Ts*(Ta + Ts));
388
389 forAll(hrad, i)
390 {
391 scalar hradTmp0 = hrad[i];
392 scalar TaLambda =
393 (htcCoeff[i] + hradTmp0)
394 /(htcCoeff[i] + hradTmp0 + 1/totalSolidRes);
395
396 scalar TsiNew = TaLambda*Ta + (1 - TaLambda)*Tp[i];
397 scalar Tsi = Ts[i];
398
399 while (mag(Tsi - TsiNew)/Tsi > 0.01)
400 {
401 Tsi = TsiNew;
402 scalar hradNew
403 (
404 eSig*((pow3(Ta) + pow3(Tsi)) + Ta*Tsi*(Ta + Tsi))
405 );
406
407 TaLambda =
408 (htcCoeff[i] + hradNew)
409 /(htcCoeff[i] + hradNew + 1/totalSolidRes);
410
411 TsiNew = TaLambda*Ta + (1 - TaLambda)*Tp[i];
412 };
413
414 hrad[i] =
415 eSig*((pow3(Ta) + pow3(Tsi)) + Ta*Tsi*(Ta + Tsi));
416 }
417 }
418 else
419 {
420 hrad = eSig*((pow3(Ta) + pow3(Tp)) + Ta*Tp*(Ta + Tp));
421 }
422 }
423
424 const scalarField hp(1/(1/(htcCoeff + hrad) + totalSolidRes));
425
426 const scalarField hpTa(hp*Ta);
427
428 const scalarField kappaDeltaCoeffs
429 (
430 this->kappa(Tp)*patch().deltaCoeffs()
431 );
432
433 refGrad() = 0;
434
435 forAll(Tp, i)
436 {
437 if (qr[i] < 0)
438 {
439 const scalar hpmqr = hp[i] - qr[i]/Tp[i];
440
441 refValue()[i] = hpTa[i]/hpmqr;
442 valueFraction()[i] = hpmqr/(hpmqr + kappaDeltaCoeffs[i]);
443 }
444 else
445 {
446 refValue()[i] = (hpTa[i] + qr[i])/hp[i];
447 valueFraction()[i] = hp[i]/(hp[i] + kappaDeltaCoeffs[i]);
448 }
449 }
450
451 break;
452 }
453 }
454
455 valueFraction() = lerp(valueFraction0, valueFraction(), relaxation_);
456 refValue() = lerp(refValue0, refValue(), relaxation_);
457
458 mixedFvPatchScalarField::updateCoeffs();
459
461 << patch().boundaryMesh().mesh().name() << ':' << patch().name() << ':'
462 << internalField().name() << " :"
463 << " heat transfer rate:" << gSum(kappa(Tp)*patch().magSf()*snGrad())
464 << " wall temperature "
465 << " min:" << gMin(*this)
466 << " max:" << gMax(*this)
467 << " avg:" << gAverage(*this) << nl;
468}
469
470
472(
473 Ostream& os
474) const
475{
477
478 os.writeEntry("mode", operationModeNames[mode_]);
480
481 if (Q_)
482 {
483 Q_->writeData(os);
484 }
485 if (q_)
486 {
487 q_->writeData(os);
488 }
489 if (h_)
490 {
491 h_->writeData(os);
492 }
493 if (Ta_)
494 {
495 Ta_->writeData(os);
496 }
497
498 switch (mode_)
499 {
500 case fixedHeatTransferCoeff:
501 {
502 if (relaxation_ < 1)
503 {
504 os.writeEntry("relaxation", relaxation_);
505 }
506
507 if (emissivity_ > 0)
508 {
509 os.writeEntry("emissivity", emissivity_);
510 }
511
512 if (thicknessLayers_.size())
513 {
514 thicknessLayers_.writeEntry("thicknessLayers", os);
515 kappaLayers_.writeEntry("kappaLayers", os);
516 }
517
518 break;
519 }
520
521 default:
522 break;
523 }
524
525 os.writeEntry("qr", qrName_);
526
527 if (qrName_ != "none")
528 {
529 os.writeEntry("qrRelaxation", qrRelaxation_);
530
531 qrPrevious_.writeEntry("qrPrevious", os);
532 }
533
534 refValue().writeEntry("refValue", os);
535 refGrad().writeEntry("refGradient", os);
536 valueFraction().writeEntry("valueFraction", os);
539
540
541// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
542
543namespace Foam
544{
546 (
548 externalWallHeatFluxTemperatureFvPatchScalarField
549 );
550}
551
552// ************************************************************************* //
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...
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition Enum.H:57
virtual label size() const =0
The size of the mapper.
@ MUST_READ
Reading required.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition Ostream.H:331
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
const Type & value() const noexcept
Return const reference to value.
This boundary condition applies a heat flux condition to temperature on an external wall in one of th...
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
externalWallHeatFluxTemperatureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
A FieldMapper for finite-volume patch fields.
virtual void write(Ostream &) const
Write.
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
Common functions used in temperature coupled boundaries.
virtual void autoMap(const fvPatchFieldMapper &)=0
Map (and resize as needed) from self given a mapping object.
void write(Ostream &os) const
Write.
virtual void rmap(const fvPatchField< scalar > &, const labelList &)=0
Reverse map the given fvPatchField onto this fvPatchField.
virtual tmp< scalarField > kappa(const scalarField &Tp) const
Given patch temperature calculate corresponding K field.
temperatureCoupledBase(const fvPatch &patch, const KMethodType method=KMethodType::mtFluidThermo)
Default construct from patch, using fluidThermo (default) or specified method.
A class for managing temporary objects.
Definition tmp.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition error.H:629
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Define a concrete fvPatchField type and add to run-time tables Example, (fvPatchScalarField,...
#define DebugInfo
Report an information message using Foam::Info.
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m2/K4].
const std::string patch
OpenFOAM patch number as a std::string.
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition fvcSnGrad.C:40
Namespace for OpenFOAM.
Type gAverage(const FieldField< Field, Type > &f, const label comm)
The global arithmetic average of a FieldField.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
Type gSum(const FieldField< Field, Type > &f)
List< label > labelList
A List of labels.
Definition List.H:62
dimensionedScalar pow3(const dimensionedScalar &ds)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
void rhs(fvMatrix< typename Expr::value_type > &m, const Expr &expression)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
Type gMin(const FieldField< Field, Type > &f)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
Type gMax(const FieldField< Field, Type > &f)
dimensioned< Type > lerp(const dimensioned< Type > &a, const dimensioned< Type > &b, const scalar t)
fvPatchField< scalar > fvPatchScalarField
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
dictionary dict
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299