coupledTemperatureFvPatchScalarField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2024 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
28 #include "volFields.H"
29 #include "fieldMapper.H"
30 #include "mappedFvPatchBaseBase.H"
32 
33 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
34 
36 (
38  tmp<scalarField>& sumKappaTByDelta,
39  tmp<scalarField>& sumKappaByDeltaNbr,
40  scalarField& sumq,
41  tmp<scalarField>& qByKappa
42 ) const
43 {
44  const thermophysicalTransportModel& ttm =
45  patch().boundaryMesh().mesh()
46  .lookupType<thermophysicalTransportModel>();
47 
48  kappa = ttm.kappaEff(patch().index());
49 
50  qByKappa = sumq/kappa();
51 
52  sumq = 0;
53 
54  tmp<scalarField> qCorr(ttm.qCorr(patch().index()));
55 
56  if (qCorr.valid())
57  {
58  sumq += qCorr;
59  }
60 }
61 
62 
64 (
65  tmp<scalarField>& sumKappaTByDeltaNbr,
66  tmp<scalarField>& sumKappaByDeltaNbr,
67  tmp<scalarField>& qNbr
68 ) const
69 {
70  const thermophysicalTransportModel& ttm =
71  patch().boundaryMesh().mesh()
72  .lookupType<thermophysicalTransportModel>();
73 
74  sumKappaByDeltaNbr = ttm.kappaEff(patch().index())*patch().deltaCoeffs();
75 
76  sumKappaTByDeltaNbr = sumKappaByDeltaNbr()*patchInternalField();
77 
78  qNbr = ttm.qCorr(patch().index());
79 }
80 
81 
83 (
84  tmp<scalarField>& TrefNbr,
85  tmp<scalarField>& qNbr
86 ) const
87 {
88  const thermophysicalTransportModel& ttm =
89  patch().boundaryMesh().mesh()
90  .lookupType<thermophysicalTransportModel>();
91 
92  const fvPatchScalarField& Tp =
93  patch().lookupPatchField<volScalarField, scalar>
94  (
95  internalField().name()
96  );
97 
98  TrefNbr = Tp;
99 
100  qNbr = ttm.qCorr(patch().index());
101 }
102 
103 
105 (
106  tmp<scalarField>& result,
107  const tmp<scalarField>& field
108 ) const
109 {
110  if (result.valid())
111  {
112  result.ref() += field;
113  }
114  else
115  {
116  if (field.isTmp())
117  {
118  result = field;
119  }
120  else
121  {
122  result = field().clone();
123  }
124  }
125 }
126 
127 
128 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
129 
132 (
133  const fvPatch& p,
135  const dictionary& dict
136 )
137 :
138  mixedFvPatchScalarField(p, iF, dict, false),
139  TnbrName_(dict.lookupOrDefault<word>("Tnbr", "T")),
140  qrNbrName_(dict.lookupOrDefault<word>("qrNbr", "none")),
141  qrName_(dict.lookupOrDefault<word>("qr", "none")),
142  thicknessLayers_(0),
143  kappaLayers_(0),
144  qs_(),
145  Qs_(0),
146  wallKappaByDelta_(0)
147 {
149  (
150  *this,
151  iF,
152  dict,
154  );
155 
156  if (dict.found("thicknessLayers"))
157  {
158  dict.lookup("thicknessLayers") >> thicknessLayers_;
159  dict.lookup("kappaLayers") >> kappaLayers_;
160 
161  if (thicknessLayers_.size() > 0)
162  {
163  // Calculate effective thermal resistance by harmonic averaging
164  forAll(thicknessLayers_, i)
165  {
166  wallKappaByDelta_ += thicknessLayers_[i]/kappaLayers_[i];
167  }
168  wallKappaByDelta_ = 1/wallKappaByDelta_;
169  }
170  }
171 
172  if (dict.found("qs"))
173  {
174  if (dict.found("Qs"))
175  {
177  << "Either qs or Qs should be specified, not both"
178  << exit(FatalIOError);
179  }
180 
181  qs_ = new scalarField("qs", dimPower/dimTime, dict, p.size());
182  }
183  else if (dict.found("Qs"))
184  {
185  Qs_ = dict.lookup<scalar>("Qs");
186  qs_ = new scalarField(p.size(), Qs_/gSum(patch().magSf()));
187  }
188 
190  (
191  scalarField("value", iF.dimensions(), dict, p.size())
192  );
193 
194  if (dict.found("refValue"))
195  {
196  // Full restart
197  refValue() = scalarField("refValue", iF.dimensions(), dict, p.size());
198  refGrad() =
200  (
201  "refGradient",
202  iF.dimensions()/dimLength,
203  dict,
204  p.size()
205  );
206  valueFraction() =
207  scalarField("valueFraction", unitFraction, dict, p.size());
208  }
209  else
210  {
211  // Start from user entered data. Assume fixedValue.
212  refValue() = *this;
213  refGrad() = 0;
214  valueFraction() = 1;
215  }
216 }
217 
218 
221 (
223  const fvPatch& p,
225  const fieldMapper& mapper
226 )
227 :
228  mixedFvPatchScalarField(psf, p, iF, mapper),
229  TnbrName_(psf.TnbrName_),
230  qrNbrName_(psf.qrNbrName_),
231  qrName_(psf.qrName_),
232  thicknessLayers_(psf.thicknessLayers_),
233  kappaLayers_(psf.kappaLayers_),
234  qs_(psf.qs_.valid() ? new scalarField(p.size()) : nullptr),
235  Qs_(psf.Qs_),
236  wallKappaByDelta_(psf.wallKappaByDelta_)
237 {
238  map(psf, mapper);
239 }
240 
241 
244 (
247 )
248 :
249  mixedFvPatchScalarField(psf, iF),
250  TnbrName_(psf.TnbrName_),
251  qrNbrName_(psf.qrNbrName_),
252  qrName_(psf.qrName_),
253  thicknessLayers_(psf.thicknessLayers_),
254  kappaLayers_(psf.kappaLayers_),
255  qs_(psf.qs_, false),
256  Qs_(psf.Qs_),
257  wallKappaByDelta_(psf.wallKappaByDelta_)
258 {}
259 
260 
261 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
262 
264 (
266  const fieldMapper& mapper
267 )
268 {
269  // Unmapped faces are considered zero-gradient/adiabatic
270  mapper(*this, ptf, [&](){ return patchInternalField(); });
271  mapper(refValue(), ptf.refValue(), [&](){ return patchInternalField(); });
272  mapper(refGrad(), ptf.refGrad(), scalar(0));
273  mapper(valueFraction(), ptf.valueFraction(), scalar(0));
274 
275  // Map the heat flux, if present
276  if (ptf.qs_.valid())
277  {
278  mapper(qs_(), ptf.qs_());
279  }
280 }
281 
282 
284 (
285  const fvPatchScalarField& ptf,
286  const fieldMapper& mapper
287 )
288 {
289  map(refCast<const coupledTemperatureFvPatchScalarField>(ptf), mapper);
290 }
291 
292 
294 (
295  const fvPatchScalarField& ptf
296 )
297 {
298  mixedFvPatchScalarField::reset(ptf);
299 
301  refCast<const coupledTemperatureFvPatchScalarField>(ptf);
302 
303  if (tiptf.qs_.valid())
304  {
305  qs_().reset(tiptf.qs_());
306  }
307 }
308 
309 
311 {
312  if (updated())
313  {
314  return;
315  }
316 
317  // Since we're inside initEvaluate/evaluate there might be processor
318  // comms underway. Change the tag we use.
319  int oldTag = UPstream::msgType();
320  UPstream::msgType() = oldTag + 1;
321 
322  // Get the mapper and the neighbouring patch
323  const mappedFvPatchBaseBase& mapper =
325  const fvPatch& patchNbr = mapper.nbrFvPatch();
326 
327  const fvPatchScalarField& TpNbr =
328  patchNbr.lookupPatchField<volScalarField, scalar>(TnbrName_);
329 
330  if (!isA<coupledTemperatureFvPatchScalarField>(TpNbr))
331  {
333  << "Patch field for " << internalField().name() << " on "
334  << this->patch().name() << " is of type "
335  << coupledTemperatureFvPatchScalarField::typeName
336  << endl << "The neighbouring patch field "
337  << internalField().name() << " on "
338  << patchNbr.name() << " is required to be the same, but is "
339  << "currently of type " << TpNbr.type() << exit(FatalError);
340  }
341 
342  const coupledTemperatureFvPatchScalarField& coupledTemperatureNbr =
343  refCast<const coupledTemperatureFvPatchScalarField>(TpNbr);
344 
345  scalarField sumq(size(), 0);
346 
347  if (qs_.valid())
348  {
349  sumq += qs_();
350  }
351 
352  if (qrName_ != "none")
353  {
354  sumq += patch().lookupPatchField<volScalarField, scalar>(qrName_);
355  }
356 
357  if (qrNbrName_ != "none")
358  {
359  sumq += mapper.fromNeighbour
360  (
361  patchNbr.lookupPatchField<volScalarField, scalar>(qrNbrName_)
362  );
363  }
364 
366  tmp<scalarField> sumKappaTByDelta;
367  tmp<scalarField> sumKappaByDelta;
368  tmp<scalarField> qByKappa;
369 
370  // q = alpha.this*sumq
371  getThis(kappa, sumKappaTByDelta, sumKappaByDelta, sumq, qByKappa);
372 
373  // Add neighbour contributions
374  {
375  tmp<scalarField> sumKappaTByDeltaNbr;
376  tmp<scalarField> sumKappaByDeltaNbr;
377  tmp<scalarField> qNbr;
378 
379  if (wallKappaByDelta_ == 0)
380  {
381  coupledTemperatureNbr.getNbr
382  (
383  sumKappaTByDeltaNbr,
384  sumKappaByDeltaNbr,
385  qNbr
386  );
387 
388  add(sumKappaTByDelta, mapper.fromNeighbour(sumKappaTByDeltaNbr));
389  add(sumKappaByDelta, mapper.fromNeighbour(sumKappaByDeltaNbr));
390  }
391  else
392  {
393  // Get the neighbour wall temperature and flux correction
394  tmp<scalarField> TwNbr;
395  coupledTemperatureNbr.getNbr(TwNbr, qNbr);
396 
397  add(sumKappaByDelta, scalarField(size(), wallKappaByDelta_));
398  add
399  (
400  sumKappaTByDelta,
401  wallKappaByDelta_*mapper.fromNeighbour(TwNbr)
402  );
403  }
404 
405  if (qNbr.valid())
406  {
407  sumq += mapper.fromNeighbour(qNbr);
408  }
409  }
410 
411  this->valueFraction() =
412  sumKappaByDelta()/(kappa()*patch().deltaCoeffs() + sumKappaByDelta());
413 
414  this->refValue() = (sumKappaTByDelta() + sumq)/sumKappaByDelta();
415 
416  this->refGrad() = qByKappa;
417 
418  mixedFvPatchScalarField::updateCoeffs();
419 
420  if (debug)
421  {
422  const scalar Q = gSum(kappa()*patch().magSf()*snGrad());
423 
424  Info<< patch().boundaryMesh().mesh().name() << ':'
425  << patch().name() << ':'
426  << this->internalField().name() << " <- "
427  << mapper.nbrMesh().name() << ':'
428  << patchNbr.name() << ':'
429  << this->internalField().name() << " :"
430  << " heat transfer rate:" << Q
431  << " walltemperature "
432  << " min:" << gMin(*this)
433  << " max:" << gMax(*this)
434  << " avg:" << gAverage(*this)
435  << endl;
436  }
437 
438  // Restore tag
439  UPstream::msgType() = oldTag;
440 }
441 
442 
444 (
445  Ostream& os
446 ) const
447 {
449 
450  writeEntryIfDifferent<word>(os, "Tnbr", "T", TnbrName_);
451  writeEntryIfDifferent<word>(os, "qrNbr", "none", qrNbrName_);
452  writeEntryIfDifferent<word>(os, "qr", "none", qrName_);
453 
454  if (Qs_ != 0)
455  {
456  writeEntry(os, "Qs", Qs_);
457  }
458  else if (qs_.valid())
459  {
460  writeEntry(os, "qs", qs_());
461  }
462 
463  if (thicknessLayers_.size())
464  {
465  writeEntry(os, "thicknessLayers", thicknessLayers_);
466  writeEntry(os, "kappaLayers", kappaLayers_);
467  }
468 }
469 
470 
471 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
472 
473 namespace Foam
474 {
476  (
479  );
480 
482  (
485  patchMapper,
486  turbulentTemperatureCoupledBaffleMixed,
487  "compressible::turbulentTemperatureCoupledBaffleMixed"
488  );
489 
491  (
494  dictionary,
495  turbulentTemperatureCoupledBaffleMixed,
496  "compressible::turbulentTemperatureCoupledBaffleMixed"
497  );
498 
499 
501  (
504  patchMapper,
505  turbulentTemperatureRadCoupledMixed,
506  "compressible::turbulentTemperatureRadCoupledMixed"
507  );
508 
510  (
513  dictionary,
514  turbulentTemperatureRadCoupledMixed,
515  "compressible::turbulentTemperatureRadCoupledMixed"
516  );
517 }
518 
519 
520 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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...
const dimensionSet & dimensions() const
Return dimensions.
Generic GeometricField class.
virtual const fileName & name() const
Return the name of the stream.
Definition: IOstream.H:294
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
static int & msgType()
Message tag of standard messages.
Definition: UPstream.H:476
Mixed boundary condition for temperature, to be used for heat-transfer with another region in a CHT c...
virtual void getThis(tmp< scalarField > &kappa, tmp< scalarField > &sumKappaTByDelta, tmp< scalarField > &sumKappaByDelta, scalarField &qTot, tmp< scalarField > &qByKappa) const
Get the patch kappa, kappa*Tc/delta and kappa/delta and also the.
virtual void reset(const fvPatchScalarField &)
Reset the fvPatchField to the given fvPatchField.
void map(const coupledTemperatureFvPatchScalarField &, const fieldMapper &)
Map the given fvPatchField onto this fvPatchField.
virtual void getNbr(tmp< scalarField > &sumKappaTByDeltaNbr, tmp< scalarField > &sumKappaByDeltaNbr, tmp< scalarField > &qNbr) const
Get the neighbour patch kappa*Tc/delta and kappa/delta.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
coupledTemperatureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
void add(tmp< scalarField > &result, const tmp< scalarField > &field) const
Add field to result which may have not been previously set.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Abstract base class for field mapping.
Definition: fieldMapper.H:48
const word & name() const
Return reference to name.
Definition: fvMesh.H:432
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:88
friend Ostream & operator(Ostream &, const fvPatchField< Type > &)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
virtual const word & name() const
Return name.
Definition: fvPatch.H:126
const GeometricField::Patch & lookupPatchField(const word &name) const
Lookup and return the patchField of the named field from the.
Base class for fv patches that provide mapping between two fv patches.
const fvPatch & nbrFvPatch() const
Get the patch to map from.
static const mappedFvPatchBaseBase & getMap(const fvPatch &patch)
Cast the given fvPatch to a mappedFvPatchBaseBase. Handle errors.
const fvMesh & nbrMesh() const
Get the mesh for the region to map from.
static void validateMapForField(const PatchField &field, const FieldType &iF, const dictionary &context, const label froms=from::any)
Validate that the map exists and is appropriate for the given.
Abstract base class for all fluid and solid thermophysical transport models.
virtual tmp< scalarField > qCorr(const label patchi) const =0
Return the patch heat flux correction [W/m^2].
virtual tmp< volScalarField > kappaEff() const =0
Effective thermal turbulent conductivity.
A class for managing temporary objects.
Definition: tmp.H:55
bool valid() const
Is this temporary object valid,.
Definition: tmpI.H:167
bool isTmp() const
Return true if this is really a temporary object.
Definition: tmpI.H:153
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
A class for handling words, derived from string.
Definition: word.H:62
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:346
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
bool valid(const PtrList< ModelType > &l)
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
tmp< SurfaceField< Type > > snGrad(const VolField< Type > &vf, const word &name)
Definition: fvcSnGrad.C:45
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Type gSum(const FieldField< Field, Type > &f)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
const dimensionSet dimPower
messageStream Info
const dimensionSet dimLength
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
addBackwardCompatibleToRunTimeSelectionTable(fvPatchScalarField, coupledTemperatureFvPatchScalarField, patchMapper, turbulentTemperatureCoupledBaffleMixed, "compressible::turbulentTemperatureCoupledBaffleMixed")
const dimensionSet dimTime
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
void add(FieldField< Field1, typename typeOfSum< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:64
IOerror FatalIOError
Type gAverage(const FieldField< Field, Type > &f)
error FatalError
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Type gMin(const FieldField< Field, Type > &f)
Type gMax(const FieldField< Field, Type > &f)
const unitConversion unitFraction
dictionary dict
volScalarField & p