polyLine.H
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-------------------------------------------------------------------------------
11License
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
27Class
28 Foam::polyLine
29
30Description
31 A series of straight line segments, which can also be interpreted as
32 a series of control points for splines, etc.
33
34 A future implementation could also handle a closed polyLine.
35
36SourceFiles
37 polyLine.C
38
39\*---------------------------------------------------------------------------*/
40
41#ifndef polyLine_H
42#define polyLine_H
43
44#include "pointField.H"
45#include "scalarList.H"
46
47// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48
49namespace Foam
50{
51
52/*---------------------------------------------------------------------------*\
53 Class polyLine Declaration
54\*---------------------------------------------------------------------------*/
56class polyLine
57{
58protected:
59
60 // Protected Data
61
62 //- The control points or ends of each segments
64
65 //- The real (total) line length
66 scalar lineLength_;
67
68 //- The rational (0-1) cumulative parameter value for each point
70
71
72 // Protected Member Functions
73
74 //- Precalculate the rational cumulative parameter value
75 //- and the line-length
76 void calcParam();
77
78 //- Return the line segment and the local parameter [0..1]
79 //- corresponding to the global lambda [0..1]
80 label localParameter(scalar& lambda) const;
81
82
83public:
84
85 // Constructors
86
87 //- Construct from components
88 explicit polyLine
89 (
90 const pointField& points,
91 const bool notImplementedClosed = false
92 );
93
94 //- Construct from begin, intermediate, end points
96 (
97 const point& start,
98 const pointField& intermediate,
99 const point& end,
100 const bool notImplementedClosed = false
101 );
102
103
104 // Static Member Functions
105
106 //- Concatenate begin, intermediate and end points
108 (
109 const point& start,
110 const pointField& intermediate,
111 const point& end
112 );
113
114
115 // Member Functions
116
117 //- Return const-access to the control-points
118 const pointField& points() const noexcept;
119
120 //- The number of line segments
121 label nSegments() const noexcept;
122
123 //- The point position corresponding to the curve parameter
124 // 0 <= lambda <= 1
125 point position(const scalar) const;
126
127 //- The point position corresponding to the local parameter
128 // 0 <= lambda <= 1 on the given segment
129 point position(const label segment, const scalar) const;
130
131 //- The length of the curve
132 scalar length() const noexcept;
133};
134
135
136// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
137
138} // End namespace Foam
139
140// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
141
142#endif
143
144// ************************************************************************* //
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: Pair.H:69
A series of straight line segments, which can also be interpreted as a series of control points for s...
Definition: polyLine.H:56
scalarList param_
The rational (0-1) cumulative parameter value for each point.
Definition: polyLine.H:68
const pointField & points() const noexcept
Return const-access to the control-points.
Definition: polyLine.C:112
label localParameter(scalar &lambda) const
Definition: polyLine.C:124
scalar length() const noexcept
The length of the curve.
Definition: polyLine.C:210
void calcParam()
Definition: polyLine.C:57
pointField points_
The control points or ends of each segments.
Definition: polyLine.H:62
static tmp< pointField > concat(const point &start, const pointField &intermediate, const point &end)
Concatenate begin, intermediate and end points.
Definition: polyLine.C:35
scalar lineLength_
The real (total) line length.
Definition: polyLine.H:65
label nSegments() const noexcept
The number of line segments.
Definition: polyLine.C:118
point position(const scalar) const
The point position corresponding to the curve parameter.
Definition: polyLine.C:158
A class for managing temporary objects.
Definition: tmp.H:65
Namespace for OpenFOAM.
const direction noexcept
Definition: Scalar.H:223
dimensionedScalar lambda("lambda", dimTime/sqr(dimLength), laminarTransport)