Loading...
Searching...
No Matches
channelIndex.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) 2022 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
29#include "channelIndex.H"
30#include "boolList.H"
31#include "syncTools.H"
32#include "OFstream.H"
33#include "meshTools.H"
34#include "Time.H"
35#include "SortableList.H"
36
37// * * * * * * * * * * * * * Static Member Data * * * * * * * * * * * * * * //
38
39const Foam::Enum
40<
42>
43Foam::channelIndex::vectorComponentsNames_
44({
45 { vector::components::X, "x" },
46 { vector::components::Y, "y" },
47 { vector::components::Z, "z" },
48});
49
50
51// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
52
53// Determines face blocking
54void Foam::channelIndex::walkOppositeFaces
55(
56 const polyMesh& mesh,
57 const labelList& startFaces,
58 boolList& blockedFace
59)
60{
61 const cellList& cells = mesh.cells();
62 const faceList& faces = mesh.faces();
63 const label nBnd = mesh.nBoundaryFaces();
64
65 DynamicList<label> frontFaces(startFaces);
66 forAll(frontFaces, i)
67 {
68 label facei = frontFaces[i];
69 blockedFace[facei] = true;
70 }
71
72 while (returnReduceOr(frontFaces.size()))
73 {
74 // Transfer across.
75 boolList isFrontBndFace(nBnd, false);
76 forAll(frontFaces, i)
77 {
78 label facei = frontFaces[i];
79
80 if (!mesh.isInternalFace(facei))
81 {
82 isFrontBndFace[facei-mesh.nInternalFaces()] = true;
83 }
84 }
86
87 // Add
88 forAll(isFrontBndFace, i)
89 {
90 label facei = mesh.nInternalFaces()+i;
91 if (isFrontBndFace[i] && !blockedFace[facei])
92 {
93 blockedFace[facei] = true;
94 frontFaces.append(facei);
95 }
96 }
97
98 // Transfer across cells
99 DynamicList<label> newFrontFaces(frontFaces.size());
100
101 forAll(frontFaces, i)
102 {
103 label facei = frontFaces[i];
104
105 {
106 const cell& ownCell = cells[mesh.faceOwner()[facei]];
107
108 label oppositeFacei = ownCell.opposingFaceLabel(facei, faces);
109
110 if (oppositeFacei == -1)
111 {
113 << "Face:" << facei << " owner cell:" << ownCell
114 << " is not a hex?" << abort(FatalError);
115 }
116 else
117 {
118 if (!blockedFace[oppositeFacei])
119 {
120 blockedFace[oppositeFacei] = true;
121 newFrontFaces.append(oppositeFacei);
122 }
123 }
124 }
125
126 if (mesh.isInternalFace(facei))
127 {
128 const cell& neiCell = mesh.cells()[mesh.faceNeighbour()[facei]];
129
130 label oppositeFacei = neiCell.opposingFaceLabel(facei, faces);
131
132 if (oppositeFacei == -1)
133 {
135 << "Face:" << facei << " neighbour cell:" << neiCell
136 << " is not a hex?" << abort(FatalError);
137 }
138 else
139 {
140 if (!blockedFace[oppositeFacei])
141 {
142 blockedFace[oppositeFacei] = true;
143 newFrontFaces.append(oppositeFacei);
144 }
145 }
146 }
147 }
148
149 frontFaces.transfer(newFrontFaces);
150 }
151}
152
153
154// Calculate regions.
155void Foam::channelIndex::calcLayeredRegions
156(
157 const polyMesh& mesh,
158 const labelList& startFaces
159)
160{
161 boolList blockedFace(mesh.nFaces(), false);
162 walkOppositeFaces
163 (
164 mesh,
165 startFaces,
166 blockedFace
167 );
168
169
170 if (false)
171 {
172 OFstream str(mesh.time().path()/"blockedFaces.obj");
173 label vertI = 0;
174 forAll(blockedFace, facei)
175 {
176 if (blockedFace[facei])
177 {
178 const face& f = mesh.faces()[facei];
179 forAll(f, fp)
180 {
181 meshTools::writeOBJ(str, mesh.points()[f[fp]]);
182 }
183 str<< 'f';
184 forAll(f, fp)
185 {
186 str << ' ' << vertI+fp+1;
187 }
188 str << nl;
189 vertI += f.size();
190 }
191 }
192 }
193
194
195 // Do analysis for connected regions
196 cellRegion_.reset(new regionSplit(mesh, blockedFace));
197
198 Info<< "Detected " << cellRegion_().nRegions() << " layers." << nl << endl;
199
200 // Sum number of entries per region
201 regionCount_ = regionSum(scalarField(mesh.nCells(), 1.0));
202
203 // Average cell centres to determine ordering.
204 pointField regionCc
205 (
207 / regionCount_
208 );
209
210 SortableList<scalar> sortComponent(regionCc.component(dir_));
211
212 sortMap_ = sortComponent.indices();
213
214 y_ = sortComponent;
215
216 if (symmetric_)
217 {
218 y_.setSize(cellRegion_().nRegions()/2);
219 }
220}
221
222
223// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
224
225Foam::channelIndex::channelIndex
226(
227 const polyMesh& mesh,
228 const dictionary& dict
229)
230:
231 symmetric_(dict.get<bool>("symmetric")),
232 dir_(vectorComponentsNames_.get("component", dict))
233{
234 const polyBoundaryMesh& patches = mesh.boundaryMesh();
235
236 const wordList patchNames(dict.get<wordList>("patches"));
237
238 label nFaces = 0;
239
241 {
242 const label patchi = patches.findPatchID(patchNames[i]);
243
244 if (patchi == -1)
245 {
247 << "Illegal patch " << patchNames[i]
248 << ". Valid patches are " << patches.name()
249 << exit(FatalError);
250 }
251
252 nFaces += patches[patchi].size();
253 }
254
255 labelList startFaces(nFaces);
256 nFaces = 0;
257
259 {
260 const polyPatch& pp = patches[patchNames[i]];
261
262 forAll(pp, j)
263 {
264 startFaces[nFaces++] = pp.start()+j;
265 }
266 }
267
268 // Calculate regions.
269 calcLayeredRegions(mesh, startFaces);
270}
271
272
273Foam::channelIndex::channelIndex
274(
275 const polyMesh& mesh,
276 const labelList& startFaces,
277 const bool symmetric,
278 const direction dir
279)
280:
281 symmetric_(symmetric),
282 dir_(dir)
283{
284 // Calculate regions.
285 calcLayeredRegions(mesh, startFaces);
286}
287
288
289// ************************************************************************* //
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition Enum.H:57
Output to file stream as an OSstream, normally using std::ofstream for the actual output.
Definition OFstream.H:75
A list that is sorted upon construction or when explicitly requested with the sort() method.
fileName path() const
The path for the case = rootPath/caseName.
Definition TimePathsI.H:102
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
A face is a list of labels corresponding to mesh vertices.
Definition face.H:71
const Time & time() const
Return the top-level database.
Definition fvMesh.H:360
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
virtual const faceList & faces() const
Return raw faces.
Definition polyMesh.C:1088
virtual const pointField & points() const
Return raw points.
Definition polyMesh.C:1063
const vectorField & cellCentres() const
label nCells() const noexcept
Number of mesh cells.
label nFaces() const noexcept
Number of mesh faces.
This class separates the mesh into distinct unconnected regions, each of which is then given a label ...
static void swapBoundaryFaceList(const polyMesh &mesh, UList< T > &faceValues, const bool parRun=UPstream::parRun())
Swap coupled boundary face values. Uses eqOp.
Definition syncTools.H:524
const polyBoundaryMesh & patches
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
const cellShapeList & cells
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of a point.
Definition meshTools.C:196
List< word > wordList
List of word.
Definition fileName.H:60
bool returnReduceOr(const bool value, const int communicator=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
List< label > labelList
A List of labels.
Definition List.H:62
messageStream Info
Information stream (stdout output on master, null elsewhere).
List< face > faceList
List of faces.
Definition faceListFwd.H:41
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
List< cell > cellList
List of cell.
Definition cellListFwd.H:41
errorManip< error > abort(error &err)
Definition errorManip.H:139
static Map< Type > regionSum(const regionSplit &regions, const Field< Type > &fld)
List< bool > boolList
A List of bools.
Definition List.H:60
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
vectorField pointField
pointField is a vectorField.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
wordList patchNames(nPatches)
labelList f(nPoints)
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299