cellLooper.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2024 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 
26 #include "cellLooper.H"
27 #include "polyMesh.H"
28 #include "ListOps.H"
29 #include "meshTools.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
36 }
37 
38 
39 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
40 
42 (
43  const label celli,
44  const label edgeI,
45  const label vertI
46 ) const
47 {
48  // Get faces connected to startEdge
49  label face0, face1;
50  meshTools::getEdgeFaces(mesh(), celli, edgeI, face0, face1);
51 
52  const labelList& pFaces = mesh().pointFaces()[vertI];
53 
54  labelList vertFaces(pFaces.size());
55  label vertFacei = 0;
56 
57  forAll(pFaces, pFacei)
58  {
59  label facei = pFaces[pFacei];
60 
61  if
62  (
63  (facei != face0)
64  && (facei != face1)
65  && (meshTools::faceOnCell(mesh(), celli, facei))
66  )
67  {
68  vertFaces[vertFacei++] = facei;
69  }
70  }
71  vertFaces.setSize(vertFacei);
72 
73  return vertFaces;
74 }
75 
76 
78 (
79  const label facei,
80  const label vertI
81 ) const
82 {
83  const labelList& fEdges = mesh().faceEdges()[facei];
84 
85  forAll(fEdges, fEdgeI)
86  {
87  label edgeI = fEdges[fEdgeI];
88 
89  const edge& e = mesh().edges()[edgeI];
90 
91  if ((e.start() == vertI) || (e.end() == vertI))
92  {
93  return edgeI;
94  }
95  }
96 
98  << "Can not find edge on face " << facei
99  << " using vertex " << vertI
100  << abort(FatalError);
101 
102  return -1;
103 }
104 
105 
107 (
108  const label celli,
109  const label facei,
110  const label vertI
111 ) const
112 {
113  const labelList& exclEdges = mesh().faceEdges()[facei];
114 
115  const labelList& pEdges = mesh().pointEdges()[vertI];
116 
117  labelList vertEdges(pEdges.size());
118  label vertEdgeI = 0;
119 
120  forAll(pEdges, pEdgeI)
121  {
122  label edgeI = pEdges[pEdgeI];
123 
124  if
125  (
126  (findIndex(exclEdges, edgeI) == -1)
127  && meshTools::edgeOnCell(mesh(), celli, edgeI)
128  )
129  {
130  vertEdges[vertEdgeI++] = edgeI;
131  }
132  }
133 
134  vertEdges.setSize(vertEdgeI);
135 
136  return vertEdges;
137 }
138 
139 
141 (
142  const vector& refDir,
143  const label celli
144 ) const
145 {
146  const labelList& cEdges = mesh().cellEdges()[celli];
147 
148  label cutEdgeI = -1;
149  scalar maxCos = -great;
150 
151  forAll(cEdges, cEdgeI)
152  {
153  label edgeI = cEdges[cEdgeI];
154 
155  scalar cosAngle = mag(refDir & meshTools::normEdgeVec(mesh(), edgeI));
156 
157  if (cosAngle > maxCos)
158  {
159  maxCos = cosAngle;
160 
161  cutEdgeI = edgeI;
162  }
163  }
164 
165  return cutEdgeI;
166 }
167 
168 
169 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
170 
172 :
173  edgeVertex(mesh)
174 {}
175 
176 
177 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
178 
180 {}
181 
182 
183 // ************************************************************************* //
Various functions to operate on Lists.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
void setSize(const label)
Reset size of List.
Definition: List.C:281
Abstract base class. Concrete implementations know how to cut a cell (i.e. determine a loop around th...
Definition: cellLooper.H:72
cellLooper(const polyMesh &mesh)
Construct from components.
Definition: cellLooper.C:171
labelList getVertEdgesNonFace(const label celli, const label facei, const label vertI) const
Get edges (on cell) connected to vertI which are not on facei.
Definition: cellLooper.C:107
label getMisAlignedEdge(const vector &refDir, const label celli) const
Return edge from cellEdges that is most perpendicular.
Definition: cellLooper.C:141
label getFirstVertEdge(const label facei, const label vertI) const
Get first edge connected to vertI and on facei.
Definition: cellLooper.C:78
labelList getVertFacesNonEdge(const label celli, const label edgeI, const label vertI) const
Get faces (on cell) connected to vertI which are not using edgeI.
Definition: cellLooper.C:42
virtual ~cellLooper()
Destructor.
Definition: cellLooper.C:179
Combines edge or vertex in single label. Used to specify cuts across cell circumference.
Definition: edgeVertex.H:53
const polyMesh & mesh() const
Definition: edgeVertex.H:93
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:61
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
const labelListList & pointFaces() const
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
void getEdgeFaces(const primitiveMesh &, const label celli, const label edgeI, label &face0, label &face1)
Get faces on cell using edgeI. Throws error if no two found.
Definition: meshTools.C:458
bool faceOnCell(const primitiveMesh &, const label celli, const label facei)
Is face used by cell.
Definition: meshTools.C:308
vector normEdgeVec(const primitiveMesh &, const label edgeI)
Normalised edge vector.
Definition: meshTools.C:189
bool edgeOnCell(const primitiveMesh &, const label celli, const label edgeI)
Is edge used by cell.
Definition: meshTools.C:285
Namespace for OpenFOAM.
const doubleScalar e
Definition: doubleScalar.H:106
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
errorManip< error > abort(error &err)
Definition: errorManip.H:131
dimensioned< scalar > mag(const dimensioned< Type > &)
defineTypeNameAndDebug(combustionModel, 0)
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.
error FatalError
Info<< "Finished reading KIVA file"<< endl;cellShapeList cellShapes(nPoints);labelList cellZoning(nPoints, -1);const cellModel &hex=*(cellModeller::lookup("hex"));labelList hexLabels(8);label activeCells=0;labelList pointMap(nPoints);forAll(pointMap, i){ pointMap[i]=i;}for(label i=0;i< nPoints;i++){ if(f[i] > 0.0) { hexLabels[0]=i;hexLabels[1]=i1tab[i];hexLabels[2]=i3tab[i1tab[i]];hexLabels[3]=i3tab[i];hexLabels[4]=i8tab[i];hexLabels[5]=i1tab[i8tab[i]];hexLabels[6]=i3tab[i1tab[i8tab[i]]];hexLabels[7]=i3tab[i8tab[i]];cellShapes[activeCells]=cellShape(hex, hexLabels);edgeList edges=cellShapes[activeCells].edges();forAll(edges, ei) { if(edges[ei].mag(points)< small) { label start=pointMap[edges[ei].start()];while(start !=pointMap[start]) { start=pointMap[start];} label end=pointMap[edges[ei].end()];while(end !=pointMap[end]) { end=pointMap[end];} label minLabel=min(start, end);pointMap[start]=pointMap[end]=minLabel;} } cellZoning[activeCells]=idreg[i];activeCells++;}}cellShapes.setSize(activeCells);cellZoning.setSize(activeCells);forAll(cellShapes, celli){ cellShape &cs=cellShapes[celli];forAll(cs, i) { cs[i]=pointMap[cs[i]];} cs.collapse();}label bcIDs[11]={-1, 0, 2, 4, -1, 5, -1, 6, 7, 8, 9};const label nBCs=12;const word *kivaPatchTypes[nBCs]={ &wallPolyPatch::typeName, &wallPolyPatch::typeName, &wallPolyPatch::typeName, &wallPolyPatch::typeName, &symmetryPolyPatch::typeName, &wedgePolyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &symmetryPolyPatch::typeName, &mergedCyclicPolyPatch::typeName};enum patchTypeNames{ PISTON, VALVE, LINER, CYLINDERHEAD, AXIS, WEDGE, INFLOW, OUTFLOW, PRESIN, PRESOUT, SYMMETRYPLANE, CYCLIC};const char *kivaPatchNames[nBCs]={ "piston", "valve", "liner", "cylinderHead", "axis", "wedge", "inflow", "outflow", "presin", "presout", "symmetryPlane", "cyclic"};List< SLList< face > > pFaces[nBCs]
Definition: readKivaGrid.H:229