tetMatcher.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 -------------------------------------------------------------------------------
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 "tetMatcher.H"
29 #include "cellMatcher.H"
30 #include "primitiveMesh.H"
31 #include "cellModel.H"
32 #include "ListOps.H"
33 
34 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35 
37 :
39  (
40  vertPerCell,
41  facePerCell,
42  maxVertPerFace,
43  "tet" // same as cellModel::modelNames[cellModel::TET]
44  )
45 {}
46 
47 
48 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
49 
51 {}
52 
53 
54 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
55 
56 
58 (
59  const bool checkOnly,
60  const faceList& faces,
61  const labelList& owner,
62  const label celli,
63  const labelList& myFaces
64 )
65 {
66  if (!faceSizeMatch(faces, myFaces))
67  {
68  return false;
69  }
70 
71  // Tet for sure now
72  if (checkOnly)
73  {
74  return true;
75  }
76 
77  // Calculate localFaces_ and mapping pointMap_, faceMap_
78  label numVert = calcLocalFaces(faces, myFaces);
79 
80  if (numVert != vertPerCell)
81  {
82  return false;
83  }
84 
85  // Set up 'edge' to face mapping.
86  calcEdgeAddressing(numVert);
87 
88  // Set up point on face to index-in-face mapping
89  calcPointFaceIndex();
90 
91  // Storage for maps -vertex to mesh and -face to mesh
92  vertLabels_.setSize(vertPerCell);
93  faceLabels_.setSize(facePerCell);
94 
95  //
96  // Try bottom face (face 3)
97  //
98 
99  label face3I = 0;
100  const face& face3 = localFaces_[face3I];
101  label face3vert0 = 0;
102 
103  //
104  // Try to follow prespecified path on faces of cell,
105  // starting at face3vert0
106  //
107 
108  vertLabels_[0] = pointMap_[face3[face3vert0]];
109  faceLabels_[3] = faceMap_[face3I];
110 
111  // Walk face 3 from vertex 0 to 1
112  label face3vert1 =
113  nextVert
114  (
115  face3vert0,
116  faceSize_[face3I],
117  !(owner[faceMap_[face3I]] == celli)
118  );
119  vertLabels_[1] = pointMap_[face3[face3vert1]];
120 
121  // Walk face 3 from vertex 1 to 2
122  label face3vert2 =
123  nextVert
124  (
125  face3vert1,
126  faceSize_[face3I],
127  !(owner[faceMap_[face3I]] == celli)
128  );
129  vertLabels_[2] = pointMap_[face3[face3vert2]];
130 
131  // Jump edge from face3 to face2
132  label face2I =
133  otherFace
134  (
135  numVert,
136  face3[face3vert0],
137  face3[face3vert1],
138  face3I
139  );
140  faceLabels_[2] = faceMap_[face2I];
141 
142  // Jump edge from face3 to face0
143  label face0I =
144  otherFace
145  (
146  numVert,
147  face3[face3vert1],
148  face3[face3vert2],
149  face3I
150  );
151  faceLabels_[0] = faceMap_[face0I];
152 
153  // Jump edge from face3 to face1
154  label face1I =
155  otherFace
156  (
157  numVert,
158  face3[face3vert2],
159  face3[face3vert0],
160  face3I
161  );
162  faceLabels_[1] = faceMap_[face1I];
163  const face& face1 = localFaces_[face1I];
164 
165  // Get index of vert0 in face 1
166  label face1vert0 = pointFaceIndex_[face3[face3vert0]][face1I];
167 
168  // Walk face 1 from vertex 0 to 3
169  label face1vert3 =
170  nextVert
171  (
172  face1vert0,
173  faceSize_[face1I],
174  (owner[faceMap_[face1I]] == celli)
175  );
176  vertLabels_[3] = pointMap_[face1[face1vert3]];
177 
178  return true;
179 }
180 
181 
183 {
184  return 4*3;
185 }
186 
187 
189 (
190  const faceList& faces,
191  const labelList& myFaces
192 ) const
193 {
194  if (myFaces.size() != 4)
195  {
196  return false;
197  }
198 
199  for (const label facei : myFaces)
200  {
201  const label size = faces[facei].size();
202 
203  if (size != 3)
204  {
205  return false;
206  }
207  }
208 
209  return true;
210 }
211 
212 
214 {
215  return matchShape
216  (
217  true,
218  mesh.faces(),
219  mesh.faceOwner(),
220  celli,
221  mesh.cells()[celli]
222  );
223 }
224 
225 
226 bool Foam::tetMatcher::isA(const faceList& faces)
227 {
228  // Do as if mesh with one cell only
229  return matchShape
230  (
231  true,
232  faces, // all faces in mesh
233  labelList(faces.size(), Zero), // cell 0 is owner of all faces
234  0, // cell label
235  identity(faces.size()) // faces of cell 0
236  );
237 }
238 
239 
241 (
242  const primitiveMesh& mesh,
243  const label celli,
244  cellShape& shape
245 )
246 {
247  if
248  (
249  matchShape
250  (
251  false,
252  mesh.faces(),
253  mesh.faceOwner(),
254  celli,
255  mesh.cells()[celli]
256  )
257  )
258  {
259  shape = cellShape(model(), vertLabels());
260 
261  return true;
262  }
263 
264  return false;
265 }
266 
267 
268 // ************************************************************************* //
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:74
Foam::tetMatcher::faceHashValue
virtual label faceHashValue() const
Hash value of all face sizes of this shape. Can be used for.
Definition: tetMatcher.C:182
Foam::meshTools::otherFace
label otherFace(const primitiveMesh &mesh, const label celli, const label facei, const label edgeI)
Return face on cell using edgeI but not facei. Throws error.
Definition: meshTools.C:555
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::primitiveMesh::cells
const cellList & cells() const
Definition: primitiveMeshCells.C:138
primitiveMesh.H
Foam::cellMatcher
Base class for cellshape matchers (hexMatch, prismMatch, etc.). These are classes which given a mesh ...
Definition: cellMatcher.H:101
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:62
cellModel.H
Foam::polyMesh::faceOwner
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1076
Foam::tetMatcher::matches
virtual bool matches(const primitiveMesh &mesh, const label celli, cellShape &shape)
Like isA but also constructs a cellShape (if shape matches)
Definition: tetMatcher.C:241
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::cellShape
An analytical geometric cellShape.
Definition: cellShape.H:71
Foam::tetMatcher::matchShape
virtual bool matchShape(const bool checkOnly, const faceList &faces, const labelList &faceOwner, const label celli, const labelList &myFaces)
Low level shape recognition. Return true if matches.
Definition: tetMatcher.C:58
Foam::tetMatcher::faceSizeMatch
virtual bool faceSizeMatch(const faceList &, const labelList &) const
Check whether number of face sizes match the shape.
Definition: tetMatcher.C:189
Foam::polyMesh::faces
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1063
Foam::tetMatcher::~tetMatcher
~tetMatcher()
Destructor.
Definition: tetMatcher.C:50
Foam::tetMatcher::isA
virtual bool isA(const primitiveMesh &mesh, const label celli)
Exact match. Uses faceSizeMatch.
Definition: tetMatcher.C:213
Foam::List< face >
cellMatcher.H
Foam::identity
labelList identity(const label len, label start=0)
Create identity map of the given length with (map[i] == i)
Definition: labelList.C:38
tetMatcher.H
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:74
ListOps.H
Various functions to operate on Lists.
Foam::tetMatcher::tetMatcher
tetMatcher()
Construct null.
Definition: tetMatcher.C:36
Foam::primitiveMesh
Cell-face mesh analysis engine.
Definition: primitiveMesh.H:78