heatExchangerModel.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 
28 #include "heatExchangerModel.H"
29 #include "coupledPolyPatch.H"
30 #include "surfaceInterpolate.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36  defineTypeNameAndDebug(heatExchangerModel, 0);
38 }
39 
40 
41 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
42 
44 (
45  const fvMesh& mesh,
46  const word& name,
47  const dictionary& coeffs
48 )
49 :
50  writeFile(mesh, name, typeName, coeffs),
51  mesh_(mesh),
52  coeffs_(coeffs),
53  name_(name),
54  UName_("U"),
55  TName_("T"),
56  phiName_("phi"),
57  faceZoneName_("unknown-faceZone"),
58  faceId_(),
59  facePatchId_(),
60  faceSign_()
61 {}
62 
63 
64 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
65 
67 {
68  const label zoneID = mesh_.faceZones().findZoneID(faceZoneName_);
69 
70  if (zoneID < 0)
71  {
73  << type() << " " << name_ << ": "
74  << " Unknown face zone name: " << faceZoneName_
75  << ". Valid face zones are: " << mesh_.faceZones().names()
76  << exit(FatalError);
77  }
78 
79  const faceZone& fZone = mesh_.faceZones()[zoneID];
80 
81  // Total number of faces selected
82  label numFaces = fZone.size();
83 
84  faceId_.resize_nocopy(numFaces);
85  facePatchId_.resize_nocopy(numFaces);
86  faceSign_.resize_nocopy(numFaces);
87 
88  numFaces = 0;
89 
90  // TDB: handle multiple zones
91  {
92  forAll(fZone, i)
93  {
94  const label meshFacei = fZone[i];
95  const label flipSign = (fZone.flipMap()[i] ? -1 : 1);
96 
97  // Internal faces
98  label faceId = meshFacei;
99  label facePatchId = -1;
100 
101  // Boundary faces
102  if (!mesh_.isInternalFace(meshFacei))
103  {
104  facePatchId = mesh_.boundaryMesh().whichPatch(meshFacei);
105  const polyPatch& pp = mesh_.boundaryMesh()[facePatchId];
106 
107  if (isA<emptyPolyPatch>(pp))
108  {
109  continue; // Ignore empty patch
110  }
111 
112  const auto* cpp = isA<coupledPolyPatch>(pp);
113 
114  if (cpp && !cpp->owner())
115  {
116  continue; // Ignore neighbour side
117  }
118 
119  faceId = pp.whichFace(meshFacei);
120  }
121 
122  if (faceId >= 0)
123  {
124  faceId_[numFaces] = faceId;
125  facePatchId_[numFaces] = facePatchId;
126  faceSign_[numFaces] = flipSign;
127 
128  ++numFaces;
129  }
130  }
131  }
132 
133  // Shrink to size used
134  faceId_.resize(numFaces);
135  facePatchId_.resize(numFaces);
136  faceSign_.resize(numFaces);
137 }
138 
139 
140 // ************************************************************************* //
labelList faceId_
Local list of face IDs.
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
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
word faceZoneName_
Name of the faceZone at the heat exchanger inlet.
virtual void initialise()
Initialise data members of the model.
void resize_nocopy(const label len)
Adjust allocated size of list without necessarily.
Definition: ListI.H:175
labelList faceSign_
List of +1/-1 representing face flip map (1 use as is, -1 negate)
bool isInternalFace(const label faceIndex) const noexcept
Return true if given face label is internal to the mesh.
Base class for heat exchanger models to handle various characteristics for the heatExchangerSource fv...
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
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 word & name_
Reference to the name of the fvOption source.
heatExchangerModel(const fvMesh &mesh, const word &name, const dictionary &coeffs)
Construct from components.
label whichPatch(const label meshFacei) const
Return patch index for a given mesh face index. Uses binary search.
const fvMesh & mesh_
Reference to the mesh.
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
const faceZoneMesh & faceZones() const noexcept
Return face zone mesh.
Definition: polyMesh.H:670
defineTypeNameAndDebug(combustionModel, 0)
labelList facePatchId_
Local list of patch IDs per face.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
wordList names() const
A list of the zone names.
Definition: ZoneMesh.C:362
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Namespace for OpenFOAM.