BSpline.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 Copyright (C) 2020 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
27\*---------------------------------------------------------------------------*/
28
29#include "BSpline.H"
30
31// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32
34(
35 const pointField& knots,
36 const bool closed
37)
38:
39 polyLine(knots, closed)
40{}
41
42
43// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
44
46{
47 // endpoints
48 if (mu < SMALL)
49 {
50 return points().first();
51 }
52 else if (mu > 1 - SMALL)
53 {
54 return points().last();
55 }
56
57 scalar lambda = mu;
58 label segment = localParameter(lambda);
59 return position(segment, lambda);
60}
61
62
64(
65 const label segment,
66 const scalar mu
67) const
68{
69 // out-of-bounds
70 if (segment < 0)
71 {
72 return points().first();
73 }
74 else if (segment > nSegments())
75 {
76 return points().last();
77 }
78
79 const point& p0 = points()[segment];
80 const point& p1 = points()[segment+1];
81
82 // Special cases - no calculation needed
83 // Note: only checks at overall endpoints
84 if ((segment == 0) && (mu <= 0.0))
85 {
86 return p0;
87 }
88 else if ((segment == nSegments() - 1) && (mu >= 1.0))
89 {
90 return p1;
91 }
92
93 // Determine the end points
94 point e0;
95 point e1;
96
97 if (segment == 0)
98 {
99 // end: simple reflection
100 e0 = 2*p0 - p1;
101 }
102 else
103 {
104 e0 = points()[segment-1];
105 }
106
107 if (segment+1 == nSegments())
108 {
109 // end: simple reflection
110 e1 = 2*p1 - p0;
111 }
112 else
113 {
114 e1 = points()[segment+2];
115 }
116
117 return 1.0/6.0 *
118 (
119 ( e0 + 4*p0 + p1 )
120 + mu *
121 (
122 ( -3*e0 + 3*p1 )
123 + mu *
124 (
125 ( 3*e0 - 6*p0 + 3*p1 )
126 + mu *
127 ( -e0 + 3*p0 - 3*p1 + e1 )
128 )
129 )
130 );
131}
132
133
134Foam::scalar Foam::BSpline::length() const
135{
137 return 1;
138}
139
140
141// ************************************************************************* //
An implementation of B-splines.
Definition: BSpline.H:79
scalar length() const
The length of the curve.
Definition: BSpline.C:134
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: Pair.H:69
T & first()
Return the first element of the list.
Definition: UListI.H:202
T & last()
Return the last element of the list.
Definition: UListI.H:216
A series of straight line segments, which can also be interpreted as a series of control points for s...
Definition: polyLine.H:56
const volScalarField & mu
const volScalarField & p0
Definition: EEqn.H:36
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:517
const pointField & points
Pair< point > segment
dimensionedScalar lambda("lambda", dimTime/sqr(dimLength), laminarTransport)