All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
polyLine.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-2018 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 "polyLine.H"
27 
28 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
29 
31 {
33 
34  if (param_.size())
35  {
36  param_[0] = 0.0;
37 
38  for (label i=1; i < param_.size(); i++)
39  {
40  param_[i] = param_[i-1] + mag(points_[i] - points_[i-1]);
41  }
42 
43  // Normalize on the interval 0-1
45  for (label i=1; i < param_.size() - 1; i++)
46  {
47  param_[i] /= lineLength_;
48  }
49  param_.last() = 1.0;
50  }
51  else
52  {
53  lineLength_ = 0.0;
54  }
55 }
56 
57 
58 
59 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
60 
61 Foam::polyLine::polyLine(const pointField& ps, const bool)
62 :
63  points_(ps),
64  lineLength_(0.0),
65  param_(0)
66 {
67  calcParam();
68 }
69 
70 
71 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
72 
74 {
75  return points_;
76 }
77 
78 
80 {
81  return points_.size()-1;
82 }
83 
84 
86 {
87  // Check endpoints
88  if (lambda < small)
89  {
90  lambda = 0;
91  return 0;
92  }
93  else if (lambda > 1 - small)
94  {
95  lambda = 1;
96  return nSegments();
97  }
98 
99  // Search table of cumulative distances to find which line-segment
100  // we are on.
101  // Check the upper bound.
102 
103  label segmentI = 1;
104  while (param_[segmentI] < lambda)
105  {
106  segmentI++;
107  }
108  segmentI--; // We want the corresponding lower bound
109 
110  // The local parameter [0-1] on this line segment
111  lambda =
112  (lambda - param_[segmentI])/(param_[segmentI+1] - param_[segmentI]);
113 
114  return segmentI;
115 }
116 
117 
118 Foam::point Foam::polyLine::position(const scalar mu) const
119 {
120  // Check end-points
121  if (mu < small)
122  {
123  return points_.first();
124  }
125  else if (mu > 1 - small)
126  {
127  return points_.last();
128  }
129 
130  scalar lambda = mu;
131  label segment = localParameter(lambda);
132  return position(segment, lambda);
133 }
134 
135 
137 (
138  const label segment,
139  const scalar mu
140 ) const
141 {
142  // Out-of-bounds
143  if (segment < 0)
144  {
145  return points_.first();
146  }
147  else if (segment > nSegments())
148  {
149  return points_.last();
150  }
151 
152  const point& p0 = points()[segment];
153  const point& p1 = points()[segment+1];
154 
155  // Special cases - no calculation needed
156  if (mu <= 0.0)
157  {
158  return p0;
159  }
160  else if (mu >= 1.0)
161  {
162  return p1;
163  }
164  else
165  {
166  // Linear interpolation
167  return points_[segment] + mu*(p1 - p0);
168  }
169 }
170 
171 
172 Foam::scalar Foam::polyLine::length() const
173 {
174  return lineLength_;
175 }
176 
177 
178 // ************************************************************************* //
scalar lineLength_
The real line length.
Definition: polyLine.H:63
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
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
Pair< point > segment
pointField points_
The control points or ends of each segments.
Definition: polyLine.H:60
scalar length() const
Return the length of the curve.
Definition: polyLine.C:172
T & first()
Return the first element of the list.
Definition: UListI.H:114
point position(const scalar) const
Return the point position corresponding to the curve parameter.
Definition: polyLine.C:118
void calcParam()
Precalculate the rational cumulative parameter value.
Definition: polyLine.C:30
const pointField & points() const
Return const-access to the control-points.
Definition: polyLine.C:73
polyLine(const pointField &, const bool notImplementedClosed=false)
Construct from components.
Definition: polyLine.C:61
scalarList param_
The rational (0-1) cumulative parameter value for each point.
Definition: polyLine.H:66
label localParameter(scalar &lambda) const
Return the line segment and the local parameter [0..1].
Definition: polyLine.C:85
void setSize(const label)
Reset size of List.
Definition: List.C:281
const dimensionedScalar & mu
Atomic mass unit.
dimensionedScalar lambda(laminarTransport.lookup("lambda"))
dimensioned< scalar > mag(const dimensioned< Type > &)
T & last()
Return the last element of the list.
Definition: UListI.H:128
label nSegments() const
Return the number of line segments.
Definition: polyLine.C:79