This repository has been archived on 2023-11-20. You can view files and clone it, but cannot push or open issues or pull requests.
foam-extend4.1-coherent-io/applications/utilities/mesh/generation/blockMesh/curvedEdges/CatmullRomSpline.C
2013-12-11 16:09:41 +00:00

140 lines
3.1 KiB
C

/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | foam-extend: Open Source CFD
\\ / O peration |
\\ / A nd | For copyright notice see file Copyright
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of foam-extend.
foam-extend 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.
foam-extend 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
along with foam-extend. If not, see <http://www.gnu.org/licenses/>.
\*---------------------------------------------------------------------------*/
#include "error.H"
#include "CatmullRomSpline.H"
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Foam::CatmullRomSpline::CatmullRomSpline
(
const pointField& knots,
const bool closed
)
:
polyLine(knots, closed)
{}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
Foam::point Foam::CatmullRomSpline::position(const scalar mu) const
{
// endpoints
if (mu < SMALL)
{
return points()[0];
}
else if (mu > 1 - SMALL)
{
return points()[points().size()-1];
}
scalar lambda = mu;
label segment = localParameter(lambda);
return position(segment, lambda);
}
Foam::point Foam::CatmullRomSpline::position
(
const label segment,
const scalar mu
) const
{
// out-of-bounds
if (segment < 0)
{
return points()[0];
}
else if (segment > nSegments())
{
return points()[points().size()-1];
}
const point& p0 = points()[segment];
const point& p1 = points()[segment+1];
// special cases - no calculation needed
if (mu <= 0.0)
{
return p0;
}
else if (mu >= 1.0)
{
return p1;
}
// 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];
}
return 0.5 *
(
( 2*p0 )
+ mu *
(
( -e0 + p1 )
+ mu *
(
( 2*e0 - 5*p0 + 4*p1 - e1 )
+ mu *
( -e0 + 3*p0 - 3*p1 + e1 )
)
)
);
}
Foam::scalar Foam::CatmullRomSpline::length() const
{
notImplemented("CatmullRomSpline::length() const");
return 1.0;
}
// ************************************************************************* //