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 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "volFields.H"
29 #include "surfaceFields.H"
31 #include "mapDistribute.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace compressible
38 {
39 
40 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
41 
42 template<class solidType>
45 (
46  const fvPatch& p,
48 )
49 :
50  mappedPatchBase(p.patch()),
51  mixedFvPatchScalarField(p, iF),
52  TName_("T"),
53  baffleActivated_(true),
54  thickness_(p.size()),
55  qs_(p.size()),
56  solidDict_(),
57  solidPtr_(nullptr),
58  qrPrevious_(p.size()),
59  qrRelaxation_(1),
60  qrName_("undefined-qr")
61 {}
62 
63 
64 template<class solidType>
67 (
69  const fvPatch& p,
71  const fvPatchFieldMapper& mapper
72 )
73 :
74  mappedPatchBase(p.patch(), ptf),
75  mixedFvPatchScalarField(ptf, p, iF, mapper),
76  TName_(ptf.TName_),
77  baffleActivated_(ptf.baffleActivated_),
78  thickness_(ptf.thickness_, mapper),
79  qs_(ptf.qs_, mapper),
80  solidDict_(ptf.solidDict_),
81  solidPtr_(ptf.solidPtr_),
82  qrPrevious_(ptf.qrPrevious_, mapper),
83  qrRelaxation_(ptf.qrRelaxation_),
84  qrName_(ptf.qrName_)
85 {}
86 
87 
88 template<class solidType>
91 (
92  const fvPatch& p,
94  const dictionary& dict
95 )
96 :
97  mappedPatchBase(p.patch(), NEARESTPATCHFACE, dict),
98  mixedFvPatchScalarField(p, iF),
99  TName_("T"),
100  baffleActivated_(dict.lookupOrDefault("baffleActivated", true)),
101  thickness_(),
102  qs_(p.size(), 0),
103  solidDict_(dict),
104  solidPtr_(),
105  qrPrevious_(p.size(), Zero),
106  qrRelaxation_
107  (
108  dict.lookupOrDefaultCompat("qrRelaxation", {{"relaxation", 1712}}, 1)
109  ),
110  qrName_(dict.lookupOrDefault<word>("qr", "none"))
111 {
112  fvPatchScalarField::operator=(scalarField("value", dict, p.size()));
113 
114  if (dict.found("thickness"))
115  {
116  thickness_ = scalarField("thickness", dict, p.size());
117  }
118 
119  if (dict.found("qs"))
120  {
121  qs_ = scalarField("qs", dict, p.size());
122  }
123 
124  if (dict.found("qrPrevious"))
125  {
126  qrPrevious_ = scalarField("qrPrevious", dict, p.size());
127  }
128 
129  if (dict.found("refValue") && baffleActivated_)
130  {
131  // Full restart
132  refValue() = scalarField("refValue", dict, p.size());
133  refGrad() = scalarField("refGradient", dict, p.size());
134  valueFraction() = scalarField("valueFraction", dict, p.size());
135  }
136  else
137  {
138  // Start from user entered data. Assume zeroGradient.
139  refValue() = *this;
140  refGrad() = 0.0;
141  valueFraction() = 0.0;
142  }
143 
144 }
145 
146 
147 template<class solidType>
148 thermalBaffle1DFvPatchScalarField<solidType>::
149 thermalBaffle1DFvPatchScalarField
150 (
152 )
153 :
154  mappedPatchBase(ptf.patch().patch(), ptf),
155  mixedFvPatchScalarField(ptf),
156  TName_(ptf.TName_),
157  baffleActivated_(ptf.baffleActivated_),
158  thickness_(ptf.thickness_),
159  qs_(ptf.qs_),
160  solidDict_(ptf.solidDict_),
161  solidPtr_(ptf.solidPtr_),
162  qrPrevious_(ptf.qrPrevious_),
163  qrRelaxation_(ptf.qrRelaxation_),
164  qrName_(ptf.qrName_)
165 {}
166 
167 
168 template<class solidType>
171 (
174 )
175 :
176  mappedPatchBase(ptf.patch().patch(), ptf),
177  mixedFvPatchScalarField(ptf, iF),
178  TName_(ptf.TName_),
179  baffleActivated_(ptf.baffleActivated_),
180  thickness_(ptf.thickness_),
181  qs_(ptf.qs_),
182  solidDict_(ptf.solidDict_),
183  solidPtr_(ptf.solidPtr_),
184  qrPrevious_(ptf.qrPrevious_),
185  qrRelaxation_(ptf.qrRelaxation_),
186  qrName_(ptf.qrName_)
187 {}
188 
189 
190 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
191 
192 template<class solidType>
194 {
195  const label patchi = patch().index();
196 
197  const label nbrPatchi = samplePolyPatch().index();
198 
199  return (patchi < nbrPatchi);
200 }
201 
202 
203 template<class solidType>
204 const solidType& thermalBaffle1DFvPatchScalarField<solidType>::solid() const
205 {
206  if (this->owner())
207  {
208  if (solidPtr_.empty())
209  {
210  solidPtr_.reset(new solidType(solidDict_));
211  }
212  return *solidPtr_;
213  }
214  else
215  {
216  const fvPatch& nbrPatch =
217  patch().boundaryMesh()[samplePolyPatch().index()];
218 
219  const thermalBaffle1DFvPatchScalarField& nbrField =
220  refCast<const thermalBaffle1DFvPatchScalarField>
221  (
222  nbrPatch.template lookupPatchField<volScalarField, scalar>(TName_)
223  );
224 
225  return nbrField.solid();
226  }
227 }
228 
229 
230 template<class solidType>
231 tmp<scalarField> thermalBaffle1DFvPatchScalarField<solidType>::
232 baffleThickness() const
233 {
234  if (this->owner())
235  {
236  if (thickness_.size() != patch().size())
237  {
238  FatalIOErrorInFunction(solidDict_)
239  << "Field thickness has not been specified"
240  " for patch " << this->patch().name()
241  << exit(FatalIOError);
242  }
243 
244  return thickness_;
245  }
246  else
247  {
248  const mapDistribute& mapDist = this->mappedPatchBase::map();
249 
250  const fvPatch& nbrPatch =
251  patch().boundaryMesh()[samplePolyPatch().index()];
252  const thermalBaffle1DFvPatchScalarField& nbrField =
253  refCast<const thermalBaffle1DFvPatchScalarField>
254  (
255  nbrPatch.template lookupPatchField<volScalarField, scalar>(TName_)
256  );
257 
258  tmp<scalarField> tthickness
259  (
260  new scalarField(nbrField.baffleThickness())
261  );
262  scalarField& thickness = tthickness.ref();
263  mapDist.distribute(thickness);
264  return tthickness;
265  }
266 }
267 
268 
269 template<class solidType>
270 tmp<scalarField> thermalBaffle1DFvPatchScalarField<solidType>::qs() const
271 {
272  if (this->owner())
273  {
274  return qs_;
275  }
276  else
277  {
278  const mapDistribute& mapDist = this->mappedPatchBase::map();
279 
280  const fvPatch& nbrPatch =
281  patch().boundaryMesh()[samplePolyPatch().index()];
282 
283  const thermalBaffle1DFvPatchScalarField& nbrField =
284  refCast<const thermalBaffle1DFvPatchScalarField>
285  (
286  nbrPatch.template lookupPatchField<volScalarField, scalar>(TName_)
287  );
288 
289  tmp<scalarField> tqs(new scalarField(nbrField.qs()));
290  scalarField& qs = tqs.ref();
291  mapDist.distribute(qs);
292  return tqs;
293  }
294 }
295 
296 
297 template<class solidType>
298 void thermalBaffle1DFvPatchScalarField<solidType>::autoMap
299 (
300  const fvPatchFieldMapper& m
301 )
302 {
304 
305  mixedFvPatchScalarField::autoMap(m);
306 
307  if (this->owner())
308  {
309  thickness_.autoMap(m);
310  qs_.autoMap(m);
311  }
312 }
313 
314 
315 template<class solidType>
317 (
318  const fvPatchScalarField& ptf,
319  const labelList& addr
320 )
321 {
322  mixedFvPatchScalarField::rmap(ptf, addr);
323 
324  const thermalBaffle1DFvPatchScalarField& tiptf =
325  refCast<const thermalBaffle1DFvPatchScalarField>(ptf);
326 
327  if (this->owner())
328  {
329  thickness_.rmap(tiptf.thickness_, addr);
330  qs_.rmap(tiptf.qs_, addr);
331  }
332 }
333 
334 
335 template<class solidType>
337 {
338  if (updated())
339  {
340  return;
341  }
342  // Since we're inside initEvaluate/evaluate there might be processor
343  // comms underway. Change the tag we use.
344  int oldTag = UPstream::msgType();
345  UPstream::msgType() = oldTag+1;
346 
347  const mapDistribute& mapDist = this->mappedPatchBase::map();
348 
349  const label patchi = patch().index();
350 
351  const label nbrPatchi = samplePolyPatch().index();
352 
353  if (baffleActivated_)
354  {
355  const fvPatch& nbrPatch = patch().boundaryMesh()[nbrPatchi];
356 
357  const compressible::turbulenceModel& turbModel =
358  db().template lookupObject<compressible::turbulenceModel>
359  (
361  );
362 
363  // local properties
364  const scalarField kappaw(turbModel.kappaEff(patchi));
365 
366  const fvPatchScalarField& Tp =
367  patch().template lookupPatchField<volScalarField, scalar>(TName_);
368 
369 
370  scalarField qr(Tp.size(), Zero);
371 
372  if (qrName_ != "none")
373  {
374  qr = patch().template lookupPatchField<volScalarField, scalar>
375  (qrName_);
376 
377  qr = qrRelaxation_*qr + (1.0 - qrRelaxation_)*qrPrevious_;
378  qrPrevious_ = qr;
379  }
380 
381  tmp<scalarField> Ti = patchInternalField();
382 
383  scalarField myKDelta(patch().deltaCoeffs()*kappaw);
384 
385  // nrb properties
386  scalarField nbrTp =
387  turbModel.transport().T().boundaryField()[nbrPatchi];
388  mapDist.distribute(nbrTp);
389 
390  // solid properties
391  scalarField kappas(patch().size(), Zero);
392  forAll(kappas, i)
393  {
394  kappas[i] = solid().kappa(0.0, (Tp[i] + nbrTp[i])/2.0);
395  }
396 
397  scalarField KDeltaSolid(kappas/baffleThickness());
398 
399  scalarField alpha(KDeltaSolid - qr/Tp);
400 
401  valueFraction() = alpha/(alpha + myKDelta);
402 
403  refValue() = (KDeltaSolid*nbrTp + qs()/2.0)/alpha;
404 
405  if (debug)
406  {
407  scalar Q = gAverage(kappaw*snGrad());
408  Info<< patch().boundaryMesh().mesh().name() << ':'
409  << patch().name() << ':'
410  << this->internalField().name() << " <- "
411  << nbrPatch.name() << ':'
412  << this->internalField().name() << " :"
413  << " heat[W]:" << Q
414  << " walltemperature "
415  << " min:" << gMin(*this)
416  << " max:" << gMax(*this)
417  << " avg:" << gAverage(*this)
418  << endl;
419  }
420  }
421 
422  // Restore tag
423  UPstream::msgType() = oldTag;
424 
425  mixedFvPatchScalarField::updateCoeffs();
426 }
427 
428 template<class solidType>
430 {
433 
434  if (this->owner())
435  {
436  baffleThickness()().writeEntry("thickness", os);
437  qs()().writeEntry("qs", os);
438  solid().write(os);
439  }
440 
441  qrPrevious_.writeEntry("qrPrevious", os);
442  os.writeEntry("qr", qrName_);
443  os.writeEntry("qrRelaxation", qrRelaxation_);
444 }
445 
446 
447 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
448 
449 } // End namespace compressible
450 } // End namespace Foam
451 
452 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::fvPatchField< scalar >
volFields.H
Foam::fvc::snGrad
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:47
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::mappedPatchBase::write
virtual void write(Ostream &) const
Write as a dictionary.
Definition: mappedPatchBase.C:1364
Foam::mappedPatchBase::clearOut
void clearOut()
Definition: mappedPatchBase.C:1200
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::dictionary::found
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Search for an entry (const access) with the given keyword.
Definition: dictionary.C:359
Foam::mappedPatchBase::map
const mapDistribute & map() const
Return reference to the parallel distribution map.
Definition: mappedPatchBaseI.H:146
Foam::gAverage
Type gAverage(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:604
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
Foam::turbulenceModel::propertiesName
static const word propertiesName
Default name of the turbulence properties dictionary.
Definition: turbulenceModel.H:100
Foam::FatalIOError
IOerror FatalIOError
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
surfaceFields.H
Foam::surfaceFields.
Foam::dictionary::lookupOrDefault
T lookupOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionary.H:1241
Foam::fvPatch::name
virtual const word & name() const
Return name.
Definition: fvPatch.H:151
Foam::mappedPatchBase
Determines a mapping between patch face centres and mesh cell or face centres and processors they're ...
Definition: mappedPatchBase.H:105
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::label
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:62
Foam::Field< scalar >
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::mapDistribute
Class containing processor-to-processor mapping information.
Definition: mapDistribute.H:163
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
Foam::mapDistribute::distribute
void distribute(List< T > &fld, const bool dummyTransform=true, const int tag=UPstream::msgType()) const
Distribute data using default commsType.
Definition: mapDistributeTemplates.C:152
compressible
bool compressible
Definition: pEqn.H:3
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::compressible::thermalBaffle1DFvPatchScalarField
This BC solves a steady 1D thermal baffle.
Definition: thermalBaffle1DFvPatchScalarField.H:111
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::ThermalDiffusivity::kappaEff
virtual tmp< volScalarField > kappaEff() const
Return the effective turbulent thermal diffusivity for temperature.
Definition: ThermalDiffusivity.H:146
Foam::UPstream::msgType
static int & msgType()
Message tag of standard messages.
Definition: UPstream.H:491
mapDistribute.H
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::ThermalDiffusivity
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
Definition: phaseCompressibleTurbulenceModelFwd.H:47
Foam::List< label >
Foam::compressible::thermalBaffle1DFvPatchScalarField::thermalBaffle1DFvPatchScalarField
thermalBaffle1DFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: thermalBaffle1DFvPatchScalarField.C:45
Foam::Ostream::writeEntry
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:219
Foam::vtk::write
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Component-wise write of a value (N times)
Definition: foamVtkOutputTemplates.C:35
FatalIOErrorInFunction
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:375
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::gMin
Type gMin(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:593
Foam::fvPatchField< scalar >::operator=
virtual void operator=(const UList< scalar > &)
Definition: fvPatchField.C:379
Foam::gMax
Type gMax(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:592
turbulentFluidThermoModel.H
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54