polyLine.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) 2020-2021 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 "polyLine.H"
30 #include "SubList.H"
31 
32 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
33 
35 (
36  const point& p0,
37  const pointField& intermediate,
38  const point& p1
39 )
40 {
41  auto tresult = tmp<pointField>::New(intermediate.size() + 2);
42  auto& result = tresult.ref();
43 
44  // Intermediate points (knots)
45  SubList<point>(result, intermediate.size(), 1) = intermediate;
46 
47  // Start/end points (knots)
48  result.first() = p0;
49  result.last() = p1;
50 
51  return tresult;
52 }
53 
54 
55 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
56 
58 {
59  lineLength_ = 0;
60  param_.resize(points_.size());
61 
62  if (param_.size())
63  {
64  param_[0] = 0;
65 
66  for (label i=1; i < param_.size(); i++)
67  {
68  param_[i] = param_[i-1] + mag(points_[i] - points_[i-1]);
69  }
70 
71  // Normalize on the interval 0-1
72  lineLength_ = param_.last();
73  for (label i=1; i < param_.size() - 1; i++)
74  {
75  param_[i] /= lineLength_;
76  }
77  param_.last() = 1;
78  }
79 }
80 
81 
82 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
83 
85 :
86  points_(p),
87  lineLength_(0),
88  param_()
89 {
90  calcParam();
91 }
92 
93 
95 (
96  const point& start,
97  const pointField& intermediate,
98  const point& end,
99  const bool
100 )
101 :
102  points_(polyLine::concat(start, intermediate, end)),
103  lineLength_(0),
104  param_()
105 {
106  calcParam();
107 }
108 
109 
110 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
111 
113 {
114  return points_;
115 }
116 
117 
118 Foam::label Foam::polyLine::nSegments() const noexcept
119 {
120  return points_.size()-1;
121 }
122 
123 
124 Foam::label Foam::polyLine::localParameter(scalar& lambda) const
125 {
126  // Check endpoints
127  if (lambda < SMALL)
128  {
129  lambda = 0;
130  return 0;
131  }
132  else if (lambda > 1 - SMALL)
133  {
134  lambda = 1;
135  return nSegments();
136  }
137 
138  // Search table of cumulative distances to find which line-segment
139  // we are on.
140  // Check the upper bound.
141  // Too small to bother with a binary search...
142 
143  label segment = 1;
144  while (param_[segment] < lambda)
145  {
146  ++segment;
147  }
148  --segment; // We want the corresponding lower bound
149 
150  // The local parameter [0-1] on this line segment
151  lambda =
152  (lambda - param_[segment])/(param_[segment+1] - param_[segment]);
153 
154  return segment;
155 }
156 
157 
159 {
160  // Check end-points
161  if (mu < SMALL)
162  {
163  return points_.first();
164  }
165  else if (mu > 1 - SMALL)
166  {
167  return points_.last();
168  }
169 
170  scalar lambda = mu;
171  label segment = localParameter(lambda);
172  return position(segment, lambda);
173 }
174 
175 
177 (
178  const label segment,
179  const scalar mu
180 ) const
181 {
182  // Out-of-bounds
183  if (segment < 0)
184  {
185  return points_.first();
186  }
187  else if (segment > nSegments())
188  {
189  return points_.last();
190  }
191 
192  const point& p0 = points_[segment];
193  const point& p1 = points_[segment+1];
194 
195  // Special cases - no calculation needed
196  if (mu <= 0.0)
197  {
198  return p0;
199  }
200  else if (mu >= 1.0)
201  {
202  return p1;
203  }
204 
205  // Linear interpolation
206  return points_[segment] + mu*(p1 - p0);
207 }
208 
209 
210 Foam::scalar Foam::polyLine::length() const noexcept
211 {
212  return lineLength_;
213 }
214 
215 
216 // ************************************************************************* //
p
volScalarField & p
Definition: createFieldRefs.H:8
SubList.H
Foam::List::resize
void resize(const label len)
Adjust allocated size of list.
Definition: ListI.H:139
Foam::constant::physicoChemical::mu
const dimensionedScalar mu
Atomic mass unit.
Definition: createFieldRefs.H:4
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::SubList
A List obtained as a section of another List.
Definition: SubList.H:54
polyLine.H
Foam::polyLine::polyLine
polyLine(const pointField &points, const bool notImplementedClosed=false)
Construct from components.
Definition: polyLine.C:84
Foam::polyLine::length
scalar length() const noexcept
The length of the curve.
Definition: polyLine.C:210
Foam::polyLine::points
const pointField & points() const noexcept
Return const-access to the control-points.
Definition: polyLine.C:112
Foam::Field< vector >
Foam::polyLine::concat
static tmp< pointField > concat(const point &start, const pointField &intermediate, const point &end)
Concatenate begin, intermediate and end points.
Definition: polyLine.C:35
Foam::polyLine::nSegments
label nSegments() const noexcept
The number of line segments.
Definition: polyLine.C:118
Foam::polyLine::param_
scalarList param_
The rational (0-1) cumulative parameter value for each point.
Definition: polyLine.H:68
Foam::polyLine::position
point position(const scalar) const
The point position corresponding to the curve parameter.
Definition: polyLine.C:158
Foam::segment
Pair< point > segment
Definition: distributedTriSurfaceMesh.H:75
Foam::polyLine::points_
pointField points_
The control points or ends of each segments.
Definition: polyLine.H:62
lambda
dimensionedScalar lambda("lambda", dimTime/sqr(dimLength), laminarTransport)
stdFoam::end
constexpr auto end(C &c) -> decltype(c.end())
Return iterator to the end of the container c.
Definition: stdFoam.H:121
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
Foam::Pair
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: Pair.H:54
Foam::Vector< scalar >
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::polyLine::lineLength_
scalar lineLength_
The real (total) line length.
Definition: polyLine.H:65
Foam::polyLine::calcParam
void calcParam()
Definition: polyLine.C:57
p0
const volScalarField & p0
Definition: EEqn.H:36
Foam::polyLine::localParameter
label localParameter(scalar &lambda) const
Definition: polyLine.C:124