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