faceZoneReferenceTemperature.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) 2022 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
29 #include "surfaceInterpolate.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace heatTransferCoeffModels
37 {
40  (
44  );
45 }
46 }
47 
48 
49 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
50 
51 void Foam::heatTransferCoeffModels::faceZoneReferenceTemperature::
52 setFaceZoneFaces(const dictionary& dict)
53 {
54  const auto& mesh =
55  mesh_.objectRegistry::db().lookupObject<fvMesh>(refRegionName_);
56 
57  const word faceZoneName(dict.get<word>("referenceFaceZone"));
58 
59  faceZonei_ = mesh.faceZones().findZoneID(faceZoneName);
60 
61  if (faceZonei_ < 0)
62  {
64  << "referenceFaceZone: " << faceZoneName
65  << " does not exist in referenceRegion: " << refRegionName_
66  << exit(FatalIOError);
67  }
68 
69  const faceZone& fZone = mesh.faceZones()[faceZonei_];
70 
71  label numFaces = fZone.size();
72 
73  if (!returnReduceOr(numFaces))
74  {
76  << "referenceFaceZone: " << faceZoneName
77  << " contains no faces."
78  << exit(FatalIOError);
79  }
80 
81  faceId_.resize_nocopy(numFaces);
82  facePatchId_.resize_nocopy(numFaces);
83 
84  numFaces = 0;
85 
86  // TDB: handle multiple zones
87  {
88  forAll(fZone, i)
89  {
90  const label meshFacei = fZone[i];
91 
92  // Internal faces
93  label faceId = meshFacei;
94  label facePatchId = -1;
95 
96  // Boundary faces
97  if (!mesh.isInternalFace(meshFacei))
98  {
99  facePatchId = mesh.boundaryMesh().whichPatch(meshFacei);
100  const polyPatch& pp = mesh.boundaryMesh()[facePatchId];
101 
102  if (isA<emptyPolyPatch>(pp))
103  {
104  continue; // Ignore empty patch
105  }
106 
107  const auto* cpp = isA<coupledPolyPatch>(pp);
108 
109  if (cpp && !cpp->owner())
110  {
111  continue; // Ignore neighbour side
112  }
113 
114  faceId = pp.whichFace(meshFacei);
115  }
116 
117  if (faceId >= 0)
118  {
119  faceId_[numFaces] = faceId;
120  facePatchId_[numFaces] = facePatchId;
121 
122  ++numFaces;
123  }
124  }
125  }
126 
127  // Shrink to size used
128  faceId_.resize(numFaces);
129  facePatchId_.resize(numFaces);
130 }
131 
132 
133 Foam::scalar Foam::heatTransferCoeffModels::faceZoneReferenceTemperature::
134 faceZoneAverageTemperature()
135 {
136  const auto& mesh =
137  mesh_.objectRegistry::db().lookupObject<fvMesh>(refRegionName_);
138 
139  const auto& T = mesh.lookupObject<volScalarField>(TName_);
141 
142  const surfaceScalarField& magSf = mesh.magSf();
143 
144  scalar Tmean = 0;
145  scalar sumMagSf = 0;
146 
147  forAll(faceId_, i)
148  {
149  const label facei = faceId_[i];
150  if (facePatchId_[i] != -1)
151  {
152  const label patchi = facePatchId_[i];
153  const scalar sf = magSf.boundaryField()[patchi][facei];
154 
155  Tmean += Tf.boundaryField()[patchi][facei]*sf;
156  sumMagSf += sf;
157  }
158  else
159  {
160  const scalar sf = magSf[facei];
161  Tmean += Tf[facei]*sf;
162  sumMagSf += sf;
163  }
164  }
165  reduce(Tmean, sumOp<scalar>());
166  reduce(sumMagSf, sumOp<scalar>());
167 
168  Tmean /= sumMagSf;
169 
170  return Tmean;
171 }
172 
173 
174 void Foam::heatTransferCoeffModels::faceZoneReferenceTemperature::htc
175 (
176  volScalarField& htc,
177  const FieldField<Field, scalar>& q
178 )
179 {
180  // Retrieve temperature boundary fields for current region
181  const auto& T = mesh_.lookupObject<volScalarField>(TName_);
182  const volScalarField::Boundary& Tbf = T.boundaryField();
183 
184  // Retrieve heat-transfer coefficient boundary fields for current region
185  volScalarField::Boundary& htcBf = htc.boundaryFieldRef();
186 
187  // Calculate area-averaged temperature field
188  // for the reference face zone and region
189  // (reference region can be different from current region)
190  const scalar Tref = faceZoneAverageTemperature();
191 
192  // Calculate heat-transfer coefficient boundary fields for current region
193  for (const label patchi : patchIDs_)
194  {
195  htcBf[patchi] = q[patchi]/(Tref - Tbf[patchi] + ROOTVSMALL);
196  }
197 }
198 
199 
200 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
201 
204 (
205  const dictionary& dict,
206  const fvMesh& mesh,
207  const word& TName
208 )
209 :
211  faceZonei_(-1),
212  refRegionName_(polyMesh::defaultRegion),
213  faceId_(),
214  facePatchId_()
215 {
217 }
218 
219 
220 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
221 
223 (
224  const dictionary& dict
225 )
226 {
228  {
229  return false;
230  }
231 
232  dict.readIfPresent("referenceRegion", refRegionName_);
233 
234  setFaceZoneFaces(dict);
235 
236  return true;
237 }
238 
239 
240 // ************************************************************************* //
dictionary dict
const Type & lookupObject(const word &name, const bool recursive=false) const
Lookup and return const reference to the object of the given Type. Fatal if not found or the wrong ty...
label faceId(-1)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
void resize(const label len)
Adjust allocated size of list.
Definition: ListI.H:160
defineTypeNameAndDebug(faceZoneReferenceTemperature, 0)
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
faceZoneReferenceTemperature(const dictionary &dict, const fvMesh &mesh, const word &TName)
Construct from components.
void resize_nocopy(const label len)
Adjust allocated size of list without necessarily.
Definition: ListI.H:175
GeometricBoundaryField< scalar, fvPatchField, volMesh > Boundary
Type of boundary fields.
virtual bool read(const dictionary &dict)
Read from dictionary.
Macros for easy insertion into run-time selection tables.
bool isInternalFace(const label faceIndex) const noexcept
Return true if given face label is internal to the mesh.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
A base class for heat transfer coefficient models.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
dynamicFvMesh & mesh
label findZoneID(const word &zoneName) const
Find zone index by name, return -1 if not found.
Definition: ZoneMesh.C:629
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition: polyMesh.H:608
A class for handling words, derived from Foam::string.
Definition: word.H:63
label size() const noexcept
The number of entries in the list.
Definition: UPtrListI.H:106
const fvMesh & mesh_
Const reference to the mesh.
label whichPatch(const label meshFacei) const
Return patch index for a given mesh face index. Uses binary search.
const fvMesh & mesh() const noexcept
Return const reference to the mesh.
const faceZoneMesh & faceZones() const noexcept
Return face zone mesh.
Definition: polyMesh.H:670
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:627
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Heat transfer coefficient calculation that employs the area-average temperature of a specified face z...
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
Reduce inplace (cf. MPI Allreduce) using specified communication schedule.
addToRunTimeSelectionTable(heatTransferCoeffModel, faceZoneReferenceTemperature, dictionary)
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
virtual bool read(const dictionary &dict)
Read from dictionary.
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
bool returnReduceOr(const bool value, const label comm=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Namespace for OpenFOAM.
IOerror FatalIOError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL IO ERROR&#39; header text and ...