geomCellLooper.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) 2019 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 "geomCellLooper.H"
30 #include "polyMesh.H"
31 #include "DynamicList.H"
32 #include "plane.H"
33 #include "meshTools.H"
34 #include "SortableList.H"
35 #include "triSurfaceTools.H"
36 #include "HashSet.H"
37 #include "ListOps.H"
38 #include "transform.H"
39 
41 
42 
43 // Extension factor of edges to make sure we catch intersections through
44 // edge endpoints
45 const Foam::scalar Foam::geomCellLooper::pointEqualTol_ = 1e-3;
46 
47 
48 // Snap cuts through edges onto edge endpoints. Fraction of edge length.
49 Foam::scalar Foam::geomCellLooper::snapTol_ = 0.1;
50 
51 
52 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
53 
54 namespace Foam
55 {
56 
57 defineTypeNameAndDebug(geomCellLooper, 0);
58 addToRunTimeSelectionTable(cellLooper, geomCellLooper, word);
59 
60 }
61 
62 
63 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
64 
65 Foam::scalar Foam::geomCellLooper::minEdgeLen(const label vertI) const
66 {
67  scalar minLen = GREAT;
68 
69  const labelList& pEdges = mesh().pointEdges()[vertI];
70 
71  forAll(pEdges, pEdgeI)
72  {
73  const edge& e = mesh().edges()[pEdges[pEdgeI]];
74 
75  minLen = min(minLen, e.mag(mesh().points()));
76  }
77  return minLen;
78 }
79 
80 
81 // Cut edge with plane. Return true and set weight to fraction between
82 // edge-start and edge-end
83 bool Foam::geomCellLooper::cutEdge
84 (
85  const plane& cutPlane,
86  const label edgeI,
87  scalar& weight
88 ) const
89 {
90  const pointField& pts = mesh().points();
91 
92  const edge& e = mesh().edges()[edgeI];
93 
94  scalar s = cutPlane.normalIntersect(pts[e.start()], e.vec(pts));
95 
96  if ((s > -pointEqualTol_) && (s < 1 + pointEqualTol_))
97  {
98  weight = s;
99 
100  return true;
101  }
102  else
103  {
104  // Make sure we don't use this value
105  weight = -GREAT;
106 
107  return false;
108  }
109 }
110 
111 
112 // If edge close enough to endpoint snap to endpoint.
113 Foam::label Foam::geomCellLooper::snapToVert
114 (
115  const scalar tol,
116  const label edgeI,
117  const scalar weight
118 ) const
119 {
120  const edge& e = mesh().edges()[edgeI];
121 
122  if (weight < tol)
123  {
124  return e.start();
125  }
126  else if (weight > (1-tol))
127  {
128  return e.end();
129  }
130 
131  return -1;
132 }
133 
134 
135 void Foam::geomCellLooper::getBase(const vector& n, vector& e0, vector& e1)
136  const
137 {
138  // Guess for vector normal to n.
139  vector base(1,0,0);
140 
141  scalar nComp = n & base;
142 
143  if (mag(nComp) > 0.8)
144  {
145  // Was bad guess. Try with different vector.
146 
147  base.x() = 0;
148  base.y() = 1;
149 
150  nComp = n & base;
151 
152  if (mag(nComp) > 0.8)
153  {
154  base.y() = 0;
155  base.z() = 1;
156 
157  nComp = n & base;
158  }
159  }
160 
161 
162  // Use component normal to n as base vector.
163  e0 = normalised(base - nComp * n);
164 
165  e1 = n ^ e0;
166 
167  //Pout<< "Coord system:" << endl
168  // << " n : " << n << ' ' << mag(n) << endl
169  // << " e0 : " << e0 << ' ' << mag(e0) << endl
170  // << " e1 : " << e1 << ' ' << mag(e1) << endl
171  // << endl;
172 }
173 
174 
175 // Return true if the cut edge at loop[index] is preceded by cuts through
176 // the edge end points.
177 bool Foam::geomCellLooper::edgeEndsCut
178 (
179  const labelList& loop,
180  const label index
181 ) const
182 {
183  label edgeI = getEdge(loop[index]);
184 
185  const edge& e = mesh().edges()[edgeI];
186 
187  const label prevCut = loop[loop.rcIndex(index)];
188  const label nextCut = loop[loop.fcIndex(index)];
189 
190  if (!isEdge(prevCut) && !isEdge(nextCut))
191  {
192  // Cuts before and after are both vertices. Check if both
193  // the edge endpoints
194  label v0 = getVertex(prevCut);
195  label v1 = getVertex(nextCut);
196 
197  if
198  (
199  (v0 == e[0] && v1 == e[1])
200  || (v0 == e[1] && v1 == e[0])
201  )
202  {
203  return true;
204  }
205  }
206  return false;
207 }
208 
209 
210 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
211 
212 // Construct from components
213 Foam::geomCellLooper::geomCellLooper(const polyMesh& mesh)
214 :
216 {}
217 
218 
219 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
220 
222 {}
223 
224 
225 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
226 
228 (
229  const vector& refDir,
230  const label celli,
231  const boolList& vertIsCut,
232  const boolList& edgeIsCut,
233  const scalarField& edgeWeight,
234 
235  labelList& loop,
236  scalarField& loopWeights
237 ) const
238 {
239  // Cut through cell centre normal to refDir.
240  return cut
241  (
242  plane(mesh().cellCentres()[celli], refDir),
243  celli,
244  vertIsCut,
245  edgeIsCut,
246  edgeWeight,
247  loop,
248  loopWeights
249  );
250 }
251 
252 
254 (
255  const plane& cutPlane,
256  const label celli,
257  const boolList&,
258  const boolList&,
259  const scalarField&,
260 
261  labelList& loop,
262  scalarField& loopWeights
263 ) const
264 {
265  const pointField& points = mesh().points();
266  const edgeList& edges = mesh().edges();
267 
268  // Find all cuts through edges.
269  // Special cases:
270  // - edge cut close to endpoint. Snap to endpoint.
271  // - edge endpoint close to plane (but edge not cut). Snap to endpoint.
272  // - both endpoints close to plane. Edge parallel to plane. No need to
273  // cut to edge.
274  // Note: any snap-to-point check must be based on all edges using a point
275  // since otherwise cut through close to point but on neighbouring edge
276  // might not be snapped.
277 
278  // Size overly big.
279  label nEstCuts = 2*mesh().cells()[celli].size();
280 
281  DynamicList<label> localLoop(nEstCuts);
282  DynamicList<scalar> localLoopWeights(nEstCuts);
283 
284  // Points checked. Used to make sure we don't cut edge and edge endpoints
285  // at the same time.
286  labelHashSet checkedPoints(nEstCuts);
287 
288  const labelList& cellEdges = mesh().cellEdges()[celli];
289 
290  forAll(cellEdges, i)
291  {
292  label edgeI = cellEdges[i];
293 
294  const edge& e = edges[edgeI];
295 
296  bool useStart = false;
297 
298  bool useEnd = false;
299 
300  //
301  // Check distance of endpoints to cutPlane
302  //
303 
304  if (!checkedPoints.found(e.start()))
305  {
306  checkedPoints.insert(e.start());
307 
308  scalar typStartLen = pointEqualTol_ * minEdgeLen(e.start());
309 
310  // Check distance of startPt to plane.
311  if (cutPlane.distance(points[e.start()]) < typStartLen)
312  {
313  // Use point.
314  localLoop.append(vertToEVert(e.start()));
315  localLoopWeights.append(-GREAT);
316 
317  useStart = true;
318  }
319  }
320  if (!checkedPoints.found(e.end()))
321  {
322  checkedPoints.insert(e.end());
323 
324  scalar typEndLen = pointEqualTol_ * minEdgeLen(e.end());
325 
326  // Check distance of endPt to plane.
327  if (cutPlane.distance(points[e.end()]) < typEndLen)
328  {
329  // Use point.
330  localLoop.append(vertToEVert(e.end()));
331  localLoopWeights.append(-GREAT);
332 
333  useEnd = true;
334  }
335  }
336 
337  //
338  // Check cut of edge
339  //
340 
341  if (!useEnd && !useStart)
342  {
343  // Edge end points not close to plane so edge not close to
344  // plane. Cut edge.
345  scalar cutWeight;
346 
347  if (cutEdge(cutPlane, edgeI, cutWeight))
348  {
349  // Snap edge cut to vertex.
350  label cutVertI = snapToVert(snapTol_, edgeI, cutWeight);
351 
352  if (cutVertI == -1)
353  {
354  // Proper cut of edge
355  localLoop.append(edgeToEVert(edgeI));
356  localLoopWeights.append(cutWeight);
357  }
358  else
359  {
360  // Cut through edge end point. Might be duplicate
361  // since connected edge might also be snapped to same
362  // endpoint so only insert if unique.
363  label cut = vertToEVert(cutVertI);
364 
365  if (!localLoop.found(cut))
366  {
367  localLoop.append(vertToEVert(cutVertI));
368  localLoopWeights.append(-GREAT);
369  }
370  }
371  }
372  }
373  }
374 
375  if (localLoop.size() <= 2)
376  {
377  return false;
378  }
379 
380  localLoop.shrink();
381  localLoopWeights.shrink();
382 
383 
384  // Get points on loop and centre of loop
385  pointField loopPoints(localLoop.size());
386  point ctr(Zero);
387 
388  forAll(localLoop, i)
389  {
390  loopPoints[i] = coord(localLoop[i], localLoopWeights[i]);
391  ctr += loopPoints[i];
392  }
393  ctr /= localLoop.size();
394 
395 
396  // Get base vectors of coordinate system normal to refDir
397  vector e0, e1;
398  getBase(cutPlane.normal(), e0, e1);
399 
400 
401  // Get sorted angles from point on loop to centre of loop.
402  SortableList<scalar> sortedAngles(localLoop.size());
403 
404  forAll(sortedAngles, i)
405  {
406  const vector toCtr = normalised(loopPoints[i] - ctr);
407 
408  sortedAngles[i] = pseudoAngle(e0, e1, toCtr);
409  }
410  sortedAngles.sort();
411 
412  loop.setSize(localLoop.size());
413  loopWeights.setSize(loop.size());
414 
415 
416  // Fill loop and loopweights according to sorted angles
417 
418  const labelList& indices = sortedAngles.indices();
419 
420  forAll(indices, i)
421  {
422  loop[i] = localLoop[indices[i]];
423  loopWeights[i] = localLoopWeights[indices[i]];
424  }
425 
426 
427  // Check for cut edges along already cut vertices.
428  bool filterLoop = false;
429 
430  forAll(loop, i)
431  {
432  label cut = loop[i];
433 
434  if (isEdge(cut) && edgeEndsCut(loop, i))
435  {
436  filterLoop = true;
437  }
438  }
439 
440  if (filterLoop)
441  {
442  // Found edges in loop where both end points are also in loop. This
443  // is superfluous information so we can remove it.
444 
445  labelList filteredLoop(loop.size());
446  scalarField filteredLoopWeights(loopWeights.size());
447  label filterI = 0;
448 
449  forAll(loop, i)
450  {
451  label cut = loop[i];
452 
453  if (isEdge(cut) && edgeEndsCut(loop, i))
454  {
455  // Superfluous cut. Edge end points cut and edge itself as well.
456  }
457  else
458  {
459  filteredLoop[filterI] = loop[i];
460  filteredLoopWeights[filterI] = loopWeights[i];
461  filterI++;
462  }
463  }
464  filteredLoop.setSize(filterI);
465  filteredLoopWeights.setSize(filterI);
466 
467  loop.transfer(filteredLoop);
468  loopWeights.transfer(filteredLoopWeights);
469  }
470 
471  if (debug&2)
472  {
473  Pout<< "cell:" << celli << endl;
474  forAll(loop, i)
475  {
476  Pout<< "At angle:" << sortedAngles[i] << endl
477  << " cut:";
478 
479  writeCut(Pout, loop[i], loopWeights[i]);
480 
481  Pout<< " coord:" << coord(loop[i], loopWeights[i]);
482 
483  Pout<< endl;
484  }
485  }
486 
487  return true;
488 }
489 
490 
491 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::plane::distance
scalar distance(const point &p) const
Return distance (magnitude) from the given point to the plane.
Definition: planeI.H:75
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:74
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
Foam::polyMesh::points
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1038
meshTools.H
Foam::plane::normal
const vector & normal() const
The plane unit normal.
Definition: planeI.H:39
s
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Definition: gmvOutputSpray.H:25
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::DynamicList< label >
Foam::geomCellLooper::cut
virtual bool cut(const vector &refDir, const label celli, const boolList &vertIsCut, const boolList &edgeIsCut, const scalarField &edgeWeight, labelList &loop, scalarField &loopWeights) const
Create cut along circumference of celli. Gets current mesh cuts.
Definition: geomCellLooper.C:228
Foam::primitiveMesh::cells
const cellList & cells() const
Definition: primitiveMeshCells.C:138
Foam::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:63
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
Foam::primitiveMesh::edges
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
Definition: primitiveMeshEdges.C:505
geomCellLooper.H
Foam::Pout
prefixOSstream Pout
An Ostream wrapper for parallel output to std::cout.
Foam::primitiveMesh::pointEdges
const labelListList & pointEdges() const
Definition: primitiveMeshEdges.C:516
polyMesh.H
Foam::HashSet< label, Hash< label > >
Foam::DynamicList::shrink
DynamicList< T, SizeMin > & shrink()
Shrink the allocated space to the number of elements used.
Definition: DynamicListI.H:376
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::plane
Geometric class that creates a 3D plane and can return the intersection point between a line and the ...
Definition: plane.H:89
n
label n
Definition: TABSMDCalcMethod2.H:31
SortableList.H
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
Foam::Field< scalar >
plane.H
Foam::primitiveMesh::cellEdges
const labelListList & cellEdges() const
Definition: primitiveMeshCellEdges.C:120
Foam::DynamicList::append
DynamicList< T, SizeMin > & append(const T &val)
Append an element to the end of this list.
Definition: DynamicListI.H:472
Foam::List::transfer
void transfer(List< T > &list)
Definition: List.C:436
Foam::cellLooper
Abstract base class. Concrete implementations know how to cut a cell (i.e. determine a loop around th...
Definition: cellLooper.H:71
cut
Patchify triangles based on orientation w.r.t other (triangulated or triangulatable) surfaces.
HashSet.H
Foam::SortableList
A list that is sorted upon construction or when explicitly requested with the sort() method.
Definition: List.H:66
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::pseudoAngle
scalar pseudoAngle(const vector &e0, const vector &e1, const vector &vec)
Estimate angle of vec in coordinate system (e0, e1, e0^e1).
Definition: transform.H:424
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:51
Foam::normalised
VectorSpace< Form, Cmpt, Ncmpts > normalised(const VectorSpace< Form, Cmpt, Ncmpts > &vs)
Definition: VectorSpaceI.H:496
Foam::geomCellLooper::~geomCellLooper
virtual ~geomCellLooper()
Destructor.
Definition: geomCellLooper.C:221
Foam::Vector< scalar >
Foam::List< bool >
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::HashSet::insert
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition: HashSet.H:182
DynamicList.H
ListOps.H
Various functions to operate on Lists.
transform.H
3D tensor transformation operations.
Foam::HashTable::found
bool found(const Key &key) const
Return true if hashed entry is found in table.
Definition: HashTableI.H:100
Foam::List::setSize
void setSize(const label newSize)
Alias for resize(const label)
Definition: ListI.H:146
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
triSurfaceTools.H
Foam::edgeVertex::mesh
const polyMesh & mesh() const
Definition: edgeVertex.H:100