MeshKernel
Public Member Functions | Public Attributes | List of all members
meshkernel::Splines Class Reference

A class describing splines. More...

#include <Splines.hpp>

Public Member Functions

 Splines ()=default
 Default constructor.
 
 Splines (Projection projection)
 Ctor, set projection. More...
 
 Splines (CurvilinearGrid const &grid)
 Ctor from grids, each gridline is converted to spline, first the first m_n horizontal lines then the m_m vertical lines. More...
 
void AddSpline (const std::vector< Point > &splines, UInt start, UInt size)
 Adds a new spline to m_splineCornerPoints. More...
 
void AddSpline (const std::vector< Point > &splines)
 Adds a new spline to m_splineCornerPoints. More...
 
void Replace (const UInt splineIndex, const std::vector< Point > &splinePoints)
 Replaces an existing spline. More...
 
void SwapSplines (const UInt firstSpline, const UInt secondSpline)
 Swap all the data for two splines. More...
 
void Reverse (const UInt splineId)
 Reverse the order of the points in the spline.
 
void SnapSpline (const size_t splineIndex, const LandBoundary &landBoundary, const int numberOfIterations=constants::numeric::defaultSnappingIterations)
 Snap the spline to the land boundary (snap_spline) More...
 
bool GetSplinesIntersection (UInt first, UInt second, double &crossProductIntersection, Point &intersectionPoint, double &firstSplineRatio, double &secondSplineRatio) const
 Computes the intersection of two splines (sect3r) More...
 
double ComputeSplineLength (UInt index, double startAdimensionalCoordinate, double endAdimensionalCoordinate, UInt numSamples=100, bool accountForCurvature=false, double height=1.0, double assignedDelta=-1.0) const
 Computes the spline length in s coordinates (GETDIS) More...
 
std::tuple< std::vector< Point >, std::vector< double > > ComputePointOnSplineFromAdimensionalDistance (UInt index, double maximumGridHeight, bool isSpacingCurvatureAdapted, const std::vector< double > &distances) const
 Compute the points on a spline lying at certain distance. More...
 
Point ComputeClosestPointOnSplineSegment (UInt index, double startSplineSegment, double endSplineSegment, Point point) const
 Computes the point on a spline segment which is the closest to another point. More...
 
Point ComputeClosestPoint (UInt index, Point point) const
 Computes the point on a spline segment which is the closest to another point. More...
 
Point Evaluate (const UInt whichSpline, const double lambda) const
 Computes the coordinate of a point on a spline, given the dimensionless distance from the first corner point (splint) More...
 
BoundingBox GetBoundingBox (const UInt splineIndex) const
 Compute the boundary box for the spline indicated by splineIndex.
 
auto GetNumSplines () const
 Get the number of splines. More...
 
UInt Size (const UInt whichSpline) const
 Get the size of a specific spline. More...
 
UInt MaxSizeIndex () const
 Get the index of the spline with the largest number of spline points.
 

Public Attributes

std::vector< std::vector< Point > > m_splineNodes
 The spline corner points.
 
std::vector< std::vector< Point > > m_splineDerivatives
 The spline derivatives at the corner points.
 
std::vector< double > m_splinesLength
 The length of each spline.
 
Projection m_projection = Projection::cartesian
 The map projection.
 

Detailed Description

A class describing splines.

This class stores the corner points of each spline. Besides the corner points, the derivatives at the corner points are also stored. The coordinates of the points between the corner points are computed in the method Splines::ComputePointOnSplineFromAdimensionalDistance.

Constructor & Destructor Documentation

◆ Splines() [1/2]

meshkernel::Splines::Splines ( Projection  projection)
explicit

Ctor, set projection.

[in] projection The map projection

◆ Splines() [2/2]

meshkernel::Splines::Splines ( CurvilinearGrid const &  grid)
explicit

Ctor from grids, each gridline is converted to spline, first the first m_n horizontal lines then the m_m vertical lines.

[in] grid The curvilinear grid

Member Function Documentation

◆ AddSpline() [1/2]

void meshkernel::Splines::AddSpline ( const std::vector< Point > &  splines)

Adds a new spline to m_splineCornerPoints.

Parameters
[in]splinesThe spline corner points

◆ AddSpline() [2/2]

void meshkernel::Splines::AddSpline ( const std::vector< Point > &  splines,
UInt  start,
UInt  size 
)

Adds a new spline to m_splineCornerPoints.

Parameters
[in]splinesThe spline corner points
[in]startThe starting index in splines
[in]sizeThe number of points making up the spline(s).

◆ ComputeClosestPoint()

Point meshkernel::Splines::ComputeClosestPoint ( UInt  index,
Point  point 
) const

