cuttingSurfaceBaseTemplates.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) 2018-2021 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
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 "volFields.H"
29#include "edgeHashes.H"
30#include "HashOps.H"
31
32// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
33
34template<class EdgeOrientIntersect, class EdgeAlphaIntersect>
36(
37 const primitiveMesh& mesh,
38 const bitSet& cellCuts,
39 const EdgeOrientIntersect& edgeOrientIntersect,
40 const EdgeAlphaIntersect& edgeAlphaIntersect,
41 const bool triangulate,
42 label nFaceCuts
43)
44{
45 // Information required from the mesh
46 const faceList& faces = mesh.faces();
47 const cellList& cells = mesh.cells();
48 const pointField& points = mesh.points();
49
50 // Dynamic lists to handle triangulation and/or missed cuts etc
51 const label nCellCuts = cellCuts.count();
52
53 DynamicList<point> dynCutPoints(4*nCellCuts);
54 DynamicList<face> dynCutFaces(4*nCellCuts);
55 DynamicList<label> dynCutCells(nCellCuts);
56
57 // No nFaceCuts provided? Use a reasonable estimate
58 if (!nFaceCuts)
59 {
60 nFaceCuts = 4*nCellCuts;
61 }
62
63 // Edge to pointId mapping
64 EdgeMap<label> handledEdges(4*nFaceCuts);
65
66 // Local scratch space for face vertices
67 DynamicList<label> localFaceLoop(16);
68
69 // Local scratch space for edge to pointId
70 EdgeMap<label> localEdges(128);
71
72 // Local scratch space for edge to faceId
74
75 // Avoid duplicates for cuts exactly through a mesh point.
76 // No other way to distinguish them, since there is no single edge
77 // that "owns" the point.
78 Map<label> endPoints;
79
80 // Hash of faces (face points) that are exactly on a cell face
81 HashSet<labelList> onCellFace;
82
83
84 // Failure handling
85
86 // Cells where walking failed (concave, degenerate, ...)
87 labelHashSet failedCells;
88
89 // To unwind insertion of end-point cutting
90 labelHashSet localEndPoints;
91
92 // Our recovery point on failure
93 label unwindPoint = 0;
94
95 // Cleanup routine for failed cell cut:
96 const auto unwindWalk =
97 [&](const label failedCellId = -1) -> void
98 {
99 // Discard points introduced
100 dynCutPoints.resize(unwindPoint);
101
102 // Discard end-point cuts
103 endPoints.erase(localEndPoints);
104
105 // Optionally record the failedCellId
106 if (failedCellId != -1)
107 {
108 failedCells.insert(failedCellId);
109 }
110 };
111
112
113 // Loop over cells that are cut
114 for (const label celli : cellCuts)
115 {
116 const cell& cFace = cells[celli];
117
118 // Reset local scratch
119 localEdges.clear();
120 localFaces.clear();
121 localFaceLoop.clear();
122 localEndPoints.clear();
123
124 unwindPoint = dynCutPoints.size();
125
126
127 // Classification for all the points cut - see intersectsFace() above
128 // for more detail
129 unsigned pointCutType = 0u;
130
131 for (const label facei : cFace)
132 {
133 const face& f = faces[facei];
134
135 forAll(f, fp)
136 {
137 edge e(f.edge(fp));
138
139 // Action #1: Orient edge (+ve gradient) and detect intersect
140 if (!edgeOrientIntersect(e))
141 {
142 continue;
143 }
144
145 // Record the edge/face combination for the edge cut.
146 // NB: the second operation is edge::insert() which places
147 // facei in the unoccupied 'slot'
148 localFaces(e).insert(facei);
149
150 // Already handled cut point in this inner-loop?
151 if (localEdges.found(e))
152 {
153 // No new edge cut required
154 continue;
155 }
156
157 // Already handled cut point in the outer-loop?
158 label cutPointId = handledEdges.lookup(e, -1);
159
160 if (cutPointId >= 0)
161 {
162 // Copy existing edge cut-point index
163 localEdges.insert(e, cutPointId);
164 continue;
165 }
166
167 // Expected id for the cut point
168 cutPointId = dynCutPoints.size();
169
170 const point& p0 = points[e.first()];
171 const point& p1 = points[e.last()];
172
173 // Action #2: edge cut alpha
174 const scalar alpha = edgeAlphaIntersect(e);
175
176 if (alpha < SMALL)
177 {
178 pointCutType |= 0x1; // Cut at 0 (first)
179
180 const label endp = e.first();
181
182 if (endPoints.insert(endp, cutPointId))
183 {
184 localEndPoints.insert(endp);
185 dynCutPoints.append(p0);
186 }
187 else
188 {
189 cutPointId = endPoints[endp];
190 }
191 }
192 else if (alpha >= (1.0 - SMALL))
193 {
194 pointCutType |= 0x2; // Cut at 1 (last)
195
196 const label endp = e.last();
197
198 if (endPoints.insert(endp, cutPointId))
199 {
200 localEndPoints.insert(endp);
201 dynCutPoints.append(p1);
202 }
203 else
204 {
205 cutPointId = endPoints[endp];
206 }
207 }
208 else
209 {
210 pointCutType |= 0x4; // Cut between
211
212 dynCutPoints.append((1-alpha)*p0 + alpha*p1);
213 }
214
215 // Introduce new edge cut point
216 localEdges.insert(e, cutPointId);
217 }
218 }
219
220
221 // The keys of localEdges, localFaces are now identical.
222 // The size of either should represent the number of points for
223 // the resulting face loop.
224
225 const label nTargetLoop = localFaces.size();
226
227 if (nTargetLoop < 3)
228 {
229 unwindWalk(celli);
230 continue;
231 }
232
233
234 // Handling cuts between two cells
235 // After the previous edgeIntersectAndOrient call, the edge is oriented
236 // according to the gradient.
237 // If we only ever cut at the same edge end we know that we have
238 // a cut coinciding with a cell face.
239
240 if (pointCutType == 1 || pointCutType == 2)
241 {
242 // Hash the face-points to avoid duplicate faces
243 if (!onCellFace.insert(HashTableOps::values(localEdges, true)))
244 {
246 <<"skip duplicate on-place cut for cell " << celli
247 << " type (" << pointCutType << ")" << endl;
248
249 // A duplicate is not failure
250 unwindWalk();
251 continue;
252 }
253 }
254
255
256 // Start somewhere.
257
258 // Since the local edges are oriented according to the gradient,
259 // they can also be used to determine the correct face orientation.
260
261 const edge refEdge = localFaces.begin().key();
262 label nextFace = localFaces.begin().val()[0];
263
265 << "search face " << nextFace << " IN " << localEdges << endl;
266
267 for
268 (
269 label loopi = 0;
270 localFaces.size() && loopi < 2*nTargetLoop;
271 ++loopi
272 )
273 {
274 bool ok = false;
275
277 {
279 << "lookup " << nextFace << " in " << iter.val() << nl;
280
281 // Find local index (0,1) or -1 on failure
282 const label got = iter.val().which(nextFace);
283
284 if (got != -1)
285 {
286 ok = true;
287
288 // The other face
289 nextFace = iter.val()[(got?0:1)];
290
291 // The edge -> cut point
292 localFaceLoop.append(localEdges[iter.key()]);
293
295 <<" faces " << iter.val()
296 << " point " << localFaceLoop.last()
297 << " edge=" << iter.key() << " nextFace=" << nextFace
298 << nl;
299
300 // Done this connection
301 localFaces.erase(iter);
302 break;
303 }
304 }
305
306 if (!ok)
307 {
308 break;
309 }
310 }
311
312
313 // Could also check if localFaces is now empty.
314
315 if (nTargetLoop != localFaceLoop.size())
316 {
318 << "Warn expected " << nTargetLoop << " but got "
319 << localFaceLoop.size() << endl;
320
321 unwindWalk(celli);
322 continue;
323 }
324
325
326 // Success
327 handledEdges += localEdges;
328
329 face f(localFaceLoop);
330
331 // Orient face to point in the same direction as the edge gradient
332 if ((f.areaNormal(dynCutPoints) & refEdge.vec(points)) < 0)
333 {
334 f.flip();
335 }
336
337 // The cut faces can be quite ugly, so optionally triangulate
338 if (triangulate)
339 {
340 label nTri = f.triangles(dynCutPoints, dynCutFaces);
341 while (nTri--)
342 {
343 dynCutCells.append(celli);
344 }
345 }
346 else
347 {
348 dynCutFaces.append(f);
349 dynCutCells.append(celli);
350 }
351 }
352
353 if (failedCells.size())
354 {
356 << "Failed cuts for " << failedCells.size() << " cells:" << nl
357 << " " << flatOutput(failedCells.sortedToc()) << nl
358 << endl;
359 }
360
361
362 // No cuts? Then no need for any of this information
363 if (dynCutCells.empty())
364 {
365 this->storedPoints().clear();
366 this->storedFaces().clear();
368 }
369 else
370 {
371 this->storedPoints().transfer(dynCutPoints);
372 this->storedFaces().transfer(dynCutFaces);
373 meshCells_.transfer(dynCutCells);
374 }
375}
376
377
378// ************************************************************************* //
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:72
void clear() noexcept
Clear the addressed list, i.e. set the size to zero.
Definition: DynamicListI.H:391
void append(const T &val)
Copy append an element to the end of this list.
Definition: DynamicListI.H:503
void resize(const label len)
Definition: DynamicListI.H:353
Map from edge (expressed as its endpoints) to value. For easier forward declaration it is currently i...
Definition: EdgeMap.H:54
A HashTable with keys but without contents that is similar to std::unordered_set.
Definition: HashSet.H:96
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition: HashSet.H:191
List< Key > sortedToc() const
The table of contents (the keys) in sorted order.
Definition: HashTable.C:137
const T & lookup(const Key &key, const T &deflt) const
Return hashed entry if it exists, or return the given default.
Definition: HashTableI.H:224
bool found(const Key &key) const
Return true if hashed entry is found in table.
Definition: HashTableI.H:100
bool insert(const Key &key, const T &obj)
Copy insert a new entry, not overwriting existing entries.
Definition: HashTableI.H:180
label size() const noexcept
The number of elements in table.
Definition: HashTableI.H:52
bool erase(const iterator &iter)
Erase an entry specified by given iterator.
Definition: HashTable.C:440
void clear()
Clear all entries from table.
Definition: HashTable.C:678
void transfer(List< T > &list)
Definition: List.C:447
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:116
A HashTable to objects of type <T> with a label key.
Definition: Map.H:60
pointField & storedPoints()
Non-const access to global points.
virtual label triangulate()
Triangulate in-place, returning the number of triangles added.
List< face > & storedFaces()
Non-const access to the faces.
const Field< point_type > & points() const noexcept
Return reference to global points.
const List< face_type > & localFaces() const
Return patch faces addressing into local point list.
T & first()
Return the first element of the list.
Definition: UListI.H:202
bool empty() const noexcept
True if the UList is empty (ie, size() is zero)
Definition: UListI.H:427
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
T & last()
Return the last element of the list.
Definition: UListI.H:216
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:66
Description of cuts across cells.
Definition: cellCuts.H:113
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:57
void walkCellCuts(const primitiveMesh &mesh, const bitSet &cellCuts, const EdgeOrientIntersect &edgeOrientIntersect, const EdgeAlphaIntersect &edgeAlphaIntersect, const bool triangulate, label nFaceCuts=0)
Walk cell cuts to create faces.
labelList meshCells_
List of the cells cut.
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:66
vector vec(const UList< point > &pts) const
Return the vector (end - start)
Definition: edgeI.H:417
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1108
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1083
Cell-face mesh analysis engine.
Definition: primitiveMesh.H:79
const cellList & cells() const
const volScalarField & p0
Definition: EEqn.H:36
dynamicFvMesh & mesh
const cellShapeList & cells
#define DebugInfo
Report an information message using Foam::Info.
#define WarningInFunction
Report a warning using Foam::Warning.
List< T > values(const HashTable< T, Key, Hash > &tbl, const bool doSort=false)
List of values from HashTable, optionally sorted.
Definition: HashOps.H:149
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Definition: FlatOutput.H:215
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
labelList f(nPoints)
volScalarField & alpha
volScalarField & e
Definition: createFields.H:11
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
#define forAllIters(container, iter)
Iterate across all elements in the container object.
Definition: stdFoam.H:260