thermalBaffle1DFvPatchScalarField.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-2022 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 
27 #include "volFields.H"
28 #include "surfaceFields.H"
30 #include "mappedPatchBase.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace compressible
37 {
38 
39 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
40 
41 template<class solidType>
44 (
45  const fvPatch& p,
47 )
48 :
49  mixedFvPatchScalarField(p, iF),
50  TName_("T"),
51  baffleActivated_(true),
52  thickness_(p.size()),
53  qs_(p.size()),
54  solidDict_(),
55  solidPtr_(nullptr),
56  qrPrevious_(p.size()),
57  qrRelaxation_(1),
58  qrName_("undefined-qr")
59 {}
60 
61 
62 template<class solidType>
65 (
66  const fvPatch& p,
68  const dictionary& dict
69 )
70 :
71  mixedFvPatchScalarField(p, iF),
72  TName_("T"),
73  baffleActivated_(dict.lookupOrDefault<bool>("baffleActivated", true)),
74  thickness_(),
75  qs_(p.size(), 0),
76  solidDict_(dict),
77  solidPtr_(),
78  qrPrevious_(p.size(), 0.0),
79  qrRelaxation_(dict.lookupOrDefault<scalar>("qrRelaxation", 1)),
80  qrName_(dict.lookupOrDefault<word>("qr", "none"))
81 {
82  if (!isA<mappedPatchBase>(this->patch().patch()))
83  {
85  << "' not type '" << mappedPatchBase::typeName << "'"
86  << "\n for patch " << p.name()
87  << " of field " << internalField().name()
88  << " in file " << internalField().objectPath()
89  << exit(FatalError);
90  }
91 
92  fvPatchScalarField::operator=(scalarField("value", dict, p.size()));
93 
94  if (dict.found("thickness"))
95  {
96  thickness_ = scalarField("thickness", dict, p.size());
97  }
98 
99  if (dict.found("qs"))
100  {
101  qs_ = scalarField("qs", dict, p.size());
102  }
103  else if (dict.found("Qs"))
104  {
105  qs_ = scalarField("Qs", dict, p.size());
106  }
107 
108  if (dict.found("qrPrevious"))
109  {
110  qrPrevious_ = scalarField("qrPrevious", dict, p.size());
111  }
112 
113  if (dict.found("refValue") && baffleActivated_)
114  {
115  // Full restart
116  refValue() = scalarField("refValue", dict, p.size());
117  refGrad() = scalarField("refGradient", dict, p.size());
118  valueFraction() = scalarField("valueFraction", dict, p.size());
119  }
120  else
121  {
122  // Start from user entered data. Assume zeroGradient.
123  refValue() = *this;
124  refGrad() = 0.0;
125  valueFraction() = 0.0;
126  }
127 }
128 
129 
130 template<class solidType>
133 (
135  const fvPatch& p,
137  const fvPatchFieldMapper& mapper
138 )
139 :
140  mixedFvPatchScalarField(ptf, p, iF, mapper),
141  TName_(ptf.TName_),
142  baffleActivated_(ptf.baffleActivated_),
143  thickness_(mapper(ptf.thickness_)),
144  qs_(mapper(ptf.qs_)),
145  solidDict_(ptf.solidDict_),
146  solidPtr_(ptf.solidPtr_),
147  qrPrevious_(mapper(ptf.qrPrevious_)),
148  qrRelaxation_(ptf.qrRelaxation_),
149  qrName_(ptf.qrName_)
150 {}
151 
152 
153 template<class solidType>
156 (
159 )
160 :
161  mixedFvPatchScalarField(ptf, iF),
162  TName_(ptf.TName_),
163  baffleActivated_(ptf.baffleActivated_),
164  thickness_(ptf.thickness_),
165  qs_(ptf.qs_),
166  solidDict_(ptf.solidDict_),
167  solidPtr_(ptf.solidPtr_),
168  qrPrevious_(ptf.qrPrevious_),
169  qrRelaxation_(ptf.qrRelaxation_),
170  qrName_(ptf.qrName_)
171 {}
172 
173 
174 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
175 
176 template<class solidType>
178 {
179  const mappedPatchBase& mpp =
180  refCast<const mappedPatchBase>(patch().patch());
181 
182  if (!mpp.sameRegion())
183  {
185  << "A" << typeName
186  << " must map to a patch field in the same region"
187  << exit(FatalError);
188  }
189 
190  return patch().patch().index() < mpp.samplePolyPatch().index();
191 }
192 
193 
194 template<class solidType>
197 {
198  const mappedPatchBase& mpp =
199  refCast<const mappedPatchBase>(patch().patch());
200  const polyMesh& nbrMesh = mpp.sampleMesh();
201  const label samplePatchi = mpp.samplePolyPatch().index();
202  const fvPatch& nbrPatch =
203  refCast<const fvMesh>(nbrMesh).boundary()[samplePatchi];
204 
205  return
206  refCast<const thermalBaffle1DFvPatchScalarField>
207  (
208  nbrPatch.template lookupPatchField<volScalarField, scalar>
209  (
210  TName_
211  )
212  );
213 }
214 
215 
216 template<class solidType>
218 {
219  if (this->owner())
220  {
221  if (solidPtr_.empty())
222  {
223  solidPtr_.reset(new solidType(solidDict_));
224  }
225 
226  return solidPtr_();
227  }
228  else
229  {
230  return nbrField().solid();
231  }
232 }
233 
234 
235 template<class solidType>
237 baffleThickness() const
238 {
239  if (this->owner())
240  {
241  if (thickness_.size() != patch().size())
242  {
243  FatalIOErrorInFunction(solidDict_)
244  << " Field thickness has not been specified "
245  << " for patch " << this->patch().name()
246  << exit(FatalIOError);
247  }
248 
249  return thickness_;
250  }
251  else
252  {
253  const mappedPatchBase& mpp =
254  refCast<const mappedPatchBase>(patch().patch());
255 
256  tmp<scalarField> nbrBaffleThickness
257  (
258  new scalarField(nbrField().baffleThickness())
259  );
260  mpp.distribute(nbrBaffleThickness.ref());
261  return nbrBaffleThickness;
262  }
263 }
264 
265 
266 template<class solidType>
268 {
269  if (this->owner())
270  {
271  return qs_;
272  }
273  else
274  {
275  const mappedPatchBase& mpp =
276  refCast<const mappedPatchBase>(patch().patch());
277 
278  tmp<scalarField> nbrQs
279  (
280  new scalarField(nbrField().qs())
281  );
282  mpp.distribute(nbrQs.ref());
283  return nbrQs;
284  }
285 }
286 
287 
288 template<class solidType>
290 (
291  const fvPatchFieldMapper& m
292 )
293 {
294  mixedFvPatchScalarField::autoMap(m);
295 
296  if (this->owner())
297  {
298  m(thickness_, thickness_);
299  m(qs_, qs_);
300  }
301 }
302 
303 
304 template<class solidType>
306 (
307  const fvPatchScalarField& ptf,
308  const labelList& addr
309 )
310 {
311  mixedFvPatchScalarField::rmap(ptf, addr);
312 
313  const thermalBaffle1DFvPatchScalarField& tiptf =
314  refCast<const thermalBaffle1DFvPatchScalarField>(ptf);
315 
316  if (this->owner())
317  {
318  thickness_.rmap(tiptf.thickness_, addr);
319  qs_.rmap(tiptf.qs_, addr);
320  }
321 }
322 
323 
324 template<class solidType>
326 (
327  const fvPatchScalarField& ptf
328 )
329 {
330  mixedFvPatchScalarField::reset(ptf);
331 
332  const thermalBaffle1DFvPatchScalarField& tiptf =
333  refCast<const thermalBaffle1DFvPatchScalarField>(ptf);
334 
335  if (this->owner())
336  {
337  thickness_.reset(tiptf.thickness_);
338  qs_.reset(tiptf.qs_);
339  }
340 }
341 
342 
343 template<class solidType>
345 {
346  if (updated())
347  {
348  return;
349  }
350 
351  // Since we're inside initEvaluate/evaluate there might be processor
352  // comms underway. Change the tag we use.
353  int oldTag = UPstream::msgType();
354  UPstream::msgType() = oldTag + 1;
355 
356  const mappedPatchBase& mpp =
357  refCast<const mappedPatchBase>(patch().patch());
358 
359  if (baffleActivated_)
360  {
361  const thermophysicalTransportModel& ttm =
362  db().objectRegistry::template lookupObject
363  <
365  >
366  (
368  (
369  thermophysicalTransportModel::typeName,
370  internalField().group()
371  )
372  );
373 
374  // Local properties
375  const fvPatchScalarField& Tp =
376  patch().template lookupPatchField<volScalarField, scalar>(TName_);
377 
378  const scalarField kappap(ttm.kappaEff(patch().index()));
379 
380  scalarField qr(Tp.size(), Zero);
381 
382  if (qrName_ != "none")
383  {
384  qr = patch().template lookupPatchField<volScalarField, scalar>
385  (qrName_);
386 
387  qr = qrRelaxation_*qr + (1.0 - qrRelaxation_)*qrPrevious_;
388  qrPrevious_ = qr;
389  }
390 
391  scalarField kappaDelta(kappap*patch().deltaCoeffs());
392 
393  // Neighbour properties
394  scalarField nbrTp(nbrField());
395  mpp.distribute(nbrTp);
396 
397  // Solid properties
398  scalarField kappas(patch().size(), 0.0);
399  forAll(kappas, i)
400  {
401  kappas[i] = solid().kappa(0.0, (Tp[i] + nbrTp[i])/2.0);
402  }
403 
404  scalarField KDeltaSolid(kappas/baffleThickness());
405 
406  scalarField alpha(KDeltaSolid - qr/Tp);
407 
408  valueFraction() = alpha/(alpha + kappaDelta);
409 
410  refValue() = (KDeltaSolid*nbrTp + qs()/2.0)/alpha;
411 
412  if (debug)
413  {
414  scalar Q = gAverage(kappap*snGrad());
415  Info<< patch().boundaryMesh().mesh().name() << ':'
416  << patch().name() << ':'
417  << this->internalField().name() << " <- "
418  << nbrField().patch().name() << ':'
419  << this->internalField().name() << " :"
420  << " heat[W]:" << Q
421  << " walltemperature "
422  << " min:" << gMin(*this)
423  << " max:" << gMax(*this)
424  << " avg:" << gAverage(*this)
425  << endl;
426  }
427  }
428 
429  // Restore tag
430  UPstream::msgType() = oldTag;
431 
432  mixedFvPatchScalarField::updateCoeffs();
433 }
434 
435 
436 template<class solidType>
438 {
440 
441  if (this->owner())
442  {
443  writeEntry(os, "thickness", baffleThickness()());
444  writeEntry(os, "qs", qs()());
445  solid().write(os);
446  }
447 
448  writeEntry(os, "qrPrevious", qrPrevious_);
449  writeEntry(os, "qr", qrName_);
450  writeEntry(os, "qrRelaxation", qrRelaxation_);
451 }
452 
453 
454 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
455 
456 } // End namespace compressible
457 } // End namespace Foam
458 
459 // ************************************************************************* //
Foam::surfaceFields.
const char *const group
Group name for atomic constants.
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:663
thermalBaffle1DFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
bool sameRegion() const
Cached sampleRegion != mesh.name()
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
Type gMin(const FieldField< Field, Type > &f)
const polyMesh & sampleMesh() const
Get the region mesh.
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
static int & msgType()
Message tag of standard messages.
Definition: UPstream.H:476
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
A class for handling words, derived from string.
Definition: word.H:59
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
static word groupName(Name name, const word &group)
Foam::fvPatchFieldMapper.
virtual void reset(const fvPatchScalarField &)
Reset the fvPatchField to the given fvPatchField.
Determines a mapping between patch face centres and mesh cell or face centres and processors they&#39;re ...
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
static const zero Zero
Definition: zero.H:97
faceListList boundary(nPatches)
virtual label size() const
Return size.
Definition: fvPatch.H:157
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
Type gMax(const FieldField< Field, Type > &f)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
Abstract base class for thermophysical transport models (RAS, LES and laminar).
void distribute(List< Type > &lst) const
Wrapper around map/interpolate data distribution.
virtual tmp< volScalarField > kappaEff() const =0
Effective thermal turbulent diffusivity for temperature.
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:318
Type gAverage(const FieldField< Field, Type > &f)
label index() const
Return the index of this patch in the boundaryMesh.
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:258
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
messageStream Info
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:76
virtual const word & name() const
Return name.
Definition: fvPatch.H:145
A class for managing temporary objects.
Definition: PtrList.H:53
const polyPatch & samplePolyPatch() const
Get the patch on the region.
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:45
Namespace for OpenFOAM.
IOerror FatalIOError