Computes the point on a spline segment which is the closest to another point.

Parameters
[in]indexThe spline index
[in]pointThe point to account for in the calculation
Returns
The point on a spline segment which is the closest to the input point

◆ ComputeClosestPointOnSplineSegment()

Point meshkernel::Splines::ComputeClosestPointOnSplineSegment ( UInt  index,
double  startSplineSegment,
double  endSplineSegment,
Point  point 
) const

Computes the point on a spline segment which is the closest to another point.

Parameters
[in]indexThe spline index
[in]startSplineSegmentThe begin of the spline segment to consider
[in]endSplineSegmentThe end of the spline segment to consider
[in]pointThe point to account for in the calculation
Returns
The point on a spline segment which is the closest to the input point

◆ ComputePointOnSplineFromAdimensionalDistance()

std::tuple<std::vector<Point>, std::vector<double> > meshkernel::Splines::ComputePointOnSplineFromAdimensionalDistance ( UInt  index,
double  maximumGridHeight,
bool  isSpacingCurvatureAdapted,
const std::vector< double > &  distances 
) const

Compute the points on a spline lying at certain distance.

Parameters
[in]indexThe spline index
[in]maximumGridHeightMaximum grid height
[in]isSpacingCurvatureAdaptedIs spacing-curvature adapted
[in]distancesThe dimensional distances of each point
Returns
The points along the spline and the adimensional distances of each point

◆ ComputeSplineLength()

double meshkernel::Splines::ComputeSplineLength ( UInt  index,
double  startAdimensionalCoordinate,
double  endAdimensionalCoordinate,
UInt  numSamples = 100,
bool  accountForCurvature = false,
double  height = 1.0,
double  assignedDelta = -1.0 
) const

Computes the spline length in s coordinates (GETDIS)

Parameters
[in]indexThe spline index
[in]startAdimensionalCoordinateAdimensional start spline
[in]endAdimensionalCoordinateAdimensional end spline
[in]numSamplesHow many intervals to use between the startAdimensionalCoordinate and endAdimensionalCoordinate
[in]accountForCurvatureAccounting for curvature
[in]heightWhen accounting for curvature, the height to use
[in]assignedDeltaWhen larger than zero, the number of intervals the spline is divided when computing the length
Returns
The computed length

◆ Evaluate()

Point meshkernel::Splines::Evaluate ( const UInt  whichSpline,
const double  lambda 
) const

Computes the coordinate of a point on a spline, given the dimensionless distance from the first corner point (splint)

Parameters
[in]whichSplineSpline index
[in]lambdaThe adimensinal coordinate where to perform the interpolation
Returns
The interpolated point

◆ GetNumSplines()

auto meshkernel::Splines::GetNumSplines ( ) const
inline

Get the number of splines.

Returns
the number of splines

◆ GetSplinesIntersection()

bool meshkernel::Splines::GetSplinesIntersection ( UInt  first,
UInt  second,
double &  crossProductIntersection,
Point intersectionPoint,
double &  firstSplineRatio,
double &  secondSplineRatio 
) const

Computes the intersection of two splines (sect3r)

Parameters
[in]firstThe index of the first spline
[in]secondThe index of the second spline
[out]crossProductIntersectionThe cross product of the intersection
[out]intersectionPointThe intersection point
[out]firstSplineRatioThe ratio of the first spline length where the intersection occurs
[out]secondSplineRatioThe ratio of the second spline length where the intersection occurs
Returns
If a valid intersection is found

◆ Replace()

void meshkernel::Splines::Replace ( const UInt  splineIndex,
const std::vector< Point > &  splinePoints 
)

Replaces an existing spline.

Parameters
[in]splineIndexThe index of the spline to be replaced
[in]splinePointsThe spline points

◆ Size()

UInt meshkernel::Splines::Size ( const UInt  whichSpline) const

Get the size of a specific spline.

Returns
the size of the desired spline

◆ SnapSpline()

void meshkernel::Splines::SnapSpline ( const size_t  splineIndex,
const LandBoundary landBoundary,
const int  numberOfIterations = constants::numeric::defaultSnappingIterations 
)

Snap the spline to the land boundary (snap_spline)

Parameters
[in]splineIndexThe index of the spline to be snapped to boundary
[in]landBoundaryThe boundary to which the spline will be snapped
[in]numberOfIterationsThe maximum number of iterations that should be performed.

◆ SwapSplines()

void meshkernel::Splines::SwapSplines ( const UInt  firstSpline,
const UInt  secondSpline 
)

Swap all the data for two splines.

On exit the contents of the two splines will be swapped


The documentation for this class was generated from the following file: