turbulentTemperatureRadCoupledMixedFvPatchScalarField.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  Copyright (C) 2017-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 #include "basicThermo.H"
35 
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40 namespace compressible
41 {
42 
43 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44 
47 (
48  const fvPatch& p,
50 )
51 :
52  mixedFvPatchScalarField(p, iF),
54  (
55  patch(),
56  "undefined",
57  "undefined",
58  "undefined-K",
59  "undefined-alpha"
60  ),
61  TnbrName_("undefined-Tnbr"),
62  qrNbrName_("undefined-qrNbr"),
63  qrName_("undefined-qr"),
64  thicknessLayers_(0),
65  kappaLayers_(0),
66  contactRes_(0),
67  thermalInertia_(false)
68 {
69  this->refValue() = 0.0;
70  this->refGrad() = 0.0;
71  this->valueFraction() = 1.0;
72 }
73 
74 
77 (
79  const fvPatch& p,
81  const fvPatchFieldMapper& mapper
82 )
83 :
84  mixedFvPatchScalarField(psf, p, iF, mapper),
86  TnbrName_(psf.TnbrName_),
87  qrNbrName_(psf.qrNbrName_),
88  qrName_(psf.qrName_),
89  thicknessLayers_(psf.thicknessLayers_),
90  kappaLayers_(psf.kappaLayers_),
91  contactRes_(psf.contactRes_),
92  thermalInertia_(psf.thermalInertia_)
93 {}
94 
95 
98 (
99  const fvPatch& p,
101  const dictionary& dict
102 )
103 :
104  mixedFvPatchScalarField(p, iF),
106  TnbrName_(dict.getOrDefault<word>("Tnbr", "T")),
107  qrNbrName_(dict.getOrDefault<word>("qrNbr", "none")),
108  qrName_(dict.getOrDefault<word>("qr", "none")),
109  thicknessLayers_(0),
110  kappaLayers_(0),
111  contactRes_(0.0),
112  thermalInertia_(dict.getOrDefault<Switch>("thermalInertia", false))
113 {
114  if (!isA<mappedPatchBase>(this->patch().patch()))
115  {
117  << "' not type '" << mappedPatchBase::typeName << "'"
118  << "\n for patch " << p.name()
119  << " of field " << internalField().name()
120  << " in file " << internalField().objectPath()
121  << exit(FatalError);
122  }
123 
124  if (dict.readIfPresent("thicknessLayers", thicknessLayers_))
125  {
126  dict.readEntry("kappaLayers", kappaLayers_);
127 
128  if (thicknessLayers_.size() > 0)
129  {
130  // Calculate effective thermal resistance by harmonic averaging
131  forAll(thicknessLayers_, iLayer)
132  {
133  contactRes_ += thicknessLayers_[iLayer]/kappaLayers_[iLayer];
134  }
135  contactRes_ = 1.0/contactRes_;
136  }
137  }
138 
139  fvPatchScalarField::operator=(scalarField("value", dict, p.size()));
140 
141  if (dict.found("refValue"))
142  {
143  // Full restart
144  refValue() = scalarField("refValue", dict, p.size());
145  refGrad() = scalarField("refGradient", dict, p.size());
146  valueFraction() = scalarField("valueFraction", dict, p.size());
147  }
148  else
149  {
150  // Start from user entered data. Assume fixedValue.
151  refValue() = *this;
152  refGrad() = 0.0;
153  valueFraction() = 1.0;
154  }
155 }
156 
157 
160 (
163 )
164 :
165  mixedFvPatchScalarField(psf, iF),
167  TnbrName_(psf.TnbrName_),
168  qrNbrName_(psf.qrNbrName_),
169  qrName_(psf.qrName_),
170  thicknessLayers_(psf.thicknessLayers_),
171  kappaLayers_(psf.kappaLayers_),
172  contactRes_(psf.contactRes_),
173  thermalInertia_(psf.thermalInertia_)
174 {}
175 
176 
177 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
178 
180 {
181  if (updated())
182  {
183  return;
184  }
185 
186  const polyMesh& mesh = patch().boundaryMesh().mesh();
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 label patchi = patch().index();
195  const mappedPatchBase& mpp =
196  refCast<const mappedPatchBase>(patch().patch());
197  const polyMesh& nbrMesh = mpp.sampleMesh();
198  const label samplePatchi = mpp.samplePolyPatch().index();
199  const fvPatch& nbrPatch =
200  refCast<const fvMesh>(nbrMesh).boundary()[samplePatchi];
201 
202 
203  scalarField Tc(patchInternalField());
204  scalarField& Tp = *this;
205 
207  nbrField = refCast
209  (
210  nbrPatch.lookupPatchField<volScalarField, scalar>(TnbrName_)
211  );
212 
213  // Swap to obtain full local values of neighbour K*delta
214  scalarField KDeltaNbr;
215  tmp<scalarField> TcNbr(new scalarField(nbrField.size(), Zero));
216  if (contactRes_ == 0.0)
217  {
218  TcNbr.ref() = nbrField.patchInternalField();
219  KDeltaNbr = nbrField.kappa(nbrField)*nbrPatch.deltaCoeffs();
220  }
221  else
222  {
223  TcNbr.ref() = nbrField;
224  KDeltaNbr.setSize(nbrField.size(), contactRes_);
225  }
226  mpp.distribute(KDeltaNbr);
227  mpp.distribute(TcNbr.ref());
228 
229  scalarField KDelta(kappa(Tp)*patch().deltaCoeffs());
230 
231  scalarField qr(Tp.size(), Zero);
232  if (qrName_ != "none")
233  {
234  qr = patch().lookupPatchField<volScalarField, scalar>(qrName_);
235  }
236 
237  scalarField qrNbr(Tp.size(), Zero);
238  if (qrNbrName_ != "none")
239  {
240  qrNbr = nbrPatch.lookupPatchField<volScalarField, scalar>(qrNbrName_);
241  mpp.distribute(qrNbr);
242  }
243 
244  // inertia therm
245  if (thermalInertia_)
246  {
247  const scalar dt = mesh.time().deltaTValue();
248  scalarField mCpDtNbr;
249 
250  {
251  const basicThermo* thermo =
253 
254  if (thermo)
255  {
256  const scalarField& ppn =
257  thermo->p().boundaryField()[samplePatchi];
258  const scalarField& Tpn =
259  thermo->T().boundaryField()[samplePatchi];
260 
261  mCpDtNbr =
262  (
263  thermo->Cp(ppn, Tpn, samplePatchi)
264  * thermo->rho(samplePatchi)
265  / nbrPatch.deltaCoeffs()/dt
266  );
267 
268  mpp.distribute(mCpDtNbr);
269  }
270  else
271  {
272  mCpDtNbr.setSize(Tp.size(), Zero);
273  }
274  }
275 
276  scalarField mCpDt;
277 
278  // Local inertia therm
279  {
280  const basicThermo* thermo =
282 
283  if (thermo)
284  {
285  const scalarField& pp = thermo->p().boundaryField()[patchi];
286 
287  mCpDt =
288  (
289  thermo->Cp(pp, Tp, patchi)
290  * thermo->rho(patchi)
291  / patch().deltaCoeffs()/dt
292  );
293  }
294  else
295  {
296  // Issue warning?
297  mCpDt.setSize(Tp.size(), Zero);
298  }
299  }
300 
301  const volScalarField& T =
302  this->db().lookupObject<volScalarField>
303  (
304  this->internalField().name()
305  );
306 
307  const fvPatchField<scalar>& TpOld = T.oldTime().boundaryField()[patchi];
308 
309  scalarField alpha(KDeltaNbr + mCpDt + mCpDtNbr);
310 
311  valueFraction() = alpha/(alpha + KDelta);
312  scalarField c(KDeltaNbr*TcNbr + (mCpDt + mCpDtNbr)*TpOld);
313  refValue() = c/alpha;
314  refGrad() = (qr + qrNbr)/kappa(Tp);
315  }
316  else
317  {
318  valueFraction() = KDeltaNbr/(KDeltaNbr + KDelta);
319  refValue() = TcNbr;
320  refGrad() = (qr + qrNbr)/kappa(Tp);
321  }
322 
323  mixedFvPatchScalarField::updateCoeffs();
324 
325  if (debug)
326  {
327  scalar Q = gSum(kappa(Tp)*patch().magSf()*snGrad());
328 
329  Info<< patch().boundaryMesh().mesh().name() << ':'
330  << patch().name() << ':'
331  << this->internalField().name() << " <- "
332  << nbrMesh.name() << ':'
333  << nbrPatch.name() << ':'
334  << this->internalField().name() << " :"
335  << " heat transfer rate:" << Q
336  << " walltemperature "
337  << " min:" << gMin(Tp)
338  << " max:" << gMax(Tp)
339  << " avg:" << gAverage(Tp)
340  << endl;
341  }
342 
343  // Restore tag
344  UPstream::msgType() = oldTag;
345 }
346 
347 
349 (
350  Ostream& os
351 ) const
352 {
354  os.writeEntry("Tnbr", TnbrName_);
355 
356  os.writeEntry("qrNbr", qrNbrName_);
357  os.writeEntry("qr", qrName_);
358  os.writeEntry("thermalInertia", thermalInertia_);
359 
360  thicknessLayers_.writeEntry("thicknessLayers", os);
361  kappaLayers_.writeEntry("kappaLayers", os);
362 
364 }
365 
366 
367 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
368 
370 (
373 );
374 
375 
376 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
377 
378 } // End namespace compressible
379 } // End namespace Foam
380 
381 
382 // ************************************************************************* //
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
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::Switch
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:77
basicThermo.H
Foam::IOobject::name
const word & name() const
Return name.
Definition: IOobjectI.H:70
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::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
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::turbulentTemperatureRadCoupledMixedFvPatchScalarField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: turbulentTemperatureRadCoupledMixedFvPatchScalarField.C:179
Foam::temperatureCoupledBase
Common functions used in temperature coupled boundaries.
Definition: temperatureCoupledBase.H:110
thermo
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
Foam::basicThermo
Abstract base-class for fluid and solid thermodynamic properties.
Definition: basicThermo.H:63
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
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:167
Foam::mappedPatchBase
Determines a mapping between patch face centres and mesh cell or face centres and processors they're ...
Definition: mappedPatchBase.H:112
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:296
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::Field< scalar >
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:65
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:58
Foam::temperatureCoupledBase::alpha
virtual tmp< scalarField > alpha(const scalarField &Tp) const
Given patch temperature calculate corresponding alphaEff field.
Definition: temperatureCoupledBase.C:284
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Foam::dictionary::dictName
word dictName() const
The local dictionary name (final part of scoped name)
Definition: dictionary.H:458
compressible
bool compressible
Definition: pEqn.H:2
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::compressible::turbulentTemperatureRadCoupledMixedFvPatchScalarField::turbulentTemperatureRadCoupledMixedFvPatchScalarField
turbulentTemperatureRadCoupledMixedFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: turbulentTemperatureRadCoupledMixedFvPatchScalarField.C:47
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
Foam::compressible::turbulentTemperatureRadCoupledMixedFvPatchScalarField::write
virtual void write(Ostream &) const
Write.
Definition: turbulentTemperatureRadCoupledMixedFvPatchScalarField.C:349
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
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
turbulentTemperatureRadCoupledMixedFvPatchScalarField.H
Foam::UPstream::msgType
static int & msgType()
Message tag of standard messages.
Definition: UPstream.H:541
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:381
Foam::objectRegistry::findObject
const Type * findObject(const word &name, const bool recursive=false) const
Return const pointer to the object of the given Type.
Definition: objectRegistryTemplates.C:401
Foam::compressible::turbulentTemperatureRadCoupledMixedFvPatchScalarField
Mixed boundary condition for temperature and radiation heat transfer to be used for in multiregion ca...
Definition: turbulentTemperatureRadCoupledMixedFvPatchScalarField.H:141
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::fvPatch::deltaCoeffs
const scalarField & deltaCoeffs() const
Definition: fvPatch.C:190
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::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
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::GeometricField< scalar, fvPatchField, volMesh >
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