turbulentTemperatureCoupledBaffleMixedFvPatchScalarField.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 
28 #include "fvPatchFieldMapper.H"
29 #include "volFields.H"
30 #include "mappedPatchBase.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace compressible
37 {
38 
39 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
40 
43 (
44  const fvPatch& p,
46 )
47 :
48  mixedFvPatchScalarField(p, iF),
49  temperatureCoupledBase(patch()),
50  TnbrName_("undefined-Tnbr"),
51  thicknessLayers_(0),
52  kappaLayers_(0),
53  qs_(p.size()),
54  contactRes_(0)
55 {
56  this->refValue() = 0.0;
57  this->refGrad() = 0.0;
58  this->valueFraction() = 1.0;
59 }
60 
61 
64 (
65  const fvPatch& p,
67  const dictionary& dict
68 )
69 :
70  mixedFvPatchScalarField(p, iF),
71  temperatureCoupledBase(patch(), dict),
72  TnbrName_(dict.lookup("Tnbr")),
73  thicknessLayers_(0),
74  kappaLayers_(0),
75  qs_(p.size(), 0),
76  contactRes_(0.0)
77 {
78  if (!isA<mappedPatchBase>(this->patch().patch()))
79  {
81  << "' not type '" << mappedPatchBase::typeName << "'"
82  << "\n for patch " << p.name()
83  << " of field " << internalField().name()
84  << " in file " << internalField().objectPath()
85  << exit(FatalError);
86  }
87 
88  if (dict.found("thicknessLayers"))
89  {
90  dict.lookup("thicknessLayers") >> thicknessLayers_;
91  dict.lookup("kappaLayers") >> kappaLayers_;
92 
93  if (thicknessLayers_.size() > 0)
94  {
95  // Calculate effective thermal resistance by harmonic averaging
96  forAll(thicknessLayers_, iLayer)
97  {
98  contactRes_ += thicknessLayers_[iLayer]/kappaLayers_[iLayer];
99  }
100  contactRes_ = 1.0/contactRes_;
101  }
102  }
103 
104  if (dict.found("qs"))
105  {
106  if (dict.found("Qs"))
107  {
109  << "Either qs or Qs should be specified, not both"
110  << exit(FatalIOError);
111  }
112 
113  qs_ = scalarField("qs", dict, p.size());
114  }
115  else if (dict.found("Qs"))
116  {
117  qs_ = scalarField
118  (
119  p.size(),
120  dict.lookup<scalar>("Qs")/gSum(patch().magSf())
121  );
122  }
123 
124  fvPatchScalarField::operator=(scalarField("value", dict, p.size()));
125 
126  if (dict.found("refValue"))
127  {
128  // Full restart
129  refValue() = scalarField("refValue", dict, p.size());
130  refGrad() = scalarField("refGradient", dict, p.size());
131  valueFraction() = scalarField("valueFraction", dict, p.size());
132  }
133  else
134  {
135  // Start from user entered data. Assume fixedValue.
136  refValue() = *this;
137  refGrad() = 0.0;
138  valueFraction() = 1.0;
139  }
140 }
141 
142 
145 (
147  const fvPatch& p,
149  const fvPatchFieldMapper& mapper
150 )
151 :
152  mixedFvPatchScalarField(ptf, p, iF, mapper),
153  temperatureCoupledBase(patch(), ptf),
154  TnbrName_(ptf.TnbrName_),
155  thicknessLayers_(ptf.thicknessLayers_),
156  kappaLayers_(ptf.kappaLayers_),
157  qs_(mapper(ptf.qs_)),
158  contactRes_(ptf.contactRes_)
159 {}
160 
161 
164 (
167 )
168 :
169  mixedFvPatchScalarField(ptf, iF),
170  temperatureCoupledBase(patch(), ptf),
171  TnbrName_(ptf.TnbrName_),
172  thicknessLayers_(ptf.thicknessLayers_),
173  kappaLayers_(ptf.kappaLayers_),
174  qs_(ptf.qs_),
175  contactRes_(ptf.contactRes_)
176 {}
177 
178 
179 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
180 
182 {
183  if (updated())
184  {
185  return;
186  }
187 
188  // Since we're inside initEvaluate/evaluate there might be processor
189  // comms underway. Change the tag we use.
190  int oldTag = UPstream::msgType();
191  UPstream::msgType() = oldTag + 1;
192 
193  // Get the coupling information from the mappedPatchBase
194  const mappedPatchBase& mpp =
195  refCast<const mappedPatchBase>(patch().patch());
196  const polyMesh& nbrMesh = mpp.sampleMesh();
197  const label samplePatchi = mpp.samplePolyPatch().index();
198  const fvPatch& nbrPatch =
199  refCast<const fvMesh>(nbrMesh).boundary()[samplePatchi];
200 
201  // Calculate the temperature by harmonic averaging
202  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
203 
205 
206  const fvPatchScalarField& nbrTp =
207  nbrPatch.lookupPatchField<volScalarField, scalar>(TnbrName_);
208 
209  if (!isA<thisType>(nbrTp))
210  {
212  << "Patch field for " << internalField().name() << " on "
213  << patch().name() << " is of type " << thisType::typeName
214  << endl << "The neighbouring patch field " << TnbrName_ << " on "
215  << nbrPatch.name() << " is required to be the same, but is "
216  << "currently of type " << nbrTp.type() << exit(FatalError);
217  }
218 
219  const thisType& nbrField = refCast<const thisType>(nbrTp);
220 
221  // Swap to obtain full local values of neighbour internal field
222  tmp<scalarField> nbrIntFld(new scalarField(nbrField.size(), 0.0));
223  tmp<scalarField> nbrKDelta(new scalarField(nbrField.size(), 0.0));
224 
225  if (contactRes_ == 0.0)
226  {
227  nbrIntFld.ref() = nbrField.patchInternalField();
228  nbrKDelta.ref() = nbrField.kappa(nbrField)*nbrPatch.deltaCoeffs();
229  }
230  else
231  {
232  nbrIntFld.ref() = nbrField;
233  nbrKDelta.ref() = contactRes_;
234  }
235 
236  mpp.distribute(nbrIntFld.ref());
237  mpp.distribute(nbrKDelta.ref());
238 
239  tmp<scalarField> myKDelta = kappa(*this)*patch().deltaCoeffs();
240 
241  // Both sides agree on
242  // - temperature : (myKDelta*fld + nbrKDelta*nbrFld)/(myKDelta+nbrKDelta)
243  // - gradient : (temperature-fld)*delta
244  // We've got a degree of freedom in how to implement this in a mixed bc.
245  // (what gradient, what fixedValue and mixing coefficient)
246  // Two reasonable choices:
247  // 1. specify above temperature on one side (preferentially the high side)
248  // and above gradient on the other. So this will switch between pure
249  // fixedvalue and pure fixedgradient
250  // 2. specify gradient and temperature such that the equations are the
251  // same on both sides. This leads to the choice of
252  // - refGradient = qs_/kappa;
253  // - refValue = neighbour value
254  // - mixFraction = nbrKDelta / (nbrKDelta + myKDelta())
255 
256  this->refValue() = nbrIntFld();
257  this->refGrad() = qs_/kappa(*this);
258  this->valueFraction() = nbrKDelta()/(nbrKDelta() + myKDelta());
259 
260  mixedFvPatchScalarField::updateCoeffs();
261 
262  if (debug)
263  {
264  scalar Q = gSum(kappa(*this)*patch().magSf()*snGrad());
265 
266  Info<< patch().boundaryMesh().mesh().name() << ':'
267  << patch().name() << ':'
268  << this->internalField().name() << " <- "
269  << nbrMesh.name() << ':'
270  << nbrPatch.name() << ':'
271  << this->internalField().name() << " :"
272  << " heat transfer rate:" << Q
273  << " walltemperature "
274  << " min:" << gMin(*this)
275  << " max:" << gMax(*this)
276  << " avg:" << gAverage(*this)
277  << endl;
278  }
279 
280  // Restore tag
281  UPstream::msgType() = oldTag;
282 }
283 
284 
286 (
287  Ostream& os
288 ) const
289 {
291  writeEntry(os, "Tnbr", TnbrName_);
292  writeEntry(os, "thicknessLayers", thicknessLayers_);
293  writeEntry(os, "kappaLayers", kappaLayers_);
294 
296 }
297 
298 
299 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
300 
302 (
305 );
306 
307 
308 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
309 
310 } // End namespace compressible
311 } // End namespace Foam
312 
313 
314 // ************************************************************************* //
Mixed boundary condition for temperature, to be used for heat-transfer on back-to-back baffles...
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:663
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
turbulentTemperatureCoupledBaffleMixedFvPatchScalarField(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:124
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
Type gMin(const FieldField< Field, Type > &f)
const polyMesh & sampleMesh() const
Get the region mesh.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
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:63
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
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
Type gSum(const FieldField< Field, Type > &f)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
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.
faceListList boundary(nPatches)
virtual label size() const
Return size.
Definition: fvPatch.H:157
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
const GeometricField::Patch & lookupPatchField(const word &name) const
Lookup and return the patchField of the named field from the.
void distribute(List< Type > &lst) const
Wrapper around map/interpolate data distribution.
Common functions used in temperature coupled boundaries.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:318
tmp< scalarField > kappa(const fvPatchScalarField &Tp) const
Given patch temperature calculate corresponding K field.
Type gAverage(const FieldField< Field, Type > &f)
const scalarField & deltaCoeffs() const
Return the face - cell distance coeffient.
Definition: fvPatch.C:170
label index() const
Return the index of this patch in the boundaryMesh.
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:258
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
messageStream Info
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:76
virtual const word & name() const
Return name.
Definition: fvPatch.H:145
A class for managing temporary objects.
Definition: PtrList.H:53
const polyPatch & samplePolyPatch() const
Get the patch on the region.
void write(Ostream &) const
Write.
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:45
makePatchTypeField(fvPatchScalarField, thermalBaffleFvPatchScalarField)
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:864
IOerror FatalIOError