BSpline.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-2021 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::BSpline
26 
27 Description
28  An implementation of B-splines.
29 
30  In this implementation, the end tangents are created automatically
31  by reflection.
32 
33  In matrix form, the \e local interpolation on the interval t=[0..1] is
34  described as follows:
35  \verbatim
36  P(t) = 1/6 * [ t^3 t^2 t 1 ] * [ -1 3 -3 1 ] * [ P-1 ]
37  [ 3 -6 3 0 ] [ P0 ]
38  [ -3 0 3 0 ] [ P1 ]
39  [ 1 4 1 0 ] [ P2 ]
40  \endverbatim
41 
42  Where P-1 and P2 represent the neighbouring points or the extrapolated
43  end points. Simple reflection is used to automatically create the end
44  points.
45 
46  The spline is discretised based on the chord length of the individual
47  segments. In rare cases (sections with very high curvatures), the
48  resulting distribution may be sub-optimal.
49 
50  A future implementation could also handle closed splines.
51 
52 See also
53  CatmullRomSpline
54 
55 SourceFiles
56  BSpline.C
57 
58 \*---------------------------------------------------------------------------*/
59 
60 #ifndef BSpline_H
61 #define BSpline_H
62 
63 #include "polyLine.H"
64 
65 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
66 
67 namespace Foam
68 {
69 
70 /*---------------------------------------------------------------------------*\
71  Class BSpline Declaration
72 \*---------------------------------------------------------------------------*/
73 
74 class BSpline
75 :
76  public polyLine
77 {
78 public:
79 
80  // Constructors
81 
82  //- Construct from components
83  BSpline
84  (
85  const pointField& knots,
86  const bool notImplementedClosed = false
87  );
88 
89  //- Disallow default bitwise copy construction
90  BSpline(const BSpline&) = delete;
91 
92 
93  // Member Functions
94 
95  //- Return the point position corresponding to the curve parameter
96  // 0 <= lambda <= 1
97  point position(const scalar lambda) const;
98 
99  //- Return the point position corresponding to the local parameter
100  // 0 <= lambda <= 1 on the given segment
101  point position(const label segment, const scalar lambda) const;
102 
103  //- Return the length of the curve
104  scalar length() const;
105 
106 
107  // Member Operators
108 
109  //- Disallow default bitwise assignment
110  void operator=(const BSpline&) = delete;
111 };
112 
113 
114 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
115 
116 } // End namespace Foam
117 
118 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
119 
120 #endif
121 
122 // ************************************************************************* //
BSpline(const pointField &knots, const bool notImplementedClosed=false)
Construct from components.
Definition: BSpline.C:31
dimensionedScalar lambda(viscosity->lookup("lambda"))
void operator=(const BSpline &)=delete
Disallow default bitwise assignment.
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: contiguous.H:49
An implementation of B-splines.
Definition: BSpline.H:73
A series of straight line segments, which can also be interpreted as a series of control points for s...
Definition: polyLine.H:53
point position(const scalar lambda) const
Return the point position corresponding to the curve parameter.
Definition: BSpline.C:42
Namespace for OpenFOAM.
scalar length() const
Return the length of the curve.
Definition: BSpline.C:132