enrichedPatchCutFaces.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) 2017-2020 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
27Description
28 Calculating cut faces of the enriched patch, together with the addressing
29 into master and slave patch.
30
31\*---------------------------------------------------------------------------*/
32
33#include "enrichedPatch.H"
34#include "boolList.H"
35#include "DynamicList.H"
36#include "labelPair.H"
37#include "primitiveMesh.H"
38#include "edgeHashes.H"
39
40// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
41
42// If the cut face gets more than this number of points, it will be checked
43const Foam::label Foam::enrichedPatch::maxFaceSizeDebug_ = 100;
44
45
46// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
47
48// Index of debug signs:
49// x - skip a point
50// l - left turn
51// r - right turn
52
53void Foam::enrichedPatch::calcCutFaces() const
54{
55 if (cutFacesPtr_ || cutFaceMasterPtr_ || cutFaceSlavePtr_)
56 {
58 << "Cut faces addressing already calculated."
59 << abort(FatalError);
60 }
61
62 const faceList& enFaces = enrichedFaces();
63 const labelList& mp = meshPoints();
64 const faceList& lf = localFaces();
65 const pointField& lp = localPoints();
66 const labelListList& pp = pointPoints();
67 // Pout<< "enFaces: " << enFaces << endl;
68 // Pout<< "lf: " << lf << endl;
69 // Pout<< "lp: " << lp << endl;
70 // Pout<< "pp: " << pp << endl;
71 const Map<labelList>& masterPointFaceAddr = masterPointFaces();
72
73 // Prepare the storage
74 DynamicList<face> cf(2*lf.size());
75 DynamicList<label> cfMaster(2*lf.size());
76 DynamicList<label> cfSlave(2*lf.size());
77
78 // Algorithm
79 // Go through all the faces
80 // 1) For each face, start at point zero and grab the edge zero.
81 // Both points are added into the cut face. Mark the first edge
82 // as used and position the current point as the end of the first
83 // edge and previous point as point zero.
84 // 2) Grab all possible points from the current point. Excluding
85 // the previous point find the point giving the biggest right
86 // turn. Add the point to the face and mark the edges as used.
87 // Continue doing this until the face is complete, i.e. the next point
88 // to add is the first point of the face.
89 // 3) Once the facet is completed, register its owner from the face
90 // it has been created from (remember that the first lot of faces
91 // inserted into the enriched faces list are the slaves, to allow
92 // matching of the other side).
93 // 4) If the facet is created from an original slave face, go
94 // through the master patch and try to identify the master face
95 // this facet belongs to. This is recognised by the fact that the
96 // faces consists exclusively of the points of the master face and
97 // the points projected onto the face.
98
99 // Create a set of edge usage parameters
100 edgeHashSet edgesUsedOnce(pp.size());
101 edgeHashSet edgesUsedTwice(pp.size()*primitiveMesh::edgesPerPoint_);
102
103
104 forAll(lf, facei)
105 {
106 const face& curLocalFace = lf[facei];
107 const face& curGlobalFace = enFaces[facei];
108
109 // Pout<< "Doing face " << facei
110 // << " local: " << curLocalFace
111 // << " or " << curGlobalFace
112 // << endl;
113
114 // if (facei < slavePatch_.size())
115 // {
116 // Pout<< "original slave: " << slavePatch_[facei]
117 // << " local: " << slavePatch_.localFaces()[facei] << endl;
118 // }
119 // else
120 // {
121 // Pout<< "original master: "
122 // << masterPatch_[facei - slavePatch_.size()] << " "
123 // << masterPatch_.localFaces()[facei - slavePatch_.size()]
124 // << endl;
125 // }
126 // {
127 // pointField facePoints = curLocalFace.points(lp);
128 // forAll(curLocalFace, pointi)
129 // {
130 // Pout<< "v " << facePoints[pointi].x() << " "
131 // << facePoints[pointi].y() << " "
132 // << facePoints[pointi].z() << endl;
133 // }
134 // }
135
136 // Track the usage of face edges. When all edges are used, the
137 // face decomposition is complete.
138 // The seed edges include all the edges of the original face + all edges
139 // of other faces that have been used in the construction of the
140 // facet. Edges from other faces can be considered as
141 // internal to the current face if used only once.
142
143 // Track the edge usage to avoid duplicate faces and reset it to unused
144 boolList usedFaceEdges(curLocalFace.size(), false);
145
146 SLList<edge> edgeSeeds;
147
148 // Insert the edges of current face into the seed list.
149 edgeList cfe = curLocalFace.edges();
150 for (const edge& e : cfe)
151 {
152 edgeSeeds.append(e);
153 }
154
155 // Grab face normal
156 const vector normal = curLocalFace.unitNormal(lp);
157
158 while (edgeSeeds.size())
159 {
160 // Pout<< "edgeSeeds.size(): " << edgeSeeds.size() << endl;
161
162 const edge curEdge = edgeSeeds.removeHead();
163
164 // Locate the edge in current face
165 const label curEdgeWhich = curLocalFace.which(curEdge.start());
166
167 // Check if the edge is in current face and if it has been
168 // used already. If so, skip it
169 if
170 (
171 curEdgeWhich > -1
172 && curLocalFace.nextLabel(curEdgeWhich) == curEdge.end()
173 )
174 {
175 // Edge is in the starting face.
176 // If unused, mark it as used; if used, skip it
177 if (usedFaceEdges[curEdgeWhich]) continue;
178
179 usedFaceEdges[curEdgeWhich] = true;
180 }
181
182 // If the edge has already been used twice, skip it
183 if (edgesUsedTwice.found(curEdge)) continue;
184
185 // Pout<< "Trying new edge (" << mp[curEdge.start()]
186 // << ", " << mp[curEdge.end()]
187 // << ") seed: " << curEdge
188 // << " used: " << edgesUsedTwice.found(curEdge)
189 // << endl;
190
191 // Estimate the size of cut face as twice the size of original face
192 DynamicList<label> cutFaceGlobalPoints(2*curLocalFace.size());
193 DynamicList<label> cutFaceLocalPoints(2*curLocalFace.size());
194
195 // Found unused edge.
196 label prevPointLabel = curEdge.start();
197 cutFaceGlobalPoints.append(mp[prevPointLabel]);
198 cutFaceLocalPoints.append(prevPointLabel);
199 // Pout<< "prevPointLabel: " << mp[prevPointLabel] << endl;
200
201 // Grab current point and append it to the list
202 label curPointLabel = curEdge.end();
203 point curPoint = lp[curPointLabel];
204 cutFaceGlobalPoints.append(mp[curPointLabel]);
205 cutFaceLocalPoints.append(curPointLabel);
206
207 bool completedCutFace = false;
208
209 label faceSizeDebug = maxFaceSizeDebug_;
210
211 do
212 {
213 // Grab the next point options
214
215 // Pout<< "curPointLabel: " << mp[curPointLabel] << endl;
216
217 const labelList& nextPoints = pp[curPointLabel];
218
219 // Pout<< "nextPoints: "
220 // << labelUIndList(mp, nextPoints)
221 // << endl;
222
223 // Get the vector along the edge and the right vector
224 vector ahead = curPoint - lp[prevPointLabel];
225 ahead.removeCollinear(normal);
226 ahead.normalise();
227
228 const vector right = normalised(normal ^ ahead);
229
230 // Pout<< "normal: " << normal
231 // << " ahead: " << ahead
232 // << " right: " << right
233 // << endl;
234
235 scalar atanTurn = -GREAT;
236 label bestAtanPoint = -1;
237
238 forAll(nextPoints, nextI)
239 {
240 // Exclude the point we are coming from; there will always
241 // be more than one edge, so this is safe
242 if (nextPoints[nextI] != prevPointLabel)
243 {
244 // Pout<< "cur point: " << curPoint
245 // << " trying for point: "
246 // << mp[nextPoints[nextI]]
247 // << " " << lp[nextPoints[nextI]];
248 vector newDir = lp[nextPoints[nextI]] - curPoint;
249 // Pout<< " newDir: " << newDir
250 // << " mag: " << mag(newDir) << flush;
251 newDir.removeCollinear(normal);
252 scalar magNewDir = mag(newDir);
253 // Pout<< " corrected: " << newDir
254 // << " mag: " << mag(newDir) << flush;
255
256 if (magNewDir < SMALL)
257 {
259 << "projection error: slave patch probably "
260 << "does not project onto master. "
261 << "Please switch on "
262 << "enriched patch debug for more info"
263 << abort(FatalError);
264 }
265
266 newDir /= magNewDir;
267
268 const scalar curAtanTurn =
269 atan2(newDir & right, newDir & ahead);
270
271 // Pout<< " atan: " << curAtanTurn << endl;
272
273 if (curAtanTurn > atanTurn)
274 {
275 bestAtanPoint = nextPoints[nextI];
276 atanTurn = curAtanTurn;
277 }
278 } // end of prev point skip
279 } // end of next point selection
280
281 // Pout<< " bestAtanPoint: " << bestAtanPoint << " or "
282 // << mp[bestAtanPoint]
283 // << endl;
284
285 // Selected next best point.
286
287 // Pout<< "cutFaceGlobalPoints: "
288 // << cutFaceGlobalPoints
289 // << endl;
290
291 // Check if the edge about to be added has been used
292 // in the current face or twice in other faces. If
293 // so, the face is bad.
294 const label whichNextPoint = curLocalFace.which(curPointLabel);
295
296 if
297 (
298 whichNextPoint > -1
299 && curLocalFace.nextLabel(whichNextPoint) == bestAtanPoint
300 && usedFaceEdges[whichNextPoint]
301 )
302 {
303 // This edge is already used in current face
304 // face cannot be good; start on a new one
305
306 // Pout<< "Double usage in current face, cannot be good"
307 // << endl;
308
309 completedCutFace = true;
310 }
311
312
313 if (edgesUsedTwice.found(edge(curPointLabel, bestAtanPoint)))
314 {
315 // This edge is already used -
316 // face cannot be good; start on a new one
317 completedCutFace = true;
318
319 // Pout<< "Double usage elsewhere, cannot be good" << endl;
320 }
321
322 if (completedCutFace) continue;
323
324 // Insert the next best point into the list
325 if (mp[bestAtanPoint] == cutFaceGlobalPoints[0])
326 {
327 // Face is completed. Save it and move on to the next
328 // available edge
329 completedCutFace = true;
330
331 if (debug)
332 {
333 Pout<< " local: " << cutFaceLocalPoints
334 << " one side: " << facei;
335 }
336
337 // Append the face
338 face cutFaceGlobal;
339 cutFaceGlobal.transfer(cutFaceGlobalPoints);
340
341 cf.append(cutFaceGlobal);
342
343 face cutFaceLocal;
344 cutFaceLocal.transfer(cutFaceLocalPoints);
345
346 // Pout<< "\ncutFaceLocal: "
347 // << cutFaceLocal.points(lp)
348 // << endl;
349
350 // Go through all edges of the cut faces.
351 // If the edge corresponds to a starting face edge,
352 // mark the starting face edge as true
353
354 forAll(cutFaceLocal, cuti)
355 {
356 const edge curCutFaceEdge
357 (
358 cutFaceLocal[cuti],
359 cutFaceLocal.nextLabel(cuti)
360 );
361
362 // Increment the usage count using two hash sets
363 auto euoIter = edgesUsedOnce.find(curCutFaceEdge);
364
365 if (!euoIter.found())
366 {
367 // Pout<< "Found edge not used before: "
368 // << curCutFaceEdge
369 // << endl;
370 edgesUsedOnce.insert(curCutFaceEdge);
371 }
372 else
373 {
374 // Pout<< "Found edge used once: "
375 // << curCutFaceEdge
376 // << endl;
377 edgesUsedOnce.erase(euoIter);
378 edgesUsedTwice.insert(curCutFaceEdge);
379 }
380
381 const label curCutFaceEdgeWhich = curLocalFace.which
382 (
383 curCutFaceEdge.start()
384 );
385
386 if
387 (
388 curCutFaceEdgeWhich > -1
389 && curLocalFace.nextLabel(curCutFaceEdgeWhich)
390 == curCutFaceEdge.end()
391 )
392 {
393 // Found edge in original face
394
395 // Pout<< "Found edge in orig face: "
396 // << curCutFaceEdge << ": "
397 // << curCutFaceEdgeWhich
398 // << endl;
399
400 usedFaceEdges[curCutFaceEdgeWhich] = true;
401 }
402 else
403 {
404 // Edge not in original face. Add it to seeds
405
406 // Pout<< "Found new edge seed: "
407 // << curCutFaceEdge
408 // << endl;
409
410 edgeSeeds.append(curCutFaceEdge.reverseEdge());
411 }
412 }
413
414 // Find out what the other side is
415
416 // Algorithm
417
418 // If the face is in the slave half of the
419 // enrichedFaces lists, it may be matched against
420 // the master face. It can be recognised by the
421 // fact that all its points belong to one master
422 // face. For this purpose:
423 // 1) Grab the master faces around the point zero
424 // of the cut face and collect all master faces
425 // around it.
426 // 2) Loop through the rest of cut face points and
427 // if the point refers to one of the faces marked
428 // by point zero, increment its count.
429 // 3) When completed, the face whose count is
430 // equal to the number of points in the cut face
431 // is the other side. If this is not the case, there is no
432 // face on the other side.
433
434 if (facei < slavePatch_.size())
435 {
436 const auto mpfAddrIter =
437 masterPointFaceAddr.cfind(cutFaceGlobal[0]);
438
439 bool otherSideFound = false;
440
441 if (mpfAddrIter.found())
442 {
443 bool miss = false;
444
445 // Create the label count using point zero
446 const labelList& masterFacesOfPZero = mpfAddrIter();
447
448 labelList hits(masterFacesOfPZero.size(), 1);
449
450 for
451 (
452 label pointi = 1;
453 pointi < cutFaceGlobal.size();
454 pointi++
455 )
456 {
457 const auto mpfAddrPointIter =
458 masterPointFaceAddr.cfind
459 (
460 cutFaceGlobal[pointi]
461 );
462
463 if (!mpfAddrPointIter.found())
464 {
465 // Point is off the master patch. Skip
466 miss = true;
467 break;
468 }
469
470 const labelList& curMasterFaces =
471 mpfAddrPointIter();
472
473 // For every current face, try to find it in the
474 // zero-list
475 for (const label facei : curMasterFaces)
476 {
477 forAll(masterFacesOfPZero, j)
478 {
479 if (facei == masterFacesOfPZero[j])
480 {
481 hits[j]++;
482 break;
483 }
484 }
485 }
486 }
487
488 // If all point are found attempt matching
489 if (!miss)
490 {
491 forAll(hits, pointi)
492 {
493 if (hits[pointi] == cutFaceGlobal.size())
494 {
495 // Found other side.
496 otherSideFound = true;
497
498 cfMaster.append
499 (
500 masterFacesOfPZero[pointi]
501 );
502
503 cfSlave.append(facei);
504
505 // Reverse the face such that it
506 // points out of the master patch
507 cf.last().flip();
508
509 if (debug)
510 {
511 Pout<< " other side: "
512 << masterFacesOfPZero[pointi]
513 << endl;
514 }
515 } // end of hits
516 } // end of for all hits
517
518 } // end of not miss
519
520 if (!otherSideFound || miss)
521 {
522 if (debug)
523 {
524 Pout<< " solo slave A" << endl;
525 }
526
527 cfMaster.append(-1);
528 cfSlave.append(facei);
529 }
530 }
531 else
532 {
533 // First point not in master patch
534 if (debug)
535 {
536 Pout<< " solo slave B" << endl;
537 }
538
539 cfMaster.append(-1);
540 cfSlave.append(facei);
541 }
542 }
543 else
544 {
545 if (debug)
546 {
547 Pout<< " master side" << endl;
548 }
549
550 cfMaster.append(facei - slavePatch_.size());
551 cfSlave.append(-1);
552 }
553 }
554 else
555 {
556 // Face not completed. Prepare for the next point search
557 prevPointLabel = curPointLabel;
558 curPointLabel = bestAtanPoint;
559 curPoint = lp[curPointLabel];
560 cutFaceGlobalPoints.append(mp[curPointLabel]);
561 cutFaceLocalPoints.append(curPointLabel);
562
563 if (debug || cutFaceGlobalPoints.size() > faceSizeDebug)
564 {
565 faceSizeDebug *= 2;
566
567 // Check for duplicate points in the face
568 forAll(cutFaceGlobalPoints, checkI)
569 {
570 for
571 (
572 label checkJ = checkI + 1;
573 checkJ < cutFaceGlobalPoints.size();
574 checkJ++
575 )
576 {
577 if
578 (
579 cutFaceGlobalPoints[checkI]
580 == cutFaceGlobalPoints[checkJ]
581 )
582 {
583 // Shrink local points for debugging
584 cutFaceLocalPoints.shrink();
585
586 face origFace;
587 face origFaceLocal;
588 if (facei < slavePatch_.size())
589 {
590 origFace = slavePatch_[facei];
591 origFaceLocal =
592 slavePatch_.localFaces()[facei];
593 }
594 else
595 {
596 origFace =
597 masterPatch_
598 [facei - slavePatch_.size()];
599
600 origFaceLocal =
601 masterPatch_.localFaces()
602 [facei - slavePatch_.size()];
603 }
604
606 << "Duplicate point found in cut face. "
607 << "Error in the face cutting "
608 << "algorithm for global face "
609 << origFace << " local face "
610 << origFaceLocal << nl
611 << "Slave size: " << slavePatch_.size()
612 << " Master size: "
613 << masterPatch_.size()
614 << " index: " << facei << ".\n"
615 << "Face: " << curGlobalFace << nl
616 << "Cut face: " << cutFaceGlobalPoints
617 << " local: " << cutFaceLocalPoints
618 << nl << "Points: "
619 << face(cutFaceLocalPoints).points(lp)
620 << abort(FatalError);
621 }
622 }
623 }
624 } // end of debug
625 }
626 } while (!completedCutFace);
627 } // end of used edges
628
629 if (debug)
630 {
631 Pout<< " Finished face " << facei << endl;
632 }
633
634 } // end of local faces
635
636 // Re-pack lists into compact storage
637 cutFacesPtr_.reset(new faceList(std::move(cf)));
638 cutFaceMasterPtr_.reset(new labelList(std::move(cfMaster)));
639 cutFaceSlavePtr_.reset(new labelList(std::move(cfSlave)));
640}
641
642
643void Foam::enrichedPatch::clearCutFaces()
644{
645 cutFacesPtr_.reset(nullptr);
646 cutFaceMasterPtr_.reset(nullptr);
647 cutFaceSlavePtr_.reset(nullptr);
648}
649
650
651// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
652
654{
655 if (!cutFacesPtr_)
656 {
657 calcCutFaces();
658 }
659
660 return *cutFacesPtr_;
661}
662
663
665{
666 if (!cutFaceMasterPtr_)
667 {
668 calcCutFaces();
669 }
670
671 return *cutFaceMasterPtr_;
672}
673
674
676{
677 if (!cutFaceSlavePtr_)
678 {
679 calcCutFaces();
680 }
681
682 return *cutFaceSlavePtr_;
683}
684
685
686// ************************************************************************* //
void append(const T &val)
Append an element at the end of the list.
Definition: ListI.H:175
const List< face_type > & localFaces() const
Return patch faces addressing into local point list.
const faceList & localFaces() const
Return local faces.
const pointField & localPoints() const
Return local points.
const Map< labelList > & masterPointFaces() const
Master point face addressing.
const labelList & cutFaceMaster() const
Return cut face master list.
const faceList & enrichedFaces() const
Return enriched faces.
const labelList & meshPoints() const
Return mesh points.
const labelListList & pointPoints() const
Return point-point addressing.
const labelList & cutFaceSlave() const
Return cut face slave list.
const faceList & cutFaces() const
Return list of cut faces.
static const unsigned edgesPerPoint_
Estimated number of edges per point.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
const dimensionedScalar mp
Proton mass.
List< label > labelList
A List of labels.
Definition: List.H:66
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
Definition: quaternionI.H:680
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
vector point
Point is a vector.
Definition: point.H:43
dimensionedScalar atan2(const dimensionedScalar &x, const dimensionedScalar &y)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:56
errorManip< error > abort(error &err)
Definition: errorManip.H:144
HashSet< edge, Hash< edge > > edgeHashSet
A HashSet with edge for its key.
Definition: edgeHashes.H:51
List< bool > boolList
A List of bools.
Definition: List.H:64
error FatalError
List< edge > edgeList
A List of edges.
Definition: edgeList.H:63
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
List< face > faceList
A List of faces.
Definition: faceListFwd.H:47
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
volScalarField & e
Definition: createFields.H:11
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333