CatmullRomSpline.C 4.65 KB
Newer Older
1
2
3
4
/*---------------------------------------------------------------------------*\
  =========                 |
  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
   \\    /   O peration     |
OpenFOAM bot's avatar
OpenFOAM bot committed
5
    \\  /    A nd           | www.openfoam.com
6
     \\/     M anipulation  |
OpenFOAM bot's avatar
OpenFOAM bot committed
7
-------------------------------------------------------------------------------
OpenFOAM bot's avatar
OpenFOAM bot committed
8
    Copyright (C) 2011-2016 OpenFOAM Foundation
9
10
11
12
-------------------------------------------------------------------------------
License
    This file is part of OpenFOAM.

13
14
15
16
    OpenFOAM is free software: you can redistribute it and/or modify it
    under the terms of the GNU General Public License as published by
    the Free Software Foundation, either version 3 of the License, or
    (at your option) any later version.
17
18
19
20
21
22
23

    OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
    ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
    FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
    for more details.

    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
27

\*---------------------------------------------------------------------------*/

28
#include "CatmullRomSpline.H"
29

30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
// * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //

Foam::scalar Foam::CatmullRomSpline::derivative
(
    const label segment,
    const scalar mu
) const
{
    const point& p0 = points()[segment];
    const point& p1 = points()[segment+1];

    // determine the end points
    point e0;
    point e1;

    if (segment == 0)
    {
        // end: simple reflection
        e0 = 2*p0 - p1;
    }
    else
    {
        e0 = points()[segment-1];
    }

    if (segment+1 == nSegments())
    {
        // end: simple reflection
        e1 = 2*p1 - p0;
    }
    else
    {
        e1 = points()[segment+2];
    }
    const point derivativePoint
    (
        0.5 *
        (
            (-e0 + p1)
          + mu *
            (
                2 * (2*e0 - 5*p0 + 4*p1 - e1)
              + mu * 3 * (-e0 + 3*p0 - 3*p1 + e1)
            )
        )
    );
    return mag(derivativePoint);
}


80
81
82
83
// * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //

Foam::CatmullRomSpline::CatmullRomSpline
(
84
85
    const pointField& knots,
    const bool closed
86
87
)
:
88
    polyLine(knots, closed)
89
{}
90

91
92
93
94

// * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //

Foam::point Foam::CatmullRomSpline::position(const scalar mu) const
95
{
96
97
    // endpoints
    if (mu < SMALL)
98
    {
99
        return points().first();
100
    }
101
    else if (mu > 1 - SMALL)
102
    {
103
        return points().last();
104
    }
105
106
107
108
109
110
111
112
113
114
115
116
117

    scalar lambda = mu;
    label segment = localParameter(lambda);
    return position(segment, lambda);
}


Foam::point Foam::CatmullRomSpline::position
(
    const label segment,
    const scalar mu
) const
{
118
119
120
121
122
123
124
125
126
127
    // out-of-bounds
    if (segment < 0)
    {
        return points().first();
    }
    else if (segment > nSegments())
    {
        return points().last();
    }

128
129
130
131
    const point& p0 = points()[segment];
    const point& p1 = points()[segment+1];

    // special cases - no calculation needed
132
    if (mu <= 0.0)
133
    {
134
        return p0;
135
    }
136
    else if (mu >= 1.0)
137
    {
138
        return p1;
139
    }
140

141

142
143
144
145
146
    // determine the end points
    point e0;
    point e1;

    if (segment == 0)
147
    {
148
        // end: simple reflection
149
        e0 = 2*p0 - p1;
150
151
152
    }
    else
    {
153
        e0 = points()[segment-1];
154
155
    }

156
    if (segment+1 == nSegments())
157
    {
158
        // end: simple reflection
159
        e1 = 2*p1 - p0;
160
    }
161
162
163
164
165
    else
    {
        e1 = points()[segment+2];
    }

166

167
168
    return
        0.5 *
169
        (
170
            (2*p0)
171
172
          + mu *
            (
173
                (-e0 + p1)
174
              + mu *
175
176
177
178
                (
                    (2*e0 - 5*p0 + 4*p1 - e1)
                  + mu*(-e0 + 3*p0 - 3*p1 + e1)
                )
179
            )
180
        );
181
182
183
}


184
Foam::scalar Foam::CatmullRomSpline::length() const
185
{
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
    const solveScalar xi[5]=
    {
        -0.9061798459386639927976,
        -0.5384693101056830910363,
        0,
        0.5384693101056830910363,
        0.9061798459386639927976
    };
    const solveScalar wi[5]=
    {
        0.2369268850561890875143,
        0.4786286704993664680413,
        0.5688888888888888888889,
        0.4786286704993664680413,
        0.2369268850561890875143
    };
    scalar sum=0;
    for (label segment=0;segment<nSegments();segment++)
    {
        for (int i=0;i<5;i++)
        {
            sum+=wi[i]*derivative(segment,(xi[i]+1.0)/2.0)/2.0;
        }
    }

    return sum;
212
213
214
215
}


// ************************************************************************* //