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 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  ),
60  TnbrName_("undefined-Tnbr"),
61  thicknessLayers_(0),
62  kappaLayers_(0),
63  contactRes_(0)
64 {
65  this->refValue() = 0.0;
66  this->refGrad() = 0.0;
67  this->valueFraction() = 1.0;
68 }
69 
70 
73 (
75  const fvPatch& p,
77  const fvPatchFieldMapper& mapper
78 )
79 :
80  mixedFvPatchScalarField(ptf, p, iF, mapper),
82  TnbrName_(ptf.TnbrName_),
83  thicknessLayers_(ptf.thicknessLayers_),
84  kappaLayers_(ptf.kappaLayers_),
85  contactRes_(ptf.contactRes_)
86 {}
87 
88 
91 (
92  const fvPatch& p,
94  const dictionary& dict
95 )
96 :
97  mixedFvPatchScalarField(p, iF),
99  TnbrName_(dict.get<word>("Tnbr")),
100  thicknessLayers_(0),
101  kappaLayers_(0),
102  contactRes_(0.0)
103 {
104  if (!isA<mappedPatchBase>(this->patch().patch()))
105  {
107  << "' not type '" << mappedPatchBase::typeName << "'"
108  << "\n for patch " << p.name()
109  << " of field " << internalField().name()
110  << " in file " << internalField().objectPath()
111  << exit(FatalError);
112  }
113 
114  if (dict.readIfPresent("thicknessLayers", thicknessLayers_))
115  {
116  dict.readEntry("kappaLayers", kappaLayers_);
117 
118  if (thicknessLayers_.size() > 0)
119  {
120  // Calculate effective thermal resistance by harmonic averaging
121  forAll(thicknessLayers_, iLayer)
122  {
123  contactRes_ += thicknessLayers_[iLayer]/kappaLayers_[iLayer];
124  }
125  contactRes_ = 1.0/contactRes_;
126  }
127  }
128 
129  fvPatchScalarField::operator=(scalarField("value", dict, p.size()));
130 
131  if (dict.found("refValue"))
132  {
133  // Full restart
134  refValue() = scalarField("refValue", dict, p.size());
135  refGrad() = scalarField("refGradient", dict, p.size());
136  valueFraction() = scalarField("valueFraction", dict, p.size());
137  }
138  else
139  {
140  // Start from user entered data. Assume fixedValue.
141  refValue() = *this;
142  refGrad() = 0.0;
143  valueFraction() = 1.0;
144  }
145 }
146 
147 
150 (
153 )
154 :
155  mixedFvPatchScalarField(wtcsf, iF),
156  temperatureCoupledBase(patch(), wtcsf),
157  TnbrName_(wtcsf.TnbrName_),
158  thicknessLayers_(wtcsf.thicknessLayers_),
159  kappaLayers_(wtcsf.kappaLayers_),
160  contactRes_(wtcsf.contactRes_)
161 {}
162 
163 
164 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
165 
167 {
168  if (updated())
169  {
170  return;
171  }
172 
173  // Since we're inside initEvaluate/evaluate there might be processor
174  // comms underway. Change the tag we use.
175  int oldTag = UPstream::msgType();
176  UPstream::msgType() = oldTag+1;
177 
178  // Get the coupling information from the mappedPatchBase
179  const mappedPatchBase& mpp =
180  refCast<const mappedPatchBase>(patch().patch());
181  const polyMesh& nbrMesh = mpp.sampleMesh();
182  const label samplePatchi = mpp.samplePolyPatch().index();
183  const fvPatch& nbrPatch =
184  refCast<const fvMesh>(nbrMesh).boundary()[samplePatchi];
185 
186  // Calculate the temperature by harmonic averaging
187  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
188 
190  refCast
191  <
193  >
194  (
195  nbrPatch.lookupPatchField<volScalarField, scalar>
196  (
197  TnbrName_
198  )
199  );
200 
201  // Swap to obtain full local values of neighbour internal field
202  tmp<scalarField> nbrIntFld(new scalarField(nbrField.size(), Zero));
203  tmp<scalarField> nbrKDelta(new scalarField(nbrField.size(), Zero));
204 
205  if (contactRes_ == 0.0)
206  {
207  nbrIntFld.ref() = nbrField.patchInternalField();
208  nbrKDelta.ref() = nbrField.kappa(nbrField)*nbrPatch.deltaCoeffs();
209  }
210  else
211  {
212  nbrIntFld.ref() = nbrField;
213  nbrKDelta.ref() = contactRes_;
214  }
215 
216  mpp.distribute(nbrIntFld.ref());
217  mpp.distribute(nbrKDelta.ref());
218 
219  tmp<scalarField> myKDelta = kappa(*this)*patch().deltaCoeffs();
220 
221 
222  // Both sides agree on
223  // - temperature : (myKDelta*fld + nbrKDelta*nbrFld)/(myKDelta+nbrKDelta)
224  // - gradient : (temperature-fld)*delta
225  // We've got a degree of freedom in how to implement this in a mixed bc.
226  // (what gradient, what fixedValue and mixing coefficient)
227  // Two reasonable choices:
228  // 1. specify above temperature on one side (preferentially the high side)
229  // and above gradient on the other. So this will switch between pure
230  // fixedvalue and pure fixedgradient
231  // 2. specify gradient and temperature such that the equations are the
232  // same on both sides. This leads to the choice of
233  // - refGradient = zero gradient
234  // - refValue = neighbour value
235  // - mixFraction = nbrKDelta / (nbrKDelta + myKDelta())
236 
237  this->refValue() = nbrIntFld();
238  this->refGrad() = 0.0;
239  this->valueFraction() = nbrKDelta()/(nbrKDelta() + myKDelta());
240 
241  mixedFvPatchScalarField::updateCoeffs();
242 
243  if (debug)
244  {
245  scalar Q = gSum(kappa(*this)*patch().magSf()*snGrad());
246 
247  Info<< patch().boundaryMesh().mesh().name() << ':'
248  << patch().name() << ':'
249  << this->internalField().name() << " <- "
250  << nbrMesh.name() << ':'
251  << nbrPatch.name() << ':'
252  << this->internalField().name() << " :"
253  << " heat transfer rate:" << Q
254  << " walltemperature "
255  << " min:" << gMin(*this)
256  << " max:" << gMax(*this)
257  << " avg:" << gAverage(*this)
258  << endl;
259  }
260 
261  // Restore tag
262  UPstream::msgType() = oldTag;
263 }
264 
265 
267 (
268  Ostream& os
269 ) const
270 {
272  os.writeEntry("Tnbr", TnbrName_);
273  thicknessLayers_.writeEntry("thicknessLayers", os);
274  kappaLayers_.writeEntry("kappaLayers", os);
275 
277 }
278 
279 
280 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
281 
283 (
286 );
287 
288 
289 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
290 
291 } // End namespace compressible
292 } // End namespace Foam
293 
294 
295 // ************************************************************************* //
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::IOobject::name
const word & name() const
Return name.
Definition: IOobjectI.H:46
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::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:1223
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:337
fvPatchFieldMapper.H
Foam::gSum
Type gSum(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:594
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
Foam::compressible::turbulentTemperatureCoupledBaffleMixedFvPatchScalarField::write
virtual void write(Ostream &) const
Write.
Definition: turbulentTemperatureCoupledBaffleMixedFvPatchScalarField.C:267
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
Foam::compressible::makePatchTypeField
makePatchTypeField(fvPatchScalarField, alphatFixedDmdtWallBoilingWallFunctionFvPatchScalarField)
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::compressible::turbulentTemperatureCoupledBaffleMixedFvPatchScalarField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: turbulentTemperatureCoupledBaffleMixedFvPatchScalarField.C:166
Foam::mappedPatchBase::sampleMesh
const polyMesh & sampleMesh() const
Get the region mesh.
Definition: mappedPatchBase.C:1210
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:258
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::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::mappedPatchBase::distribute
void distribute(List< Type > &lst) const
Wrapper around map/interpolate data distribution.
Definition: mappedPatchBaseTemplates.C:29
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
Foam::compressible::turbulentTemperatureCoupledBaffleMixedFvPatchScalarField
Mixed boundary condition for temperature, to be used for heat-transfer on back-to-back baffles....
Definition: turbulentTemperatureCoupledBaffleMixedFvPatchScalarField.H:126
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
compressible
bool compressible
Definition: pEqn.H:3
Foam::fvPatch::lookupPatchField
const GeometricField::Patch & lookupPatchField(const word &name, const GeometricField *=nullptr, const Type *=nullptr) const
Definition: fvPatchFvMeshTemplates.C:34
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
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::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
turbulentTemperatureCoupledBaffleMixedFvPatchScalarField.H
Foam::UPstream::msgType
static int & msgType()
Message tag of standard messages.
Definition: UPstream.H:491
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::fvPatch::deltaCoeffs
const scalarField & deltaCoeffs() const
Definition: fvPatch.C:160
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: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
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:133
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