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-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 
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  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", 1)),
59  qrName_(dict.lookupOrDefault<word>("qr", "none"))
60 {
62  (
63  *this,
64  iF,
65  dict,
68  );
69 
71 
72  if (dict.found("thickness"))
73  {
74  thickness_ = scalarField("thickness", dict, p.size());
75  }
76 
77  if (dict.found("qs"))
78  {
79  qs_ = scalarField("qs", dict, p.size());
80  }
81  else if (dict.found("Qs"))
82  {
83  qs_ = scalarField("Qs", dict, p.size());
84  }
85 
86  if (dict.found("qrPrevious"))
87  {
88  qrPrevious_ = scalarField("qrPrevious", dict, p.size());
89  }
90 
91  if (dict.found("refValue") && baffleActivated_)
92  {
93  // Full restart
94  refValue() = scalarField("refValue", dict, p.size());
95  refGrad() = scalarField("refGradient", dict, p.size());
96  valueFraction() = scalarField("valueFraction", dict, p.size());
97  }
98  else
99  {
100  // Start from user entered data. Assume zeroGradient.
101  refValue() = *this;
102  refGrad() = 0.0;
103  valueFraction() = 0.0;
104  }
105 }
106 
107 
108 template<class solidType>
111 (
113  const fvPatch& p,
115  const fvPatchFieldMapper& mapper
116 )
117 :
118  mixedFvPatchScalarField(ptf, p, iF, mapper),
119  TName_(ptf.TName_),
120  baffleActivated_(ptf.baffleActivated_),
121  thickness_(mapper(ptf.thickness_)),
122  qs_(mapper(ptf.qs_)),
123  solidDict_(ptf.solidDict_),
124  solidPtr_(ptf.solidPtr_),
125  qrPrevious_(mapper(ptf.qrPrevious_)),
126  qrRelaxation_(ptf.qrRelaxation_),
127  qrName_(ptf.qrName_)
128 {}
129 
130 
131 template<class solidType>
134 (
137 )
138 :
139  mixedFvPatchScalarField(ptf, iF),
140  TName_(ptf.TName_),
141  baffleActivated_(ptf.baffleActivated_),
142  thickness_(ptf.thickness_),
143  qs_(ptf.qs_),
144  solidDict_(ptf.solidDict_),
145  solidPtr_(ptf.solidPtr_),
146  qrPrevious_(ptf.qrPrevious_),
147  qrRelaxation_(ptf.qrRelaxation_),
148  qrName_(ptf.qrName_)
149 {}
150 
151 
152 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
153 
154 template<class solidType>
156 {
157  const mappedPatchBase& mpp = mappedPatchBase::getMap(patch().patch());
158  return patch().patch().index() < mpp.nbrPolyPatch().index();
159 }
160 
161 
162 template<class solidType>
163 const thermalBaffle1DFvPatchScalarField<solidType>&
164 thermalBaffle1DFvPatchScalarField<solidType>::nbrField() const
165 {
166  const mappedPatchBase& mpp = mappedPatchBase::getMap(patch().patch());
167  const polyMesh& nbrMesh = mpp.nbrMesh();
168  const label nbrPatchi = mpp.nbrPolyPatch().index();
169  const fvPatch& nbrPatch =
170  refCast<const fvMesh>(nbrMesh).boundary()[nbrPatchi];
171 
172  return
173  refCast<const thermalBaffle1DFvPatchScalarField>
174  (
175  nbrPatch.template lookupPatchField<volScalarField, scalar>
176  (
177  TName_
178  )
179  );
180 }
181 
182 
183 template<class solidType>
184 const solidType& thermalBaffle1DFvPatchScalarField<solidType>::solid() const
185 {
186  if (this->owner())
187  {
188  if (solidPtr_.empty())
189  {
190  solidPtr_.reset(new solidType("solid", solidDict_));
191  }
192 
193  return solidPtr_();
194  }
195  else
196  {
197  return nbrField().solid();
198  }
199 }
200 
201 
202 template<class solidType>
203 tmp<scalarField> thermalBaffle1DFvPatchScalarField<solidType>::
204 baffleThickness() const
205 {
206  if (this->owner())
207  {
208  if (thickness_.size() != patch().size())
209  {
210  FatalIOErrorInFunction(solidDict_)
211  << " Field thickness has not been specified "
212  << " for patch " << this->patch().name()
213  << exit(FatalIOError);
214  }
215 
216  return thickness_;
217  }
218  else
219  {
220  const mappedPatchBase& mpp = mappedPatchBase::getMap(patch().patch());
221  return mpp.fromNeighbour(nbrField().baffleThickness());
222  }
223 }
224 
225 
226 template<class solidType>
227 tmp<scalarField> thermalBaffle1DFvPatchScalarField<solidType>::qs() const
228 {
229  if (this->owner())
230  {
231  return qs_;
232  }
233  else
234  {
235  const mappedPatchBase& mpp = mappedPatchBase::getMap(patch().patch());
236  return mpp.fromNeighbour(nbrField().qs());
237  }
238 }
239 
240 
241 template<class solidType>
243 (
244  const fvPatchScalarField& ptf,
245  const fvPatchFieldMapper& mapper
246 )
247 {
248  mixedFvPatchScalarField::map(ptf, mapper);
249 
250  const thermalBaffle1DFvPatchScalarField& tiptf =
251  refCast<const thermalBaffle1DFvPatchScalarField>(ptf);
252 
253  if (this->owner())
254  {
255  mapper(thickness_, tiptf.thickness_);
256  mapper(qs_, tiptf.qs_);
257  }
258 }
259 
260 
261 template<class solidType>
263 (
264  const fvPatchScalarField& ptf
265 )
266 {
267  mixedFvPatchScalarField::reset(ptf);
268 
269  const thermalBaffle1DFvPatchScalarField& tiptf =
270  refCast<const thermalBaffle1DFvPatchScalarField>(ptf);
271 
272  if (this->owner())
273  {
274  thickness_.reset(tiptf.thickness_);
275  qs_.reset(tiptf.qs_);
276  }
277 }
278 
279 
280 template<class solidType>
282 {
283  if (updated())
284  {
285  return;
286  }
287 
288  // Since we're inside initEvaluate/evaluate there might be processor
289  // comms underway. Change the tag we use.
290  int oldTag = UPstream::msgType();
291  UPstream::msgType() = oldTag + 1;
292 
293  const mappedPatchBase& mpp = mappedPatchBase::getMap(patch().patch());
294 
295  if (baffleActivated_)
296  {
297  const thermophysicalTransportModel& ttm =
298  db().objectRegistry::template lookupObject
299  <
301  >
302  (
304  (
305  thermophysicalTransportModel::typeName,
306  internalField().group()
307  )
308  );
309 
310  // Local properties
311  const fvPatchScalarField& Tp =
312  patch().template lookupPatchField<volScalarField, scalar>(TName_);
313 
314  const scalarField kappap(ttm.kappaEff(patch().index()));
315 
316  scalarField qr(Tp.size(), Zero);
317 
318  if (qrName_ != "none")
319  {
320  qr = patch().template lookupPatchField<volScalarField, scalar>
321  (qrName_);
322 
323  qr = qrRelaxation_*qr + (1.0 - qrRelaxation_)*qrPrevious_;
324  qrPrevious_ = qr;
325  }
326 
327  const scalarField kappaDelta(kappap*patch().deltaCoeffs());
328 
329  // Neighbour properties
330  const scalarField nbrTp(mpp.fromNeighbour(nbrField()));
331 
332  // Solid properties
333  scalarField kappas(patch().size(), 0.0);
334  forAll(kappas, i)
335  {
336  kappas[i] = solid().kappa(0.0, (Tp[i] + nbrTp[i])/2.0);
337  }
338 
339  const scalarField kappaByDeltaSolid(kappas/baffleThickness());
340 
341  const scalarField alpha(kappaByDeltaSolid - qr/Tp);
342 
343  valueFraction() = alpha/(alpha + kappaDelta);
344 
345  refValue() = (kappaByDeltaSolid*nbrTp + qs()/2.0)/alpha;
346 
347  if (debug)
348  {
349  scalar Q = gAverage(kappap*snGrad());
350  Info<< patch().boundaryMesh().mesh().name() << ':'
351  << patch().name() << ':'
352  << this->internalField().name() << " <- "
353  << nbrField().patch().name() << ':'
354  << this->internalField().name() << " :"
355  << " heat[W]:" << Q
356  << " walltemperature "
357  << " min:" << gMin(*this)
358  << " max:" << gMax(*this)
359  << " avg:" << gAverage(*this)
360  << endl;
361  }
362  }
363 
364  // Restore tag
365  UPstream::msgType() = oldTag;
366 
367  mixedFvPatchScalarField::updateCoeffs();
368 }
369 
370 
371 template<class solidType>
373 {
375 
376  if (this->owner())
377  {
378  writeEntry(os, "thickness", baffleThickness()());
379  writeEntry(os, "qs", qs()());
380  solid().write(os);
381  }
382 
383  writeEntry(os, "qrPrevious", qrPrevious_);
384  writeEntry(os, "qr", qrName_);
385  writeEntry(os, "qrRelaxation", qrRelaxation_);
386 }
387 
388 
389 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
390 
391 } // End namespace compressible
392 } // End namespace Foam
393 
394 // ************************************************************************* //
#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...
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 fvPatchFieldMapper &)
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: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
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.
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:318
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:251
messageStream Info
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)
Type gMin(const FieldField< Field, Type > &f)
Type gMax(const FieldField< Field, Type > &f)
dictionary dict
volScalarField & p
static const label differentPatch
static const label sameRegion
Foam::surfaceFields.