FreeFOAM The Cross-Platform CFD Toolkit
Hosted by SourceForge:
Get FreeFOAM at SourceForge.net.
            Fast, secure and Free Open Source software downloads

CatmullRomSpline.H

Go to the documentation of this file.
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::CatmullRomSpline
00026 
00027 Description
00028     An implementation of Catmull-Rom splines
00029     (sometimes known as Overhauser splines).
00030 
00031     In this implementation, the end tangents are created automatically
00032     by reflection.
00033 
00034     In matrix form, the @e local interpolation on the interval t=[0..1] is
00035     described as follows:
00036     @verbatim
00037     P(t) = 1/2 * [ t^3 t^2 t 1 ] * [ -1  3 -3  1 ] * [ P-1 ]
00038                                    [  2 -5  4 -1 ]   [ P0 ]
00039                                    [ -1  0  1  0 ]   [ P1 ]
00040                                    [  0  2  0  0 ]   [ P2 ]
00041     @endverbatim
00042 
00043     Where P-1 and P2 represent the neighbouring points or the extrapolated
00044     end points. Simple reflection is used to automatically create the end
00045     points.
00046 
00047     The spline is discretized based on the chord length of the individual
00048     segments. In rare cases (sections with very high curvatures), the
00049     resulting distribution may be sub-optimal.
00050 
00051 SeeAlso
00052     http://www.algorithmist.net/catmullrom.html provides a nice
00053     introduction
00054 
00055 ToDo
00056     A future implementation could also handle closed splines.
00057 
00058 SourceFiles
00059     CatmullRomSpline.C
00060 
00061 \*---------------------------------------------------------------------------*/
00062 
00063 #ifndef CatmullRomSpline_H
00064 #define CatmullRomSpline_H
00065 
00066 #include "polyLine.H"
00067 
00068 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00069 
00070 namespace Foam
00071 {
00072 
00073 /*---------------------------------------------------------------------------*\
00074                       Class CatmullRomSpline Declaration
00075 \*---------------------------------------------------------------------------*/
00076 
00077 class CatmullRomSpline
00078 :
00079     public polyLine
00080 {
00081     // Private Member Functions
00082 
00083         //- Disallow default bitwise copy construct
00084         CatmullRomSpline(const CatmullRomSpline&);
00085 
00086         //- Disallow default bitwise assignment
00087         void operator=(const CatmullRomSpline&);
00088 
00089 
00090 public:
00091 
00092     // Constructors
00093 
00094         //- Construct from components
00095         CatmullRomSpline
00096         (
00097             const pointField& knots,
00098             const bool notImplementedClosed = false
00099         );
00100 
00101 
00102     // Member Functions
00103 
00104         //- Return the point position corresponding to the curve parameter
00105         //  0 <= lambda <= 1
00106         point position(const scalar lambda) const;
00107 
00108         //- Return the point position corresponding to the local parameter
00109         //  0 <= lambda <= 1 on the given segment
00110         point position(const label segment, const scalar lambda) const;
00111 
00112         //- Return the length of the curve - not implemented
00113         scalar length() const;
00114 };
00115 
00116 
00117 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00118 
00119 } // End namespace Foam
00120 
00121 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00122 
00123 #endif
00124 
00125 // ************************************************************************* //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines