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 -------------------------------------------------------------------------------
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"
33 #include "mappedPatchBase.H"
34 
35 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39 namespace 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")
59 {
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),
77  mappedPatchFieldBase<scalar>
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),
101  mappedPatchFieldBase<scalar>
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
132  thicknessLayer_ = PatchFunction1<scalar>::NewIfPresent
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 //
173 }
174 
175 
178 (
179  const turbulentTemperatureCoupledBaffleMixedFvPatchScalarField& wtcsf,
180  const DimensionedField<scalar, volMesh>& iF
181 )
182 :
183  mixedFvPatchScalarField(wtcsf, iF),
184  temperatureCoupledBase(patch(), wtcsf),
185  mappedPatchFieldBase<scalar>
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),
207  mappedPatchFieldBase<scalar>
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_)
232  {
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  Info<< patch().boundaryMesh().mesh().name() << ':'
389  << patch().name() << ':'
390  << this->internalField().name() << " <- "
391  << mpp.sampleRegion() << ':'
392  << mpp.samplePatch() << ':'
393  << this->internalField().name() << " :"
394  << " heat transfer rate:" << Q
395  << " walltemperature "
396  << " min:" << gMin(*this)
397  << " max:" << gMax(*this)
398  << " avg:" << gAverage(*this)
399  << endl;
400  }
401 
402  UPstream::msgType(oldTag); // Restore tag
403 }
404 
405 
407 (
408  fvMatrix<scalar>& m,
409  const label iMatrix,
410  const direction cmpt
411 )
412 {
414  << "This BC does not support energy coupling "
415  << "Use compressible::turbulentTemperatureRadCoupledMixed "
416  << "which has more functionalities and it can handle "
417  << "the assemble coupled option for energy. "
418  << abort(FatalError);
419 }
420 
421 
422 tmp<Field<scalar>>
423 turbulentTemperatureCoupledBaffleMixedFvPatchScalarField::coeffs
424 (
425  fvMatrix<scalar>& matrix,
426  const Field<scalar>& coeffs,
427  const label mat
428 ) const
429 {
431  << "This BC does not support energy coupling "
432  << "Use compressible::turbulentTemperatureRadCoupledMixed "
433  << "which has more functionalities and it can handle "
434  << "the assemble coupled option for energy. "
435  << abort(FatalError);
436  /*
437  const label index(this->patch().index());
438 
439  const label nSubFaces(matrix.lduMesh().cellBoundMap()[mat][index].size());
440 
441  Field<scalar> mapCoeffs(nSubFaces, Zero);
442 
443  label subFaceI = 0;
444  forAll(*this, faceI)
445  {
446 
447  }
448  */
449  return tmp<Field<scalar>>(new Field<scalar>());
450 }
451 
452 
454 (
455  Ostream& os
456 ) const
457 {
459  os.writeEntry("Tnbr", TnbrName_);
460  if (thicknessLayer_)
461  {
462  thicknessLayer_().writeData(os);
463  kappaLayer_().writeData(os);
464  }
465  if (thicknessLayers_.size())
466  {
467  thicknessLayers_.writeEntry("thicknessLayers", os);
468  kappaLayers_.writeEntry("kappaLayers", os);
469  }
472 }
473 
474 
475 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
476 
478 (
480  turbulentTemperatureCoupledBaffleMixedFvPatchScalarField
481 );
482 
483 
484 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
485 
486 } // End namespace compressible
487 } // End namespace Foam
488 
489 
490 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
dictionary dict
Mixed boundary condition for temperature, to be used for heat-transfer on back-to-back baffles...
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
uint8_t direction
Definition: direction.H:46
turbulentTemperatureCoupledBaffleMixedFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
static int incrMsgType(int val=1) noexcept
Increment the message tag for standard messages.
Definition: UPstream.H:1251
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
static const mappedPatchBase & mapper(const fvPatch &p, const DimensionedField< Type, volMesh > &iF)
Check that patch is of correct type.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
Type gMin(const FieldField< Field, Type > &f)
virtual void write(Ostream &os) const
Write.
Type & refCast(U &obj)
A dynamic_cast (for references). Generates a FatalError on failed casts and uses the virtual type() m...
Definition: typeInfo.H:159
virtual tmp< scalarField > kappa(const scalarField &Tp) const
Given patch temperature calculate corresponding K field.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
virtual void rmap(const fvPatchField< scalar > &, const labelList &)=0
Reverse map the given fvPatchField onto this fvPatchField.
static int & msgType() noexcept
Message tag of standard messages.
Definition: UPstream.H:1229
virtual void autoMap(const fvPatchFieldMapper &)=0
Map (and resize as needed) from self given a mapping object.
void write(Ostream &os) const
Write.
Macros for easy insertion into run-time selection tables.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
makePatchTypeField(fvPatchScalarField, alphatFixedDmdtWallBoilingWallFunctionFvPatchScalarField)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
Type gSum(const FieldField< Field, Type > &f)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
fvPatchField< scalar > fvPatchScalarField
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
A FieldMapper for finite-volume patch fields.
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:64
errorManip< error > abort(error &err)
Definition: errorManip.H:139
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
bool compressible
Definition: pEqn.H:2
int debug
Static debugging option.
Type gMax(const FieldField< Field, Type > &f)
OBJstream os(runTime.globalPath()/outputName)
Common functions used in temperature coupled boundaries.
#define WarningInFunction
Report a warning using Foam::Warning.
Type gAverage(const FieldField< Field, Type > &f)
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
const std::string patch
OpenFOAM patch number as a std::string.
static autoPtr< PatchFunction1< Type > > NewIfPresent(const polyPatch &pp, const word &entryName, const dictionary &dict, const bool faceValues=true)
An optional selector.
messageStream Info
Information stream (stdout output on master, null elsewhere)
void writeEntry(Ostream &os) const
Write the UList with its compound type.
Definition: UListIO.C:30
virtual void write(Ostream &) const
Write.
virtual void rmap(const fvPatchField< scalar > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
virtual tmp< scalarField > kappa(const scalarField &Tp) const
Given patch temperature calculate corresponding K field. Override temperatureCoupledBase::kappa to in...
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
void distribute(const word &fieldName, Field< T > &newValues) const
Wrapper for mapDistribute::distribute that knows about dabase mapping.
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:40
Namespace for OpenFOAM.
virtual void manipulateMatrix(fvMatrix< scalar > &m, const label iMatrix, const direction cmpt)
Manipulate matrix.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127