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() = 0.0;
61  this->refGrad() = 0.0;
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  fvPatchScalarField::operator=(scalarField("value", dict, p.size()));
147 
148  if (dict.found("refValue"))
149  {
150  // Full restart
151  refValue() = scalarField("refValue", dict, p.size());
152  refGrad() = scalarField("refGradient", dict, p.size());
153  valueFraction() = scalarField("valueFraction", dict, p.size());
154  }
155  else
156  {
157  // Start from user entered data. Assume fixedValue.
158  refValue() = *this;
159  refGrad() = 0.0;
160  valueFraction() = 1.0;
161  }
162 
163 // This blocks (crashes) with more than two worlds!
164 //
176 }
177 
178 
181 (
182  const turbulentTemperatureCoupledBaffleMixedFvPatchScalarField& wtcsf,
183  const DimensionedField<scalar, volMesh>& iF
184 )
185 :
186  mixedFvPatchScalarField(wtcsf, iF),
187  temperatureCoupledBase(patch(), wtcsf),
188  mappedPatchFieldBase<scalar>
189  (
190  mappedPatchFieldBase<scalar>::mapper(patch(), iF),
191  *this,
192  wtcsf
193  ),
194  TnbrName_(wtcsf.TnbrName_),
195  thicknessLayers_(wtcsf.thicknessLayers_),
196  thicknessLayer_(wtcsf.thicknessLayer_.clone(patch().patch())),
197  kappaLayers_(wtcsf.kappaLayers_),
198  kappaLayer_(wtcsf.kappaLayer_.clone(patch().patch()))
199 {}
200 
201 
204 (
206 )
207 :
208  mixedFvPatchScalarField(wtcsf),
209  temperatureCoupledBase(patch(), wtcsf),
210  mappedPatchFieldBase<scalar>
211  (
212  mappedPatchFieldBase<scalar>::mapper(patch(), wtcsf.internalField()),
213  *this,
214  wtcsf
215  ),
216  TnbrName_(wtcsf.TnbrName_),
217  thicknessLayers_(wtcsf.thicknessLayers_),
218  thicknessLayer_(wtcsf.thicknessLayer_.clone(patch().patch())),
219  kappaLayers_(wtcsf.kappaLayers_),
220  kappaLayer_(wtcsf.kappaLayer_.clone(patch().patch()))
221 {}
222 
223 
224 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
225 
227 (
228  const fvPatchFieldMapper& mapper
229 )
230 {
231  mixedFvPatchScalarField::autoMap(mapper);
233  //mappedPatchFieldBase<scalar>::autoMap(mapper);
234  if (thicknessLayer_)
235  {
236  thicknessLayer_().autoMap(mapper);
237  kappaLayer_().autoMap(mapper);
238  }
239 }
240 
241 
243 (
244  const fvPatchField<scalar>& ptf,
245  const labelList& addr
246 )
247 {
248  mixedFvPatchScalarField::rmap(ptf, addr);
249 
251  refCast
252  <
254  >(ptf);
255 
256  temperatureCoupledBase::rmap(tiptf, addr);
257  //mappedPatchFieldBase<scalar>::rmap(ptf, addr);
258  if (thicknessLayer_)
259  {
260  thicknessLayer_().rmap(tiptf.thicknessLayer_(), addr);
261  kappaLayer_().rmap(tiptf.kappaLayer_(), addr);
262  }
263 }
264 
265 
268 (
269  const scalarField& Tp
270 ) const
271 {
272  // Get kappa from relevant thermo
274 
275  // Optionally modify with explicit resistance
276  if (thicknessLayer_ || thicknessLayers_.size())
277  {
278  scalarField KDelta(tk*patch().deltaCoeffs());
279 
280  // Harmonic averaging of kappa*deltaCoeffs
281  {
282  KDelta = 1.0/KDelta;
283  if (thicknessLayer_)
284  {
285  const scalar t = db().time().timeOutputValue();
286  KDelta +=
287  thicknessLayer_().value(t)
288  /kappaLayer_().value(t);
289  }
290  if (thicknessLayers_.size())
291  {
292  forAll(thicknessLayers_, iLayer)
293  {
294  KDelta += thicknessLayers_[iLayer]/kappaLayers_[iLayer];
295  }
296  }
297  KDelta = 1.0/KDelta;
298  }
299 
300  // Update kappa from KDelta
301  tk = KDelta/patch().deltaCoeffs();
302  }
303 
304  return tk;
305 }
306 
307 
309 {
310  if (updated())
311  {
312  return;
313  }
314 
315  // Since we're inside initEvaluate/evaluate there might be processor
316  // comms underway. Change the tag we use.
317  const int oldTag = UPstream::msgType();
318  UPstream::msgType() = oldTag+1;
319 
320  // Get the coupling information from the mappedPatchBase
321  const mappedPatchBase& mpp =
323  (
324  patch(),
325  this->internalField()
326  );
327 
328  const scalarField& Tp = *this;
329  const scalarField kappaTp(kappa(Tp));
330  const tmp<scalarField> myKDelta = kappaTp*patch().deltaCoeffs();
331 
332 
333  scalarField nbrIntFld;
334  scalarField nbrKDelta;
335  if (mpp.sameWorld())
336  {
337  // Same world so lookup
338  const auto& nbrMesh = refCast<const fvMesh>(mpp.sampleMesh());
339  const label nbrPatchID = mpp.samplePolyPatch().index();
340  const auto& nbrPatch = nbrMesh.boundary()[nbrPatchID];
341 
343  nbrField =
344  refCast
345  <
347  >
348  (
349  nbrPatch.lookupPatchField<volScalarField, scalar>
350  (
351  TnbrName_
352  )
353  );
354 
355  // Swap to obtain full local values of neighbour K*delta
356  nbrIntFld = nbrField.patchInternalField();
357  nbrKDelta = nbrField.kappa(nbrField)*nbrPatch.deltaCoeffs();
358  }
359  else
360  {
361  // Different world so use my region,patch. Distribution below will
362  // do the reordering.
363  nbrIntFld = patchInternalField();
364  nbrKDelta = myKDelta.ref();
365  }
366  distribute(this->internalField().name() + "_value", nbrIntFld);
367  distribute(this->internalField().name() + "_weights", nbrKDelta);
368 
369 
370  // Both sides agree on
371  // - temperature : (myKDelta*fld + nbrKDelta*nbrFld)/(myKDelta+nbrKDelta)
372  // - gradient : (temperature-fld)*delta
373  // We've got a degree of freedom in how to implement this in a mixed bc.
374  // (what gradient, what fixedValue and mixing coefficient)
375  // Two reasonable choices:
376  // 1. specify above temperature on one side (preferentially the high side)
377  // and above gradient on the other. So this will switch between pure
378  // fixedvalue and pure fixedgradient
379  // 2. specify gradient and temperature such that the equations are the
380  // same on both sides. This leads to the choice of
381  // - refGradient = zero gradient
382  // - refValue = neighbour value
383  // - mixFraction = nbrKDelta / (nbrKDelta + myKDelta())
384 
385  this->refValue() = nbrIntFld;
386  this->refGrad() = 0.0;
387  this->valueFraction() = nbrKDelta/(nbrKDelta + myKDelta());
388 
389  mixedFvPatchScalarField::updateCoeffs();
390 
391  if (debug)
392  {
393  scalar Q = gSum(kappaTp*patch().magSf()*snGrad());
394 
395  Info<< patch().boundaryMesh().mesh().name() << ':'
396  << patch().name() << ':'
397  << this->internalField().name() << " <- "
398  << mpp.sampleRegion() << ':'
399  << mpp.samplePatch() << ':'
400  << this->internalField().name() << " :"
401  << " heat transfer rate:" << Q
402  << " walltemperature "
403  << " min:" << gMin(*this)
404  << " max:" << gMax(*this)
405  << " avg:" << gAverage(*this)
406  << endl;
407  }
409  // Restore tag
410  UPstream::msgType() = oldTag;
411 }
412 
413 
415 (
416  fvMatrix<scalar>& m,
417  const label iMatrix,
418  const direction cmpt
419 )
420 {
422  << "This BC does not support energy coupling "
423  << "Use compressible::turbulentTemperatureRadCoupledMixed "
424  << "which has more functionalities and it can handle "
425  << "the assemble coupled option for energy. "
426  << abort(FatalError);
427 }
428 
429 
431 turbulentTemperatureCoupledBaffleMixedFvPatchScalarField::coeffs
432 (
433  fvMatrix<scalar>& matrix,
434  const Field<scalar>& coeffs,
435  const label mat
436 ) const
437 {
439  << "This BC does not support energy coupling "
440  << "Use compressible::turbulentTemperatureRadCoupledMixed "
441  << "which has more functionalities and it can handle "
442  << "the assemble coupled option for energy. "
443  << abort(FatalError);
444  /*
445  const label index(this->patch().index());
446 
447  const label nSubFaces(matrix.lduMesh().cellBoundMap()[mat][index].size());
448 
449  Field<scalar> mapCoeffs(nSubFaces, Zero);
450 
451  label subFaceI = 0;
452  forAll(*this, faceI)
453  {
454 
455  }
456  */
457  return tmp<Field<scalar>>(new Field<scalar>());
458 }
459 
460 
462 (
463  Ostream& os
464 ) const
465 {
467  os.writeEntry("Tnbr", TnbrName_);
468  if (thicknessLayer_)
469  {
470  thicknessLayer_().writeData(os);
471  kappaLayer_().writeData(os);
472  }
473  if (thicknessLayers_.size())
474  {
475  thicknessLayers_.writeEntry("thicknessLayers", os);
476  kappaLayers_.writeEntry("kappaLayers", os);
477  }
480 }
481 
482 
483 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
484 
486 (
488  turbulentTemperatureCoupledBaffleMixedFvPatchScalarField
489 );
490 
491 
492 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
493 
494 } // End namespace compressible
495 } // End namespace Foam
496 
497 
498 // ************************************************************************* //
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:118
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
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:120
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:578
Type gMin(const FieldField< Field, Type > &f)
virtual void write(Ostream &os) const
Write.
Type & refCast(U &obj)
A dynamic_cast (for references) that generates FatalError on failed casts, uses the virtual type() me...
Definition: typeInfo.H:151
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:487
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:68
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:806
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:413
makePatchTypeField(fvPatchScalarField, alphatFixedDmdtWallBoilingWallFunctionFvPatchScalarField)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:84
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Component-wise write of a value (N times)
Type gSum(const FieldField< Field, Type > &f)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:52
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:55
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)
virtual void operator=(const UList< scalar > &)
Definition: fvPatchField.C:352
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:31
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.