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 | Copyright (C) 2011-2015 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 {
35  defineTypeNameAndDebug(cellLooper, 0);
36  defineRunTimeSelectionTable(cellLooper, word);
37 }
38 
39 
40 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41 
43 (
44  const word& type,
45  const polyMesh& mesh
46 )
47 {
48  wordConstructorTable::iterator cstrIter =
49  wordConstructorTablePtr_->find(type);
50 
51  if (cstrIter == wordConstructorTablePtr_->end())
52  {
54  << "Unknown set type "
55  << type << nl << nl
56  << "Valid cellLooper types : " << endl
57  << wordConstructorTablePtr_->sortedToc()
58  << exit(FatalError);
59  }
60 
61  return autoPtr<cellLooper>(cstrIter()(mesh));
62 }
63 
64 
65 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
66 
67 // Get faces (on cell) connected to vertI which are not using edgeI
69 (
70  const label cellI,
71  const label edgeI,
72  const label vertI
73 ) const
74 {
75  // Get faces connected to startEdge
76  label face0, face1;
77  meshTools::getEdgeFaces(mesh(), cellI, edgeI, face0, face1);
78 
79  const labelList& pFaces = mesh().pointFaces()[vertI];
80 
81  labelList vertFaces(pFaces.size());
82  label vertFaceI = 0;
83 
84  forAll(pFaces, pFaceI)
85  {
86  label faceI = pFaces[pFaceI];
87 
88  if
89  (
90  (faceI != face0)
91  && (faceI != face1)
92  && (meshTools::faceOnCell(mesh(), cellI, faceI))
93  )
94  {
95  vertFaces[vertFaceI++] = faceI;
96  }
97  }
98  vertFaces.setSize(vertFaceI);
99 
100  return vertFaces;
101 }
102 
103 
104 // Get first edge connected to vertI and on faceI
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 
134 // Get edges (on cell) connected to vertI which are not on faceI
136 (
137  const label cellI,
138  const label faceI,
139  const label vertI
140 ) const
141 {
142  const labelList& exclEdges = mesh().faceEdges()[faceI];
143 
144  const labelList& pEdges = mesh().pointEdges()[vertI];
145 
146  labelList vertEdges(pEdges.size());
147  label vertEdgeI = 0;
148 
149  forAll(pEdges, pEdgeI)
150  {
151  label edgeI = pEdges[pEdgeI];
152 
153  if
154  (
155  (findIndex(exclEdges, edgeI) == -1)
156  && meshTools::edgeOnCell(mesh(), cellI, edgeI)
157  )
158  {
159  vertEdges[vertEdgeI++] = edgeI;
160  }
161  }
162 
163  vertEdges.setSize(vertEdgeI);
164 
165  return vertEdges;
166 }
167 
168 
169 // Return edge from cellEdges that is most perpendicular
170 // to refinement direction.
172 (
173  const vector& refDir,
174  const label cellI
175 ) const
176 {
177  const labelList& cEdges = mesh().cellEdges()[cellI];
178 
179  label cutEdgeI = -1;
180  scalar maxCos = -GREAT;
181 
182  forAll(cEdges, cEdgeI)
183  {
184  label edgeI = cEdges[cEdgeI];
185 
186  scalar cosAngle = mag(refDir & meshTools::normEdgeVec(mesh(), edgeI));
187 
188  if (cosAngle > maxCos)
189  {
190  maxCos = cosAngle;
191 
192  cutEdgeI = edgeI;
193  }
194  }
195 
196  return cutEdgeI;
197 }
198 
199 
200 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
201 
203 :
205 {}
206 
207 
208 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
209 
211 {}
212 
213 
214 // ************************************************************************* //
meshTools.H
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
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:136
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::findIndex
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurence of given element and return index,.
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::cellLooper::~cellLooper
virtual ~cellLooper()
Destructor.
Definition: cellLooper.C:210
Foam::defineRunTimeSelectionTable
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
cellLooper.H
Foam::meshTools::edgeOnCell
bool edgeOnCell(const primitiveMesh &, const label cellI, const label edgeI)
Is edge used by cell.
Definition: meshTools.C:285
Foam::cellLooper::getMisAlignedEdge
label getMisAlignedEdge(const vector &refDir, const label cellI) const
Return edge from cellEdges that is most perpendicular.
Definition: cellLooper.C:172
Foam::cellLooper::cellLooper
cellLooper(const cellLooper &)
Disallow default bitwise copy construct.
Foam::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
Foam::edgeVertex
Combines edge or vertex in single label. Used to specify cuts across cell circumference.
Definition: edgeVertex.H:52
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
polyMesh.H
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
Foam::meshTools::faceOnCell
bool faceOnCell(const primitiveMesh &, const label cellI, const label faceI)
Is face used by cell.
Definition: meshTools.C:307
pFaces
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, &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
Foam::label
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
Foam::nl
static const char nl
Definition: Ostream.H:260
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:69
Foam::FatalError
error FatalError
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::meshTools::getEdgeFaces
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:456
Foam::List::setSize
void setSize(const label)
Reset size of List.
Foam::autoPtr
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:117
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::meshTools::normEdgeVec
vector normEdgeVec(const primitiveMesh &, const label edgeI)
Normalized edge vector.
Definition: meshTools.C:189
Foam::Vector< scalar >
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: doubleFloat.H:94
Foam::cellLooper::New
static autoPtr< cellLooper > New(const word &type, const polyMesh &mesh)
Return a reference to the selected cellLooper.
Definition: cellLooper.C:43
ListOps.H
Various functions to operate on Lists.
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::type
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:588
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)