cellLooper.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) 2019-2021 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 "cellLooper.H"
30 #include "polyMesh.H"
31 #include "ListOps.H"
32 #include "meshTools.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38  defineTypeNameAndDebug(cellLooper, 0);
39  defineRunTimeSelectionTable(cellLooper, word);
40 }
41 
42 
43 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44 
46 (
47  const word& type,
48  const polyMesh& mesh
49 )
50 {
51  auto* ctorPtr = wordConstructorTable(type);
52 
53  if (!ctorPtr)
54  {
56  (
57  "cellLooper",
58  type,
59  *wordConstructorTablePtr_
60  ) << exit(FatalError);
61  }
62 
63  return autoPtr<cellLooper>(ctorPtr(mesh));
64 }
65 
66 
67 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
68 
70 (
71  const label celli,
72  const label edgeI,
73  const label vertI
74 ) const
75 {
76  // Get faces connected to startEdge
77  label face0, face1;
78  meshTools::getEdgeFaces(mesh(), celli, edgeI, face0, face1);
79 
80  const labelList& pFaces = mesh().pointFaces()[vertI];
81 
82  labelList vertFaces(pFaces.size());
83  label vertFacei = 0;
84 
85  forAll(pFaces, pFacei)
86  {
87  label facei = pFaces[pFacei];
88 
89  if
90  (
91  (facei != face0)
92  && (facei != face1)
93  && (meshTools::faceOnCell(mesh(), celli, facei))
94  )
95  {
96  vertFaces[vertFacei++] = facei;
97  }
98  }
99  vertFaces.setSize(vertFacei);
100 
101  return vertFaces;
102 }
103 
104 
106 (
107  const label facei,
108  const label vertI
109 ) const
110 {
111  const labelList& fEdges = mesh().faceEdges()[facei];
112 
113  forAll(fEdges, fEdgeI)
114  {
115  label edgeI = fEdges[fEdgeI];
116 
117  const edge& e = mesh().edges()[edgeI];
118 
119  if ((e.start() == vertI) || (e.end() == vertI))
120  {
121  return edgeI;
122  }
123  }
124 
126  << "Can not find edge on face " << facei
127  << " using vertex " << vertI
128  << abort(FatalError);
129 
130  return -1;
131 }
132 
133 
135 (
136  const label celli,
137  const label facei,
138  const label vertI
139 ) const
140 {
141  const labelList& exclEdges = mesh().faceEdges()[facei];
142 
143  const labelList& pEdges = mesh().pointEdges()[vertI];
144 
145  labelList vertEdges(pEdges.size());
146  label vertEdgeI = 0;
147 
148  forAll(pEdges, pEdgeI)
149  {
150  label edgeI = pEdges[pEdgeI];
151 
152  if
153  (
154  !exclEdges.found(edgeI)
155  && meshTools::edgeOnCell(mesh(), celli, edgeI)
156  )
157  {
158  vertEdges[vertEdgeI++] = edgeI;
159  }
160  }
161 
162  vertEdges.setSize(vertEdgeI);
163 
164  return vertEdges;
165 }
166 
167 
169 (
170  const vector& refDir,
171  const label celli
172 ) const
173 {
174  const labelList& cEdges = mesh().cellEdges()[celli];
175 
176  label cutEdgeI = -1;
177  scalar maxCos = -GREAT;
178 
179  forAll(cEdges, cEdgeI)
180  {
181  label edgeI = cEdges[cEdgeI];
182 
183  scalar cosAngle = mag(refDir & meshTools::normEdgeVec(mesh(), edgeI));
184 
185  if (cosAngle > maxCos)
186  {
187  maxCos = cosAngle;
188 
189  cutEdgeI = edgeI;
190  }
191  }
192 
193  return cutEdgeI;
194 }
195 
196 
197 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
198 
200 :
202 {}
203 
204 
205 // ************************************************************************* //
meshTools.H
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::cellLooper::getVertEdgesNonFace
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:135
Foam::defineRunTimeSelectionTable
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
cellLooper.H
Foam::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:63
Foam::edgeVertex
Combines edge or vertex in single label. Used to specify cuts across cell circumference.
Definition: edgeVertex.H:55
Foam::cellLooper::getFirstVertEdge
label getFirstVertEdge(const label facei, const label vertI) const
Get first edge connected to vertI and on facei.
Definition: cellLooper.C:106
Foam::meshTools::edgeOnCell
bool edgeOnCell(const primitiveMesh &mesh, const label celli, const label edgeI)
Is edge used by cell.
Definition: meshTools.C:308
Foam::cellLooper::cellLooper
cellLooper(const cellLooper &)=delete
No copy construct.
polyMesh.H
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::meshTools::faceOnCell
bool faceOnCell(const primitiveMesh &mesh, const label celli, const label facei)
Is face used by cell.
Definition: meshTools.C:330
Foam::List::setSize
void setSize(const label n)
Alias for resize()
Definition: List.H:222
Foam::meshTools::getEdgeFaces
void getEdgeFaces(const primitiveMesh &mesh, 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:479
Foam::FatalError
error FatalError
FatalErrorInLookup
#define FatalErrorInLookup(lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalError.
Definition: error.H:457
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::meshTools::normEdgeVec
vector normEdgeVec(const primitiveMesh &, const label edgeI)
Normalized edge vector.
Definition: meshTools.C:193
Foam::Vector< scalar >
Foam::List< label >
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::cellLooper::getVertFacesNonEdge
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:70
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::cellLooper::New
static autoPtr< cellLooper > New(const word &type, const polyMesh &mesh)
Return a reference to the selected cellLooper.
Definition: cellLooper.C:46
Foam::roots::type
type
Types of root.
Definition: Roots.H:54
ListOps.H
Various functions to operate on Lists.
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::cellLooper::getMisAlignedEdge
label getMisAlignedEdge(const vector &refDir, const label celli) const
Return edge from cellEdges that is most perpendicular.
Definition: cellLooper.C:169
pFaces
Info<< "Finished reading KIVA file"<< endl;cellShapeList cellShapes(nPoints);labelList cellZoning(nPoints, -1);const cellModel &hex=cellModel::ref(cellModel::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].reset(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, &oldCyclicPolyPatch::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:235