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