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  Copyright (C) 2020 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 "tetMatcher.H"
30 #include "cellMatcher.H"
31 #include "primitiveMesh.H"
32 #include "cellModel.H"
33 #include "ListOps.H"
34 
35 // * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39 
40 // Check (4 tri)
41 static inline bool checkFaceSizeMatch(const UList<face>& faces)
42 {
43  if (faces.size() != 4) // facePerCell
44  {
45  return false;
46  }
47 
48  for (const face& f : faces)
49  {
50  if (f.size() != 3) // tri
51  {
52  return false;
53  }
54  }
55 
56  return true;
57 }
58 
59 
60 // Check (4 tri)
61 static inline bool checkFaceSizeMatch
62 (
63  const UList<face>& meshFaces,
64  const labelUList& cellFaces
65 )
66 {
67  if (cellFaces.size() != 4) // facePerCell
68  {
69  return false;
70  }
71 
72  for (const label facei : cellFaces)
73  {
74  if (meshFaces[facei].size() != 3) // tri
75  {
76  return false;
77  }
78  }
79 
80  return true;
81 }
82 
83 
84 } // End namespace Foam
85 
86 
87 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
88 
90 {
91  return checkFaceSizeMatch(faces);
92 }
93 
94 
95 bool Foam::tetMatcher::test(const primitiveMesh& mesh, const label celli)
96 {
97  return checkFaceSizeMatch(mesh.faces(), mesh.cells()[celli]);
98 }
99 
100 
101 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
102 
104 :
106  (
107  vertPerCell,
108  facePerCell,
109  maxVertPerFace,
110  "tet" // == cellModel::modelNames[cellModel::TET]
111  )
112 {}
113 
114 
115 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
116 
118 (
119  const bool checkOnly,
120  const faceList& faces,
121  const labelList& owner,
122  const label celli,
123  const labelList& myFaces
124 )
125 {
126  if (!faceSizeMatch(faces, myFaces))
127  {
128  return false;
129  }
130 
131  // Tet for sure now
132  if (checkOnly)
133  {
134  return true;
135  }
136 
137  // Calculate localFaces_ and mapping pointMap_, faceMap_
138  label numVert = calcLocalFaces(faces, myFaces);
139 
140  if (numVert != vertPerCell)
141  {
142  return false;
143  }
144 
145  // Set up 'edge' to face mapping.
146  calcEdgeAddressing(numVert);
147 
148  // Set up point on face to index-in-face mapping
149  calcPointFaceIndex();
150 
151  // Storage for maps -vertex to mesh and -face to mesh
152  vertLabels_.setSize(vertPerCell);
153  faceLabels_.setSize(facePerCell);
154 
155  //
156  // Try bottom face (face 3)
157  //
158 
159  label face3I = 0;
160  const face& face3 = localFaces_[face3I];
161  label face3vert0 = 0;
162 
163  //
164  // Try to follow prespecified path on faces of cell,
165  // starting at face3vert0
166  //
167 
168  vertLabels_[0] = pointMap_[face3[face3vert0]];
169  faceLabels_[3] = faceMap_[face3I];
170 
171  // Walk face 3 from vertex 0 to 1
172  label face3vert1 =
173  nextVert
174  (
175  face3vert0,
176  faceSize_[face3I],
177  !(owner[faceMap_[face3I]] == celli)
178  );
179  vertLabels_[1] = pointMap_[face3[face3vert1]];
180 
181  // Walk face 3 from vertex 1 to 2
182  label face3vert2 =
183  nextVert
184  (
185  face3vert1,
186  faceSize_[face3I],
187  !(owner[faceMap_[face3I]] == celli)
188  );
189  vertLabels_[2] = pointMap_[face3[face3vert2]];
190 
191  // Jump edge from face3 to face2
192  label face2I =
193  otherFace
194  (
195  numVert,
196  face3[face3vert0],
197  face3[face3vert1],
198  face3I
199  );
200  faceLabels_[2] = faceMap_[face2I];
201 
202  // Jump edge from face3 to face0
203  label face0I =
204  otherFace
205  (
206  numVert,
207  face3[face3vert1],
208  face3[face3vert2],
209  face3I
210  );
211  faceLabels_[0] = faceMap_[face0I];
212 
213  // Jump edge from face3 to face1
214  label face1I =
215  otherFace
216  (
217  numVert,
218  face3[face3vert2],
219  face3[face3vert0],
220  face3I
221  );
222  faceLabels_[1] = faceMap_[face1I];
223  const face& face1 = localFaces_[face1I];
224 
225  // Get index of vert0 in face 1
226  label face1vert0 = pointFaceIndex_[face3[face3vert0]][face1I];
227 
228  // Walk face 1 from vertex 0 to 3
229  label face1vert3 =
230  nextVert
231  (
232  face1vert0,
233  faceSize_[face1I],
234  (owner[faceMap_[face1I]] == celli)
235  );
236  vertLabels_[3] = pointMap_[face1[face1vert3]];
237 
238  return true;
239 }
240 
241 
243 {
244  return 4*3;
245 }
246 
247 
249 (
250  const faceList& meshFaces,
251  const labelList& cellFaces
252 ) const
253 {
254  return checkFaceSizeMatch(meshFaces, cellFaces);
255 }
256 
257 
259 (
260  const primitiveMesh& mesh,
261  const label celli,
262  cellShape& shape
263 )
264 {
265  if
266  (
267  matchShape
268  (
269  false,
270  mesh.faces(),
271  mesh.faceOwner(),
272  celli,
273  mesh.cells()[celli]
274  )
275  )
276  {
277  shape.reset(model(), vertLabels());
278  return true;
279  }
280 
281  return false;
282 }
283 
284 
285 // ************************************************************************* //
Foam::tetMatcher::faceHashValue
virtual label faceHashValue() const
Hash value of all face sizes of this shape. Can be used for.
Definition: tetMatcher.C:242
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::checkFaceSizeMatch
static bool checkFaceSizeMatch(const UList< face > &faces)
Definition: hexMatcher.C:39
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:99
cellModel.H
Foam::polyMesh::faceOwner
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1107
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:259
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::cellShape
An analytical geometric cellShape.
Definition: cellShape.H:69
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
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:118
Foam::cellShape::reset
void reset(const cellModel &model, const labelUList &labels, const bool doCollapse=false)
Reset from components.
Definition: cellShapeI.H:300
Foam::tetMatcher::faceSizeMatch
virtual bool faceSizeMatch(const faceList &, const labelList &) const
Check whether number of face sizes match the shape.
Definition: tetMatcher.C:249
Foam::polyMesh::faces
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1094
f
labelList f(nPoints)
Foam::List< face >
cellMatcher.H
Foam::UList< face >
tetMatcher.H
Foam::tetMatcher::test
static bool test(const UList< face > &faces)
Definition: tetMatcher.C:89
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:72
ListOps.H
Various functions to operate on Lists.
Foam::UList::size
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
Foam::tetMatcher::tetMatcher
tetMatcher()
Default construct.
Definition: tetMatcher.C:103
Foam::primitiveMesh
Cell-face mesh analysis engine.
Definition: primitiveMesh.H:78