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