00001 /*---------------------------------------------------------------------------*\ 00002 ========= | 00003 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox 00004 \\ / O peration | 00005 \\ / A nd | Copyright (C) 2009-2010 OpenCFD Ltd. 00006 \\/ M anipulation | 00007 ------------------------------------------------------------------------------- 00008 License 00009 This file is part of OpenFOAM. 00010 00011 OpenFOAM is free software: you can redistribute it and/or modify it 00012 under the terms of the GNU General Public License as published by 00013 the Free Software Foundation, either version 3 of the License, or 00014 (at your option) any later version. 00015 00016 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT 00017 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or 00018 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License 00019 for more details. 00020 00021 You should have received a copy of the GNU General Public License 00022 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>. 00023 00024 Class 00025 Foam::BSpline 00026 00027 Description 00028 An implementation of B-splines. 00029 00030 In this implementation, the end tangents are created automatically 00031 by reflection. 00032 00033 In matrix form, the @e local interpolation on the interval t=[0..1] is 00034 described as follows: 00035 @verbatim 00036 P(t) = 1/6 * [ t^3 t^2 t 1 ] * [ -1 3 -3 1 ] * [ P-1 ] 00037 [ 3 -6 3 0 ] [ P0 ] 00038 [ -3 0 3 0 ] [ P1 ] 00039 [ 1 4 1 0 ] [ P2 ] 00040 @endverbatim 00041 00042 Where P-1 and P2 represent the neighbouring points or the extrapolated 00043 end points. Simple reflection is used to automatically create the end 00044 points. 00045 00046 The spline is discretized based on the chord length of the individual 00047 segments. In rare cases (sections with very high curvatures), the 00048 resulting distribution may be sub-optimal. 00049 00050 SeeAlso 00051 CatmullRomSpline 00052 00053 ToDo 00054 A future implementation could also handle closed splines. 00055 00056 SourceFiles 00057 BSpline.C 00058 00059 \*---------------------------------------------------------------------------*/ 00060 00061 #ifndef BSpline_H 00062 #define BSpline_H 00063 00064 #include "polyLine.H" 00065 00066 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // 00067 00068 namespace Foam 00069 { 00070 00071 /*---------------------------------------------------------------------------*\ 00072 Class BSpline Declaration 00073 \*---------------------------------------------------------------------------*/ 00074 00075 class BSpline 00076 : 00077 public polyLine 00078 { 00079 // Private Member Functions 00080 00081 //- Disallow default bitwise copy construct 00082 BSpline(const BSpline&); 00083 00084 //- Disallow default bitwise assignment 00085 void operator=(const BSpline&); 00086 00087 00088 public: 00089 00090 // Constructors 00091 00092 //- Construct from components 00093 BSpline 00094 ( 00095 const pointField& knots, 00096 const bool notImplementedClosed = false 00097 ); 00098 00099 00100 // Member Functions 00101 00102 //- Return the point position corresponding to the curve parameter 00103 // 0 <= lambda <= 1 00104 point position(const scalar lambda) const; 00105 00106 //- Return the point position corresponding to the local parameter 00107 // 0 <= lambda <= 1 on the given segment 00108 point position(const label segment, const scalar lambda) const; 00109 00110 //- Return the length of the curve 00111 scalar length() const; 00112 }; 00113 00114 00115 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // 00116 00117 } // End namespace Foam 00118 00119 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // 00120 00121 #endif 00122 00123 // ************************ vim: set sw=4 sts=4 et: ************************ //