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 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "polyLine.H"
29 
30 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
31 
33 {
34  param_.setSize(points_.size());
35 
36  if (param_.size())
37  {
38  param_[0] = 0.0;
39 
40  for (label i=1; i < param_.size(); i++)
41  {
42  param_[i] = param_[i-1] + mag(points_[i] - points_[i-1]);
43  }
44 
45  // Normalize on the interval 0-1
46  lineLength_ = param_.last();
47  for (label i=1; i < param_.size() - 1; i++)
48  {
49  param_[i] /= lineLength_;
50  }
51  param_.last() = 1.0;
52  }
53  else
54  {
55  lineLength_ = 0.0;
56  }
57 }
58 
59 
60 
61 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
62 
63 Foam::polyLine::polyLine(const pointField& ps, const bool)
64 :
65  points_(ps),
66  lineLength_(0.0),
67  param_(0)
68 {
69  calcParam();
70 }
71 
72 
73 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
74 
76 {
77  return points_;
78 }
79 
80 
82 {
83  return points_.size()-1;
84 }
85 
86 
88 {
89  // Check endpoints
90  if (lambda < SMALL)
91  {
92  lambda = 0;
93  return 0;
94  }
95  else if (lambda > 1 - SMALL)
96  {
97  lambda = 1;
98  return nSegments();
99  }
100 
101  // Search table of cumulative distances to find which line-segment
102  // we are on.
103  // Check the upper bound.
104 
105  label segmentI = 1;
106  while (param_[segmentI] < lambda)
107  {
108  segmentI++;
109  }
110  segmentI--; // We want the corresponding lower bound
111 
112  // The local parameter [0-1] on this line segment
113  lambda =
114  (lambda - param_[segmentI])/(param_[segmentI+1] - param_[segmentI]);
115 
116  return segmentI;
117 }
118 
119 
121 {
122  // Check end-points
123  if (mu < SMALL)
124  {
125  return points_.first();
126  }
127  else if (mu > 1 - SMALL)
128  {
129  return points_.last();
130  }
131 
132  scalar lambda = mu;
133  label segment = localParameter(lambda);
134  return position(segment, lambda);
135 }
136 
137 
139 (
140  const label segment,
141  const scalar mu
142 ) const
143 {
144  // Out-of-bounds
145  if (segment < 0)
146  {
147  return points_.first();
148  }
149  else if (segment > nSegments())
150  {
151  return points_.last();
152  }
153 
154  const point& p0 = points()[segment];
155  const point& p1 = points()[segment+1];
156 
157  // Special cases - no calculation needed
158  if (mu <= 0.0)
159  {
160  return p0;
161  }
162  else if (mu >= 1.0)
163  {
164  return p1;
165  }
166  else
167  {
168  // Linear interpolation
169  return points_[segment] + mu*(p1 - p0);
170  }
171 }
172 
173 
174 Foam::scalar Foam::polyLine::length() const
175 {
176  return lineLength_;
177 }
178 
179 
180 // ************************************************************************* //
Foam::polyLine::points
const pointField & points() const
Return const-access to the control-points.
Definition: polyLine.C:75
Foam::constant::physicoChemical::mu
const dimensionedScalar mu
Atomic mass unit.
Definition: createFieldRefs.H:4
polyLine.H
Foam::polyLine::nSegments
label nSegments() const
Return the number of line segments.
Definition: polyLine.C:81
Foam::label
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:62
Foam::Field< vector >
Foam::polyLine::param_
scalarList param_
The rational (0-1) cumulative parameter value for each point.
Definition: polyLine.H:77
Foam::polyLine::position
point position(const scalar) const
Return the point position corresponding to the curve parameter.
Definition: polyLine.C:120
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:71
lambda
dimensionedScalar lambda("lambda", dimTime/sqr(dimLength), laminarTransport)
Foam::Pair
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: Pair.H:54
Foam::Vector< scalar >
Foam::polyLine::lineLength_
scalar lineLength_
The real line length.
Definition: polyLine.H:74
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::polyLine::calcParam
void calcParam()
Precalculate the rational cumulative parameter value.
Definition: polyLine.C:32
p0
const volScalarField & p0
Definition: EEqn.H:36
Foam::polyLine::localParameter
label localParameter(scalar &lambda) const
Return the line segment and the local parameter [0..1].
Definition: polyLine.C:87
Foam::List::setSize
void setSize(const label newSize)
Alias for resize(const label)
Definition: ListI.H:146
Foam::polyLine::length
scalar length() const
Return the length of the curve.
Definition: polyLine.C:174