extrudedCircle_searchableSurface.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) 2016-2025 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 
27 #include "Time.H"
28 #include "edgeMesh.H"
29 #include "indexedOctree.H"
30 #include "treeDataEdge.H"
32 #include "quaternion.H"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39  namespace searchableSurfaces
40  {
42 
44  (
48  );
49 
51  (
54  dictionary,
55  searchableExtrudedCircle,
56  "searchableExtrudedCircle"
57  );
58  }
59 }
60 
61 
62 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
63 
65 (
66  const IOobject& io,
67  const dictionary& dict
68 )
69 :
71  eMeshPtr_
72  (
73  edgeMesh::New
74  (
75  IOobject
76  (
77  dict.lookup("file"),
78  io.time().constant(),
79  searchableSurface::geometryDir(),
80  io.time(),
81  IOobject::MUST_READ,
82  IOobject::NO_WRITE,
83  false
84  ).objectPath(true)
85  )
86  ),
87  radius_(dict.lookup<scalar>("radius", dimLength))
88 {
89  const edgeMesh& eMesh = eMeshPtr_();
90 
91  const pointField& points = eMesh.points();
92  const edgeList& edges = eMesh.edges();
93  bounds() = boundBox(points, false);
94 
95  vector halfSpan(0.5*bounds().span());
96  point ctr(bounds().midpoint());
97 
98  bounds().min() = ctr - mag(halfSpan)*vector(1, 1, 1);
99  bounds().max() = ctr + mag(halfSpan)*vector(1, 1, 1);
100 
101  // Calculate bb of all points
102  treeBoundBox bb(bounds());
103 
104  // Slightly extended bb. Slightly off-centred just so on symmetric
105  // geometry there are less face/edge aligned items.
106  bb.min() -= point(rootVSmall, rootVSmall, rootVSmall);
107  bb.max() += point(rootVSmall, rootVSmall, rootVSmall);
108 
109  edgeTree_.reset
110  (
112  (
114  (
115  false, // do not cache bb
116  edges,
117  points,
118  identityMap(edges.size())
119  ),
120  bb, // overall search domain
121  8, // maxLevel
122  10, // leafsize
123  3.0 // duplicity
124  )
125  );
126 }
127 
128 
129 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
130 
132 {}
133 
134 
135 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
136 
138 {
139  if (regions_.empty())
140  {
141  regions_.setSize(1);
142  regions_[0] = "region0";
143  }
144  return regions_;
145 }
146 
147 
149 {
150  return eMeshPtr_().points().size();
151 }
152 
153 
156 {
157  return eMeshPtr_().points();
158 }
159 
160 
162 (
163  pointField& centres,
164  scalarField& radiusSqr
165 ) const
166 {
167  centres = eMeshPtr_().points();
168  radiusSqr.setSize(centres.size());
169  radiusSqr = Foam::sqr(radius_);
170  // Add a bit to make sure all points are tested inside
171  radiusSqr += Foam::sqr(small);
172 }
173 
174 
176 (
177  const pointField& samples,
178  const scalarField& nearestDistSqr,
179  List<pointIndexHit>& info
180 ) const
181 {
182  const indexedOctree<treeDataEdge>& tree = edgeTree_();
183 
184  info.setSize(samples.size());
185 
186  forAll(samples, i)
187  {
188  info[i] = tree.findNearest(samples[i], nearestDistSqr[i]);
189 
190  if (info[i].hit())
191  {
192  vector d(samples[i]-info[i].hitPoint());
193  info[i].setPoint(info[i].hitPoint() + d/mag(d)*radius_);
194  }
195  }
196 }
197 
198 
200 (
201  const point& start,
202  const point& end,
203  const scalarField& rawLambdas,
204  const scalarField& nearestDistSqr,
205  List<pointIndexHit>& info
206 ) const
207 {
208  const edgeMesh& mesh = eMeshPtr_();
209  const indexedOctree<treeDataEdge>& tree = edgeTree_();
210  const edgeList& edges = mesh.edges();
211  const pointField& points = mesh.points();
212  const labelListList& pointEdges = mesh.pointEdges();
213 
214  const scalar maxDistSqr(Foam::magSqr(bounds().span()));
215 
216  // Normalise lambdas
217  const scalarField lambdas
218  (
219  (rawLambdas-rawLambdas[0])
220  /(rawLambdas.last()-rawLambdas[0])
221  );
222 
223 
224  // Nearest point on curve and local axis direction
225  pointField curvePoints(lambdas.size());
226  vectorField axialVecs(lambdas.size());
227 
228  const pointIndexHit startInfo = tree.findNearest(start, maxDistSqr);
229  curvePoints[0] = startInfo.hitPoint();
230  axialVecs[0] = edges[startInfo.index()].vec(points);
231 
232  const pointIndexHit endInfo = tree.findNearest(end, maxDistSqr);
233  curvePoints.last() = endInfo.hitPoint();
234  axialVecs.last() = edges[endInfo.index()].vec(points);
235 
236 
237 
238  scalarField curveLambdas(points.size(), -1.0);
239 
240  {
241  scalar endDistance = -1.0;
242 
243  // Determine edge lengths walking from start to end.
244 
245  const point& start = curvePoints[0];
246  const point& end = curvePoints.last();
247 
248  label edgei = startInfo.index();
249  const edge& startE = edges[edgei];
250 
251  label pointi = startE[0];
252  if ((startE.vec(points)&(end-start)) > 0)
253  {
254  pointi = startE[1];
255  }
256 
257  curveLambdas[pointi] = mag(points[pointi]-curvePoints[0]);
258  label otherPointi = startE.otherVertex(pointi);
259  curveLambdas[otherPointi] = -mag(points[otherPointi]-curvePoints[0]);
260 
261  // Pout<< "for point:" << points[pointi] << " have distance "
262  // << curveLambdas[pointi] << endl;
263 
264 
265  while (true)
266  {
267  const labelList& pEdges = pointEdges[pointi];
268  if (pEdges.size() == 1)
269  {
270  break;
271  }
272  else if (pEdges.size() != 2)
273  {
274  FatalErrorInFunction << "Curve " << name()
275  << " is not a single path. This is not supported"
276  << exit(FatalError);
277  break;
278  }
279 
280  label oldEdgei = edgei;
281  if (pEdges[0] == oldEdgei)
282  {
283  edgei = pEdges[1];
284  }
285  else
286  {
287  edgei = pEdges[0];
288  }
289 
290  if (edgei == endInfo.index())
291  {
292  endDistance = curveLambdas[pointi] + mag(end-points[pointi]);
293 
294  // Pout<< "Found end edge:" << edges[edgei].centre(points)
295  // << " endPt:" << end
296  // << " point before:" << points[pointi]
297  // << " accumulated length:" << endDistance << endl;
298  }
299 
300 
301  label oldPointi = pointi;
302  pointi = edges[edgei].otherVertex(oldPointi);
303 
304  if (curveLambdas[pointi] >= 0)
305  {
306  break;
307  }
308 
309  curveLambdas[pointi] =
310  curveLambdas[oldPointi] + edges[edgei].mag(points);
311  }
312 
313  // Normalise curveLambdas
314  forAll(curveLambdas, i)
315  {
316  if (curveLambdas[i] >= 0)
317  {
318  curveLambdas[i] /= endDistance;
319  }
320  }
321  }
322 
323 
324 
325  // Interpolation engine
326  linearInterpolationWeights interpolator(curveLambdas);
327 
328  // Find wanted location along curve
329  labelList indices;
330  scalarField weights;
331  for (label i = 1; i < curvePoints.size()-1; i++)
332  {
333  interpolator.valueWeights(lambdas[i], indices, weights);
334 
335  if (indices.size() == 1)
336  {
337  // On outside of curve. Choose one of the connected edges.
338  label pointi = indices[0];
339  const point& p0 = points[pointi];
340  label edge0 = pointEdges[pointi][0];
341  const edge& e0 = edges[edge0];
342  axialVecs[i] = e0.vec(points);
343  curvePoints[i] = weights[0]*p0;
344  }
345  else if (indices.size() == 2)
346  {
347  const point& p0 = points[indices[0]];
348  const point& p1 = points[indices[1]];
349  axialVecs[i] = p1-p0;
350  curvePoints[i] = weights[0]*p0+weights[1]*p1;
351  }
352  }
353  axialVecs /= mag(axialVecs);
354 
355 
356  // Now we have the lambdas, curvePoints and axialVecs.
357 
358 
359 
360  info.setSize(lambdas.size());
361  info = pointIndexHit();
362 
363  // Given the current lambdas interpolate radial direction in between
364  // endpoints (all projected onto the starting coordinate system)
365  quaternion qStart;
366  vector radialStart;
367  {
368  radialStart = start-curvePoints[0];
369  radialStart -= (radialStart&axialVecs[0])*axialVecs[0];
370  radialStart /= mag(radialStart);
371  qStart = quaternion(radialStart, 0.0);
372 
373  info[0] = pointIndexHit(true, start, 0);
374  }
375 
376  quaternion qProjectedEnd;
377  {
378  vector radialEnd(end-curvePoints.last());
379  radialEnd -= (radialEnd&axialVecs.last())*axialVecs.last();
380  radialEnd /= mag(radialEnd);
381 
382  vector projectedEnd = radialEnd;
383  projectedEnd -= (projectedEnd&axialVecs[0])*axialVecs[0];
384  projectedEnd /= mag(projectedEnd);
385  qProjectedEnd = quaternion(projectedEnd, 0.0);
386 
387  info.last() = pointIndexHit(true, end, 0);
388  }
389 
390  for (label i = 1; i < lambdas.size()-1; i++)
391  {
392  quaternion q(slerp(qStart, qProjectedEnd, lambdas[i]));
393  vector radialDir(q.transform(radialStart));
394 
395  radialDir -= (radialDir&axialVecs[i])*axialVecs.last();
396  radialDir /= mag(radialDir);
397 
398  info[i] = pointIndexHit(true, curvePoints[i]+radius_*radialDir, 0);
399  }
400 }
401 
402 
404 (
405  const List<pointIndexHit>& info,
406  labelList& region
407 ) const
408 {
409  region.setSize(info.size());
410  region = 0;
411 }
412 
413 
415 (
416  const List<pointIndexHit>& info,
417  vectorField& normal
418 ) const
419 {
420  const edgeMesh& mesh = eMeshPtr_();
421  const indexedOctree<treeDataEdge>& tree = edgeTree_();
422  const edgeList& edges = mesh.edges();
423  const pointField& points = mesh.points();
424 
425  normal.setSize(info.size());
426  normal = Zero;
427 
428  forAll(info, i)
429  {
430  if (info[i].hit())
431  {
432  // Find nearest on curve
433  pointIndexHit curvePt = tree.findNearest
434  (
435  info[i].hitPoint(),
436  Foam::magSqr(bounds().span())
437  );
438 
439  normal[i] = info[i].hitPoint()-curvePt.hitPoint();
440 
441  // Subtract axial direction
442  vector axialVec = edges[curvePt.index()].vec(points);
443  axialVec /= mag(axialVec);
444  normal[i] -= (normal[i]&axialVec)*axialVec;
445  normal[i] /= mag(normal[i]);
446  }
447  }
448 }
449 
450 
451 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
Macros for easy insertion into run-time selection tables.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
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
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
Definition: PointIndexHit.H:54
const Point & hitPoint() const
Return hit point.
label index() const
Return index.
T & last()
Return the last element of the list.
Definition: UListI.H:128
A bounding box defined in terms of the points at its extremities.
Definition: boundBox.H:59
const point & min() const
Minimum point defining the bounding box.
Definition: boundBoxI.H:60
const point & max() const
Maximum point defining the bounding box.
Definition: boundBoxI.H:66
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
Points connected by edges.
Definition: edgeMesh.H:72
const edgeList & edges() const
Return edges.
Definition: edgeMeshI.H:68
const pointField & points() const
Return points.
Definition: edgeMeshI.H:62
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:61
vector vec(const pointField &) const
Return the vector (end - start)
Definition: edgeI.H:175
label otherVertex(const label a) const
Given one vertex, return the other.
Definition: edgeI.H:140
Non-pointer based hierarchical recursive searching.
Definition: indexedOctree.H:72
virtual bool valueWeights(const scalar t, labelList &indices, scalarField &weights) const
Calculate weights and indices to calculate t from samples.
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1331
const labelListList & pointEdges() const
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
Quaternion class used to perform rotations in 3D space.
Definition: quaternion.H:61
vector transform(const vector &v) const
Rotate the given vector.
Definition: quaternionI.H:296
Base class of (analytical or triangulated) surface. Encapsulates all the search routines....
const boundBox & bounds() const
Return const reference to boundBox.
Surface geometry with a tube shape, which can be used with snappyHexMesh. The geometry is formed from...
virtual label size() const
Range of local indices that can be returned.
virtual void boundingSpheres(pointField &centres, scalarField &radiusSqr) const
Get bounding spheres (centre and radius squared), one per element.
virtual tmp< pointField > points() const
Get the points that define the surface.
virtual void getNormal(const List< pointIndexHit > &, vectorField &normal) const
From a set of points and indices get the normal.
virtual void getRegion(const List< pointIndexHit > &, labelList &region) const
From a set of points and indices get the region.
virtual const wordList & regions() const
Names of regions.
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 void findNearest(const pointField &sample, const scalarField &nearestDistSqr, List< pointIndexHit > &) const
extrudedCircle(const IOobject &io, const dictionary &dict)
Construct from dictionary (used by searchableSurface)
virtual tmp< pointField > coordinates() const
Get representative set of element coordinates.
A class for managing temporary objects.
Definition: tmp.H:55
Standard boundBox + extra functionality for use in octree.
Definition: treeBoundBox.H:90
Holds data for octree to work on an edges subset.
Definition: treeDataEdge.H:54
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
const pointField & points
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
addToRunTimeSelectionTable(searchableSurface, box, dictionary)
addBackwardCompatibleToRunTimeSelectionTable(searchableSurface, box, dictionary, searchableBox, "searchableBox")
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static const zero Zero
Definition: zero.H:97
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
PointIndexHit< point > pointIndexHit
Definition: pointIndexHit.H:42
const dimensionSet dimLength
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
void sqr(LagrangianPatchField< typename outerProduct< Type, Type >::type > &f, const LagrangianPatchField< Type > &f1)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
error FatalError
void magSqr(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
labelList identityMap(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
quaternion slerp(const quaternion &qa, const quaternion &qb, const scalar t)
Spherical linear interpolation of quaternions.
Definition: quaternion.C:55
dictionary dict
scalarField samples(nIntervals, 0)