thermalBaffle1DFvPatchScalarField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "volFields.H"
30 #include "surfaceFields.H"
32 #include "mapDistribute.H"
33 
34 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 namespace compressible
39 {
40 
41 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
42 
43 template<class solidType>
46 (
47  const fvPatch& p,
49 )
50 :
51  mappedPatchBase(p.patch()),
52  mixedFvPatchScalarField(p, iF),
53  TName_("T"),
54  baffleActivated_(true),
55  thickness_(p.size()),
56  qs_(p.size()),
57  solidDict_(),
58  solidPtr_(nullptr),
59  qrPrevious_(p.size()),
60  qrRelaxation_(1),
61  qrName_("undefined-qr")
62 {}
63 
64 
65 template<class solidType>
68 (
70  const fvPatch& p,
72  const fvPatchFieldMapper& mapper
73 )
74 :
75  mappedPatchBase(p.patch(), ptf),
76  mixedFvPatchScalarField(ptf, p, iF, mapper),
77  TName_(ptf.TName_),
78  baffleActivated_(ptf.baffleActivated_),
79  thickness_(ptf.thickness_, mapper),
80  qs_(ptf.qs_, mapper),
81  solidDict_(ptf.solidDict_),
82  solidPtr_(ptf.solidPtr_),
83  qrPrevious_(ptf.qrPrevious_, mapper),
84  qrRelaxation_(ptf.qrRelaxation_),
85  qrName_(ptf.qrName_)
86 {}
87 
88 
89 template<class solidType>
92 (
93  const fvPatch& p,
95  const dictionary& dict
96 )
97 :
98  mappedPatchBase(p.patch(), NEARESTPATCHFACE, dict),
99  mixedFvPatchScalarField(p, iF),
100  TName_("T"),
101  baffleActivated_(dict.getOrDefault("baffleActivated", true)),
102  thickness_(),
103  qs_(p.size(), Zero),
104  solidDict_(dict),
105  solidPtr_(),
106  qrPrevious_(p.size(), Zero),
107  qrRelaxation_
108  (
109  dict.getOrDefaultCompat("qrRelaxation", {{"relaxation", 1712}}, 1)
110  ),
111  qrName_(dict.getOrDefault<word>("qr", "none"))
112 {
113  this->readValueEntry(dict, IOobjectOption::MUST_READ);
114 
115  if (dict.found("thickness"))
116  {
117  thickness_ = scalarField("thickness", dict, p.size());
118  }
119 
120  if (dict.found("qs"))
121  {
122  qs_ = scalarField("qs", dict, p.size());
123  }
124 
125  if (dict.found("qrPrevious"))
126  {
127  qrPrevious_ = scalarField("qrPrevious", dict, p.size());
128  }
129 
130  if (baffleActivated_ && this->readMixedEntries(dict))
131  {
132  // Full restart
133  }
134  else
135  {
136  // Start from user entered data. Assume zeroGradient.
137  refValue() = *this;
138  refGrad() = 0.0;
139  valueFraction() = 0.0;
140  }
142 }
143 
144 
145 template<class solidType>
148 (
150 )
151 :
152  mappedPatchBase(ptf.patch().patch(), ptf),
153  mixedFvPatchScalarField(ptf),
154  TName_(ptf.TName_),
155  baffleActivated_(ptf.baffleActivated_),
156  thickness_(ptf.thickness_),
157  qs_(ptf.qs_),
158  solidDict_(ptf.solidDict_),
159  solidPtr_(ptf.solidPtr_),
160  qrPrevious_(ptf.qrPrevious_),
161  qrRelaxation_(ptf.qrRelaxation_),
162  qrName_(ptf.qrName_)
163 {}
164 
165 
166 template<class solidType>
169 (
172 )
173 :
174  mappedPatchBase(ptf.patch().patch(), ptf),
175  mixedFvPatchScalarField(ptf, iF),
176  TName_(ptf.TName_),
177  baffleActivated_(ptf.baffleActivated_),
178  thickness_(ptf.thickness_),
179  qs_(ptf.qs_),
180  solidDict_(ptf.solidDict_),
181  solidPtr_(ptf.solidPtr_),
182  qrPrevious_(ptf.qrPrevious_),
183  qrRelaxation_(ptf.qrRelaxation_),
184  qrName_(ptf.qrName_)
185 {}
186 
187 
188 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
189 
190 template<class solidType>
191 bool thermalBaffle1DFvPatchScalarField<solidType>::owner() const
192 {
193  const label patchi = patch().index();
194 
195  const label nbrPatchi = samplePolyPatch().index();
196 
197  return (patchi < nbrPatchi);
198 }
199 
200 
201 template<class solidType>
202 const solidType& thermalBaffle1DFvPatchScalarField<solidType>::solid() const
203 {
204  if (this->owner())
205  {
206  if (!solidPtr_)
207  {
208  solidPtr_.reset(new solidType(solidDict_));
209  }
210  return *solidPtr_;
211  }
212  else
213  {
214  const fvPatch& nbrPatch =
215  patch().boundaryMesh()[samplePolyPatch().index()];
216 
217  const thermalBaffle1DFvPatchScalarField& nbrField =
218  refCast<const thermalBaffle1DFvPatchScalarField>
219  (
220  nbrPatch.template lookupPatchField<volScalarField>(TName_)
221  );
222 
223  return nbrField.solid();
224  }
225 }
226 
227 
228 template<class solidType>
229 tmp<scalarField> thermalBaffle1DFvPatchScalarField<solidType>::
230 baffleThickness() const
231 {
232  if (this->owner())
233  {
234  if (thickness_.size() != patch().size())
235  {
236  FatalIOErrorInFunction(solidDict_)
237  << "Field thickness has not been specified"
238  " for patch " << this->patch().name()
239  << exit(FatalIOError);
240  }
241 
242  return thickness_;
243  }
244  else
245  {
246  const mapDistribute& mapDist = this->mappedPatchBase::map();
247 
248  const fvPatch& nbrPatch =
249  patch().boundaryMesh()[samplePolyPatch().index()];
250  const thermalBaffle1DFvPatchScalarField& nbrField =
251  refCast<const thermalBaffle1DFvPatchScalarField>
252  (
253  nbrPatch.template lookupPatchField<volScalarField>(TName_)
254  );
255 
256  tmp<scalarField> tthickness
257  (
258  new scalarField(nbrField.baffleThickness())
259  );
260  scalarField& thickness = tthickness.ref();
261  mapDist.distribute(thickness);
262  return tthickness;
263  }
264 }
265 
266 
267 template<class solidType>
268 tmp<scalarField> thermalBaffle1DFvPatchScalarField<solidType>::qs() const
269 {
270  if (this->owner())
271  {
272  return qs_;
273  }
274  else
275  {
276  const mapDistribute& mapDist = this->mappedPatchBase::map();
277 
278  const fvPatch& nbrPatch =
279  patch().boundaryMesh()[samplePolyPatch().index()];
280 
281  const thermalBaffle1DFvPatchScalarField& nbrField =
282  refCast<const thermalBaffle1DFvPatchScalarField>
283  (
284  nbrPatch.template lookupPatchField<volScalarField>(TName_)
285  );
286 
287  tmp<scalarField> tqs(new scalarField(nbrField.qs()));
288  scalarField& qs = tqs.ref();
289  mapDist.distribute(qs);
290  return tqs;
291  }
292 }
293 
294 
295 template<class solidType>
297 (
298  const fvPatchFieldMapper& m
299 )
300 {
302 
303  mixedFvPatchScalarField::autoMap(m);
304 
305  if (this->owner())
306  {
307  thickness_.autoMap(m);
308  qs_.autoMap(m);
309  }
310 }
311 
312 
313 template<class solidType>
315 (
316  const fvPatchScalarField& ptf,
317  const labelList& addr
318 )
319 {
320  mixedFvPatchScalarField::rmap(ptf, addr);
321 
322  const thermalBaffle1DFvPatchScalarField& tiptf =
323  refCast<const thermalBaffle1DFvPatchScalarField>(ptf);
324 
325  if (this->owner())
326  {
327  thickness_.rmap(tiptf.thickness_, addr);
328  qs_.rmap(tiptf.qs_, addr);
329  }
330 }
331 
332 
333 template<class solidType>
335 {
336  if (updated())
337  {
338  return;
339  }
340  // Since we're inside initEvaluate/evaluate there might be processor
341  // comms underway. Change the tag we use.
342  const int oldTag = UPstream::incrMsgType();
343 
344  const mapDistribute& mapDist = this->mappedPatchBase::map();
345 
346  const label patchi = patch().index();
347 
348  const label nbrPatchi = samplePolyPatch().index();
349 
350  if (baffleActivated_)
351  {
352  const fvPatch& nbrPatch = patch().boundaryMesh()[nbrPatchi];
353 
354  const compressible::turbulenceModel& turbModel =
355  db().template lookupObject<compressible::turbulenceModel>
356  (
358  );
359 
360  // local properties
361  const scalarField kappaw(turbModel.kappaEff(patchi));
362 
363  const fvPatchScalarField& Tp =
364  patch().template lookupPatchField<volScalarField>(TName_);
365 
366 
367  scalarField qr(Tp.size(), Zero);
368 
369  if (qrName_ != "none")
370  {
371  qr = patch().template lookupPatchField<volScalarField>(qrName_);
372 
373  qr = qrRelaxation_*qr + (1.0 - qrRelaxation_)*qrPrevious_;
374  qrPrevious_ = qr;
375  }
376 
377  tmp<scalarField> Ti = patchInternalField();
378 
379  scalarField myKDelta(patch().deltaCoeffs()*kappaw);
380 
381  // nrb properties
382  scalarField nbrTp =
383  turbModel.transport().T().boundaryField()[nbrPatchi];
384  mapDist.distribute(nbrTp);
385 
386  // solid properties
387  scalarField kappas(patch().size(), Zero);
388  forAll(kappas, i)
389  {
390  kappas[i] = solid().kappa(0.0, (Tp[i] + nbrTp[i])/2.0);
391  }
392 
393  scalarField KDeltaSolid(kappas/baffleThickness());
394 
395  scalarField alpha(KDeltaSolid - qr/Tp);
396 
397  valueFraction() = alpha/(alpha + myKDelta);
398 
399  refValue() = (KDeltaSolid*nbrTp + qs()/2.0)/alpha;
400 
401  if (debug)
402  {
403  scalar Q = gAverage(kappaw*snGrad());
404  Info<< patch().boundaryMesh().mesh().name() << ':'
405  << patch().name() << ':'
406  << this->internalField().name() << " <- "
407  << nbrPatch.name() << ':'
408  << this->internalField().name() << " :"
409  << " heat[W]:" << Q
410  << " walltemperature "
411  << " min:" << gMin(*this)
412  << " max:" << gMax(*this)
413  << " avg:" << gAverage(*this)
414  << endl;
415  }
416  }
418  UPstream::msgType(oldTag); // Restore tag
419 
420  mixedFvPatchScalarField::updateCoeffs();
421 }
422 
423 template<class solidType>
425 {
428 
429  if (this->owner())
430  {
431  baffleThickness()().writeEntry("thickness", os);
432  qs()().writeEntry("qs", os);
433  solid().write(os);
434  }
435 
436  qrPrevious_.writeEntry("qrPrevious", os);
437  os.writeEntry("qr", qrName_);
438  os.writeEntry("qrRelaxation", qrRelaxation_);
439 }
440 
441 
442 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
443 
444 } // End namespace compressible
445 } // End namespace Foam
446 
447 // ************************************************************************* //
Foam::surfaceFields.
dictionary dict
thermalBaffle1DFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
static int incrMsgType(int val=1) noexcept
Increment the message tag for standard messages.
Definition: UPstream.H:1251
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
Type gMin(const FieldField< Field, Type > &f)
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
static int & msgType() noexcept
Message tag of standard messages.
Definition: UPstream.H:1229
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.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
fvPatchField< scalar > fvPatchScalarField
static const word propertiesName
Default name of the turbulence properties dictionary.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
A FieldMapper for finite-volume patch fields.
Determines a mapping between patch face centres and mesh cell or face centres and processors they&#39;re ...
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
bool compressible
Definition: pEqn.H:2
int debug
Static debugging option.
Type gMax(const FieldField< Field, Type > &f)
OBJstream os(runTime.globalPath()/outputName)
virtual void write(Ostream &os) const
Write as a dictionary.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:627
Type gAverage(const FieldField< Field, Type > &f)
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
const std::string patch
OpenFOAM patch number as a std::string.
messageStream Info
Information stream (stdout output on master, null elsewhere)
tmp< Field< Type > > T() const
Return the field transpose (only defined for second rank tensors)
Definition: Field.C:711
virtual void write(Ostream &) const
Write.
volScalarField & p
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:40
Namespace for OpenFOAM.
IOerror FatalIOError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL IO ERROR&#39; header text and ...
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127