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 | Copyright (C) 2011-2015 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 \*---------------------------------------------------------------------------*/
25 
26 #include "error.H"
27 #include "BSpline.H"
28 
29 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
30 
31 Foam::BSpline::BSpline
32 (
33  const pointField& knots,
34  const bool closed
35 )
36 :
37  polyLine(knots, closed)
38 {}
39 
40 
41 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
42 
43 Foam::point Foam::BSpline::position(const scalar mu) const
44 {
45  // endpoints
46  if (mu < SMALL)
47  {
48  return points().first();
49  }
50  else if (mu > 1 - SMALL)
51  {
52  return points().last();
53  }
54 
55  scalar lambda = mu;
56  label segment = localParameter(lambda);
57  return position(segment, lambda);
58 }
59 
60 
62 (
63  const label segment,
64  const scalar mu
65 ) const
66 {
67  // out-of-bounds
68  if (segment < 0)
69  {
70  return points().first();
71  }
72  else if (segment > nSegments())
73  {
74  return points().last();
75  }
76 
77  const point& p0 = points()[segment];
78  const point& p1 = points()[segment+1];
79 
80  // special cases - no calculation needed
81  if (mu <= 0.0)
82  {
83  return p0;
84  }
85  else if (mu >= 1.0)
86  {
87  return p1;
88  }
89 
90 
91  // determine the end points
92  point e0;
93  point e1;
94 
95  if (segment == 0)
96  {
97  // end: simple reflection
98  e0 = 2*p0 - p1;
99  }
100  else
101  {
102  e0 = points()[segment-1];
103  }
104 
105  if (segment+1 == nSegments())
106  {
107  // end: simple reflection
108  e1 = 2*p1 - p0;
109  }
110  else
111  {
112  e1 = points()[segment+2];
113  }
114 
115 
116  return 1.0/6.0 *
117  (
118  ( e0 + 4*p0 + p1 )
119  + mu *
120  (
121  ( -3*e0 + 3*p1 )
122  + mu *
123  (
124  ( 3*e0 - 6*p0 + 3*p1 )
125  + mu *
126  ( -e0 + 3*p0 - 3*p1 + e1 )
127  )
128  )
129  );
130 }
131 
132 
133 Foam::scalar Foam::BSpline::length() const
134 {
136  return 1.0;
137 }
138 
139 
140 // ************************************************************************* //
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
label nSegments() const
Return the number of line segments.
Definition: polyLine.C:79
const pointField & points() const
Return const-access to the control-points.
Definition: polyLine.C:73
scalar length() const
Return the length of the curve.
Definition: BSpline.C:133
Pair< point > segment
T & first()
Return the first element of the list.
Definition: UListI.H:114
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:43
const dimensionedScalar mu
Atomic mass unit.
label localParameter(scalar &lambda) const
Return the line segment and the local parameter [0..1].
Definition: polyLine.C:85
dimensionedScalar lambda(laminarTransport.lookup("lambda"))
T & last()
Return the last element of the list.
Definition: UListI.H:128
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:366