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),
53  (
54  patch(),
55  "undefined",
56  "undefined",
57  "undefined-K",
58  "undefined-alpha"
59  ),
61  (
63  *this
64  ),
65  TnbrName_("undefined-Tnbr")
66 {
67  this->refValue() = 0.0;
68  this->refGrad() = 0.0;
69  this->valueFraction() = 1.0;
70 }
71 
72 
75 (
77  const fvPatch& p,
79  const fvPatchFieldMapper& mapper
80 )
81 :
82  mixedFvPatchScalarField(ptf, p, iF, mapper),
85  (
87  *this,
88  ptf
89  ),
90  TnbrName_(ptf.TnbrName_),
91  thicknessLayers_(ptf.thicknessLayers_),
92  thicknessLayer_(ptf.thicknessLayer_.clone(p.patch())),
93  kappaLayers_(ptf.kappaLayers_),
94  kappaLayer_(ptf.kappaLayer_.clone(p.patch()))
95 {}
96 
97 
100 (
101  const fvPatch& p,
103  const dictionary& dict
104 )
105 :
106  mixedFvPatchScalarField(p, iF),
109  (
111  *this,
112  dict
113  ),
114  TnbrName_(dict.get<word>("Tnbr"))
115 {
116  if (!isA<mappedPatchBase>(this->patch().patch()))
117  {
119  << "' not type '" << mappedPatchBase::typeName << "'"
120  << "\n for patch " << p.name()
121  << " of field " << internalField().name()
122  << " in file " << internalField().objectPath()
123  << exit(FatalError);
124  }
125 
127  << "This BC has been superseded by "
128  << "compressible::turbulentTemperatureRadCoupledMixed "
129  << "which has more functionalities and it can handle "
130  << "the assemble coupled option for energy. "
131  << endl;
132 
133  // Read list of layers
134  if (dict.readIfPresent("thicknessLayers", thicknessLayers_))
135  {
136  dict.readEntry("kappaLayers", kappaLayers_);
137  }
138  // Read single additional PatchFunction1
139  thicknessLayer_ = PatchFunction1<scalar>::NewIfPresent
140  (
141  p.patch(),
142  "thicknessLayer",
143  dict
144  );
146  (
147  p.patch(),
148  "kappaLayer",
149  dict
150  );
151 
152 
153  fvPatchScalarField::operator=(scalarField("value", dict, p.size()));
154 
155  if (dict.found("refValue"))
156  {
157  // Full restart
158  refValue() = scalarField("refValue", dict, p.size());
159  refGrad() = scalarField("refGradient", dict, p.size());
160  valueFraction() = scalarField("valueFraction", dict, p.size());
161  }
162  else
163  {
164  // Start from user entered data. Assume fixedValue.
165  refValue() = *this;
166  refGrad() = 0.0;
167  valueFraction() = 1.0;
168  }
169 
170 // This blocks (crashes) with more than two worlds!
171 //
183 }
184 
185 
188 (
191 )
192 :
193  mixedFvPatchScalarField(wtcsf, iF),
194  temperatureCoupledBase(patch(), wtcsf),
196  (
198  *this,
199  wtcsf
200  ),
201  TnbrName_(wtcsf.TnbrName_),
202  thicknessLayers_(wtcsf.thicknessLayers_),
203  thicknessLayer_(wtcsf.thicknessLayer_.clone(patch().patch())),
204  kappaLayers_(wtcsf.kappaLayers_),
205  kappaLayer_(wtcsf.kappaLayer_.clone(patch().patch()))
206 {}
207 
208 
211 (
213 )
214 :
215  mixedFvPatchScalarField(wtcsf),
216  temperatureCoupledBase(patch(), wtcsf),
218  (
219  mappedPatchFieldBase<scalar>::mapper(patch(), wtcsf.internalField()),
220  *this,
221  wtcsf
222  ),
223  TnbrName_(wtcsf.TnbrName_),
224  thicknessLayers_(wtcsf.thicknessLayers_),
225  thicknessLayer_(wtcsf.thicknessLayer_.clone(patch().patch())),
226  kappaLayers_(wtcsf.kappaLayers_),
227  kappaLayer_(wtcsf.kappaLayer_.clone(patch().patch()))
228 {}
229 
230 
231 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
232 
234 (
235  const fvPatchFieldMapper& mapper
236 )
237 {
238  mixedFvPatchScalarField::autoMap(mapper);
240  //mappedPatchFieldBase<scalar>::autoMap(mapper);
241  if (thicknessLayer_)
242  {
243  thicknessLayer_().autoMap(mapper);
244  kappaLayer_().autoMap(mapper);
245  }
246 }
247 
248 
250 (
251  const fvPatchField<scalar>& ptf,
252  const labelList& addr
253 )
254 {
255  mixedFvPatchScalarField::rmap(ptf, addr);
256 
258  refCast
259  <
261  >(ptf);
262 
263  temperatureCoupledBase::rmap(tiptf, addr);
264  //mappedPatchFieldBase<scalar>::rmap(ptf, addr);
265  if (thicknessLayer_)
266  {
267  thicknessLayer_().rmap(tiptf.thicknessLayer_(), addr);
268  kappaLayer_().rmap(tiptf.kappaLayer_(), addr);
269  }
270 }
271 
272 
275 (
276  const scalarField& Tp
277 ) const
278 {
279  // Get kappa from relevant thermo
281 
282  // Optionally modify with explicit resistance
283  if (thicknessLayer_ || thicknessLayers_.size())
284  {
285  scalarField KDelta(tk*patch().deltaCoeffs());
286 
287  // Harmonic averaging of kappa*deltaCoeffs
288  {
289  KDelta = 1.0/KDelta;
290  if (thicknessLayer_)
291  {
292  const scalar t = db().time().timeOutputValue();
293  KDelta +=
294  thicknessLayer_().value(t)
295  /kappaLayer_().value(t);
296  }
297  if (thicknessLayers_.size())
298  {
299  forAll(thicknessLayers_, iLayer)
300  {
301  KDelta += thicknessLayers_[iLayer]/kappaLayers_[iLayer];
302  }
303  }
304  KDelta = 1.0/KDelta;
305  }
306 
307  // Update kappa from KDelta
308  tk = KDelta/patch().deltaCoeffs();
309  }
310 
311  return tk;
312 }
313 
314 
316 {
317  if (updated())
318  {
319  return;
320  }
321 
322  // Since we're inside initEvaluate/evaluate there might be processor
323  // comms underway. Change the tag we use.
324  const int oldTag = UPstream::msgType();
325  UPstream::msgType() = oldTag+1;
326 
327  // Get the coupling information from the mappedPatchBase
328  const mappedPatchBase& mpp =
330  (
331  patch(),
332  this->internalField()
333  );
334 
335  const scalarField& Tp = *this;
336  const scalarField kappaTp(kappa(Tp));
337  const tmp<scalarField> myKDelta = kappaTp*patch().deltaCoeffs();
338 
339 
340  scalarField nbrIntFld;
341  scalarField nbrKDelta;
342  if (mpp.sameWorld())
343  {
344  // Same world so lookup
345  const auto& nbrMesh = refCast<const fvMesh>(mpp.sampleMesh());
346  const label nbrPatchID = mpp.samplePolyPatch().index();
347  const auto& nbrPatch = nbrMesh.boundary()[nbrPatchID];
348 
350  nbrField =
351  refCast
352  <
354  >
355  (
356  nbrPatch.lookupPatchField<volScalarField, scalar>
357  (
358  TnbrName_
359  )
360  );
361 
362  // Swap to obtain full local values of neighbour K*delta
363  nbrIntFld = nbrField.patchInternalField();
364  nbrKDelta = nbrField.kappa(nbrField)*nbrPatch.deltaCoeffs();
365  }
366  else
367  {
368  // Different world so use my region,patch. Distribution below will
369  // do the reordering.
370  nbrIntFld = patchInternalField();
371  nbrKDelta = myKDelta.ref();
372  }
373  distribute(this->internalField().name() + "_value", nbrIntFld);
374  distribute(this->internalField().name() + "_weights", nbrKDelta);
375 
376 
377  // Both sides agree on
378  // - temperature : (myKDelta*fld + nbrKDelta*nbrFld)/(myKDelta+nbrKDelta)
379  // - gradient : (temperature-fld)*delta
380  // We've got a degree of freedom in how to implement this in a mixed bc.
381  // (what gradient, what fixedValue and mixing coefficient)
382  // Two reasonable choices:
383  // 1. specify above temperature on one side (preferentially the high side)
384  // and above gradient on the other. So this will switch between pure
385  // fixedvalue and pure fixedgradient
386  // 2. specify gradient and temperature such that the equations are the
387  // same on both sides. This leads to the choice of
388  // - refGradient = zero gradient
389  // - refValue = neighbour value
390  // - mixFraction = nbrKDelta / (nbrKDelta + myKDelta())
391 
392  this->refValue() = nbrIntFld;
393  this->refGrad() = 0.0;
394  this->valueFraction() = nbrKDelta/(nbrKDelta + myKDelta());
395 
396  mixedFvPatchScalarField::updateCoeffs();
397 
398  if (debug)
399  {
400  scalar Q = gSum(kappaTp*patch().magSf()*snGrad());
401 
402  Info<< patch().boundaryMesh().mesh().name() << ':'
403  << patch().name() << ':'
404  << this->internalField().name() << " <- "
405  << mpp.sampleRegion() << ':'
406  << mpp.samplePatch() << ':'
407  << this->internalField().name() << " :"
408  << " heat transfer rate:" << Q
409  << " walltemperature "
410  << " min:" << gMin(*this)
411  << " max:" << gMax(*this)
412  << " avg:" << gAverage(*this)
413  << endl;
414  }
415 
416  // Restore tag
417  UPstream::msgType() = oldTag;
418 }
419 
420 
422 (
423  fvMatrix<scalar>& m,
424  const label iMatrix,
425  const direction cmpt
426 )
427 {
429  << "This BC does not support energy coupling "
430  << "Use compressible::turbulentTemperatureRadCoupledMixed "
431  << "which has more functionalities and it can handle "
432  << "the assemble coupled option for energy. "
433  << abort(FatalError);
434 }
435 
436 
438 turbulentTemperatureCoupledBaffleMixedFvPatchScalarField::coeffs
439 (
440  fvMatrix<scalar>& matrix,
441  const Field<scalar>& coeffs,
442  const label mat
443 ) const
444 {
446  << "This BC does not support energy coupling "
447  << "Use compressible::turbulentTemperatureRadCoupledMixed "
448  << "which has more functionalities and it can handle "
449  << "the assemble coupled option for energy. "
450  << abort(FatalError);
451  /*
452  const label index(this->patch().index());
453 
454  const label nSubFaces(matrix.lduMesh().cellBoundMap()[mat][index].size());
455 
456  Field<scalar> mapCoeffs(nSubFaces, Zero);
457 
458  label subFaceI = 0;
459  forAll(*this, faceI)
460  {
461 
462  }
463  */
464  return tmp<Field<scalar>>(new Field<scalar>());
465 }
466 
467 
469 (
470  Ostream& os
471 ) const
472 {
474  os.writeEntry("Tnbr", TnbrName_);
475  if (thicknessLayer_)
476  {
477  thicknessLayer_().writeData(os);
478  kappaLayer_().writeData(os);
479  }
480  if (thicknessLayers_.size())
481  {
482  thicknessLayers_.writeEntry("thicknessLayers", os);
483  kappaLayers_.writeEntry("kappaLayers", os);
484  }
487 }
488 
489 
490 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
491 
493 (
496 );
497 
498 
499 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
500 
501 } // End namespace compressible
502 } // End namespace Foam
503 
504 
505 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::fvPatchField< scalar >
volFields.H
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
Foam::compressible::turbulentTemperatureCoupledBaffleMixedFvPatchScalarField::manipulateMatrix
virtual void manipulateMatrix(fvMatrix< scalar > &m, const label iMatrix, const direction cmpt)
Manipulate matrix.
Definition: turbulentTemperatureCoupledBaffleMixedFvPatchScalarField.C:422
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::gAverage
Type gAverage(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:604
Foam::mappedPatchBase::samplePolyPatch
const polyPatch & samplePolyPatch() const
Get the patch on the region.
Definition: mappedPatchBase.C:1662
Foam::compressible::turbulentTemperatureCoupledBaffleMixedFvPatchScalarField::turbulentTemperatureCoupledBaffleMixedFvPatchScalarField
turbulentTemperatureCoupledBaffleMixedFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: turbulentTemperatureCoupledBaffleMixedFvPatchScalarField.C:46
Foam::temperatureCoupledBase
Common functions used in temperature coupled boundaries.
Definition: temperatureCoupledBase.H:135
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
fvPatchFieldMapper.H
Foam::mappedPatchBase::sameWorld
bool sameWorld() const
Is sample world the local world?
Definition: mappedPatchBaseI.H:169
Foam::gSum
Type gSum(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:594
Foam::mappedPatchBase
Determines a mapping between patch face centres and mesh cell or face centres and processors they're ...
Definition: mappedPatchBase.H:112
Foam::temperatureCoupledBase::autoMap
virtual void autoMap(const fvPatchFieldMapper &)=0
Map (and resize as needed) from self given a mapping object.
Definition: temperatureCoupledBase.C:172
Foam::compressible::makePatchTypeField
makePatchTypeField(fvPatchScalarField, alphatFixedDmdtWallBoilingWallFunctionFvPatchScalarField)
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::compressible::turbulentTemperatureCoupledBaffleMixedFvPatchScalarField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: turbulentTemperatureCoupledBaffleMixedFvPatchScalarField.C:315
Foam::mappedPatchBase::sampleMesh
const polyMesh & sampleMesh() const
Get the region mesh.
Definition: mappedPatchBase.C:1649
Foam::compressible::turbulentTemperatureCoupledBaffleMixedFvPatchScalarField::kappa
virtual tmp< scalarField > kappa(const scalarField &Tp) const
Definition: turbulentTemperatureCoupledBaffleMixedFvPatchScalarField.C:275
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:227
Foam::Field< scalar >
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::mappedPatchFieldBase< scalar >::distribute
void distribute(const word &fieldName, Field< T > &newValues) const
Definition: mappedPatchFieldBase.C:524
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
Foam::compressible::turbulentTemperatureCoupledBaffleMixedFvPatchScalarField
Mixed boundary condition for temperature, to be used for heat-transfer on back-to-back baffles....
Definition: turbulentTemperatureCoupledBaffleMixedFvPatchScalarField.H:140
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Foam::compressible::turbulentTemperatureCoupledBaffleMixedFvPatchScalarField::write
virtual void write(Ostream &os) const
Write.
Definition: turbulentTemperatureCoupledBaffleMixedFvPatchScalarField.C:469
Foam::compressible::turbulentTemperatureCoupledBaffleMixedFvPatchScalarField::rmap
virtual void rmap(const fvPatchField< scalar > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Definition: turbulentTemperatureCoupledBaffleMixedFvPatchScalarField.C:250
compressible
bool compressible
Definition: pEqn.H:2
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::mappedPatchFieldBase< scalar >
Foam::mappedPatchFieldBase::mapper
static const mappedPatchBase & mapper(const fvPatch &p, const DimensionedField< Type, volMesh > &iF)
Check that patch is of correct type.
Definition: mappedPatchFieldBase.C:959
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
os
OBJstream os(runTime.globalPath()/outputName)
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::temperatureCoupledBase::rmap
virtual void rmap(const fvPatchField< scalar > &, const labelList &)=0
Reverse map the given fvPatchField onto this fvPatchField.
Definition: temperatureCoupledBase.C:188
Foam::patchIdentifier::index
label index() const noexcept
The index of this patch in the boundaryMesh.
Definition: patchIdentifier.H:147
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::refCast
To & refCast(From &r)
Reference type cast template function.
Definition: typeInfo.H:131
Foam::mappedPatchBase::sampleRegion
const word & sampleRegion() const
Region to sample.
Definition: mappedPatchBaseI.H:42
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::compressible::turbulentTemperatureCoupledBaffleMixedFvPatchScalarField::autoMap
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Definition: turbulentTemperatureCoupledBaffleMixedFvPatchScalarField.C:234
turbulentTemperatureCoupledBaffleMixedFvPatchScalarField.H
Foam::mappedPatchFieldBase::write
virtual void write(Ostream &os) const
Write.
Definition: mappedPatchFieldBase.C:1007
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::UPstream::msgType
static int & msgType() noexcept
Message tag of standard messages.
Definition: UPstream.H:540
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::List< label >
Foam::temperatureCoupledBase::kappa
virtual tmp< scalarField > kappa(const scalarField &Tp) const
Given patch temperature calculate corresponding K field.
Definition: temperatureCoupledBase.C:210
Foam::direction
uint8_t direction
Definition: direction.H:52
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:68
Foam::vtk::write
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Component-wise write of a value (N times)
Definition: foamVtkOutputTemplates.C:36
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
Foam::temperatureCoupledBase::write
void write(Ostream &os) const
Write.
Definition: temperatureCoupledBase.C:512
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::gMin
Type gMin(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:593
Foam::fvPatchField< scalar >::operator=
virtual void operator=(const UList< scalar > &)
Definition: fvPatchField.C:404
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:328
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
Foam::PatchFunction1::NewIfPresent
static autoPtr< PatchFunction1< Type > > NewIfPresent(const polyPatch &pp, const word &entryName, const dictionary &dict, const bool faceValues=true)
An optional selector.
Definition: PatchFunction1New.C:214
mappedPatchBase.H
Foam::mappedPatchBase::samplePatch
const word & samplePatch() const
Patch (only if NEARESTPATCHFACE)
Definition: mappedPatchBaseI.H:68