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  solidPtr_
56  (
57  this->owner()
58  ? new solidType("solid", dict)
59  : nullptr
60  ),
61  qrPrevious_(p.size(), 0.0),
62  qrRelaxation_(dict.lookupOrDefault<scalar>("qrRelaxation", dimless, 1)),
63  qrName_(dict.lookupOrDefault<word>("qr", "none"))
64 {
66  (
67  *this,
68  iF,
69  dict,
72  );
73 
75  (
76  scalarField("value", iF.dimensions(), dict, p.size())
77  );
78 
79  if (dict.found("thickness"))
80  {
81  thickness_ = scalarField("thickness", dimLength, dict, p.size());
82  }
83 
84  if (dict.found("qs"))
85  {
86  qs_ = scalarField("qs", dimPower/dimArea, dict, p.size());
87  }
88 
89  if (dict.found("qrPrevious"))
90  {
91  qrPrevious_ =
92  scalarField("qrPrevious", dimPower/dimArea, dict, p.size());
93  }
94 
95  if (dict.found("refValue") && baffleActivated_)
96  {
97  // Full restart
98  refValue() = scalarField("refValue", iF.dimensions(), dict, p.size());
99  refGrad() =
101  (
102  "refGradient",
103  iF.dimensions()/dimLength,
104  dict,
105  p.size()
106  );
107  valueFraction() =
108  scalarField("valueFraction", unitFraction, dict, p.size());
109  }
110  else
111  {
112  // Start from user entered data. Assume zeroGradient.
113  refValue() = *this;
114  refGrad() = 0.0;
115  valueFraction() = 0.0;
116  }
117 }
118 
119 
120 template<class solidType>
123 (
125  const fvPatch& p,
127  const fieldMapper& mapper
128 )
129 :
130  mixedFvPatchScalarField(ptf, p, iF, mapper),
131  TName_(ptf.TName_),
132  baffleActivated_(ptf.baffleActivated_),
133  thickness_(mapper(ptf.thickness_)),
134  qs_(mapper(ptf.qs_)),
135  solidPtr_(ptf.solidPtr_),
136  qrPrevious_(mapper(ptf.qrPrevious_)),
137  qrRelaxation_(ptf.qrRelaxation_),
138  qrName_(ptf.qrName_)
139 {}
140 
141 
142 template<class solidType>
145 (
148 )
149 :
150  mixedFvPatchScalarField(ptf, iF),
151  TName_(ptf.TName_),
152  baffleActivated_(ptf.baffleActivated_),
153  thickness_(ptf.thickness_),
154  qs_(ptf.qs_),
155  solidPtr_(ptf.solidPtr_),
156  qrPrevious_(ptf.qrPrevious_),
157  qrRelaxation_(ptf.qrRelaxation_),
158  qrName_(ptf.qrName_)
159 {}
160 
161 
162 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
163 
164 template<class solidType>
166 {
167  const mappedFvPatchBaseBase& mapper =
169  return patch().index() < mapper.nbrFvPatch().index();
170 }
171 
172 
173 template<class solidType>
174 const thermalBaffle1DFvPatchScalarField<solidType>&
175 thermalBaffle1DFvPatchScalarField<solidType>::nbrField() const
176 {
177  const mappedFvPatchBaseBase& mapper =
179  const polyMesh& nbrMesh = mapper.nbrMesh();
180  const label nbrPatchi = mapper.nbrFvPatch().index();
181  const fvPatch& nbrPatch =
182  refCast<const fvMesh>(nbrMesh).boundary()[nbrPatchi];
183 
184  return
185  refCast<const thermalBaffle1DFvPatchScalarField>
186  (
187  nbrPatch.template lookupPatchField<volScalarField, scalar>
188  (
189  TName_
190  )
191  );
192 }
193 
194 
195 template<class solidType>
196 const solidType& thermalBaffle1DFvPatchScalarField<solidType>::solid() const
197 {
198  if (this->owner())
199  {
200  if (solidPtr_.empty())
201  {
203  << "solid not allocated" << exit(FatalError);
204  }
205 
206  return solidPtr_();
207  }
208  else
209  {
210  return nbrField().solid();
211  }
212 }
213 
214 
215 template<class solidType>
216 tmp<scalarField> thermalBaffle1DFvPatchScalarField<solidType>::
217 baffleThickness() const
218 {
219  if (this->owner())
220  {
221  if (thickness_.size() != patch().size())
222  {
224  << " Field thickness has not been specified "
225  << " for patch " << this->patch().name()
226  << exit(FatalError);
227  }
228 
229  return thickness_;
230  }
231  else
232  {
233  const mappedFvPatchBaseBase& mapper =
235  return mapper.fromNeighbour(nbrField().baffleThickness());
236  }
237 }
238 
239 
240 template<class solidType>
241 tmp<scalarField> thermalBaffle1DFvPatchScalarField<solidType>::qs() const
242 {
243  if (this->owner())
244  {
245  return qs_;
246  }
247  else
248  {
249  const mappedFvPatchBaseBase& mapper =
251  return mapper.fromNeighbour(nbrField().qs());
252  }
253 }
254 
255 
256 template<class solidType>
258 (
259  const fvPatchScalarField& ptf,
260  const fieldMapper& mapper
261 )
262 {
263  mixedFvPatchScalarField::map(ptf, mapper);
264 
265  const thermalBaffle1DFvPatchScalarField& tiptf =
266  refCast<const thermalBaffle1DFvPatchScalarField>(ptf);
267 
268  if (this->owner())
269  {
270  mapper(thickness_, tiptf.thickness_);
271  mapper(qs_, tiptf.qs_);
272  }
273 }
274 
275 
276 template<class solidType>
278 (
279  const fvPatchScalarField& ptf
280 )
281 {
282  mixedFvPatchScalarField::reset(ptf);
283 
284  const thermalBaffle1DFvPatchScalarField& tiptf =
285  refCast<const thermalBaffle1DFvPatchScalarField>(ptf);
286 
287  if (this->owner())
288  {
289  thickness_.reset(tiptf.thickness_);
290  qs_.reset(tiptf.qs_);
291  }
292 }
293 
294 
295 template<class solidType>
297 {
298  if (updated())
299  {
300  return;
301  }
302 
303  // Since we're inside initEvaluate/evaluate there might be processor
304  // comms underway. Change the tag we use.
305  int oldTag = UPstream::msgType();
306  UPstream::msgType() = oldTag + 1;
307 
308  const mappedFvPatchBaseBase& mapper =
310 
311  if (baffleActivated_)
312  {
313  const thermophysicalTransportModel& ttm =
314  db().objectRegistry::template lookupObject
315  <
317  >
318  (
320  (
321  thermophysicalTransportModel::typeName,
322  internalField().group()
323  )
324  );
325 
326  // Local properties
327  const fvPatchScalarField& Tp =
328  patch().template lookupPatchField<volScalarField, scalar>(TName_);
329 
330  const scalarField kappap(ttm.kappaEff(patch().index()));
331 
332  scalarField qr(Tp.size(), Zero);
333 
334  if (qrName_ != "none")
335  {
336  qr = patch().template lookupPatchField<volScalarField, scalar>
337  (qrName_);
338 
339  qr = qrRelaxation_*qr + (1.0 - qrRelaxation_)*qrPrevious_;
340  qrPrevious_ = qr;
341  }
342 
343  const scalarField kappaDelta(kappap*patch().deltaCoeffs());
344 
345  // Neighbour properties
346  const scalarField nbrTp(mapper.fromNeighbour(nbrField()));
347 
348  // Solid properties
349  scalarField kappas(patch().size(), 0.0);
350  forAll(kappas, i)
351  {
352  kappas[i] = solid().kappa(0.0, (Tp[i] + nbrTp[i])/2.0);
353  }
354 
355  const scalarField kappaByDeltaSolid(kappas/baffleThickness());
356 
357  const scalarField alpha(kappaByDeltaSolid - qr/Tp);
358 
359  valueFraction() = alpha/(alpha + kappaDelta);
360 
361  refValue() = (kappaByDeltaSolid*nbrTp + qs()/2.0)/alpha;
362 
363  if (debug)
364  {
365  scalar Q = gAverage(kappap*snGrad());
366  Info<< patch().boundaryMesh().mesh().name() << ':'
367  << patch().name() << ':'
368  << this->internalField().name() << " <- "
369  << nbrField().patch().name() << ':'
370  << this->internalField().name() << " :"
371  << " heat[W]:" << Q
372  << " walltemperature "
373  << " min:" << gMin(*this)
374  << " max:" << gMax(*this)
375  << " avg:" << gAverage(*this)
376  << endl;
377  }
378  }
379 
380  // Restore tag
381  UPstream::msgType() = oldTag;
382 
383  mixedFvPatchScalarField::updateCoeffs();
384 }
385 
386 
387 template<class solidType>
389 {
391 
392  if (this->owner())
393  {
394  writeEntry(os, "thickness", baffleThickness()());
395  writeEntry(os, "qs", qs()());
396  solid().write(os);
397  }
398 
399  writeEntry(os, "qrPrevious", qrPrevious_);
400  writeEntry(os, "qr", qrName_);
401  writeEntry(os, "qrRelaxation", qrRelaxation_);
402 }
403 
404 
405 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
406 
407 } // End namespace compressible
408 } // End namespace Foam
409 
410 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
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 keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
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:91
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 FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
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:258
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
Type gAverage(const FieldField< Field, Type > &f)
error FatalError
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.