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 | Copyright (C) 2011-2015 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"
28 #include "mappedPatchBase.H"
30 #include "mapDistribute.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 )
48 :
50  mixedFvPatchScalarField(p, iF),
51  TName_("T"),
52  baffleActivated_(true),
53  thickness_(p.size()),
54  Qs_(p.size()),
55  solidDict_(),
56  solidPtr_(NULL),
57  QrPrevious_(p.size()),
58  QrRelaxation_(1),
59  QrName_("undefined-Qr")
60 {}
61 
62 
63 template<class solidType>
66 (
68  const fvPatch& p,
70  const fvPatchFieldMapper& mapper
71 )
72 :
73  mappedPatchBase(p.patch(), ptf),
74  mixedFvPatchScalarField(ptf, p, iF, mapper),
75  TName_(ptf.TName_),
76  baffleActivated_(ptf.baffleActivated_),
77  thickness_(ptf.thickness_, mapper),
78  Qs_(ptf.Qs_, mapper),
79  solidDict_(ptf.solidDict_),
80  solidPtr_(ptf.solidPtr_),
81  QrPrevious_(ptf.QrPrevious_, mapper),
82  QrRelaxation_(ptf.QrRelaxation_),
83  QrName_(ptf.QrName_)
84 {}
85 
86 
87 template<class solidType>
90 (
91  const fvPatch& p,
93  const dictionary& dict
94 )
95 :
96  mappedPatchBase(p.patch(), NEARESTPATCHFACE, dict),
97  mixedFvPatchScalarField(p, iF),
98  TName_("T"),
99  baffleActivated_(dict.lookupOrDefault<bool>("baffleActivated", true)),
100  thickness_(),
101  Qs_(p.size(), 0),
102  solidDict_(dict),
103  solidPtr_(),
104  QrPrevious_(p.size(), 0.0),
105  QrRelaxation_(dict.lookupOrDefault<scalar>("relaxation", 1)),
106  QrName_(dict.lookupOrDefault<word>("Qr", "none"))
107 {
108  fvPatchScalarField::operator=(scalarField("value", dict, p.size()));
109 
110  if (dict.found("thickness"))
111  {
112  thickness_ = scalarField("thickness", dict, p.size());
113  }
114 
115  if (dict.found("Qs"))
116  {
117  Qs_ = scalarField("Qs", dict, p.size());
118  }
119 
120  if (dict.found("QrPrevious"))
121  {
122  QrPrevious_ = scalarField("QrPrevious", dict, p.size());
123  }
124 
125  if (dict.found("refValue") && baffleActivated_)
126  {
127  // Full restart
128  refValue() = scalarField("refValue", dict, p.size());
129  refGrad() = scalarField("refGradient", dict, p.size());
130  valueFraction() = scalarField("valueFraction", dict, p.size());
131  }
132  else
133  {
134  // Start from user entered data. Assume zeroGradient.
135  refValue() = *this;
136  refGrad() = 0.0;
137  valueFraction() = 0.0;
138  }
139 
140 }
141 
142 
143 template<class solidType>
146 (
148 )
149 :
150  mappedPatchBase(ptf.patch().patch(), ptf),
151  mixedFvPatchScalarField(ptf),
152  TName_(ptf.TName_),
153  baffleActivated_(ptf.baffleActivated_),
154  thickness_(ptf.thickness_),
155  Qs_(ptf.Qs_),
156  solidDict_(ptf.solidDict_),
157  solidPtr_(ptf.solidPtr_),
158  QrPrevious_(ptf.QrPrevious_),
159  QrRelaxation_(ptf.QrRelaxation_),
160  QrName_(ptf.QrName_)
161 {}
162 
163 
164 template<class solidType>
167 (
170 )
171 :
172  mappedPatchBase(ptf.patch().patch(), ptf),
173  mixedFvPatchScalarField(ptf, iF),
174  TName_(ptf.TName_),
175  baffleActivated_(ptf.baffleActivated_),
176  thickness_(ptf.thickness_),
177  Qs_(ptf.Qs_),
178  solidDict_(ptf.solidDict_),
179  solidPtr_(ptf.solidPtr_),
180  QrPrevious_(ptf.QrPrevious_),
181  QrRelaxation_(ptf.QrRelaxation_),
182  QrName_(ptf.QrName_)
183 {}
184 
185 
186 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
187 
188 template<class solidType>
190 {
191  const label patchi = patch().index();
192 
193  const label nbrPatchi = samplePolyPatch().index();
194 
195  return (patchi < nbrPatchi);
196 }
197 
198 
199 template<class solidType>
201 {
202  if (this->owner())
203  {
204  if (solidPtr_.empty())
205  {
206  solidPtr_.reset(new solidType(solidDict_));
207  }
208  return solidPtr_();
209  }
210  else
211  {
212  const fvPatch& nbrPatch =
213  patch().boundaryMesh()[samplePolyPatch().index()];
214 
215  const thermalBaffle1DFvPatchScalarField& nbrField =
216  refCast<const thermalBaffle1DFvPatchScalarField>
217  (
218  nbrPatch.template lookupPatchField<volScalarField, scalar>(TName_)
219  );
220 
221  return nbrField.solid();
222  }
223 }
224 
225 
226 template<class solidType>
228 baffleThickness() const
229 {
230  if (this->owner())
231  {
232  if (thickness_.size() != patch().size())
233  {
235  (
236  " template<class solidType>"
237  " tmp<scalarField> thermalBaffle1DFvPatchScalarField<solidType>"
238  " baffleThickness() const",
239  solidDict_
240  )<< " Field thickness has not been specified "
241  << " for patch " << this->patch().name()
242  << exit(FatalIOError);
243  }
244 
245  return thickness_;
246  }
247  else
248  {
249  const mapDistribute& mapDist = this->mappedPatchBase::map();
250 
251  const fvPatch& nbrPatch =
252  patch().boundaryMesh()[samplePolyPatch().index()];
253  const thermalBaffle1DFvPatchScalarField& nbrField =
254  refCast<const thermalBaffle1DFvPatchScalarField>
255  (
256  nbrPatch.template lookupPatchField<volScalarField, scalar>(TName_)
257  );
258 
259  tmp<scalarField> tthickness
260  (
261  new scalarField(nbrField.baffleThickness())
262  );
263  scalarField& thickness = tthickness();
264  mapDist.distribute(thickness);
265  return tthickness;
266  }
267 }
268 
269 
270 template<class solidType>
272 {
273  if (this->owner())
274  {
275  return Qs_;
276  }
277  else
278  {
279  const mapDistribute& mapDist = this->mappedPatchBase::map();
280 
281  const fvPatch& nbrPatch =
282  patch().boundaryMesh()[samplePolyPatch().index()];
283 
284  const thermalBaffle1DFvPatchScalarField& nbrField =
285  refCast<const thermalBaffle1DFvPatchScalarField>
286  (
287  nbrPatch.template lookupPatchField<volScalarField, scalar>(TName_)
288  );
289 
290  tmp<scalarField> tQs(new scalarField(nbrField.Qs()));
291  scalarField& Qs = tQs();
292  mapDist.distribute(Qs);
293  return tQs;
294  }
295 }
296 
297 
298 template<class solidType>
300 (
301  const fvPatchFieldMapper& m
302 )
303 {
304  mixedFvPatchScalarField::autoMap(m);
305 
306  if (this->owner())
307  {
308  thickness_.autoMap(m);
309  Qs_.autoMap(m);
310  }
311 }
312 
313 
314 template<class solidType>
316 (
317  const fvPatchScalarField& ptf,
318  const labelList& addr
319 )
320 {
321  mixedFvPatchScalarField::rmap(ptf, addr);
322 
323  const thermalBaffle1DFvPatchScalarField& tiptf =
324  refCast<const thermalBaffle1DFvPatchScalarField>(ptf);
325 
326  if (this->owner())
327  {
328  thickness_.rmap(tiptf.thickness_, addr);
329  Qs_.rmap(tiptf.Qs_, addr);
330  }
331 }
332 
333 
334 template<class solidType>
336 {
337  if (updated())
338  {
339  return;
340  }
341  // Since we're inside initEvaluate/evaluate there might be processor
342  // comms underway. Change the tag we use.
343  int oldTag = UPstream::msgType();
344  UPstream::msgType() = oldTag+1;
345 
346  const mapDistribute& mapDist = this->mappedPatchBase::map();
347 
348  const label patchi = patch().index();
349 
350  const label nbrPatchi = samplePolyPatch().index();
351 
352  if (baffleActivated_)
353  {
354  const fvPatch& nbrPatch = patch().boundaryMesh()[nbrPatchi];
355 
356  const compressible::turbulenceModel& turbModel =
357  db().template lookupObject<compressible::turbulenceModel>
358  (
360  );
361 
362  // local properties
363  const scalarField kappaw(turbModel.kappaEff(patchi));
364 
365  const fvPatchScalarField& Tp =
366  patch().template lookupPatchField<volScalarField, scalar>(TName_);
367 
368 
369  scalarField Qr(Tp.size(), 0.0);
370 
371  if (QrName_ != "none")
372  {
373  Qr = patch().template lookupPatchField<volScalarField, scalar>
374  (QrName_);
375 
376  Qr = QrRelaxation_*Qr + (1.0 - QrRelaxation_)*QrPrevious_;
377  QrPrevious_ = Qr;
378  }
379 
380  tmp<scalarField> Ti = patchInternalField();
381 
382  scalarField myKDelta(patch().deltaCoeffs()*kappaw);
383 
384  // nrb properties
385  scalarField nbrTp =
386  turbModel.transport().T().boundaryField()[nbrPatchi];
387  mapDist.distribute(nbrTp);
388 
389  // solid properties
390  scalarField kappas(patch().size(), 0.0);
391  forAll(kappas, i)
392  {
393  kappas[i] = solid().kappa(0.0, (Tp[i] + nbrTp[i])/2.0);
394  }
395 
396  scalarField KDeltaSolid(kappas/baffleThickness());
397 
398  scalarField alpha(KDeltaSolid - Qr/Tp);
399 
400  valueFraction() = alpha/(alpha + myKDelta);
401 
402  refValue() = (KDeltaSolid*nbrTp + Qs()/2.0)/alpha;
403 
404  if (debug)
405  {
406  scalar Q = gAverage(kappaw*snGrad());
407  Info<< patch().boundaryMesh().mesh().name() << ':'
408  << patch().name() << ':'
409  << this->dimensionedInternalField().name() << " <- "
410  << nbrPatch.name() << ':'
411  << this->dimensionedInternalField().name() << " :"
412  << " heat[W]:" << Q
413  << " walltemperature "
414  << " min:" << gMin(*this)
415  << " max:" << gMax(*this)
416  << " avg:" << gAverage(*this)
417  << endl;
418  }
419  }
420 
421  // Restore tag
422  UPstream::msgType() = oldTag;
423 
424  mixedFvPatchScalarField::updateCoeffs();
425 }
426 
427 template<class solidType>
429 {
432 
433  if (this->owner())
434  {
435  baffleThickness()().writeEntry("thickness", os);
436  Qs()().writeEntry("Qs", os);
437  solid().write(os);
438  }
439 
440  QrPrevious_.writeEntry("QrPrevious", os);
441  os.writeKeyword("Qr")<< QrName_ << token::END_STATEMENT << nl;
442  os.writeKeyword("relaxation")<< QrRelaxation_
443  << token::END_STATEMENT << nl;
444 }
445 
446 
447 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
448 
449 } // End namespace compressible
450 } // End namespace Foam
451 
452 // ************************************************************************* //
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:306
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:387
Foam::surfaceFields.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
virtual tmp< volScalarField > kappaEff() const
Return the effective turbulent thermal diffusivity for temperature.
static void distribute(const Pstream::commsTypes commsType, const List< labelPair > &schedule, const label constructSize, const labelListList &subMap, const labelListList &constructMap, List< T > &, const int tag=UPstream::msgType())
Distribute data. Note:schedule only used for Pstream::scheduled.
Class containing processor-to-processor mapping information.
A class for handling words, derived from string.
Definition: word.H:59
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
Foam::fvPatchFieldMapper.
messageStream Info
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
Namespace for OpenFOAM.
const word & name() const
Return name.
Definition: fvPatch.H:149
runTime write()
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const mapDistribute & map() const
Return reference to the parallel distribution map.
const polyPatch & patch() const
Return the polyPatch.
Definition: fvPatch.H:143
dictionary dict
static const char nl
Definition: Ostream.H:260
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
IOerror FatalIOError
thermalBaffle1DFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
const fvBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: fvPatch.H:185
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
#define forAll(list, i)
Definition: UList.H:421
label patchi
static int & msgType()
Message tag of standard messages.
Definition: UPstream.H:451
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:45
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
Type gMin(const FieldField< Field, Type > &f)
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:59
rDeltaT dimensionedInternalField()
bool compressible
Definition: pEqn.H:40
Type gAverage(const FieldField< Field, Type > &f)
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
static const word propertiesName
Default name of the turbulence properties dictionary.
Determines a mapping between patch face centres and mesh cell or face centres and processors they&#39;re ...
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
virtual label size() const
Return size.
Definition: fvPatch.H:161
This function object calculates and outputs the second invariant of the velocity gradient tensor [1/s...
Definition: Q.H:66
#define FatalIOErrorIn(functionName, ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:325
Type gMax(const FieldField< Field, Type > &f)
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
A class for managing temporary objects.
Definition: PtrList.H:118
virtual void write(Ostream &) const
Write as a dictionary.