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-2023 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 "fvPatchFieldMapper.H"
30 #include "mappedPatchBase.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", 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 
189  fvPatchScalarField::operator=(scalarField("value", dict, p.size()));
190 
191  if (dict.found("refValue"))
192  {
193  // Full restart
194  refValue() = scalarField("refValue", dict, p.size());
195  refGrad() = scalarField("refGradient", dict, p.size());
196  valueFraction() = scalarField("valueFraction", dict, p.size());
197  }
198  else
199  {
200  // Start from user entered data. Assume fixedValue.
201  refValue() = *this;
202  refGrad() = 0;
203  valueFraction() = 1;
204  }
205 }
206 
207 
210 (
212  const fvPatch& p,
214  const fvPatchFieldMapper& mapper
215 )
216 :
217  mixedFvPatchScalarField(psf, p, iF, mapper),
218  TnbrName_(psf.TnbrName_),
219  qrNbrName_(psf.qrNbrName_),
220  qrName_(psf.qrName_),
221  thicknessLayers_(psf.thicknessLayers_),
222  kappaLayers_(psf.kappaLayers_),
223  qs_
224  (
225  psf.qs_.valid()
226  ? mapper(psf.qs_()).ptr()
227  : nullptr
228  ),
229  Qs_(psf.Qs_),
230  wallKappaByDelta_(psf.wallKappaByDelta_)
231 {}
232 
233 
236 (
239 )
240 :
241  mixedFvPatchScalarField(psf, iF),
242  TnbrName_(psf.TnbrName_),
243  qrNbrName_(psf.qrNbrName_),
244  qrName_(psf.qrName_),
245  thicknessLayers_(psf.thicknessLayers_),
246  kappaLayers_(psf.kappaLayers_),
247  qs_(psf.qs_, false),
248  Qs_(psf.Qs_),
249  wallKappaByDelta_(psf.wallKappaByDelta_)
250 {}
251 
252 
253 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
254 
256 {
257  if (updated())
258  {
259  return;
260  }
261 
262  // Since we're inside initEvaluate/evaluate there might be processor
263  // comms underway. Change the tag we use.
264  int oldTag = UPstream::msgType();
265  UPstream::msgType() = oldTag + 1;
266 
267  // Get the coupling information from the mappedPatchBase
268  const mappedPatchBase& mpp = mappedPatchBase::getMap(patch().patch());
269  const label patchiNbr = mpp.nbrPolyPatch().index();
270  const fvPatch& patchNbr =
271  refCast<const fvMesh>(mpp.nbrMesh()).boundary()[patchiNbr];
272 
273  const fvPatchScalarField& TpNbr =
274  patchNbr.lookupPatchField<volScalarField, scalar>(TnbrName_);
275 
276  if (!isA<coupledTemperatureFvPatchScalarField>(TpNbr))
277  {
279  << "Patch field for " << internalField().name() << " on "
280  << this->patch().name() << " is of type "
281  << coupledTemperatureFvPatchScalarField::typeName
282  << endl << "The neighbouring patch field "
283  << internalField().name() << " on "
284  << patchNbr.name() << " is required to be the same, but is "
285  << "currently of type " << TpNbr.type() << exit(FatalError);
286  }
287 
288  const coupledTemperatureFvPatchScalarField& coupledTemperatureNbr =
289  refCast<const coupledTemperatureFvPatchScalarField>(TpNbr);
290 
291  scalarField sumq(size(), 0);
292 
293  if (qs_.valid())
294  {
295  sumq += qs_();
296  }
297 
298  if (qrName_ != "none")
299  {
300  sumq += patch().lookupPatchField<volScalarField, scalar>(qrName_);
301  }
302 
303  if (qrNbrName_ != "none")
304  {
305  sumq += mpp.fromNeighbour
306  (
307  patchNbr.lookupPatchField<volScalarField, scalar>(qrNbrName_)
308  );
309  }
310 
312  tmp<scalarField> sumKappaTByDelta;
313  tmp<scalarField> sumKappaByDelta;
314  tmp<scalarField> qByKappa;
315 
316  // q = alpha.this*sumq
317  getThis(kappa, sumKappaTByDelta, sumKappaByDelta, sumq, qByKappa);
318 
319  // Add neighbour contributions
320  {
321  tmp<scalarField> sumKappaTByDeltaNbr;
322  tmp<scalarField> sumKappaByDeltaNbr;
323  tmp<scalarField> qNbr;
324 
325  if (wallKappaByDelta_ == 0)
326  {
327  coupledTemperatureNbr.getNbr
328  (
329  sumKappaTByDeltaNbr,
330  sumKappaByDeltaNbr,
331  qNbr
332  );
333 
334  add(sumKappaTByDelta, mpp.fromNeighbour(sumKappaTByDeltaNbr));
335  add(sumKappaByDelta, mpp.fromNeighbour(sumKappaByDeltaNbr));
336  }
337  else
338  {
339  // Get the neighbour wall temperature and flux correction
340  tmp<scalarField> TwNbr;
341  coupledTemperatureNbr.getNbr(TwNbr, qNbr);
342 
343  add(sumKappaByDelta, scalarField(size(), wallKappaByDelta_));
344  add(sumKappaTByDelta, wallKappaByDelta_*mpp.fromNeighbour(TwNbr));
345  }
346 
347  if (qNbr.valid())
348  {
349  sumq += mpp.fromNeighbour(qNbr);
350  }
351  }
352 
353  this->valueFraction() =
354  sumKappaByDelta()/(kappa()*patch().deltaCoeffs() + sumKappaByDelta());
355 
356  this->refValue() = (sumKappaTByDelta() + sumq)/sumKappaByDelta();
357 
358  this->refGrad() = qByKappa;
359 
360  mixedFvPatchScalarField::updateCoeffs();
361 
362  if (debug)
363  {
364  const scalar Q = gSum(kappa()*patch().magSf()*snGrad());
365 
366  Info<< patch().boundaryMesh().mesh().name() << ':'
367  << patch().name() << ':'
368  << this->internalField().name() << " <- "
369  << mpp.nbrMesh().name() << ':'
370  << patchNbr.name() << ':'
371  << this->internalField().name() << " :"
372  << " heat transfer rate:" << Q
373  << " walltemperature "
374  << " min:" << gMin(*this)
375  << " max:" << gMax(*this)
376  << " avg:" << gAverage(*this)
377  << endl;
378  }
379 
380  // Restore tag
381  UPstream::msgType() = oldTag;
382 }
383 
384 
386 (
387  Ostream& os
388 ) const
389 {
391 
392  writeEntryIfDifferent<word>(os, "Tnbr", "T", TnbrName_);
393  writeEntryIfDifferent<word>(os, "qrNbr", "none", qrNbrName_);
394  writeEntryIfDifferent<word>(os, "qr", "none", qrName_);
395 
396  if (Qs_ != 0)
397  {
398  writeEntry(os, "Qs", Qs_);
399  }
400  else if (qs_.valid())
401  {
402  writeEntry(os, "qs", qs_());
403  }
404 
405  if (thicknessLayers_.size())
406  {
407  writeEntry(os, "thicknessLayers", thicknessLayers_);
408  writeEntry(os, "kappaLayers", kappaLayers_);
409  }
410 }
411 
412 
413 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
414 
415 namespace Foam
416 {
418  (
421  );
422 
424  (
427  patchMapper,
428  turbulentTemperatureCoupledBaffleMixed,
429  "compressible::turbulentTemperatureCoupledBaffleMixed"
430  );
431 
433  (
436  dictionary,
437  turbulentTemperatureCoupledBaffleMixed,
438  "compressible::turbulentTemperatureCoupledBaffleMixed"
439  );
440 
441 
443  (
446  patchMapper,
447  turbulentTemperatureRadCoupledMixed,
448  "compressible::turbulentTemperatureRadCoupledMixed"
449  );
450 
452  (
455  dictionary,
456  turbulentTemperatureRadCoupledMixed,
457  "compressible::turbulentTemperatureRadCoupledMixed"
458  );
459 }
460 
461 
462 // ************************************************************************* //
#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...
Generic GeometricField class.
const word & name() const
Return name.
Definition: IOobject.H:310
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 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:160
Foam::fvPatchFieldMapper.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:87
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:251
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
virtual const word & name() const
Return name.
Definition: fvPatch.H:145
const GeometricField::Patch & lookupPatchField(const word &name) const
Lookup and return the patchField of the named field from the.
Engine which provides mapping between two patches.
const polyPatch & nbrPolyPatch() const
Get the patch to map from.
static const mappedPatchBase & getMap(const polyPatch &patch)
Cast the given polyPatch to a mappedPatchBase. Handle errors.
const polyMesh & nbrMesh() const
Get the mesh for the region to map from.
static void validateMapForField(const PatchFieldType &field, const FieldType &iF, const dictionary &context, const label froms=from::any)
Validate that the map exists and is appropriate for the given.
tmp< Field< Type > > fromNeighbour(const Field< Type > &nbrFld) const
Map/interpolate the neighbour patch field to this patch.
label index() const
Return the index of this patch in the boundaryMesh.
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:318
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
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)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
messageStream Info
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
addBackwardCompatibleToRunTimeSelectionTable(fvPatchScalarField, coupledTemperatureFvPatchScalarField, patchMapper, turbulentTemperatureCoupledBaffleMixed, "compressible::turbulentTemperatureCoupledBaffleMixed")
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:61
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)
dictionary dict
volScalarField & p
static const label differentPatch