BUG: UListIO: byteSize overflowing on really big faceLists
[OpenFOAM-2.0.x.git] / src / mesh / blockMesh / curvedEdges / BSpline.C
blobe101bf22e58571252e1660f6ae7950d7992e02d2
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
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.
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.
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/>.
24 \*---------------------------------------------------------------------------*/
26 #include "error.H"
27 #include "BSpline.H"
29 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
31 Foam::BSpline::BSpline
33     const pointField& knots,
34     const bool closed
37     polyLine(knots, closed)
41 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
43 Foam::point Foam::BSpline::position(const scalar mu) const
45     // endpoints
46     if (mu < SMALL)
47     {
48         return points().first();
49     }
50     else if (mu > 1 - SMALL)
51     {
52         return points().last();
53     }
55     scalar lambda = mu;
56     label segment = localParameter(lambda);
57     return position(segment, lambda);
61 Foam::point Foam::BSpline::position
63     const label segment,
64     const scalar mu
65 ) const
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     }
77     const point& p0 = points()[segment];
78     const point& p1 = points()[segment+1];
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     }
91     // determine the end points
92     point e0;
93     point e1;
95     if (segment == 0)
96     {
97         // end: simple reflection
98         e0 = 2*p0 - p1;
99     }
100     else
101     {
102         e0 = points()[segment-1];
103     }
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     }
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     );
133 Foam::scalar Foam::BSpline::length() const
135     notImplemented("BSpline::length() const");
136     return 1.0;
140 // ************************************************************************* //