turbulentTemperatureCoupledBaffleMixedFvPatchScalarField.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-2016 OpenFOAM Foundation
9  Copyright (C) 2019-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 
31 #include "fvPatchFieldMapper.H"
32 #include "volFields.H"
33 #include "mappedPatchBase.H"
34 
35 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39 namespace compressible
40 {
41 
42 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
43 
46 (
47  const fvPatch& p,
49 )
50 :
51  mixedFvPatchScalarField(p, iF),
53  (
54  patch(),
55  "undefined",
56  "undefined",
57  "undefined-K",
58  "undefined-alpha"
59  ),
61  (
63  *this
64  ),
65  TnbrName_("undefined-Tnbr"),
66  thicknessLayers_(0),
67  kappaLayers_(0),
68  contactRes_(0)
69 {
70  this->refValue() = 0.0;
71  this->refGrad() = 0.0;
72  this->valueFraction() = 1.0;
73 }
74 
75 
78 (
80  const fvPatch& p,
82  const fvPatchFieldMapper& mapper
83 )
84 :
85  mixedFvPatchScalarField(ptf, p, iF, mapper),
88  (
90  *this,
91  ptf
92  ),
93  TnbrName_(ptf.TnbrName_),
94  thicknessLayers_(ptf.thicknessLayers_),
95  kappaLayers_(ptf.kappaLayers_),
96  contactRes_(ptf.contactRes_)
97 {}
98 
99 
102 (
103  const fvPatch& p,
105  const dictionary& dict
106 )
107 :
108  mixedFvPatchScalarField(p, iF),
111  (
113  *this,
114  dict
115  ),
116  TnbrName_(dict.get<word>("Tnbr")),
117  thicknessLayers_(0),
118  kappaLayers_(0),
119  contactRes_(0.0)
120 {
121  if (dict.readIfPresent("thicknessLayers", thicknessLayers_))
122  {
123  dict.readEntry("kappaLayers", kappaLayers_);
124 
125  if (thicknessLayers_.size() > 0)
126  {
127  // Calculate effective thermal resistance by harmonic averaging
128  forAll(thicknessLayers_, iLayer)
129  {
130  contactRes_ += thicknessLayers_[iLayer]/kappaLayers_[iLayer];
131  }
132  contactRes_ = 1.0/contactRes_;
133  }
134  }
135 
136  fvPatchScalarField::operator=(scalarField("value", dict, p.size()));
137 
138  if (dict.found("refValue"))
139  {
140  // Full restart
141  refValue() = scalarField("refValue", dict, p.size());
142  refGrad() = scalarField("refGradient", dict, p.size());
143  valueFraction() = scalarField("valueFraction", dict, p.size());
144  }
145  else
146  {
147  // Start from user entered data. Assume fixedValue.
148  refValue() = *this;
149  refGrad() = 0.0;
150  valueFraction() = 1.0;
151  }
152 
153  // Store patch value as initial guess when running in database mode
155  (
156  this->internalField().name(),
157  *this
158  );
160  (
161  this->internalField().name() + "_weights",
162  this->patch().deltaCoeffs()
163  );
164 }
165 
166 
169 (
172 )
173 :
174  mixedFvPatchScalarField(wtcsf, iF),
175  temperatureCoupledBase(patch(), wtcsf),
177  (
179  *this,
180  wtcsf
181  ),
182  TnbrName_(wtcsf.TnbrName_),
183  thicknessLayers_(wtcsf.thicknessLayers_),
184  kappaLayers_(wtcsf.kappaLayers_),
185  contactRes_(wtcsf.contactRes_)
186 {}
187 
188 
191 (
193 )
194 :
195  mixedFvPatchScalarField(wtcsf),
196  temperatureCoupledBase(patch(), wtcsf),
198  (
199  mappedPatchFieldBase<scalar>::mapper(patch(), wtcsf.internalField()),
200  *this,
201  wtcsf
202  ),
203  TnbrName_(wtcsf.TnbrName_),
204  thicknessLayers_(wtcsf.thicknessLayers_),
205  kappaLayers_(wtcsf.kappaLayers_),
206  contactRes_(wtcsf.contactRes_)
207 {}
208 
209 
210 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
211 
213 {
214  if (updated())
215  {
216  return;
217  }
218 
219  // Since we're inside initEvaluate/evaluate there might be processor
220  // comms underway. Change the tag we use.
221  const int oldTag = UPstream::msgType();
222  UPstream::msgType() = oldTag+1;
223 
224  // Get the coupling information from the mappedPatchBase
225  const mappedPatchBase& mpp =
227  (
228  patch(),
229  this->internalField()
230  );
231 
232  const tmp<scalarField> myKDelta = kappa(*this)*patch().deltaCoeffs();
233  tmp<scalarField> nbrIntFld;
234  tmp<scalarField> nbrKDelta;
235 
236  //Pout<< "updateCoeffs() : mpp.sameWorld():" << mpp.sameWorld() << endl;
237 
238 
239  if (mpp.sameWorld())
240  {
241  // Same world so lookup
242  const auto& nbrMesh = refCast<const fvMesh>(mpp.sampleMesh());
243  const label nbrPatchID = mpp.samplePolyPatch().index();
244  const auto& nbrPatch = nbrMesh.boundary()[nbrPatchID];
245 
247  nbrField =
248  refCast
249  <
251  >
252  (
253  nbrPatch.lookupPatchField<volScalarField, scalar>
254  (
255  TnbrName_
256  )
257  );
258 
259  if (contactRes_ == 0.0)
260  {
261  // Get neighbour internal data in local order. Does all
262  // comms/reordering already
263  nbrIntFld = this->mappedInternalField();
264  // Calculate neighbour weighting (in neighbouring ordering)
265  nbrKDelta = nbrField.kappa(nbrField)*nbrPatch.deltaCoeffs();
266  }
267  else
268  {
269  // Get neighbour patch values
270  nbrIntFld = new scalarField(nbrField);
271  // Comms/reorder to local
272  distribute(this->internalField().name(), nbrIntFld.ref());
273 
274  // Constant neighbour weighting. Reorder/comms below
275  nbrKDelta = new scalarField(nbrField.size(), contactRes_);
276  }
277  }
278  else
279  {
280  // Different world so use my region,patch. Distribution below will
281  // do the reordering
282  if (contactRes_ == 0.0)
283  {
284  nbrIntFld = this->mappedInternalField();
285  nbrKDelta = new scalarField(myKDelta());
286  }
287  else
288  {
289  nbrIntFld = *this;
290  nbrKDelta = new scalarField(this->size(), contactRes_);
291  }
292  }
293 
294 
295  //Pout<< "updateCoeffs() : nbrIntFld:" << flatOutput(nbrIntFld()) << endl;
296  //Pout<< "updateCoeffs() : nbrKDelta BEFORE:" << flatOutput(nbrKDelta())
297  // << endl;
298 
299  distribute(this->internalField().name() + "_weights", nbrKDelta.ref());
300  //Pout<< "updateCoeffs() : nbrKDelta AFTER:" << flatOutput(nbrKDelta())
301  // << endl;
302  //Pout<< "updateCoeffs() : myKDelta:" << flatOutput(myKDelta())
303  // << endl;
304 
305 
306  // Both sides agree on
307  // - temperature : (myKDelta*fld + nbrKDelta*nbrFld)/(myKDelta+nbrKDelta)
308  // - gradient : (temperature-fld)*delta
309  // We've got a degree of freedom in how to implement this in a mixed bc.
310  // (what gradient, what fixedValue and mixing coefficient)
311  // Two reasonable choices:
312  // 1. specify above temperature on one side (preferentially the high side)
313  // and above gradient on the other. So this will switch between pure
314  // fixedvalue and pure fixedgradient
315  // 2. specify gradient and temperature such that the equations are the
316  // same on both sides. This leads to the choice of
317  // - refGradient = zero gradient
318  // - refValue = neighbour value
319  // - mixFraction = nbrKDelta / (nbrKDelta + myKDelta())
320 
321  this->refValue() = nbrIntFld();
322  this->refGrad() = 0.0;
323  this->valueFraction() = nbrKDelta()/(nbrKDelta() + myKDelta());
324 
325  mixedFvPatchScalarField::updateCoeffs();
326 
327  if (debug)
328  {
329  scalar Q = gSum(kappa(*this)*patch().magSf()*snGrad());
330 
331  Info<< patch().boundaryMesh().mesh().name() << ':'
332  << patch().name() << ':'
333  << this->internalField().name() << " <- "
334  << mpp.sampleRegion() << ':'
335  << mpp.samplePatch() << ':'
336  << this->internalField().name() << " :"
337  << " heat transfer rate:" << Q
338  << " walltemperature "
339  << " min:" << gMin(*this)
340  << " max:" << gMax(*this)
341  << " avg:" << gAverage(*this)
342  << endl;
343  }
344 
345  // Restore tag
346  UPstream::msgType() = oldTag;
347 }
348 
349 
351 (
352  Ostream& os
353 ) const
354 {
356  os.writeEntry("Tnbr", TnbrName_);
357  thicknessLayers_.writeEntry("thicknessLayers", os);
358  kappaLayers_.writeEntry("kappaLayers", os);
359 
362 }
363 
364 
365 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
366 
368 (
371 );
372 
373 
374 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
375 
376 } // End namespace compressible
377 } // End namespace Foam
378 
379 
380 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: volSurfaceMapping.H:50
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
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:61
Foam::gAverage
Type gAverage(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:604
Foam::mappedPatchBase::samplePolyPatch
const polyPatch & samplePolyPatch() const
Get the patch on the region.
Definition: mappedPatchBase.C:1619
Foam::compressible::turbulentTemperatureCoupledBaffleMixedFvPatchScalarField::turbulentTemperatureCoupledBaffleMixedFvPatchScalarField
turbulentTemperatureCoupledBaffleMixedFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: turbulentTemperatureCoupledBaffleMixedFvPatchScalarField.C:46
Foam::temperatureCoupledBase
Common functions used in temperature coupled boundaries.
Definition: temperatureCoupledBase.H:110
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
fvPatchFieldMapper.H
Foam::mappedPatchBase::sameWorld
bool sameWorld() const
Is world the local world.
Definition: mappedPatchBaseI.H:153
Foam::gSum
Type gSum(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:594
Foam::mappedPatchBase
Determines a mapping between patch face centres and mesh cell or face centres and processors they're ...
Definition: mappedPatchBase.H:112
Foam::compressible::turbulentTemperatureCoupledBaffleMixedFvPatchScalarField::write
virtual void write(Ostream &) const
Write.
Definition: turbulentTemperatureCoupledBaffleMixedFvPatchScalarField.C:351
Foam::compressible::makePatchTypeField
makePatchTypeField(fvPatchScalarField, alphatFixedDmdtWallBoilingWallFunctionFvPatchScalarField)
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::compressible::turbulentTemperatureCoupledBaffleMixedFvPatchScalarField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: turbulentTemperatureCoupledBaffleMixedFvPatchScalarField.C:212
Foam::mappedPatchBase::sampleMesh
const polyMesh & sampleMesh() const
Get the region mesh.
Definition: mappedPatchBase.C:1606
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:228
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::mappedPatchFieldBase< scalar >::distribute
void distribute(const word &fieldName, Field< T > &newValues) const
Definition: mappedPatchFieldBase.C:411
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
Foam::compressible::turbulentTemperatureCoupledBaffleMixedFvPatchScalarField
Mixed boundary condition for temperature, to be used for heat-transfer on back-to-back baffles....
Definition: turbulentTemperatureCoupledBaffleMixedFvPatchScalarField.H:127
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
compressible
bool compressible
Definition: pEqn.H:2
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::mappedPatchFieldBase< scalar >
Foam::mappedPatchFieldBase::mapper
static const mappedPatchBase & mapper(const fvPatch &p, const DimensionedField< Type, volMesh > &iF)
Check that patch is of correct type.
Definition: mappedPatchFieldBase.C:797
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::refCast
To & refCast(From &r)
Reference type cast template function.
Definition: typeInfo.H:131
Foam::mappedPatchBase::sampleRegion
const word & sampleRegion() const
Region to sample.
Definition: mappedPatchBaseI.H:42
turbulentTemperatureCoupledBaffleMixedFvPatchScalarField.H
Foam::mappedPatchFieldBase::write
virtual void write(Ostream &os) const
Write.
Definition: mappedPatchFieldBase.C:842
Foam::UPstream::msgType
static int & msgType()
Message tag of standard messages.
Definition: UPstream.H:541
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::mappedPatchFieldBase< scalar >::mappedInternalField
virtual tmp< Field< scalar > > mappedInternalField() const
Map internal of sampleField onto *this patch.
Definition: mappedPatchFieldBase.C:661
Foam::temperatureCoupledBase::kappa
virtual tmp< scalarField > kappa(const scalarField &Tp) const
Given patch temperature calculate corresponding K field.
Definition: temperatureCoupledBase.C:138
Foam::Ostream::writeEntry
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:232
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
Foam::mappedPatchFieldBase::initRetrieveField
void initRetrieveField(const objectRegistry &obr, const word &region, const word &patch, const mapDistribute &map, const word &fieldName, const Field< T > &fld) const
Construct field from registered elements.
Definition: mappedPatchFieldBase.C:170
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
Foam::temperatureCoupledBase::write
void write(Ostream &os) const
Write.
Definition: temperatureCoupledBase.C:426
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::patchIdentifier::index
label index() const
The index of this patch in the boundaryMesh.
Definition: patchIdentifier.H:158
Foam::gMax
Type gMax(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:592
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54
mappedPatchBase.H
Foam::mappedPatchBase::samplePatch
const word & samplePatch() const
Patch (only if NEARESTPATCHFACE)
Definition: mappedPatchBaseI.H:68