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 | 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 
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(), "undefined", "undefined", "undefined-K"),
50  TnbrName_("undefined-Tnbr"),
51  QrNbrName_("undefined-QrNbr"),
52  QrName_("undefined-Qr"),
53  thicknessLayers_(0),
54  kappaLayers_(0),
55  contactRes_(0)
56 {
57  this->refValue() = 0.0;
58  this->refGrad() = 0.0;
59  this->valueFraction() = 1.0;
60 }
61 
62 
65 (
67  const fvPatch& p,
69  const fvPatchFieldMapper& mapper
70 )
71 :
72  mixedFvPatchScalarField(psf, p, iF, mapper),
73  temperatureCoupledBase(patch(), psf),
74  TnbrName_(psf.TnbrName_),
75  QrNbrName_(psf.QrNbrName_),
76  QrName_(psf.QrName_),
77  thicknessLayers_(psf.thicknessLayers_),
78  kappaLayers_(psf.kappaLayers_),
79  contactRes_(psf.contactRes_)
80 {}
81 
82 
85 (
86  const fvPatch& p,
88  const dictionary& dict
89 )
90 :
91  mixedFvPatchScalarField(p, iF),
92  temperatureCoupledBase(patch(), dict),
93  TnbrName_(dict.lookupOrDefault<word>("Tnbr", "T")),
94  QrNbrName_(dict.lookupOrDefault<word>("QrNbr", "none")),
95  QrName_(dict.lookupOrDefault<word>("Qr", "none")),
96  thicknessLayers_(0),
97  kappaLayers_(0),
98  contactRes_(0.0)
99 {
100  if (!isA<mappedPatchBase>(this->patch().patch()))
101  {
103  << "' not type '" << mappedPatchBase::typeName << "'"
104  << "\n for patch " << p.name()
105  << " of field " << dimensionedInternalField().name()
106  << " in file " << dimensionedInternalField().objectPath()
107  << exit(FatalError);
108  }
109 
110  if (dict.found("thicknessLayers"))
111  {
112  dict.lookup("thicknessLayers") >> thicknessLayers_;
113  dict.lookup("kappaLayers") >> kappaLayers_;
114 
115  if (thicknessLayers_.size() > 0)
116  {
117  // Calculate effective thermal resistance by harmonic averaging
118  forAll (thicknessLayers_, iLayer)
119  {
120  contactRes_ += thicknessLayers_[iLayer]/kappaLayers_[iLayer];
121  }
122  contactRes_ = 1.0/contactRes_;
123  }
124  }
125 
126  fvPatchScalarField::operator=(scalarField("value", dict, p.size()));
127 
128  if (dict.found("refValue"))
129  {
130  // Full restart
131  refValue() = scalarField("refValue", dict, p.size());
132  refGrad() = scalarField("refGradient", dict, p.size());
133  valueFraction() = scalarField("valueFraction", dict, p.size());
134  }
135  else
136  {
137  // Start from user entered data. Assume fixedValue.
138  refValue() = *this;
139  refGrad() = 0.0;
140  valueFraction() = 1.0;
141  }
142 }
143 
144 
147 (
150 )
151 :
152  mixedFvPatchScalarField(psf, iF),
153  temperatureCoupledBase(patch(), psf),
154  TnbrName_(psf.TnbrName_),
155  QrNbrName_(psf.QrNbrName_),
156  QrName_(psf.QrName_),
157  thicknessLayers_(psf.thicknessLayers_),
158  kappaLayers_(psf.kappaLayers_),
159  contactRes_(psf.contactRes_)
160 {}
161 
162 
163 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
164 
166 {
167  if (updated())
168  {
169  return;
170  }
171 
172  // Since we're inside initEvaluate/evaluate there might be processor
173  // comms underway. Change the tag we use.
174  int oldTag = UPstream::msgType();
175  UPstream::msgType() = oldTag+1;
176 
177  // Get the coupling information from the mappedPatchBase
178  const mappedPatchBase& mpp =
179  refCast<const mappedPatchBase>(patch().patch());
180  const polyMesh& nbrMesh = mpp.sampleMesh();
181  const label samplePatchI = mpp.samplePolyPatch().index();
182  const fvPatch& nbrPatch =
183  refCast<const fvMesh>(nbrMesh).boundary()[samplePatchI];
184 
185 
186  scalarField Tc(patchInternalField());
187  scalarField& Tp = *this;
188 
190  nbrField = refCast
192  (
193  nbrPatch.lookupPatchField<volScalarField, scalar>(TnbrName_)
194  );
195 
196  // Swap to obtain full local values of neighbour internal field
197  scalarField TcNbr(nbrField.patchInternalField());
198  mpp.distribute(TcNbr);
199 
200 
201  // Swap to obtain full local values of neighbour K*delta
202  scalarField KDeltaNbr;
203  if (contactRes_ == 0.0)
204  {
205  KDeltaNbr = nbrField.kappa(nbrField)*nbrPatch.deltaCoeffs();
206  }
207  else
208  {
209  KDeltaNbr.setSize(nbrField.size(), contactRes_);
210  }
211  mpp.distribute(KDeltaNbr);
212 
213  scalarField KDelta(kappa(Tp)*patch().deltaCoeffs());
214 
215  scalarField Qr(Tp.size(), 0.0);
216  if (QrName_ != "none")
217  {
218  Qr = patch().lookupPatchField<volScalarField, scalar>(QrName_);
219  }
220 
221  scalarField QrNbr(Tp.size(), 0.0);
222  if (QrNbrName_ != "none")
223  {
224  QrNbr = nbrPatch.lookupPatchField<volScalarField, scalar>(QrNbrName_);
225  mpp.distribute(QrNbr);
226  }
227 
228  valueFraction() = KDeltaNbr/(KDeltaNbr + KDelta);
229  refValue() = TcNbr;
230  refGrad() = (Qr + QrNbr)/kappa(Tp);
231 
232  mixedFvPatchScalarField::updateCoeffs();
233 
234  if (debug)
235  {
236  scalar Q = gSum(kappa(Tp)*patch().magSf()*snGrad());
237 
238  Info<< patch().boundaryMesh().mesh().name() << ':'
239  << patch().name() << ':'
240  << this->dimensionedInternalField().name() << " <- "
241  << nbrMesh.name() << ':'
242  << nbrPatch.name() << ':'
243  << this->dimensionedInternalField().name() << " :"
244  << " heat transfer rate:" << Q
245  << " walltemperature "
246  << " min:" << gMin(Tp)
247  << " max:" << gMax(Tp)
248  << " avg:" << gAverage(Tp)
249  << endl;
250  }
251 
252  // Restore tag
253  UPstream::msgType() = oldTag;
254 }
255 
256 
258 (
259  Ostream& os
260 ) const
261 {
263  os.writeKeyword("Tnbr")<< TnbrName_ << token::END_STATEMENT << nl;
264  os.writeKeyword("QrNbr")<< QrNbrName_ << token::END_STATEMENT << nl;
265  os.writeKeyword("Qr")<< QrName_ << token::END_STATEMENT << nl;
266  thicknessLayers_.writeEntry("thicknessLayers", os);
267  kappaLayers_.writeEntry("kappaLayers", os);
268 
270 }
271 
272 
273 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
274 
276 (
279 );
280 
281 
282 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
283 
284 } // End namespace compressible
285 } // End namespace Foam
286 
287 
288 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
volFields.H
Foam::fvPatchField::operator=
virtual void operator=(const UList< Type > &)
Foam::fvc::snGrad
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:45
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:48
Foam::token::END_STATEMENT
@ END_STATEMENT
Definition: token.H:99
Foam::compressible::makePatchTypeField
makePatchTypeField(fvPatchScalarField, thermalBaffleFvPatchScalarField)
p
p
Definition: pEqn.H:62
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::compressible::turbulentTemperatureRadCoupledMixedFvPatchScalarField::kappaLayers_
scalarList kappaLayers_
Conductivity of layers.
Definition: turbulentTemperatureRadCoupledMixedFvPatchScalarField.H:158
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::gAverage
Type gAverage(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:571
dimensionedInternalField
rDeltaT dimensionedInternalField()
Foam::mappedPatchBase::samplePolyPatch
const polyPatch & samplePolyPatch() const
Get the patch on the region.
Definition: mappedPatchBase.C:1246
Foam::compressible::turbulentTemperatureRadCoupledMixedFvPatchScalarField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: turbulentTemperatureRadCoupledMixedFvPatchScalarField.C:165
Foam::temperatureCoupledBase
Common functions for use in temperature coupled boundaries.
Definition: temperatureCoupledBase.H:104
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
fvPatchFieldMapper.H
Foam::gSum
Type gSum(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:564
Foam::mappedPatchBase
Determines a mapping between patch face centres and mesh cell or face centres and processors they're ...
Definition: mappedPatchBase.H:101
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
Foam::fvPatch::name
const word & name() const
Return name.
Definition: fvPatch.H:149
Foam::Q
This function object calculates and outputs the second invariant of the velocity gradient tensor [1/s...
Definition: Q.H:119
Foam::compressible::turbulentTemperatureRadCoupledMixedFvPatchScalarField::contactRes_
scalar contactRes_
Total contact resistance.
Definition: turbulentTemperatureRadCoupledMixedFvPatchScalarField.H:161
Foam::mappedPatchBase::sampleMesh
const polyMesh & sampleMesh() const
Get the region mesh.
Definition: mappedPatchBase.C:1237
Foam::compressible::turbulentTemperatureRadCoupledMixedFvPatchScalarField::thicknessLayers_
scalarList thicknessLayers_
Thickness of layers.
Definition: turbulentTemperatureRadCoupledMixedFvPatchScalarField.H:155
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:59
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
Foam::mappedPatchBase::distribute
void distribute(List< Type > &lst) const
Wrapper around map/interpolate data distribution.
Definition: mappedPatchBaseTemplates.C:27
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
Qr
volScalarField Qr(IOobject("Qr", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE), mesh, dimensionedScalar("Qr", dimMass/pow3(dimTime), 0.0))
Foam::compressible::turbulentTemperatureRadCoupledMixedFvPatchScalarField::TnbrName_
const word TnbrName_
Name of field on the neighbour region.
Definition: turbulentTemperatureRadCoupledMixedFvPatchScalarField.H:146
Foam::IOobject::name
const word & name() const
Return name.
Definition: IOobject.H:273
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
Foam::compressible::turbulentTemperatureRadCoupledMixedFvPatchScalarField::QrName_
const word QrName_
Name of the radiative heat flux in local region.
Definition: turbulentTemperatureRadCoupledMixedFvPatchScalarField.H:152
Foam::compressible::turbulentTemperatureRadCoupledMixedFvPatchScalarField::QrNbrName_
const word QrNbrName_
Name of the radiative heat flux in the neighbout region.
Definition: turbulentTemperatureRadCoupledMixedFvPatchScalarField.H:149
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:43
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::compressible::turbulentTemperatureRadCoupledMixedFvPatchScalarField::write
virtual void write(Ostream &) const
Write.
Definition: turbulentTemperatureRadCoupledMixedFvPatchScalarField.C:258
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::refCast
To & refCast(From &r)
Reference type cast template function.
Definition: typeInfo.H:106
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
turbulentTemperatureRadCoupledMixedFvPatchScalarField.H
Foam::UPstream::msgType
static int & msgType()
Message tag of standard messages.
Definition: UPstream.H:452
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::compressible::turbulentTemperatureRadCoupledMixedFvPatchScalarField
Mixed boundary condition for temperature and radiation heat transfer to be used for in multiregion ca...
Definition: turbulentTemperatureRadCoupledMixedFvPatchScalarField.H:138
Foam::fvPatch::deltaCoeffs
const scalarField & deltaCoeffs() const
Return the face - cell distance coeffient.
Definition: fvPatch.C:164
Foam::temperatureCoupledBase::kappa
tmp< scalarField > kappa(const scalarField &Tp) const
Given patch temperature calculate corresponding K field.
Definition: temperatureCoupledBase.C:101
compressible
bool compressible
Definition: pEqn.H:40
Foam::Ostream::writeKeyword
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:59
Foam::fvPatch::lookupPatchField
const GeometricField::PatchFieldType & lookupPatchField(const word &name, const GeometricField *=NULL, const Type *=NULL) const
Lookup and return the patchField of the named field from the.
Definition: fvPatchFvMeshTemplates.C:32
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:45
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
Foam::gMin
Type gMin(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:563
Foam::temperatureCoupledBase::write
void write(Ostream &) const
Write.
Definition: temperatureCoupledBase.C:225
Foam::patchIdentifier::index
label index() const
Return the index of this patch in the boundaryMesh.
Definition: patchIdentifier.H:133
write
Tcoeff write()
Foam::gMax
Type gMax(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:562
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:51
mappedPatchBase.H