Loading...
Searching...
No Matches
turbulentTemperatureCoupledBaffleMixedFvPatchScalarField.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-2016 OpenFOAM Foundation
9 Copyright (C) 2019-2021 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"
33#include "mappedPatchBase.H"
34
35// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36
37namespace Foam
38{
39namespace compressible
40{
41
42// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
43
46(
47 const fvPatch& p,
49)
50:
51 mixedFvPatchScalarField(p, iF),
52 temperatureCoupledBase(patch()), // default method (fluidThermo)
54 (
56 *this
57 ),
58 TnbrName_("undefined-Tnbr")
60 this->refValue() = Zero;
61 this->refGrad() = Zero;
62 this->valueFraction() = 1.0;
63}
64
65
68(
70 const fvPatch& p,
72 const fvPatchFieldMapper& mapper
73)
74:
75 mixedFvPatchScalarField(ptf, p, iF, mapper),
76 temperatureCoupledBase(patch(), ptf),
78 (
79 mappedPatchFieldBase<scalar>::mapper(p, iF),
80 *this,
81 ptf
82 ),
83 TnbrName_(ptf.TnbrName_),
84 thicknessLayers_(ptf.thicknessLayers_),
85 thicknessLayer_(ptf.thicknessLayer_.clone(p.patch())),
86 kappaLayers_(ptf.kappaLayers_),
87 kappaLayer_(ptf.kappaLayer_.clone(p.patch()))
88{}
89
90
93(
94 const fvPatch& p,
96 const dictionary& dict
97)
98:
99 mixedFvPatchScalarField(p, iF),
102 (
103 mappedPatchFieldBase<scalar>::mapper(p, iF),
104 *this,
105 dict
106 ),
107 TnbrName_(dict.get<word>("Tnbr"))
108{
109 if (!isA<mappedPatchBase>(this->patch().patch()))
110 {
112 << "' not type '" << mappedPatchBase::typeName << "'"
113 << "\n for patch " << p.name()
114 << " of field " << internalField().name()
115 << " in file " << internalField().objectPath()
116 << exit(FatalError);
117 }
118
120 << "This BC has been superseded by "
121 << "compressible::turbulentTemperatureRadCoupledMixed "
122 << "which has more functionalities and it can handle "
123 << "the assemble coupled option for energy. "
124 << endl;
125
126 // Read list of layers
127 if (dict.readIfPresent("thicknessLayers", thicknessLayers_))
128 {
129 dict.readEntry("kappaLayers", kappaLayers_);
130 }
131 // Read single additional PatchFunction1
133 (
134 p.patch(),
135 "thicknessLayer",
136 dict
137 );
139 (
140 p.patch(),
141 "kappaLayer",
142 dict
143 );
144
145
146 this->readValueEntry(dict, IOobjectOption::MUST_READ);
147
148 if (this->readMixedEntries(dict))
149 {
150 // Full restart
151 }
152 else
153 {
154 // Start from user entered data. Assume fixedValue.
155 refValue() = *this;
156 refGrad() = Zero;
157 valueFraction() = 1.0;
158 }
159
160// This blocks (crashes) with more than two worlds!
161//
167
169/// (
170/// this->internalField().name() + "_weights",
171/// this->patch().deltaCoeffs()
172/// );
173}
174
175
178(
181)
182:
183 mixedFvPatchScalarField(wtcsf, iF),
184 temperatureCoupledBase(patch(), wtcsf),
186 (
187 mappedPatchFieldBase<scalar>::mapper(patch(), iF),
188 *this,
189 wtcsf
190 ),
191 TnbrName_(wtcsf.TnbrName_),
192 thicknessLayers_(wtcsf.thicknessLayers_),
193 thicknessLayer_(wtcsf.thicknessLayer_.clone(patch().patch())),
194 kappaLayers_(wtcsf.kappaLayers_),
195 kappaLayer_(wtcsf.kappaLayer_.clone(patch().patch()))
196{}
197
198
201(
203)
204:
205 mixedFvPatchScalarField(wtcsf),
206 temperatureCoupledBase(patch(), wtcsf),
208 (
209 mappedPatchFieldBase<scalar>::mapper(patch(), wtcsf.internalField()),
210 *this,
211 wtcsf
212 ),
213 TnbrName_(wtcsf.TnbrName_),
214 thicknessLayers_(wtcsf.thicknessLayers_),
215 thicknessLayer_(wtcsf.thicknessLayer_.clone(patch().patch())),
216 kappaLayers_(wtcsf.kappaLayers_),
217 kappaLayer_(wtcsf.kappaLayer_.clone(patch().patch()))
218{}
219
220
221// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
222
224(
225 const fvPatchFieldMapper& mapper
226)
227{
228 mixedFvPatchScalarField::autoMap(mapper);
230 //mappedPatchFieldBase<scalar>::autoMap(mapper);
231 if (thicknessLayer_)
233 thicknessLayer_().autoMap(mapper);
234 kappaLayer_().autoMap(mapper);
235 }
236}
237
238
240(
241 const fvPatchField<scalar>& ptf,
242 const labelList& addr
243)
244{
245 mixedFvPatchScalarField::rmap(ptf, addr);
246
248 refCast
249 <
251 >(ptf);
252
253 temperatureCoupledBase::rmap(tiptf, addr);
254 //mappedPatchFieldBase<scalar>::rmap(ptf, addr);
255 if (thicknessLayer_)
256 {
257 thicknessLayer_().rmap(tiptf.thicknessLayer_(), addr);
258 kappaLayer_().rmap(tiptf.kappaLayer_(), addr);
259 }
260}
261
262
265(
266 const scalarField& Tp
267) const
268{
269 // Get kappa from relevant thermo
271
272 // Optionally modify with explicit resistance
273 if (thicknessLayer_ || thicknessLayers_.size())
274 {
275 scalarField KDelta(tk*patch().deltaCoeffs());
276
277 // Harmonic averaging of kappa*deltaCoeffs
278 {
279 KDelta = 1.0/KDelta;
280 if (thicknessLayer_)
281 {
282 const scalar t = db().time().timeOutputValue();
283 KDelta +=
284 thicknessLayer_().value(t)
285 /kappaLayer_().value(t);
286 }
287 if (thicknessLayers_.size())
288 {
289 forAll(thicknessLayers_, iLayer)
290 {
291 KDelta += thicknessLayers_[iLayer]/kappaLayers_[iLayer];
292 }
293 }
294 KDelta = 1.0/KDelta;
295 }
296
297 // Update kappa from KDelta
298 tk = KDelta/patch().deltaCoeffs();
299 }
300
301 return tk;
302}
303
304
306{
307 if (updated())
308 {
309 return;
310 }
311
312 // Since we're inside initEvaluate/evaluate there might be processor
313 // comms underway. Change the tag we use.
314 const int oldTag = UPstream::incrMsgType();
315
316 // Get the coupling information from the mappedPatchBase
317 const mappedPatchBase& mpp =
319 (
320 patch(),
321 this->internalField()
322 );
323
324 const scalarField& Tp = *this;
325 const scalarField kappaTp(kappa(Tp));
326 const tmp<scalarField> myKDelta = kappaTp*patch().deltaCoeffs();
327
328
329 scalarField nbrIntFld;
330 scalarField nbrKDelta;
331 if (mpp.sameWorld())
332 {
333 // Same world so lookup
334 const auto& nbrMesh = refCast<const fvMesh>(mpp.sampleMesh());
335 const label nbrPatchID = mpp.samplePolyPatch().index();
336 const auto& nbrPatch = nbrMesh.boundary()[nbrPatchID];
337
339 nbrField =
340 refCast
341 <
343 >
344 (
345 nbrPatch.lookupPatchField<volScalarField>(TnbrName_)
346 );
347
348 // Swap to obtain full local values of neighbour K*delta
349 nbrIntFld = nbrField.patchInternalField();
350 nbrKDelta = nbrField.kappa(nbrField)*nbrPatch.deltaCoeffs();
351 }
352 else
353 {
354 // Different world so use my region,patch. Distribution below will
355 // do the reordering.
356 nbrIntFld = patchInternalField();
357 nbrKDelta = myKDelta.ref();
358 }
359 distribute(this->internalField().name() + "_value", nbrIntFld);
360 distribute(this->internalField().name() + "_weights", nbrKDelta);
361
362
363 // Both sides agree on
364 // - temperature : (myKDelta*fld + nbrKDelta*nbrFld)/(myKDelta+nbrKDelta)
365 // - gradient : (temperature-fld)*delta
366 // We've got a degree of freedom in how to implement this in a mixed bc.
367 // (what gradient, what fixedValue and mixing coefficient)
368 // Two reasonable choices:
369 // 1. specify above temperature on one side (preferentially the high side)
370 // and above gradient on the other. So this will switch between pure
371 // fixedvalue and pure fixedgradient
372 // 2. specify gradient and temperature such that the equations are the
373 // same on both sides. This leads to the choice of
374 // - refGradient = zero gradient
375 // - refValue = neighbour value
376 // - mixFraction = nbrKDelta / (nbrKDelta + myKDelta())
377
378 this->refValue() = nbrIntFld;
379 this->refGrad() = Zero;
380 this->valueFraction() = nbrKDelta/(nbrKDelta + myKDelta());
381
382 mixedFvPatchScalarField::updateCoeffs();
383
384 if (debug)
385 {
386 scalar Q = gSum(kappaTp*patch().magSf()*snGrad());
387
388 auto limits = gMinMax(*this);
389 auto avg = gAverage(*this);
390
391 Info<< patch().boundaryMesh().mesh().name() << ':'
392 << patch().name() << ':'
393 << this->internalField().name() << " <- "
394 << mpp.sampleRegion() << ':'
395 << mpp.samplePatch() << ':'
396 << this->internalField().name() << " :"
397 << " heat transfer rate:" << Q
398 << " walltemperature "
399 << " min:" << limits.min()
400 << " max:" << limits.max()
401 << " avg:" << avg
403 }
404
405 UPstream::msgType(oldTag); // Restore tag
406}
407
408
410(
412 const label iMatrix,
413 const direction cmpt
414)
415{
417 << "This BC does not support energy coupling "
418 << "Use compressible::turbulentTemperatureRadCoupledMixed "
419 << "which has more functionalities and it can handle "
420 << "the assemble coupled option for energy. "
421 << abort(FatalError);
422}
423
424
425tmp<Field<scalar>>
426turbulentTemperatureCoupledBaffleMixedFvPatchScalarField::coeffs
427(
428 fvMatrix<scalar>& matrix,
429 const Field<scalar>& coeffs,
430 const label mat
431) const
432{
434 << "This BC does not support energy coupling "
435 << "Use compressible::turbulentTemperatureRadCoupledMixed "
436 << "which has more functionalities and it can handle "
437 << "the assemble coupled option for energy. "
438 << abort(FatalError);
439
440 return nullptr;
441}
442
443
445(
446 Ostream& os
447) const
448{
450 os.writeEntry("Tnbr", TnbrName_);
451 if (thicknessLayer_)
452 {
453 thicknessLayer_().writeData(os);
454 kappaLayer_().writeData(os);
455 }
456 if (thicknessLayers_.size())
457 {
458 thicknessLayers_.writeEntry("thicknessLayers", os);
459 kappaLayers_.writeEntry("kappaLayers", os);
460 }
463}
464
465
466// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
467
469(
471 turbulentTemperatureCoupledBaffleMixedFvPatchScalarField
472);
473
474
475// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
476
477} // End namespace compressible
478} // End namespace Foam
479
480
481// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
if(patchID !=-1)
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
@ MUST_READ
Reading required.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
static autoPtr< PatchFunction1< Type > > NewIfPresent(const polyPatch &pp, const word &entryName, const dictionary &dict, const bool faceValues=true)
An optional selector.
void writeEntry(Ostream &os) const
Write the UList with its compound type.
Definition UListIO.C:29
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
static int & msgType() noexcept
Message tag of standard messages.
Definition UPstream.H:1926
static int incrMsgType(int val=1) noexcept
Increment the message tag for standard messages.
Definition UPstream.H:1948
Mixed boundary condition for temperature, to be used for heat-transfer on back-to-back baffles....
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
virtual tmp< scalarField > kappa(const scalarField &Tp) const
Given patch temperature calculate corresponding K field. Override temperatureCoupledBase::kappa to in...
virtual void rmap(const fvPatchField< scalar > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
virtual void manipulateMatrix(fvMatrix< scalar > &m, const label iMatrix, const direction cmpt)
Manipulate matrix.
turbulentTemperatureCoupledBaffleMixedFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition fvMatrix.H:118
A FieldMapper for finite-volume patch fields.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
Determines a mapping between patch face centres and mesh cell or face centres and processors they're ...
const polyMesh & sampleMesh() const
Get the region mesh.
const polyPatch & samplePolyPatch() const
Get the patch on the region.
const word & samplePatch() const
Patch (only if NEARESTPATCHFACE).
const word & sampleRegion() const
Region to sample.
bool sameWorld() const
Is sample world the local world?
Functionality for sampling fields using mappedPatchBase. Every call to mappedField() returns a sample...
virtual void write(Ostream &os) const
Write.
static const mappedPatchBase & mapper(const fvPatch &p, const DimensionedField< scalar, volMesh > &iF)
mappedPatchFieldBase(const mappedPatchBase &mapper, const fvPatchField< scalar > &patchField, const word &fieldName, const bool setAverage, const scalar average, const word &interpolationScheme)
void distribute(const word &fieldName, Field< T > &newValues) const
virtual void write(Ostream &) const
Write.
label index() const noexcept
The index of this patch in the boundaryMesh.
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
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition tmpI.H:235
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
auto limits
Definition setRDeltaT.H:186
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Define a concrete fvPatchField type and add to run-time tables Example, (fvPatchScalarField,...
auto & name
#define WarningInFunction
Report a warning using Foam::Warning.
Namespace for handling debugging switches.
Definition debug.C:45
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
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
errorManip< error > abort(error &err)
Definition errorManip.H:139
uint8_t direction
Definition direction.H:49
MinMax< Type > gMinMax(const FieldField< Field, Type > &f)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
fvPatchField< scalar > fvPatchScalarField
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299