All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
polyLine.H
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-2019 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 Class
25  Foam::polyLine
26 
27 Description
28  A series of straight line segments, which can also be interpreted as
29  a series of control points for splines, etc.
30 
31  A future implementation could also handle a closed polyLine.
32 
33 SourceFiles
34  polyLine.C
35 
36 \*---------------------------------------------------------------------------*/
37 
38 #ifndef polyLine_H
39 #define polyLine_H
40 
41 #include "pointField.H"
42 #include "scalarList.H"
43 
44 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45 
46 namespace Foam
47 {
48 
49 /*---------------------------------------------------------------------------*\
50  Class polyLine Declaration
51 \*---------------------------------------------------------------------------*/
52 
53 
54 class polyLine
55 {
56 protected:
57 
58  // Protected data
59 
60  //- The control points or ends of each segments
62 
63  //- The real line length
64  scalar lineLength_;
65 
66  //- The rational (0-1) cumulative parameter value for each point
68 
69 
70  // Protected Member Functions
71 
72  //- Precalculate the rational cumulative parameter value
73  // and the line-length
74  void calcParam();
75 
76  //- Return the line segment and the local parameter [0..1]
77  // corresponding to the global lambda [0..1]
78  label localParameter(scalar& lambda) const;
79 
80 
81 public:
82 
83  // Constructors
84 
85  //- Construct from components
86  polyLine
87  (
88  const pointField&,
89  const bool notImplementedClosed = false
90  );
91 
92  //- Disallow default bitwise copy construction
93  polyLine(const polyLine&) = delete;
94 
95 
96  // Member Functions
97 
98  //- Return const-access to the control-points
99  const pointField& points() const;
100 
101  //- Return the number of line segments
102  label nSegments() const;
103 
104  //- Return the point position corresponding to the curve parameter
105  // 0 <= lambda <= 1
106  point position(const scalar) const;
107 
108  //- Return the point position corresponding to the local parameter
109  // 0 <= lambda <= 1 on the given segment
110  point position(const label segment, const scalar) const;
111 
112  //- Return the length of the curve
113  scalar length() const;
114 
115 
116  // Member Operators
117 
118  //- Disallow default bitwise assignment
119  void operator=(const polyLine&) = delete;
120 };
121 
122 
123 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
124 
125 } // End namespace Foam
126 
127 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
128 
129 #endif
130 
131 // ************************************************************************* //
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
dimensionedScalar lambda(laminarTransport.lookup("lambda"))
pointField points_
The control points or ends of each segments.
Definition: polyLine.H:60
void operator=(const polyLine &)=delete
Disallow default bitwise assignment.
scalar length() const
Return the length of the curve.
Definition: polyLine.C:172
point position(const scalar) const
Return the point position corresponding to the curve parameter.
Definition: polyLine.C:118
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: contiguous.H:49
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
A series of straight line segments, which can also be interpreted as a series of control points for s...
Definition: polyLine.H:53
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
label nSegments() const
Return the number of line segments.
Definition: polyLine.C:79
Namespace for OpenFOAM.