searchableExtrudedCircle.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) 2016-2017 OpenFOAM Foundation
9  Copyright (C) 2019-2022 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 
31 #include "Time.H"
32 #include "edgeMesh.H"
33 #include "indexedOctree.H"
34 #include "treeDataEdge.H"
36 #include "quaternion.H"
37 
38 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39 
40 namespace Foam
41 {
42  defineTypeNameAndDebug(searchableExtrudedCircle, 0);
44  (
45  searchableSurface,
46  searchableExtrudedCircle,
47  dict
48  );
50  (
51  searchableSurface,
52  searchableExtrudedCircle,
53  dict,
54  extrudedCircle
55  );
56 }
57 
58 
59 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
60 
61 Foam::searchableExtrudedCircle::searchableExtrudedCircle
62 (
63  const IOobject& io,
64  const dictionary& dict
65 )
66 :
68  eMeshPtr_
69  (
70  edgeMesh::New
71  (
72  IOobject
73  (
74  dict.get<word>("file"), // name
75  io.time().constant(), // instance
76  "geometry", // local
77  io.time(), // registry
78  IOobject::MUST_READ,
79  IOobject::NO_WRITE,
80  IOobject::NO_REGISTER
81  ).objectPath()
82  )
83  ),
84  radius_(dict.get<scalar>("radius"))
85 {
86  const edgeMesh& eMesh = eMeshPtr_();
87 
88  const pointField& points = eMesh.points();
89  const edgeList& edges = eMesh.edges();
90  bounds() = boundBox(points, false);
91 
92  // Make the boundBox into a perfect cube around its centre
93  const scalar halfWidth = mag(0.5*bounds().span());
94 
95  bounds().reset(bounds().centre());
96  bounds().grow(halfWidth);
97 
98  // Slightly extended bb. Slightly off-centred just so on symmetric
99  // geometry there are less face/edge aligned items.
100  treeBoundBox bb(bounds());
101  bb.grow(ROOTVSMALL);
102 
103  edgeTree_.reset
104  (
106  (
107  treeDataEdge(edges, points), // All edges
108 
109  bb, // overall search domain
110  8, // maxLevel
111  10, // leafsize
112  3.0 // duplicity
113  )
114  );
115 }
116 
117 
118 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
121 {}
122 
123 
124 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
125 
127 {
128  if (regions_.empty())
129  {
130  regions_.setSize(1);
131  regions_[0] = "region0";
132  }
133  return regions_;
134 }
135 
137 Foam::label Foam::searchableExtrudedCircle::size() const
138 {
139  return eMeshPtr_().points().size();
140 }
141 
142 
144 {
145  return eMeshPtr_().points();
146 }
147 
148 
150 (
151  pointField& centres,
152  scalarField& radiusSqr
153 ) const
154 {
155  centres = eMeshPtr_().points();
156  radiusSqr.setSize(centres.size());
157  radiusSqr = Foam::sqr(radius_);
158  // Add a bit to make sure all points are tested inside
159  radiusSqr += Foam::sqr(SMALL);
160 }
161 
162 
164 (
165  const pointField& samples,
166  const scalarField& nearestDistSqr,
167  List<pointIndexHit>& info
168 ) const
169 {
170  const indexedOctree<treeDataEdge>& tree = edgeTree_();
171 
172  info.setSize(samples.size());
173 
174  forAll(samples, i)
175  {
176  const scalar nearestDist = Foam::sqrt(nearestDistSqr[i]);
177  const scalar searchDistSqr = Foam::sqr(nearestDist+radius_);
178 
179  // Find nearest on central edge
180  info[i] = tree.findNearest(samples[i], searchDistSqr);
181 
182  if (info[i].hit())
183  {
184  // Derive distance to nearest surface from distance to nearest edge
185  const vector d(samples[i] - info[i].point());
186  const scalar s(mag(d));
187 
188  if (s < ROOTVSMALL)
189  {
190  // Point is on edge. TBD.
191  info[i].setMiss();
192  }
193  else
194  {
195  const scalar distToSurface = radius_-s;
196  if (mag(distToSurface) > nearestDist)
197  {
198  info[i].setMiss();
199  }
200  else
201  {
202  info[i].setPoint(info[i].point() + d/s*radius_);
203  }
204  }
205  }
206  }
207 }
208 
209 
211 (
212  const point& start,
213  const point& end,
214  const scalarField& rawLambdas,
215  const scalarField& nearestDistSqr,
216  List<pointIndexHit>& info
217 ) const
218 {
219  const edgeMesh& mesh = eMeshPtr_();
220  const indexedOctree<treeDataEdge>& tree = edgeTree_();
221  const auto& treeData = tree.shapes();
222  const edgeList& edges = mesh.edges();
223  const pointField& points = mesh.points();
224  const labelListList& pointEdges = mesh.pointEdges();
225 
226  const scalar maxDistSqr = bounds().magSqr();
227 
228  // Normalise lambdas
229  const scalarField lambdas
230  (
231  (rawLambdas-rawLambdas[0])
232  /(rawLambdas.last()-rawLambdas[0])
233  );
234 
235 
236  // Nearest point on curve and local axis direction
237  pointField curvePoints(lambdas.size());
238  vectorField axialVecs(lambdas.size());
239 
240  const pointIndexHit startInfo = tree.findNearest(start, maxDistSqr);
241  curvePoints[0] = startInfo.hitPoint();
242  axialVecs[0] = treeData.line(startInfo.index()).vec();
243 
244  const pointIndexHit endInfo = tree.findNearest(end, maxDistSqr);
245  curvePoints.last() = endInfo.hitPoint();
246  axialVecs.last() = treeData.line(endInfo.index()).vec();
247 
248 
249 
250  scalarField curveLambdas(points.size(), -1.0);
251 
252  {
253  scalar endDistance = -1.0;
254 
255  // Determine edge lengths walking from start to end.
256 
257  const point& start = curvePoints[0];
258  const point& end = curvePoints.last();
259 
260  label edgei = startInfo.index();
261  const edge& startE = edges[edgei];
262 
263  label pointi = startE[0];
264  if ((startE.vec(points)&(end-start)) > 0)
265  {
266  pointi = startE[1];
267  }
268 
269  curveLambdas[pointi] = mag(points[pointi]-curvePoints[0]);
270  label otherPointi = startE.otherVertex(pointi);
271  curveLambdas[otherPointi] = -mag(points[otherPointi]-curvePoints[0]);
272 
273  //Pout<< "for point:" << points[pointi] << " have distance "
274  // << curveLambdas[pointi] << endl;
275 
276 
277  while (true)
278  {
279  const labelList& pEdges = pointEdges[pointi];
280  if (pEdges.size() == 1)
281  {
282  break;
283  }
284  else if (pEdges.size() != 2)
285  {
286  FatalErrorInFunction << "Curve " << name()
287  << " is not a single path. This is not supported"
288  << exit(FatalError);
289  break;
290  }
291 
292  label oldEdgei = edgei;
293  if (pEdges[0] == oldEdgei)
294  {
295  edgei = pEdges[1];
296  }
297  else
298  {
299  edgei = pEdges[0];
300  }
301 
302  if (edgei == endInfo.index())
303  {
304  endDistance = curveLambdas[pointi] + mag(end-points[pointi]);
305 
306  //Pout<< "Found end edge:" << edges[edgei].centre(points)
307  // << " endPt:" << end
308  // << " point before:" << points[pointi]
309  // << " accumulated length:" << endDistance << endl;
310  }
311 
312 
313  label oldPointi = pointi;
314  pointi = edges[edgei].otherVertex(oldPointi);
315 
316  if (curveLambdas[pointi] >= 0)
317  {
318  break;
319  }
320 
321  curveLambdas[pointi] =
322  curveLambdas[oldPointi] + edges[edgei].mag(points);
323  }
324 
325  // Normalise curveLambdas
326  forAll(curveLambdas, i)
327  {
328  if (curveLambdas[i] >= 0)
329  {
330  curveLambdas[i] /= endDistance;
331  }
332  }
333  }
334 
335 
336 
337  // Interpolation engine
338  linearInterpolationWeights interpolator(curveLambdas);
339 
340  // Find wanted location along curve
341  labelList indices;
342  scalarField weights;
343  for (label i = 1; i < curvePoints.size()-1; i++)
344  {
345  interpolator.valueWeights(lambdas[i], indices, weights);
346 
347  if (indices.size() == 1)
348  {
349  // On outside of curve. Choose one of the connected edges.
350  label pointi = indices[0];
351  const point& p0 = points[pointi];
352  label edge0 = pointEdges[pointi][0];
353  const edge& e0 = edges[edge0];
354  axialVecs[i] = e0.vec(points);
355  curvePoints[i] = weights[0]*p0;
356  }
357  else if (indices.size() == 2)
358  {
359  const point& p0 = points[indices[0]];
360  const point& p1 = points[indices[1]];
361  axialVecs[i] = p1-p0;
362  curvePoints[i] = weights[0]*p0+weights[1]*p1;
363  }
364  }
365  axialVecs /= mag(axialVecs);
366 
367 
368  // Now we have the lambdas, curvePoints and axialVecs.
369 
370 
371 
372  info.setSize(lambdas.size());
373  info = pointIndexHit();
374 
375  // Given the current lambdas interpolate radial direction inbetween
376  // endpoints (all projected onto the starting coordinate system)
377  quaternion qStart;
378  vector radialStart;
379  {
380  radialStart = start-curvePoints[0];
381  radialStart.removeCollinear(axialVecs[0]);
382  radialStart.normalise();
383 
384  qStart = quaternion(radialStart, 0.0);
385 
386  info[0] = pointIndexHit(true, start, 0);
387  }
388 
389  quaternion qProjectedEnd;
390  {
391  vector radialEnd(end-curvePoints.last());
392  radialEnd.removeCollinear(axialVecs.last());
393  radialEnd.normalise();
394 
395  vector projectedEnd = radialEnd;
396  projectedEnd.removeCollinear(axialVecs[0]);
397  projectedEnd.normalise();
398 
399  qProjectedEnd = quaternion(projectedEnd, 0.0);
400 
401  info.last() = pointIndexHit(true, end, 0);
402  }
403 
404  for (label i = 1; i < lambdas.size()-1; i++)
405  {
406  quaternion q(slerp(qStart, qProjectedEnd, lambdas[i]));
407  vector radialDir(q.transform(radialStart));
408 
409  radialDir.removeCollinear(axialVecs[i]);
410  radialDir.normalise();
412  info[i] = pointIndexHit(true, curvePoints[i]+radius_*radialDir, 0);
413  }
414 }
415 
416 
418 (
419  const List<pointIndexHit>& info,
420  labelList& region
421 ) const
422 {
423  region.setSize(info.size());
424  region = 0;
425 }
426 
427 
429 (
430  const List<pointIndexHit>& info,
431  vectorField& normal
432 ) const
433 {
434  const edgeMesh& mesh = eMeshPtr_();
435  const indexedOctree<treeDataEdge>& tree = edgeTree_();
436  const edgeList& edges = mesh.edges();
437  const pointField& points = mesh.points();
438 
439  normal.setSize(info.size());
440  normal = Zero;
441 
442  const scalar distSqr = bounds().magSqr();
443 
444  forAll(info, i)
445  {
446  if (info[i].hit())
447  {
448  // Find nearest on curve
449  pointIndexHit curvePt = tree.findNearest(info[i].point(), distSqr);
450 
451  normal[i] = info[i].hitPoint()-curvePt.hitPoint();
452 
453  // Subtract axial direction
454  const vector axialVec = edges[curvePt.index()].unitVec(points);
455 
456  normal[i].removeCollinear(axialVec);
457  normal[i].normalise();
458  }
459  }
460 }
461 
462 
463 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
dictionary dict
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
virtual ~searchableExtrudedCircle()
Destructor.
const pointField & points() const noexcept
Return points.
Definition: edgeMeshI.H:92
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
virtual label size() const
Range of local indices that can be returned.
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
List< edge > edgeList
List of edge.
Definition: edgeList.H:32
const labelListList & pointEdges() const
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const point_type & hitPoint() const
Return hit point. Fatal if not hit.
scalarField samples(nIntervals, Zero)
dimensionedScalar sqrt(const dimensionedScalar &ds)
PointIndexHit< point > pointIndexHit
A PointIndexHit with a 3D point.
Definition: pointIndexHit.H:58
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
Holds data for octree to work on an edges subset.
Definition: treeDataEdge.H:53
virtual void getRegion(const List< pointIndexHit > &, labelList &region) const
From a set of points and indices get the region.
A bounding box defined in terms of min/max extrema points.
Definition: boundBox.H:63
This class describes the interaction of an object (often a face) and a point. It carries the info of ...
Definition: pointIndexHit.H:44
Vector< Cmpt > & normalise(const scalar tol=ROOTVSMALL)
Inplace normalise the vector by its magnitude.
Definition: VectorI.H:114
Base class of (analytical or triangulated) surface. Encapsulates all the search routines. WIP.
virtual const boundBox & bounds() const
Return const reference to boundBox.
List< labelList > labelListList
List of labelList.
Definition: labelList.H:38
Macros for easy insertion into run-time selection tables.
addNamedToRunTimeSelectionTable(topoSetCellSource, badQualityToCell, word, badQuality)
virtual void findParametricNearest(const point &start, const point &end, const scalarField &lambdas, const scalarField &nearestDistSqr, List< pointIndexHit > &) const
Unique to parametric geometry: given points find.
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1078
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
virtual void getNormal(const List< pointIndexHit > &, vectorField &normal) const
From a set of points and indices get the normal.
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
void setSize(const label n)
Alias for resize()
Definition: List.H:316
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
const pointField & points
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
virtual void boundingSpheres(pointField &centres, scalarField &radiusSqr) const
Get bounding spheres (centre and radius squared), one per element.
Tree tree(triangles.begin(), triangles.end())
Vector< scalar > vector
Definition: vector.H:57
label index() const noexcept
Return the hit index.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
constexpr auto end(C &c) -> decltype(c.end())
Return iterator to the end of the container c.
Definition: stdFoam.H:201
virtual const wordList & regions() const
Names of regions.
defineTypeNameAndDebug(combustionModel, 0)
Mesh data needed to do the Finite Area discretisation.
Definition: edgeFaMesh.H:47
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
virtual tmp< pointField > coordinates() const
Get representative set of element coordinates.
vector point
Point is a vector.
Definition: point.H:37
Non-pointer based hierarchical recursive searching.
const edgeList & edges() const noexcept
Return edges.
Definition: edgeMeshI.H:98
quaternion slerp(const quaternion &qa, const quaternion &qb, const scalar t)
Spherical linear interpolation of quaternions.
Definition: quaternion.C:75
Standard boundBox with extra functionality for use in octree.
Definition: treeBoundBox.H:90
virtual void findNearest(const pointField &sample, const scalarField &nearestDistSqr, List< pointIndexHit > &) const
Field< vector > vectorField
Specialisation of Field<T> for vector.
void reset()
Reset to an inverted box.
Definition: boundBoxI.H:295
List< label > labelList
A List of labels.
Definition: List.H:62
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, IOobject::NO_REGISTER)
A class for managing temporary objects.
Definition: HashPtrTable.H:50
void normalise()
Inplace normalise this field. Generally a no-op except for vector fields.
Definition: Field.C:601
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))
virtual tmp< pointField > points() const
Get the points that define the surface.
Vector< Cmpt > & removeCollinear(const Vector< Cmpt > &unitVec)
Inplace removal of components that are collinear to the given unit vector.
Definition: VectorI.H:136
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
void grow(const scalar delta)
Expand box by adjusting min/max by specified amount in each dimension.
Definition: boundBoxI.H:367
const volScalarField & p0
Definition: EEqn.H:36
Namespace for OpenFOAM.
addToRunTimeSelectionTable(functionObject, pointHistory, dictionary)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127