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-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 
27 #include "volFields.H"
28 #include "surfaceFields.H"
30 #include "mappedFvPatchBaseBase.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  const dictionary& dict
48 )
49 :
50  mixedFvPatchScalarField(p, iF, dict, false),
51  TName_("T"),
52  baffleActivated_(dict.lookupOrDefault<bool>("baffleActivated", true)),
53  thickness_(),
54  qs_(p.size(), 0),
55  solidDict_(dict),
56  solidPtr_(),
57  qrPrevious_(p.size(), 0.0),
58  qrRelaxation_(dict.lookupOrDefault<scalar>("qrRelaxation", dimless, 1)),
59  qrName_(dict.lookupOrDefault<word>("qr", "none"))
60 {
62  (
63  *this,
64  iF,
65  dict,
68  );
69 
71  (
72  scalarField("value", iF.dimensions(), dict, p.size())
73  );
74 
75  if (dict.found("thickness"))
76  {
77  thickness_ = scalarField("thickness", dimLength, dict, p.size());
78  }
79 
80  if (dict.found("qs"))
81  {
82  qs_ = scalarField("qs", dimPower/dimArea, dict, p.size());
83  }
84 
85  if (dict.found("qrPrevious"))
86  {
87  qrPrevious_ =
88  scalarField("qrPrevious", dimPower/dimArea, dict, p.size());
89  }
90 
91  if (dict.found("refValue") && baffleActivated_)
92  {
93  // Full restart
94  refValue() = scalarField("refValue", iF.dimensions(), dict, p.size());
95  refGrad() =
97  (
98  "refGradient",
99  iF.dimensions()/dimLength,
100  dict,
101  p.size()
102  );
103  valueFraction() =
104  scalarField("valueFraction", unitFraction, dict, p.size());
105  }
106  else
107  {
108  // Start from user entered data. Assume zeroGradient.
109  refValue() = *this;
110  refGrad() = 0.0;
111  valueFraction() = 0.0;
112  }
113 }
114 
115 
116 template<class solidType>
119 (
121  const fvPatch& p,
123  const fieldMapper& mapper
124 )
125 :
126  mixedFvPatchScalarField(ptf, p, iF, mapper),
127  TName_(ptf.TName_),
128  baffleActivated_(ptf.baffleActivated_),
129  thickness_(mapper(ptf.thickness_)),
130  qs_(mapper(ptf.qs_)),
131  solidDict_(ptf.solidDict_),
132  solidPtr_(ptf.solidPtr_),
133  qrPrevious_(mapper(ptf.qrPrevious_)),
134  qrRelaxation_(ptf.qrRelaxation_),
135  qrName_(ptf.qrName_)
136 {}
137 
138 
139 template<class solidType>
142 (
145 )
146 :
147  mixedFvPatchScalarField(ptf, iF),
148  TName_(ptf.TName_),
149  baffleActivated_(ptf.baffleActivated_),
150  thickness_(ptf.thickness_),
151  qs_(ptf.qs_),
152  solidDict_(ptf.solidDict_),
153  solidPtr_(ptf.solidPtr_),
154  qrPrevious_(ptf.qrPrevious_),
155  qrRelaxation_(ptf.qrRelaxation_),
156  qrName_(ptf.qrName_)
157 {}
158 
159 
160 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
161 
162 template<class solidType>
164 {
165  const mappedFvPatchBaseBase& mapper =
167  return patch().index() < mapper.nbrFvPatch().index();
168 }
169 
170 
171 template<class solidType>
172 const thermalBaffle1DFvPatchScalarField<solidType>&
173 thermalBaffle1DFvPatchScalarField<solidType>::nbrField() const
174 {
175  const mappedFvPatchBaseBase& mapper =
177  const polyMesh& nbrMesh = mapper.nbrMesh();
178  const label nbrPatchi = mapper.nbrFvPatch().index();
179  const fvPatch& nbrPatch =
180  refCast<const fvMesh>(nbrMesh).boundary()[nbrPatchi];
181 
182  return
183  refCast<const thermalBaffle1DFvPatchScalarField>
184  (
185  nbrPatch.template lookupPatchField<volScalarField, scalar>
186  (
187  TName_
188  )
189  );
190 }
191 
192 
193 template<class solidType>
194 const solidType& thermalBaffle1DFvPatchScalarField<solidType>::solid() const
195 {
196  if (this->owner())
197  {
198  if (solidPtr_.empty())
199  {
200  solidPtr_.reset(new solidType("solid", solidDict_));
201  }
202 
203  return solidPtr_();
204  }
205  else
206  {
207  return nbrField().solid();
208  }
209 }
210 
211 
212 template<class solidType>
213 tmp<scalarField> thermalBaffle1DFvPatchScalarField<solidType>::
214 baffleThickness() const
215 {
216  if (this->owner())
217  {
218  if (thickness_.size() != patch().size())
219  {
220  FatalIOErrorInFunction(solidDict_)
221  << " Field thickness has not been specified "
222  << " for patch " << this->patch().name()
223  << exit(FatalIOError);
224  }
225 
226  return thickness_;
227  }
228  else
229  {
230  const mappedFvPatchBaseBase& mapper =
232  return mapper.fromNeighbour(nbrField().baffleThickness());
233  }
234 }
235 
236 
237 template<class solidType>
238 tmp<scalarField> thermalBaffle1DFvPatchScalarField<solidType>::qs() const
239 {
240  if (this->owner())
241  {
242  return qs_;
243  }
244  else
245  {
246  const mappedFvPatchBaseBase& mapper =
248  return mapper.fromNeighbour(nbrField().qs());
249  }
250 }
251 
252 
253 template<class solidType>
255 (
256  const fvPatchScalarField& ptf,
257  const fieldMapper& mapper
258 )
259 {
260  mixedFvPatchScalarField::map(ptf, mapper);
261 
262  const thermalBaffle1DFvPatchScalarField& tiptf =
263  refCast<const thermalBaffle1DFvPatchScalarField>(ptf);
264 
265  if (this->owner())
266  {
267  mapper(thickness_, tiptf.thickness_);
268  mapper(qs_, tiptf.qs_);
269  }
270 }
271 
272 
273 template<class solidType>
275 (
276  const fvPatchScalarField& ptf
277 )
278 {
279  mixedFvPatchScalarField::reset(ptf);
280 
281  const thermalBaffle1DFvPatchScalarField& tiptf =
282  refCast<const thermalBaffle1DFvPatchScalarField>(ptf);
283 
284  if (this->owner())
285  {
286  thickness_.reset(tiptf.thickness_);
287  qs_.reset(tiptf.qs_);
288  }
289 }
290 
291 
292 template<class solidType>
294 {
295  if (updated())
296  {
297  return;
298  }
299 
300  // Since we're inside initEvaluate/evaluate there might be processor
301  // comms underway. Change the tag we use.
302  int oldTag = UPstream::msgType();
303  UPstream::msgType() = oldTag + 1;
304 
305  const mappedFvPatchBaseBase& mapper =
307 
308  if (baffleActivated_)
309  {
310  const thermophysicalTransportModel& ttm =
311  db().objectRegistry::template lookupObject
312  <
314  >
315  (
317  (
318  thermophysicalTransportModel::typeName,
319  internalField().group()
320  )
321  );
322 
323  // Local properties
324  const fvPatchScalarField& Tp =
325  patch().template lookupPatchField<volScalarField, scalar>(TName_);
326 
327  const scalarField kappap(ttm.kappaEff(patch().index()));
328 
329  scalarField qr(Tp.size(), Zero);
330 
331  if (qrName_ != "none")
332  {
333  qr = patch().template lookupPatchField<volScalarField, scalar>
334  (qrName_);
335 
336  qr = qrRelaxation_*qr + (1.0 - qrRelaxation_)*qrPrevious_;
337  qrPrevious_ = qr;
338  }
339 
340  const scalarField kappaDelta(kappap*patch().deltaCoeffs());
341 
342  // Neighbour properties
343  const scalarField nbrTp(mapper.fromNeighbour(nbrField()));
344 
345  // Solid properties
346  scalarField kappas(patch().size(), 0.0);
347  forAll(kappas, i)
348  {
349  kappas[i] = solid().kappa(0.0, (Tp[i] + nbrTp[i])/2.0);
350  }
351 
352  const scalarField kappaByDeltaSolid(kappas/baffleThickness());
353 
354  const scalarField alpha(kappaByDeltaSolid - qr/Tp);
355 
356  valueFraction() = alpha/(alpha + kappaDelta);
357 
358  refValue() = (kappaByDeltaSolid*nbrTp + qs()/2.0)/alpha;
359 
360  if (debug)
361  {
362  scalar Q = gAverage(kappap*snGrad());
363  Info<< patch().boundaryMesh().mesh().name() << ':'
364  << patch().name() << ':'
365  << this->internalField().name() << " <- "
366  << nbrField().patch().name() << ':'
367  << this->internalField().name() << " :"
368  << " heat[W]:" << Q
369  << " walltemperature "
370  << " min:" << gMin(*this)
371  << " max:" << gMax(*this)
372  << " avg:" << gAverage(*this)
373  << endl;
374  }
375  }
376 
377  // Restore tag
378  UPstream::msgType() = oldTag;
379 
380  mixedFvPatchScalarField::updateCoeffs();
381 }
382 
383 
384 template<class solidType>
386 {
388 
389  if (this->owner())
390  {
391  writeEntry(os, "thickness", baffleThickness()());
392  writeEntry(os, "qs", qs()());
393  solid().write(os);
394  }
395 
396  writeEntry(os, "qrPrevious", qrPrevious_);
397  writeEntry(os, "qr", qrName_);
398  writeEntry(os, "qrRelaxation", qrRelaxation_);
399 }
400 
401 
402 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
403 
404 } // End namespace compressible
405 } // End namespace Foam
406 
407 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const dimensionSet & dimensions() const
Return dimensions.
static word groupName(Name name, const word &group)
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
virtual void reset(const fvPatchScalarField &)
Reset the fvPatchField to the given fvPatchField.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void map(const fvPatchScalarField &, const fieldMapper &)
Map the given fvPatchField onto this fvPatchField.
thermalBaffle1DFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
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
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
label index() const
Return the index of this patch in the fvBoundaryMesh.
Definition: fvPatch.H:156
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.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
Abstract base class for all fluid and solid thermophysical transport models.
virtual tmp< volScalarField > kappaEff() const =0
Effective thermal turbulent conductivity.
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
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
const char *const group
Group name for atomic constants.
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
static const zero Zero
Definition: zero.H:97
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:257
const dimensionSet dimPower
const dimensionSet dimless
messageStream Info
const dimensionSet dimLength
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
IOerror FatalIOError
Type gAverage(const FieldField< Field, Type > &f)
const dimensionSet dimArea
Type gMin(const FieldField< Field, Type > &f)
Type gMax(const FieldField< Field, Type > &f)
const unitConversion unitFraction
dictionary dict
volScalarField & p
Foam::surfaceFields.