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-2021 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 
26 #include "volFields.H"
27 #include "surfaceFields.H"
29 #include "mapDistribute.H"
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace compressible
36 {
37 
38 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
39 
40 template<class solidType>
43 (
44  const fvPatch& p,
46 )
47 :
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  mappedPatchBase(p.patch(), NEARESTPATCHFACE, dict),
72  mixedFvPatchScalarField(p, iF),
73  TName_("T"),
74  baffleActivated_(dict.lookupOrDefault<bool>("baffleActivated", true)),
75  thickness_(),
76  Qs_(p.size(), 0),
77  solidDict_(dict),
78  solidPtr_(),
79  qrPrevious_(p.size(), 0.0),
80  qrRelaxation_(dict.lookupOrDefault<scalar>("qrRelaxation", 1)),
81  qrName_(dict.lookupOrDefault<word>("qr", "none"))
82 {
83  fvPatchScalarField::operator=(scalarField("value", dict, p.size()));
84 
85  if (dict.found("thickness"))
86  {
87  thickness_ = scalarField("thickness", dict, p.size());
88  }
89 
90  if (dict.found("Qs"))
91  {
92  Qs_ = scalarField("Qs", dict, p.size());
93  }
94 
95  if (dict.found("qrPrevious"))
96  {
97  qrPrevious_ = scalarField("qrPrevious", dict, p.size());
98  }
99 
100  if (dict.found("refValue") && baffleActivated_)
101  {
102  // Full restart
103  refValue() = scalarField("refValue", dict, p.size());
104  refGrad() = scalarField("refGradient", dict, p.size());
105  valueFraction() = scalarField("valueFraction", dict, p.size());
106  }
107  else
108  {
109  // Start from user entered data. Assume zeroGradient.
110  refValue() = *this;
111  refGrad() = 0.0;
112  valueFraction() = 0.0;
113  }
114 }
115 
116 
117 template<class solidType>
120 (
122  const fvPatch& p,
124  const fvPatchFieldMapper& mapper
125 )
126 :
127  mappedPatchBase(p.patch(), ptf),
128  mixedFvPatchScalarField(ptf, p, iF, mapper),
129  TName_(ptf.TName_),
130  baffleActivated_(ptf.baffleActivated_),
131  thickness_(mapper(ptf.thickness_)),
132  Qs_(mapper(ptf.Qs_)),
133  solidDict_(ptf.solidDict_),
134  solidPtr_(ptf.solidPtr_),
135  qrPrevious_(mapper(ptf.qrPrevious_)),
136  qrRelaxation_(ptf.qrRelaxation_),
137  qrName_(ptf.qrName_)
138 {}
139 
140 
141 template<class solidType>
144 (
147 )
148 :
149  mappedPatchBase(ptf.patch().patch(), ptf),
150  mixedFvPatchScalarField(ptf, iF),
151  TName_(ptf.TName_),
152  baffleActivated_(ptf.baffleActivated_),
153  thickness_(ptf.thickness_),
154  Qs_(ptf.Qs_),
155  solidDict_(ptf.solidDict_),
156  solidPtr_(ptf.solidPtr_),
157  qrPrevious_(ptf.qrPrevious_),
158  qrRelaxation_(ptf.qrRelaxation_),
159  qrName_(ptf.qrName_)
160 {}
161 
162 
163 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
164 
165 template<class solidType>
167 {
168  const label patchi = patch().index();
169 
170  const label nbrPatchi = samplePolyPatch().index();
171 
172  return (patchi < nbrPatchi);
173 }
174 
175 
176 template<class solidType>
178 {
179  if (this->owner())
180  {
181  if (solidPtr_.empty())
182  {
183  solidPtr_.reset(new solidType(solidDict_));
184  }
185  return solidPtr_();
186  }
187  else
188  {
189  const fvPatch& nbrPatch =
190  patch().boundaryMesh()[samplePolyPatch().index()];
191 
192  const thermalBaffle1DFvPatchScalarField& nbrField =
193  refCast<const thermalBaffle1DFvPatchScalarField>
194  (
195  nbrPatch.template lookupPatchField<volScalarField, scalar>(TName_)
196  );
197 
198  return nbrField.solid();
199  }
200 }
201 
202 
203 template<class solidType>
205 baffleThickness() const
206 {
207  if (this->owner())
208  {
209  if (thickness_.size() != patch().size())
210  {
212  (
213  solidDict_
214  )<< " Field thickness has not been specified "
215  << " for patch " << this->patch().name()
216  << exit(FatalIOError);
217  }
218 
219  return thickness_;
220  }
221  else
222  {
223  const mapDistribute& mapDist = this->mappedPatchBase::map();
224 
225  const fvPatch& nbrPatch =
226  patch().boundaryMesh()[samplePolyPatch().index()];
227  const thermalBaffle1DFvPatchScalarField& nbrField =
228  refCast<const thermalBaffle1DFvPatchScalarField>
229  (
230  nbrPatch.template lookupPatchField<volScalarField, scalar>(TName_)
231  );
232 
233  tmp<scalarField> tthickness
234  (
235  new scalarField(nbrField.baffleThickness())
236  );
237  scalarField& thickness = tthickness.ref();
238  mapDist.distribute(thickness);
239  return tthickness;
240  }
241 }
242 
243 
244 template<class solidType>
246 {
247  if (this->owner())
248  {
249  return Qs_;
250  }
251  else
252  {
253  const mapDistribute& mapDist = this->mappedPatchBase::map();
254 
255  const fvPatch& nbrPatch =
256  patch().boundaryMesh()[samplePolyPatch().index()];
257 
258  const thermalBaffle1DFvPatchScalarField& nbrField =
259  refCast<const thermalBaffle1DFvPatchScalarField>
260  (
261  nbrPatch.template lookupPatchField<volScalarField, scalar>(TName_)
262  );
263 
264  tmp<scalarField> tQs(new scalarField(nbrField.Qs()));
265  scalarField& Qs = tQs.ref();
266  mapDist.distribute(Qs);
267  return tQs;
268  }
269 }
270 
271 
272 template<class solidType>
274 (
275  const fvPatchFieldMapper& m
276 )
277 {
279 
280  mixedFvPatchScalarField::autoMap(m);
281 
282  if (this->owner())
283  {
284  m(thickness_, thickness_);
285  m(Qs_, Qs_);
286  }
287 }
288 
289 
290 template<class solidType>
292 (
293  const fvPatchScalarField& ptf,
294  const labelList& addr
295 )
296 {
297  mixedFvPatchScalarField::rmap(ptf, addr);
298 
299  const thermalBaffle1DFvPatchScalarField& tiptf =
300  refCast<const thermalBaffle1DFvPatchScalarField>(ptf);
301 
302  if (this->owner())
303  {
304  thickness_.rmap(tiptf.thickness_, addr);
305  Qs_.rmap(tiptf.Qs_, addr);
306  }
307 }
308 
309 
310 template<class solidType>
312 {
313  if (updated())
314  {
315  return;
316  }
317  // Since we're inside initEvaluate/evaluate there might be processor
318  // comms underway. Change the tag we use.
319  int oldTag = UPstream::msgType();
320  UPstream::msgType() = oldTag+1;
321 
322  const mapDistribute& mapDist = this->mappedPatchBase::map();
323 
324  const label patchi = patch().index();
325 
326  const label nbrPatchi = samplePolyPatch().index();
327 
328  if (baffleActivated_)
329  {
330  const fvPatch& nbrPatch = patch().boundaryMesh()[nbrPatchi];
331 
332  const thermophysicalTransportModel& ttm =
333  db().objectRegistry::template lookupObject
334  <
336  >
337  (
339  (
340  thermophysicalTransportModel::typeName,
341  internalField().group()
342  )
343  );
344 
345  // local properties
346  const scalarField kappaw(ttm.kappaEff(patchi));
347 
348  const fvPatchScalarField& Tp =
349  patch().template lookupPatchField<volScalarField, scalar>(TName_);
350 
351 
352  scalarField qr(Tp.size(), 0.0);
353 
354  if (qrName_ != "none")
355  {
356  qr = patch().template lookupPatchField<volScalarField, scalar>
357  (qrName_);
358 
359  qr = qrRelaxation_*qr + (1.0 - qrRelaxation_)*qrPrevious_;
360  qrPrevious_ = qr;
361  }
362 
363  tmp<scalarField> Ti = patchInternalField();
364 
365  scalarField myKDelta(patch().deltaCoeffs()*kappaw);
366 
367  // nrb properties
368  scalarField nbrTp = ttm.thermo().T().boundaryField()[nbrPatchi];
369  mapDist.distribute(nbrTp);
370 
371  // solid properties
372  scalarField kappas(patch().size(), 0.0);
373  forAll(kappas, i)
374  {
375  kappas[i] = solid().kappa(0.0, (Tp[i] + nbrTp[i])/2.0);
376  }
377 
378  scalarField KDeltaSolid(kappas/baffleThickness());
379 
380  scalarField alpha(KDeltaSolid - qr/Tp);
381 
382  valueFraction() = alpha/(alpha + myKDelta);
383 
384  refValue() = (KDeltaSolid*nbrTp + Qs()/2.0)/alpha;
385 
386  if (debug)
387  {
388  scalar Q = gAverage(kappaw*snGrad());
389  Info<< patch().boundaryMesh().mesh().name() << ':'
390  << patch().name() << ':'
391  << this->internalField().name() << " <- "
392  << nbrPatch.name() << ':'
393  << this->internalField().name() << " :"
394  << " heat[W]:" << Q
395  << " walltemperature "
396  << " min:" << gMin(*this)
397  << " max:" << gMax(*this)
398  << " avg:" << gAverage(*this)
399  << endl;
400  }
401  }
402 
403  // Restore tag
404  UPstream::msgType() = oldTag;
405 
406  mixedFvPatchScalarField::updateCoeffs();
407 }
408 
409 template<class solidType>
411 {
414 
415  if (this->owner())
416  {
417  writeEntry(os, "thickness", baffleThickness()());
418  writeEntry(os, "Qs", Qs()());
419  solid().write(os);
420  }
421 
422  writeEntry(os, "qrPrevious", qrPrevious_);
423  writeEntry(os, "qr", qrName_);
424  writeEntry(os, "qrRelaxation", qrRelaxation_);
425 }
426 
427 
428 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
429 
430 } // End namespace compressible
431 } // End namespace Foam
432 
433 // ************************************************************************* //
Foam::surfaceFields.
const fvBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: fvPatch.H:180
const char *const group
Group name for atomic constants.
dictionary dict
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:643
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
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
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
Type gMin(const FieldField< Field, Type > &f)
virtual const fluidThermo & thermo() const =0
Access function to incompressible transport model.
void distribute(List< T > &fld, const bool dummyTransform=true, const int tag=UPstream::msgType()) const
Distribute data using default commsType.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
virtual void write(Ostream &) const
Write as a dictionary.
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:62
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.
const mapDistribute & map() const
Return reference to the parallel distribution map.
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)
const polyPatch & patch() const
Return the polyPatch.
Definition: fvPatch.H:138
Foam::fvPatchFieldMapper.
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.
virtual label size() const
Return size.
Definition: fvPatch.H:156
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).
label patchi
Class containing processor-to-processor mapping information.
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,.
virtual const volScalarField & T() const =0
Temperature [K].
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:335
Type gAverage(const FieldField< Field, Type > &f)
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:275
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
messageStream Info
virtual const word & name() const
Return name.
Definition: fvPatch.H:144
A class for managing temporary objects.
Definition: PtrList.H:53
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:45
Namespace for OpenFOAM.
IOerror FatalIOError