MeshKernel
Classes | Typedefs | Enumerations | Functions | Variables
meshkernel Namespace Reference

Contains the logic of the C++ static library. More...

Classes

class  AlgorithmError
 A class for throwing algorithm exceptions. More...
 
class  AveragingInterpolation
 The class used to interpolate based on averaging. More...
 
class  BilinearInterpolationOnGriddedSamples
 A class for performing bilinear interpolation on gridded samples. More...
 
class  BoundingBox
 A class defining a bounding box. More...
 
struct  Cartesian3DPoint
 A struct describing the three coordinates in a cartesian projection. More...
 
class  CasulliDeRefinement
 Compute the Casulli de-refinement for a mesh. More...
 
class  CasulliRefinement
 Compute the Casulli refinement for a mesh. More...
 
class  ConnectMeshes
 Connects grids across element faces. More...
 
class  ConstraintError
 An exception class thrown when an attempt is made that violates a range constraint. More...
 
class  Contacts
 A class describing an 1d-2d contacts. More...
 
class  ConvertCartesianToSpherical
 Converts points from spherical to Cartesian coordinate system. More...
 
class  ConvertCartesianToSphericalBase
 Converts points from spherical to Cartesian coordinate system. More...
 
class  ConvertCartesianToSphericalEPSG
 Converts points from spherical to Cartesian coordinate system using an EPSG code. More...
 
class  ConvertSphericalToCartesian
 Converts points from spherical to Cartesian coordinate system. More...
 
class  ConvertSphericalToCartesianBase
 Namespace alias for boost::geometry. More...
 
class  ConvertSphericalToCartesianEPSG
 Converts points from spherical to Cartesian coordinate system using an ESPG code. More...
 
struct  CurvilinearParameters
 A struct used to describe parameters for generating a curvilinear grid in a C-compatible manner. More...
 
struct  EdgeMeshPolyLineIntersection
 An intersection with a mesh edge. More...
 
class  ErrorCategory
 Contains error category information. More...
 
struct  FaceMeshPolyLineIntersection
 An intersection with a mesh face. More...
 
class  FlipEdges
 A class used to improve mesh connectivity. More...
 
class  FormatString
 Manages the format string and source location. More...
 
struct  FuncAdimensionalToDimensionalDistanceOnSpline
 This structure is used to create a function for converting an adimensional distance on a spline to a dimensional one. More...
 
struct  FuncDistanceFromAPoint
 This structure is used to compute the point on a spline closest to another point. More...
 
class  Hessian
 The hessian values. More...
 
class  LandBoundaries
 A class describing land boundaries. These are used to visualise the land-water interface. More...
 
class  LandBoundary
 A class containing the land boundary polylines. More...
 
class  LinearAlgebraError
 A class for throwing linear algebra exceptions. More...
 
struct  MakeGridParameters
 This struct describes the necessary parameters to create a new curvilinear grid in a C-compatible manner. More...
 
class  Mesh
 A class describing an unstructured mesh. This class contains the shared functionality between 1d or 2d meshes. More...
 
class  Mesh1D
 A class derived from Mesh, which describes 1d meshes. More...
 
class  Mesh2D
 A class derived from Mesh, which describes unstructures 2d meshes. More...
 
class  Mesh2DGenerateGlobal
 Construct a global grid in spherical coordinates, as a base for later mesh refinements. More...
 
class  Mesh2DIntersections
 Compute the intersections of polygon inner and outer perimeters. More...
 
class  Mesh2DToCurvilinear
 Construct a curvilinear grid from an unstructured mesh. More...
 
class  MeshConversion
 Apply a conversion to nodes of a mesh. More...
 
class  MeshGeometryError
 A class for throwing mesh geometry errors. More...
 
class  MeshInterpolation
 Interface for interpolation methods. More...
 
class  MeshKernelError
 A class for throwing general MeshKernel exceptions. More...
 
class  MeshRefinement
 A class used to refine a Mesh2D instance. More...
 
struct  MeshRefinementParameters
 A struct used to describe the mesh refinement parameters in a C-compatible manner. More...
 
class  MeshTransformation
 Apply a transformation to a mesh. More...
 
class  Network1D
 A class describing a network 1d. More...
 
class  NotImplementedError
 A class for throwing not implemented exceptions. More...
 
class  OrthogonalizationAndSmoothing
 Orthogonalizion (optimize the aspect ratios) and mesh smoothing (optimize internal face angles or area). More...
 
struct  OrthogonalizationParameters
 A struct used to describe the orthogonalization parameters in a C-compatible manner. More...
 
class  Orthogonalizer
 Orthogonalizion (optimize the aspect ratios) and mesh smoothing (optimize internal face angles or area). More...
 
class  Point
 A struct describing a point in a two-dimensional space. More...
 
class  Polygon
 A closed polygon. More...
 
class  PolygonalEnclosure
 A region enclosed by a polygonal permieter. More...
 
class  Polygons
 A class containing a list of polygonaly enclosed regions. More...
 
class  RangeError
 A class for throwing range error exceptions. More...
 
class  RemoveDisconnectedRegions
 Removes any small regions that are disconnected from the main part of the domain. More...
 
class  RigidBodyTransformation
 A composition of translation and rotation transformations. More...
 
class  Rotation
 Apply a rotation transformation to a point or a vector. More...
 
class  Sample
 A struct describing a sample with two coordinates and a value. More...
 
class  SamplesHessianCalculator
 A class implementing the computation of the real component the local hessian eigenvalues. More...
 
class  Smoother
 Orthogonalizion (optimize the aspect ratios) and mesh smoothing (optimize internal face angles or area). More...
 
class  SplineAlgorithms
 Provide algorithms operating on splines. More...
 
class  Splines
 A class describing splines. More...
 
struct  SplinesToCurvilinearParameters
 A struct used to describe the spline to curvilinear grid parameters in a C-compatible manner. More...
 
class  SplitRowColumnOfMesh
 Split the row or column connected to the given edge. More...
 
class  Translation
 Apply a translation transformation to a point or a vector. More...
 
class  TriangulationInterpolation
 A class used for triangulation interpolation. More...
 
struct  TriangulationWrapper
 Wrapper around the Triangle library. More...
 
class  Vector
 A class defining a vector. More...
 

Typedefs

using UInt = std::uint32_t
 Integer type used when indexing mesh graph entities.
 
using Edge = std::pair< UInt, UInt >
 Describes an edge with two indices.
 
using EdgeFaces = std::array< UInt, 2 >
 Contains the ID's of elements either side of an edge.
 
using HessianDimension = std::array< UInt, 3 >
 Array containing dimensions of the hessian.
 
using MatrixColMajor = lin_alg::Matrix< double, Eigen::ColMajor >
 Define column major orientation.
 

Enumerations

enum  Projection { cartesian = 0, spherical = 1, sphericalAccurate = 2 }
 Enumerator describing the supported projections.
 
enum  TraversalDirection { TraversalDirection::Clockwise, TraversalDirection::AntiClockwise }
 Indicator for traversal direction of the points specifying a polygon. More...
 
enum  Location { Location::Faces = 0, Location::Nodes = 1, Location::Edges = 2, Location::Unknown = 3 }
 Mesh locations enumeration. More...
 
enum  CurvilinearDirection { CurvilinearDirection::M, CurvilinearDirection::N }
 Direction to use in curvilinear grid algorithms. More...
 
enum  InterpolationDataTypes { InterpolationDataTypes::Short = 0, InterpolationDataTypes::Float = 1, InterpolationDataTypes::Int = 2, InterpolationDataTypes::Double = 3 }
 The possible types of the values to be interpolated in the gridded sample. More...
 
enum  ExitCode {
  Success = 0, MeshKernelErrorCode = 1, NotImplementedErrorCode = 2, AlgorithmErrorCode = 3,
  ConstraintErrorCode = 4, MeshGeometryErrorCode = 5, LinearAlgebraErrorCode = 6, RangeErrorCode = 7,
  StdLibExceptionCode = 8, UnknownExceptionCode = 9
}
 Enumeration of exit codes. More...
 
enum  Interpolants { AveragingInterpolation = 0, BilinearInterpolationOnGriddedSamples = 1 }
 The interpolant types.
 

Functions

const std::vector< int > & GetValidProjections ()
 Gets the valid projectionbs as vector of integers.
 
const std::vector< int > & GetValidDeletionOptions ()
 Gets the valid deletion options as vector of integers.
 
Projection GetProjectionValue (int projection)
 Convert an integer value to the Projection enumeration type. More...
 
const std::string & ProjectionToString (Projection projection)
 Get the string representation of the Projection enumeration values.
 
CurvilinearDirection GetCurvilinearDirectionValue (int direction)
 Convert an integer value to the CurvilinearDirection enumeration type. More...
 
const std::string & CurvilinearDirectionToString (CurvilinearDirection direction)
 Get the string representation of the CurvilinearDirection enumeration values.
 
template<typename T >
void ResizeAndFill2DVector (std::vector< std::vector< T >> &v, UInt const &firstDimension, UInt const &secondDimension, bool fill=false, const T &fillValue={})
 Resizes and fills a two dimensional vector. More...
 
template<typename T >
void ResizeAndFill3DVector (std::vector< std::vector< std::vector< T >>> &v, UInt const &firstDimension, UInt const &secondDimension, UInt const &thirdDim, bool fill=false, const T &fillValue={})
 Resizes and fills a three dimensional vector. More...
 
Cartesian3DPoint VectorProduct (const Cartesian3DPoint &a, const Cartesian3DPoint &b)
 Defines vector product for cartesian 3D-space. More...
 
double InnerProduct (const Cartesian3DPoint &a, const Cartesian3DPoint &b)
 Defines inner product in cartesian 3D-space. More...
 
Cartesian3DPoint operator* (const Cartesian3DPoint &p, const double value)
 Multiply Cartesian point p by a scalar value.
 
Cartesian3DPoint operator* (const double value, const Cartesian3DPoint &p)
 Multiply Cartesian point p by a scalar value.
 
Cartesian3DPoint operator/ (const Cartesian3DPoint &p, const double value)
 Divide Cartesian point p by a scalar value.
 
Cartesian3DPoint operator+ (const Cartesian3DPoint &p1, const Cartesian3DPoint &p2)
 Add Cartesian point p2 to p1.
 
Cartesian3DPoint operator- (const Cartesian3DPoint &p1, const Cartesian3DPoint &p2)
 Subtract Cartesian point p2 from p1.
 
template<typename T >
UInt FindIndex (const std::vector< T > &vec, T el)
 Find index of a certain element. More...
 
template<typename T >
UInt FindNextIndex (const std::vector< T > &vec, UInt current)
 Find the next index in the vector, wraps around when current is the last index.
 
template<typename T >
UInt FindPreviousIndex (const std::vector< T > &vec, UInt current)
 Find the previous index in the vector, wraps around when current is the first index.
 
std::vector< std::pair< UInt, UInt > > FindIndices (const std::vector< Point > &vec, size_t start, size_t end, double separator)
 Find all start-end positions in a vector separated by a separator. More...
 
UInt InvalidPointCount (const std::vector< Point > &points)
 Determine if the number of invalid points in the point array.
 
UInt InvalidPointCount (const std::vector< Point > &points, size_t start, size_t end)
 Determine if the number of invalid points in the section of the point array. More...
 
template<typename T >
std::vector< UIntSortedIndices (const std::vector< T > &v)
 Sort a vector and return the sorted indices. More...
 
template<typename T >
auto ReorderVector (const std::vector< T > &v, const std::vector< UInt > &order)
 Reorder vector accordingly to a specific order. More...
 
template<typename F >
double FindFunctionRootWithGoldenSectionSearch (F func, double min, double max)
 Algorithm performing the zero's search using the golden section algorithm's. More...
 
UInt NextCircularForwardIndex (UInt currentIndex, UInt size)
 Get the next forward index, for a range in 0 .. size-1. More...
 
UInt NextCircularBackwardIndex (UInt currentIndex, UInt size)
 Get the next backward index, for a range in 0 .. size-1. More...
 
bool IsPointOnPole (const Point &point)
 Determines if a point is close to the poles (latitude close to 90 degrees). More...
 
Cartesian3DPoint SphericalToCartesian3D (const Point &sphericalPoint)
 Transforms 2D point in spherical coordinates to 3D cartesian coordinates. More...
 
Point Cartesian3DToSpherical (const Cartesian3DPoint &cartesianPoint, double referenceLongitude)
 Transforms 3D cartesian coordinates to 2D point in spherical coordinates. More...
 
double crossProduct (const Point &firstSegmentFirstPoint, const Point &firstSegmentSecondPoint, const Point &secondSegmentFirstPoint, const Point &secondSegmentSecondPoint, const Projection &projection)
 Computes the cross product between two segments (duitpl) More...
 
bool IsPointInPolygonNodes (const Point &point, const std::vector< Point > &polygonNodes, const Projection &projection, Point polygonCenter={constants::missing::doubleValue, constants::missing::doubleValue}, UInt startNode=constants::missing::uintValue, UInt endNode=constants::missing::uintValue)
 Checks if a point is in polygonNodes using the winding number method. More...
 
void ComputeThreeBaseComponents (const Point &point, std::array< double, 3 > &exxp, std::array< double, 3 > &eyyp, std::array< double, 3 > &ezzp)
 Computes three base components.
 
void ComputeTwoBaseComponents (const Point &point, double(&elambda)[3], double(&ephi)[3])
 Computes two base components.
 
double GetDx (const Point &firstPoint, const Point &secondPoint, const Projection &projection)
 Gets dx for the given projection. More...
 
Vector GetDelta (const Point &firstPoint, const Point &secondPoint, const Projection &projection)
 Get the delta (dx, dy) for the given projection. More...
 
Vector ComputeNormalToline (const Point &start, const Point &end, const Projection projection)
 Get the normal to the line described by the two points. More...
 
double GetDy (const Point &firstPoint, const Point &secondPoint, const Projection &projection)
 Gets dy for the given projection. More...
 
double OuterProductTwoSegments (const Point &firstPointFirstSegment, const Point &secondPointFirstSegment, const Point &firstPointSecondSegment, const Point &secondPointSecondSegment, const Projection &projection)
 Outer product of two segments (dprodout)
 
Point ComputeMiddlePoint (const Point &firstPoint, const Point &secondPoint, const Projection &projection)
 Computes the middle point. More...
 
Point ComputeMiddlePointAccountingForPoles (const Point &firstPoint, const Point &secondPoint, const Projection &projection)
 Computes the middle point (account for poles, latitudes close to 90 degrees) More...
 
Point NormalVector (const Point &firstPoint, const Point &secondPoint, const Point &insidePoint, const Projection &projection)
 Normalized vector of a segment in direction 1 -> 2 with the insidePoint orientation. More...
 
void TransformGlobalVectorToLocal (const Point &reference, const Point &globalCoordinates, const Point &globalComponents, const Projection &projection, Point &localComponents)
 Transforms vector with components in global spherical coordinate directions(xglob, yglob) to local coordinate directions(xloc, yloc) around reference point(xref, yref)
 
Point NormalVectorOutside (const Point &firstPoint, const Point &secondPoint, const Projection &projection)
 Computes the normal vector outside (normalout) More...
 
void NormalVectorInside (const Point &firstPoint, const Point &secondPoint, const Point &insidePoint, Point &normal, bool &flippedNormal, const Projection &projection)
 Computes the normal vector to a line 1-2, which is outward w.r.t. an 'inside' point 3. More...
 
void AddIncrementToPoint (const Point &normal, double increment, const Point &referencePoint, const Projection &projection, Point &point)
 Moves a point by adding an increment vector to it. More...
 
void TranslateSphericalCoordinates (std::vector< Point > &polygon)
 For a given polygon the function may shift the input coordinates. More...
 
Point ReferencePoint (const std::vector< Point > &polygon, const Projection &projection)
 For a given polygon compute a reference point. More...
 
Point ReferencePoint (const std::vector< Point > &nodes, const std::vector< UInt > &polygonIndices, const Projection &projection)
 For a given polygon compute a reference point. More...
 
double ComputeSquaredDistance (const Point &firstPoint, const Point &secondPoint, const Projection &projection)
 Computes the squared distance between two points This is faster than ComputeDistance because it does not take the square root. More...
 
double ComputeDistance (const Point &firstPoint, const Point &secondPoint, const Projection &projection)
 Computes the distance between two points (dbdistance) More...
 
std::tuple< double, Point, double > DistanceFromLine (const Point &point, const Point &firstNode, const Point &secondNode, const Projection &projection)
 Computes the perpendicular distance of a point from a segment firstNode - secondNode (dlinedis3) More...
 
double InnerProductTwoSegments (const Point &firstPointFirstSegment, const Point &secondPointFirstSegment, const Point &firstPointSecondSegment, const Point &secondPointSecondSegment, const Projection &projection)
 Inner product of two segments (dprodin) More...
 
double NormalizedInnerProductTwoSegments (const Point &firstPointFirstSegment, const Point &secondPointFirstSegment, const Point &firstPointSecondSegment, const Point &secondPointSecondSegment, const Projection &projection)
 The normalized inner product of two segments (dcosphi) More...
 
Point CircumcenterOfTriangle (const Point &firstNode, const Point &secondNode, const Point &thirdNode, const Projection &projection)
 Computes the circumcenter of a triangle. More...
 
std::tuple< bool, Point, double, double, double > AreSegmentsCrossing (const Point &firstSegmentFirstPoint, const Point &firstSegmentSecondPoint, const Point &secondSegmentFirstPoint, const Point &secondSegmentSecondPoint, bool adimensionalCrossProduct, const Projection &projection)
 Determines if two segments are crossing (cross, cross3D) More...
 
template<typename T >
ComputePointOnSplineAtAdimensionalDistance (const std::vector< T > &coordinates, const std::vector< T > &coordinatesDerivatives, double pointAdimensionalCoordinate)
 Computes the coordinate of a point on a spline, given the dimensionless distance from the first corner point (splint) More...
 
std::tuple< std::vector< double >, double > ComputeAdimensionalDistancesFromPointSerie (const std::vector< Point > &v, const Projection &projection)
 Computes dimensionless distances of a vector of points such as the first entry has distance 0 and the last entry has distance 1. More...
 
template<typename T >
int sgn (T val)
 Computes the sign of a type. More...
 
lin_alg::Matrix< PointDiscretizeTransfinite (const std::vector< Point > &leftDiscretization, const std::vector< Point > &rightDiscretization, const std::vector< Point > &bottomDiscretization, const std::vector< Point > &upperDiscretization, const Projection &projection, UInt numM, UInt numN)
 Computes the transfinite discretization inside the area defined by 4 sides, each one discretized with a series of points (tranfn2). More...
 
std::vector< PointComputeEdgeCenters (const std::vector< Point > &nodes, const std::vector< Edge > &edges)
 Computes the edge centers. More...
 
double LinearInterpolationInTriangle (const Point &interpolationPoint, const std::vector< Point > &polygon, const std::vector< double > &values, const Projection &projection)
 Given a triangles with values on each node, computes the interpolated value inside the triangle, using linear interpolation. More...
 
Point ComputeAverageCoordinate (const std::vector< Point > &points, const Projection &projection)
 Given a series of point computes the average coordinate. More...
 
std::tuple< Point, double, bool > OrthogonalProjectionOnSegment (Point const &firstNode, Point const &secondNode, Point const &pointToProject)
 Cartesian projection of a point on a segment defined by other two points. More...
 
UInt AbsoluteDifference (UInt number_1, UInt number_2)
 Calculates the absolute difference between to Index numbers. More...
 
std::vector< PointComputePolyLineDiscretization (std::vector< Point > const &polyline, std::vector< double > &chainages, Projection projection)
 Computes the discretization points along a polyline. More...
 
std::vector< double > ComputePolyLineNodalChainages (std::vector< Point > const &polyline, Projection projection)
 Computes the chainages of each polyline node. More...
 
std::vector< double > ComputePolyLineEdgesLengths (std::vector< Point > const &polyline, Projection projection)
 Computes the lengths of each polyline segment. More...
 
double MatrixNorm (const std::vector< double > &x, const std::vector< double > &y, const std::vector< double > &matCoefficients)
 Compute the matrix norm. More...
 
void IncrementValidValue (UInt &value, const UInt increment)
 Increment a valid value by an increment.
 
double lengthSquared (const Point &p1)
 Compute the dot product of a point with itself. More...
 
double dot (const Point &p1, const Point &p2)
 Compute the dot product of two points. More...
 
double dot (const Point &p, const Vector &v)
 Compute the dot product of a point and a vector. More...
 
double dot (const Vector &v, const Point &p)
 Compute the dot product of a point and a vector. More...
 
double cross (const Point &p1, const Point &p2)
 Compute the cross product of two points (as vectors). More...
 
Point operator- (const Point &pnt)
 Unary minus. More...
 
Point operator+ (const Point &p1, const Point &p2)
 Add two points. More...
 
Point operator+ (const Point &pnt, const Vector &vec)
 Add vector to a point. More...
 
Point operator+ (const Point &p, const double x)
 Add points and scalar. More...
 
Point operator+ (const double x, const Point &p)
 Add points and scalar. More...
 
Point operator- (const Point &p1, const Point &p2)
 Subtract two points. More...
 
Point operator- (const Point &pnt, const Vector &vec)
 Subtract vector from a point. More...
 
Point operator- (const Point &p, const double x)
 Subtract scalar from point. More...
 
Point operator- (const double x, const Point &p)
 Subtract point from scalar. More...
 
Point operator* (const Point &p1, const Point &p2)
 Multiply point by a point. More...
 
Point operator* (const double x, const Point &p)
 Multiply point by a scalar. More...
 
Point operator* (const Point &p, const double x)
 Multiply point by a scalar. More...
 
Point operator/ (const Point &p1, const Point &p2)
 Divide point by a point. More...
 
Point operator/ (const Point &p, const double x)
 Divide point by a scalar. More...
 
bool operator== (const Point &p1, const Point &p2)
 Test points for equality using a default tolerance. More...
 
bool operator!= (const Point &p1, const Point &p2)
 Test points for equality using a default tolerance. More...
 
bool IsEqual (const Point &p1, const Point &p2, const double epsilon)
 Test points for equality upto a tolerance. More...
 
Point PointAlongLine (const Point &startPoint, const Point &endPoint, const double lambda)
 Compute the point at some position along the line connecting start- and end-point. More...
 
double GetDeltaXCartesian (const Point &p1, const Point &p2)
 Get the delta-x in Cartesian coordinate system.
 
double GetDeltaYCartesian (const Point &p1, const Point &p2)
 Get the delta-y in Cartesian coordinate system.
 
Vector GetDeltaCartesian (const Point &p1, const Point &p2)
 Get the delta-x and -y in Cartesian coordinate system.
 
void Triangulation (int jatri, double const *const xs, double const *const ys, int ns, int *const indx, int *const numtri, int *const edgeidx, int *const numedge, int *const triedge, double *const xs3, double *const ys3, int *const ns3, double trisize)
 Function of the Triangle library. More...
 
Vector normalise (const Vector &vec)
 Return the normalised vector. More...
 
double dot (const Vector &v1, const Vector &v2)
 Compute the dot product of two vectors.
 
double angleBetween (const Vector &v1, const Vector &v2)
 Compute the angle between two vector.
 
Vector operator- (const Vector &vec)
 Unary minus. More...
 
Vector operator+ (const Vector &v1, const Vector &v2)
 Add two vectors.
 
Vector operator- (const Vector &v1, const Vector &v2)
 Subtract vector from another.
 
Vector operator* (const double alpha, const Vector &vec)
 Multiply vector by a scalar.
 
Vector operator* (const Vector &vec, const double alpha)
 Multiply vector by a scalar.
 
Vector operator/ (const Vector &vec, const double alpha)
 Divide vector by a scalar.
 

Variables

template<typename T >
concept InterpolatableType
 Defines the iterpolatable data types. More...
 
template<typename Operation >
concept HasTransformationOperation = requires(Operation op, Point p) {{ op(p)} -> std::same_as<Point>; }
 Ensure any instantiation of the MeshTransformation Compute function is with the correct operation.
 
template<typename Operation >
concept HasConversionProjection
 Ensure any instantiation of the MeshConversion Compute function is able to determine source and target projection types. More...
 
template<typename Operation >
concept ConversionFunctor = HasTransformationOperation<Operation> && HasConversionProjection<Operation>
 Ensure the MeshConversion template parameter has a valid interface.
 
template<typename Function >
concept HasTransformationFunction = requires(Function op, Point p) {{ op(p)} -> std::same_as<Point>; }
 Ensure any instantiation of the MeshTransformation Compute function is with the correct operation.
 
template<typename Function >
concept HasTransformationProjection = requires(Function op) {{ op.TransformationProjection()} -> std::same_as<Projection>; }
 Ensure any instantiation of the MeshTransformation Compute function is able to determine the projection.
 
template<typename Function >
concept TransformationFunction = HasTransformationFunction<Function> && HasTransformationProjection<Function>
 To ensure the MeshTransformation Compute template parameter has a valid interface.
 

Detailed Description

Contains the logic of the C++ static library.

Enumeration Type Documentation

◆ CurvilinearDirection

Direction to use in curvilinear grid algorithms.

Enumerator

M-direction.

N-direction.

◆ ExitCode

Enumeration of exit codes.

Enumerator
Success 

Success.

MeshKernelErrorCode 

MehKernel error.

NotImplementedErrorCode 

Not implemented error.

AlgorithmErrorCode 

Algorithm error.

ConstraintErrorCode 

Constraint error.

MeshGeometryErrorCode 

Geometry error.

LinearAlgebraErrorCode 

Linear algebra error.

RangeErrorCode 

Range error.

StdLibExceptionCode 

Standrad library exception.

UnknownExceptionCode 

Unknown exception.

◆ InterpolationDataTypes

The possible types of the values to be interpolated in the gridded sample.

Enumerator
Short 

short type

Float 

float type

Int 

int type

Double 

double type

◆ Location

enum meshkernel::Location
strong

Mesh locations enumeration.

Enumerator
Faces 

Faces.

Nodes 

Nodes.

Edges 

Edges.

Unknown 

Unknown.

◆ TraversalDirection

Indicator for traversal direction of the points specifying a polygon.

Enumerator
Clockwise 

Points define a clockwise traversal of the polygon.

AntiClockwise 

Points define a anti-clockwise (counter-clockwise) traversal of the polygon.

Function Documentation

◆ AbsoluteDifference()

UInt meshkernel::AbsoluteDifference ( UInt  number_1,
UInt  number_2 
)

Calculates the absolute difference between to Index numbers.

Parameters
[in]number_1The first number
[in]number_2The second number

◆ AddIncrementToPoint()

void meshkernel::AddIncrementToPoint ( const Point normal,
double  increment,
const Point referencePoint,
const Projection projection,
Point point 
)

Moves a point by adding an increment vector to it.

Parameters
[in]normalThe increment direction.
[in]incrementThe increment to use.
[in]referencePointThe reference point containing the reference latitude to use.
[in]projectionThe coordinate system projection.
[in,out]pointThe point to be incremented.

◆ AreSegmentsCrossing()

std::tuple<bool, Point, double, double, double> meshkernel::AreSegmentsCrossing ( const Point firstSegmentFirstPoint,
const Point firstSegmentSecondPoint,
const Point secondSegmentFirstPoint,
const Point secondSegmentSecondPoint,
bool  adimensionalCrossProduct,
const Projection projection 
)

Determines if two segments are crossing (cross, cross3D)

Parameters
[in]firstSegmentFirstPointThe first point of the first segment
[in]firstSegmentSecondPointThe second point of the first segment
[in]secondSegmentFirstPointThe first point of the second segment
[in]secondSegmentSecondPointThe second point of the second segment
[in]adimensionalCrossProductWhether to compute the dimensionless cross product
[in]projectionThe coordinate system projection
Returns
A tuple with: If the two segments are crossing The intersection point The cross product of the intersection The distance of the intersection from the first node of the first segment, expressed as a ratio of the segment length The distance of the intersection from the first node of the second segment, expressed as a ratio of the segment length

◆ Cartesian3DToSpherical()

Point meshkernel::Cartesian3DToSpherical ( const Cartesian3DPoint cartesianPoint,
double  referenceLongitude 
)

Transforms 3D cartesian coordinates to 2D point in spherical coordinates.

Parameters
[in]cartesianPointThe 3d cartesian point
[in]referenceLongitudeThe reference longitude
Returns
The spherical coordinate

◆ CircumcenterOfTriangle()

Point meshkernel::CircumcenterOfTriangle ( const Point firstNode,
const Point secondNode,
const Point thirdNode,
const Projection projection 
)

Computes the circumcenter of a triangle.

Parameters
[in]firstNodeThe first triangle node
[in]secondNodeThe second triangle node
[in]thirdNodeThe third triangle node
[in]projectionThe coordinate system projection
Returns
The resulting circumcenter

◆ ComputeAdimensionalDistancesFromPointSerie()

std::tuple<std::vector<double>, double> meshkernel::ComputeAdimensionalDistancesFromPointSerie ( const std::vector< Point > &  v,
const Projection projection 
)

Computes dimensionless distances of a vector of points such as the first entry has distance 0 and the last entry has distance 1.

Parameters
[in]vThe vector of points
[in]projectionThe projection to use.
Returns
The dimensionless distances and the dimensional total distance (used for normalization)

◆ ComputeAverageCoordinate()

Point meshkernel::ComputeAverageCoordinate ( const std::vector< Point > &  points,
const Projection projection 
)

Given a series of point computes the average coordinate.

Parameters
[in]pointsThe point series.
[in]projectionThe projection to use.
Returns
The average coordinate.

◆ ComputeDistance()

double meshkernel::ComputeDistance ( const Point firstPoint,
const Point secondPoint,
const Projection projection 
)

Computes the distance between two points (dbdistance)

Parameters
[in]firstPointThe first point.
[in]secondPointThe second point.
[in]projectionThe coordinate system projection.
Returns
The distance

◆ ComputeEdgeCenters()

std::vector<Point> meshkernel::ComputeEdgeCenters ( const std::vector< Point > &  nodes,
const std::vector< Edge > &  edges 
)

Computes the edge centers.

Parameters
[in]nodesThe vector of edge nodes.
[in]edgesThe vector of edge indices.
Returns
The vector containing the edge centers.

◆ ComputeMiddlePoint()

Point meshkernel::ComputeMiddlePoint ( const Point firstPoint,
const Point secondPoint,
const Projection projection 
)

Computes the middle point.

Parameters
[in]firstPointThe first point of the segment.
[in]secondPointThe second point of the segment.
[in]projectionThe coordinate system projection.
Returns
The middle point.

◆ ComputeMiddlePointAccountingForPoles()

Point meshkernel::ComputeMiddlePointAccountingForPoles ( const Point firstPoint,
const Point secondPoint,
const Projection projection 
)

Computes the middle point (account for poles, latitudes close to 90 degrees)

Parameters
[in]firstPointThe first point of the segment.
[in]secondPointThe second point of the segment.
[in]projectionThe coordinate system projection.
Returns
The middle point.

◆ ComputeNormalToline()

Vector meshkernel::ComputeNormalToline ( const Point start,
const Point end,
const Projection  projection 
)

Get the normal to the line described by the two points.

The vector is normalised.

Parameters
[in]startPoint at the start of the line
[in]endPoint at the end of the line
[in]projectionThe coordinate system projection.

◆ ComputePointOnSplineAtAdimensionalDistance()

template<typename T >
T meshkernel::ComputePointOnSplineAtAdimensionalDistance ( const std::vector< T > &  coordinates,
const std::vector< T > &  coordinatesDerivatives,
double  pointAdimensionalCoordinate 
)

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

Parameters
[in]coordinatesThe spline node coordinates
[in]coordinatesDerivativesThe spline nodal derivatives
[in]pointAdimensionalCoordinateThe adimensinal coordinate where to perform the interpolation
Returns
The interpolated point

◆ ComputePolyLineDiscretization()

std::vector<Point> meshkernel::ComputePolyLineDiscretization ( std::vector< Point > const &  polyline,
std::vector< double > &  chainages,
Projection  projection 
)

Computes the discretization points along a polyline.

Parameters
polylineA polyline described by its nodes
chainagesThe chainages used for dicretizing the current polyline
projectionThe projection to use
Returns
The discretized polyline

◆ ComputePolyLineEdgesLengths()

std::vector<double> meshkernel::ComputePolyLineEdgesLengths ( std::vector< Point > const &  polyline,
Projection  projection 
)

Computes the lengths of each polyline segment.

Parameters
polylineA polyline described by its nodes
projectionThe projection to use
Returns
A vector containing the lengths of each polyline segment

◆ ComputePolyLineNodalChainages()

std::vector<double> meshkernel::ComputePolyLineNodalChainages ( std::vector< Point > const &  polyline,
Projection  projection 
)

Computes the chainages of each polyline node.

Parameters
polylineA polyline described by its nodes
projectionThe projection to use
Returns
A vector containing the chainage volau of the polyline nodes

◆ ComputeSquaredDistance()

double meshkernel::ComputeSquaredDistance ( const Point firstPoint,
const Point secondPoint,
const Projection projection 
)

Computes the squared distance between two points This is faster than ComputeDistance because it does not take the square root.

Parameters
[in]firstPointThe first point.
[in]secondPointThe second point.
[in]projectionThe coordinate system projection.
Returns
The squared distance

◆ cross()

double meshkernel::cross ( const Point p1,
const Point p2 
)
inline

Compute the cross product of two points (as vectors).

This is mainly a convenience function

◆ crossProduct()

double meshkernel::crossProduct ( const Point firstSegmentFirstPoint,
const Point firstSegmentSecondPoint,
const Point secondSegmentFirstPoint,
const Point secondSegmentSecondPoint,
const Projection projection 
)

Computes the cross product between two segments (duitpl)

Parameters
[in]firstSegmentFirstPointThe first point of the first segment
[in]firstSegmentSecondPointThe second point of the first segment
[in]secondSegmentFirstPointThe second point of the second segment
[in]secondSegmentSecondPointThe second point of the second segment
[in]projectionThe coordinate system projection
Returns
The cross product value

◆ DiscretizeTransfinite()

lin_alg::Matrix<Point> meshkernel::DiscretizeTransfinite ( const std::vector< Point > &  leftDiscretization,
const std::vector< Point > &  rightDiscretization,
const std::vector< Point > &  bottomDiscretization,
const std::vector< Point > &  upperDiscretization,
const Projection projection,
UInt  numM,
UInt  numN 
)

Computes the transfinite discretization inside the area defined by 4 sides, each one discretized with a series of points (tranfn2).

Parameters
[in]leftDiscretizationThe discretization of the left side.
[in]rightDiscretizationThe discretization of the right side.
[in]bottomDiscretizationThe discretization of the bottom.
[in]upperDiscretizationThe discretization of the upper side.
[in]projectionThe projection to use.
[in]numMThe number of columns to generate.
[in]numNThe number of rows to generate.
Returns
The resulting dicretization (expressed as number of points).

◆ DistanceFromLine()

std::tuple<double, Point, double> meshkernel::DistanceFromLine ( const Point point,
const Point firstNode,
const Point secondNode,
const Projection projection 
)

Computes the perpendicular distance of a point from a segment firstNode - secondNode (dlinedis3)

Parameters
[in]pointThe point to consider in the distance calculation.
[in]firstNodeThe first point of the segment.
[in]secondNodeThe second point of the segment.
[in]projectionThe coordinate system projection.
Returns
The normal distance from the segment, the intersection of the normal projection on the segment, the distance from the first node, expressed as ratio of the segment length

◆ dot() [1/3]

double meshkernel::dot ( const Point p,
const Vector v 
)
inline

Compute the dot product of a point and a vector.

This is mainly a convenience function

◆ dot() [2/3]

double meshkernel::dot ( const Point p1,
const Point p2 
)
inline

Compute the dot product of two points.

This is mainly a convenience function

◆ dot() [3/3]

double meshkernel::dot ( const Vector v,
const Point p 
)
inline

Compute the dot product of a point and a vector.

This is mainly a convenience function

◆ FindFunctionRootWithGoldenSectionSearch()

template<typename F >
double meshkernel::FindFunctionRootWithGoldenSectionSearch ( func,
double  min,
double  max 
)

Algorithm performing the zero's search using the golden section algorithm's.

Parameters
[in]funcFunction to search for a root
[in]minThe minimum value of the interval
[in]maxThe maximum value of the interval
Returns
The value where the function is approximately 0

◆ FindIndex()

template<typename T >
UInt meshkernel::FindIndex ( const std::vector< T > &  vec,
el 
)

Find index of a certain element.

Parameters
[in]vecThe vector to search in
[in]elThe element to search for
Returns
The index of element

◆ FindIndices()

std::vector<std::pair<UInt, UInt> > meshkernel::FindIndices ( const std::vector< Point > &  vec,
size_t  start,
size_t  end,
double  separator 
)

Find all start-end positions in a vector separated by a separator.

Parameters
[in]vecThe vector with separator
[in]startThe start of the range to search for
[in]endThe end of the range to search for
[in]separatorThe value of the separator
Returns
Indices of elements

◆ GetCurvilinearDirectionValue()

CurvilinearDirection meshkernel::GetCurvilinearDirectionValue ( int  direction)

Convert an integer value to the CurvilinearDirection enumeration type.

If the integer direction value does not correspond to an enumeration value then a ConstraintError will be thrown

◆ GetDelta()

Vector meshkernel::GetDelta ( const Point firstPoint,
const Point secondPoint,
const Projection projection 
)

Get the delta (dx, dy) for the given projection.

Parameters
[in]firstPoint
[in]secondPoint
[in]projectionThe coordinate system projection.

◆ GetDx()

double meshkernel::GetDx ( const Point firstPoint,
const Point secondPoint,
const Projection projection 
)

Gets dx for the given projection.

Parameters
[in]firstPoint
[in]secondPoint
[in]projectionThe coordinate system projection.

◆ GetDy()

double meshkernel::GetDy ( const Point firstPoint,
const Point secondPoint,
const Projection projection 
)

Gets dy for the given projection.

Parameters
[in]firstPoint
[in]secondPoint
[in]projectionThe coordinate system projection.

◆ GetProjectionValue()

Projection meshkernel::GetProjectionValue ( int  projection)

Convert an integer value to the Projection enumeration type.

If the integer projection value does not correspond to an enumeration value then a ConstraintError will be thrown

◆ InnerProduct()

double meshkernel::InnerProduct ( const Cartesian3DPoint a,
const Cartesian3DPoint b 
)

Defines inner product in cartesian 3D-space.

Parameters
[in]aThe first cartesian 3D point
[in]bThe second cartesian 3D point
Returns
The resulting inner product

◆ InnerProductTwoSegments()

double meshkernel::InnerProductTwoSegments ( const Point firstPointFirstSegment,
const Point secondPointFirstSegment,
const Point firstPointSecondSegment,
const Point secondPointSecondSegment,
const Projection projection 
)

Inner product of two segments (dprodin)

Parameters
[in]firstPointFirstSegmentThe first point of the first segment
[in]secondPointFirstSegmentThe second point of the first segment
[in]firstPointSecondSegmentThe first point of the second segment
[in]secondPointSecondSegmentThe second point of the second segment
[in]projectionThe coordinate system projection
Returns
The resulting inner product

◆ InvalidPointCount()

UInt meshkernel::InvalidPointCount ( const std::vector< Point > &  points,
size_t  start,
size_t  end 
)

Determine if the number of invalid points in the section of the point array.

Parameters
[in]pointsThe array of points
[in]startThe start of the range in which to search
[in]endThe end of the range in which to search

◆ IsEqual()

bool meshkernel::IsEqual ( const Point p1,
const Point p2,
const double  epsilon 
)
inline

Test points for equality upto a tolerance.

Parameters
[in]p1First point to compare
[in]p2Second point to compare
[in]epsilonRelative tolerance to which the values are compared.
Returns
Boolean value indicating where the points p1 and p2 are equal upto a relative tolerance.

◆ IsPointInPolygonNodes()

bool meshkernel::IsPointInPolygonNodes ( const Point point,
const std::vector< Point > &  polygonNodes,
const Projection projection,
Point  polygonCenter = {constants::missing::doubleValue, constants::missing::doubleValue},
UInt  startNode = constants::missing::uintValue,
UInt  endNode = constants::missing::uintValue 
)

Checks if a point is in polygonNodes using the winding number method.

Parameters
[in]pointThe point to check
[in]polygonNodesA series of closed polygons
[in]startNodeThe start index in polygonNodes
[in]endNodeThe end index in polygonNodes
[in]projectionThe coordinate system projection.
[in]polygonCenterA coordinate needed in case of sphericalAccurate projection
Returns
If point is inside the designed polygon

◆ IsPointOnPole()

bool meshkernel::IsPointOnPole ( const Point point)

Determines if a point is close to the poles (latitude close to 90 degrees).

Parameters
[in]pointThe current point.
Returns
If the point is on the pole.

◆ lengthSquared()

double meshkernel::lengthSquared ( const Point p1)
inline

Compute the dot product of a point with itself.

This is mainly a convenience function

◆ LinearInterpolationInTriangle()

double meshkernel::LinearInterpolationInTriangle ( const Point interpolationPoint,
const std::vector< Point > &  polygon,
const std::vector< double > &  values,
const Projection projection 
)

Given a triangles with values on each node, computes the interpolated value inside the triangle, using linear interpolation.

Parameters
[in]interpolationPointThe point where to interpolate.
[in]polygonThe polygon containing the triangle nodes.
[in]valuesThe values at each node.
[in]projectionThe projection to use.
Returns
The interpolated value.

◆ MatrixNorm()

double meshkernel::MatrixNorm ( const std::vector< double > &  x,
const std::vector< double > &  y,
const std::vector< double > &  matCoefficients 
)

Compute the matrix norm.

x [in] To be detailed

y [in] To be detailed

matCoefficients [in] To be detailed

Returns
The computed matrix norm

◆ NextCircularBackwardIndex()

UInt meshkernel::NextCircularBackwardIndex ( UInt  currentIndex,
UInt  size 
)

Get the next backward index, for a range in 0 .. size-1.

Parameters
[in]currentIndexThe current index.
[in]sizeThe size of the vector.
Returns
The next backward index.

◆ NextCircularForwardIndex()

UInt meshkernel::NextCircularForwardIndex ( UInt  currentIndex,
UInt  size 
)

Get the next forward index, for a range in 0 .. size-1.

Parameters
[in]currentIndexThe current index.
[in]sizeThe size of the vector.
Returns
The next forward index.

◆ normalise()

meshkernel::Vector meshkernel::normalise ( const Vector vec)
inline

Return the normalised vector.

Note
No check is made that the length is 0

◆ NormalizedInnerProductTwoSegments()

double meshkernel::NormalizedInnerProductTwoSegments ( const Point firstPointFirstSegment,
const Point secondPointFirstSegment,
const Point firstPointSecondSegment,
const Point secondPointSecondSegment,
const Projection projection 
)

The normalized inner product of two segments (dcosphi)

Parameters
[in]firstPointFirstSegmentThe first point of the first segment
[in]secondPointFirstSegmentThe second point of the first segment
[in]firstPointSecondSegmentThe first point of the second segment
[in]secondPointSecondSegmentThe second point of the second segment
[in]projectionThe coordinate system projection
Returns
The resulting normalized inner product

◆ NormalVector()

Point meshkernel::NormalVector ( const Point firstPoint,
const Point secondPoint,
const Point insidePoint,
const Projection projection 
)

Normalized vector of a segment in direction 1 -> 2 with the insidePoint orientation.

Parameters
[in]firstPointThe first point of the segment.
[in]secondPointThe second point of the segment.
[in]insidePointThe inside point of the segment
[in]projectionThe coordinate system projection.
Returns
The Normal vector

◆ NormalVectorInside()

void meshkernel::NormalVectorInside ( const Point firstPoint,
const Point secondPoint,
const Point insidePoint,
Point normal,
bool &  flippedNormal,
const Projection projection 
)

Computes the normal vector to a line 1-2, which is outward w.r.t. an 'inside' point 3.

Similar to NormalVectorOutside, except that the normal vector may be flipped based on the 'inside' point.

◆ NormalVectorOutside()

Point meshkernel::NormalVectorOutside ( const Point firstPoint,
const Point secondPoint,
const Projection projection 
)

Computes the normal vector outside (normalout)

See also
NormalVectorInside

◆ operator!=()

bool meshkernel::operator!= ( const Point p1,
const Point p2 
)
inline

Test points for equality using a default tolerance.

Returns
\( p1.x = p2.x \wedge p1.y = p2.y)\)

◆ operator*() [1/3]

meshkernel::Point meshkernel::operator* ( const double  x,
const Point p 
)
inline

Multiply point by a scalar.

Returns
\( (p.x * x, p.y * x)\)

◆ operator*() [2/3]

meshkernel::Point meshkernel::operator* ( const Point p,
const double  x 
)
inline

Multiply point by a scalar.

Returns
\( (p.x * x, p.y * x)\)

◆ operator*() [3/3]

meshkernel::Point meshkernel::operator* ( const Point p1,
const Point p2 
)
inline

Multiply point by a point.

Computes a point-wise product.

Returns
\( (p1.x * p2.x, p1.y * p2.y)\)

◆ operator+() [1/4]

meshkernel::Point meshkernel::operator+ ( const double  x,
const Point p 
)
inline

Add points and scalar.

Returns
\( (p.x + x, p.y + x)\)

◆ operator+() [2/4]

meshkernel::Point meshkernel::operator+ ( const Point p,
const double  x 
)
inline

Add points and scalar.

Returns
\( (p.x + x, p.y + x)\)

◆ operator+() [3/4]

meshkernel::Point meshkernel::operator+ ( const Point p1,
const Point p2 
)
inline

Add two points.

Returns
\( (p1.x + p2.x, p1.y + p2.y)\)

◆ operator+() [4/4]

meshkernel::Point meshkernel::operator+ ( const Point pnt,
const Vector vec 
)
inline

Add vector to a point.

Returns
\( (p.x + v.x, p.y + v.y)\)

◆ operator-() [1/6]

meshkernel::Point meshkernel::operator- ( const double  x,
const Point p 
)
inline

Subtract point from scalar.

Returns
\( (x - p.x, x - p.y)\)

◆ operator-() [2/6]

meshkernel::Point meshkernel::operator- ( const Point p,
const double  x 
)
inline

Subtract scalar from point.

Returns
\( (p.x - x, p.y - x)\)

◆ operator-() [3/6]

meshkernel::Point meshkernel::operator- ( const Point p1,
const Point p2 
)
inline

Subtract two points.

Returns
\( (p1.x - p2.x, p1.y - p2.y)\)

◆ operator-() [4/6]

meshkernel::Point meshkernel::operator- ( const Point pnt)
inline

Unary minus.

Returns
\( (-p1.x, -p1.y)\)

◆ operator-() [5/6]

meshkernel::Point meshkernel::operator- ( const Point pnt,
const Vector vec 
)
inline

Subtract vector from a point.

Returns
\( (p.x - v.x, p.y - v.y)\)

◆ operator-() [6/6]

meshkernel::Vector meshkernel::operator- ( const Vector vec)
inline

Unary minus.

Returns
\( (-vec.x, -vec.y)\)

◆ operator/() [1/2]

meshkernel::Point meshkernel::operator/ ( const Point p,
const double  x 
)
inline

Divide point by a scalar.

Note
No check is made to determine is divisor is zero.
Returns
\( (p.x / x, p.y / x)\)

◆ operator/() [2/2]

meshkernel::Point meshkernel::operator/ ( const Point p1,
const Point p2 
)
inline

Divide point by a point.

Computes a point-wise division.

Returns
\( (p1.x / p2.x, p1.y / p2.y)\)
Note
No check is made to determine is divisors are zero.

◆ operator==()

bool meshkernel::operator== ( const Point p1,
const Point p2 
)
inline

Test points for equality using a default tolerance.

Returns
\( p1.x = p2.x \wedge p1.y = p2.y)\)

◆ OrthogonalProjectionOnSegment()

std::tuple<Point, double, bool> meshkernel::OrthogonalProjectionOnSegment ( Point const &  firstNode,
Point const &  secondNode,
Point const &  pointToProject 
)

Cartesian projection of a point on a segment defined by other two points.

Parameters
firstNodeThe first node of the segment
secondNodeThe second node of the segment
pointToProjectThe point to project
Returns
The projected point, the distance of the projection from the first point (expressed as a ratio), a boolean indicating if the projected point is within the first and second point.

◆ PointAlongLine()

meshkernel::Point meshkernel::PointAlongLine ( const Point startPoint,
const Point endPoint,
const double  lambda 
)
inline

Compute the point at some position along the line connecting start- and end-point.

\( r = (1-\lambda) s + \lambda e \)

◆ ReferencePoint() [1/2]

Point meshkernel::ReferencePoint ( const std::vector< Point > &  nodes,
const std::vector< UInt > &  polygonIndices,
const Projection projection 
)

For a given polygon compute a reference point.

Parameters
[in]nodesThe input nodes
[in]polygonIndicesThe polygon node indices.
[in]projectionThe coordinate system projection.
Returns
The reference point

◆ ReferencePoint() [2/2]

Point meshkernel::ReferencePoint ( const std::vector< Point > &  polygon,
const Projection projection 
)

For a given polygon compute a reference point.

Parameters
[in]polygonThe input polygon.
[in]projectionThe coordinate system projection.
Returns
The reference point

◆ ReorderVector()

template<typename T >
auto meshkernel::ReorderVector ( const std::vector< T > &  v,
const std::vector< UInt > &  order 
)

Reorder vector accordingly to a specific order.

Parameters
[in]vThe vector to reorder
[in]orderThe order to use
Returns
The reordered vector

◆ ResizeAndFill2DVector()

template<typename T >
void meshkernel::ResizeAndFill2DVector ( std::vector< std::vector< T >> &  v,
UInt const &  firstDimension,
UInt const &  secondDimension,
bool  fill = false,
const T &  fillValue = {} 
)

Resizes and fills a two dimensional vector.

Template Parameters
TThe type of the vector elements
Parameters
[in]vThe input two dimensional vector
[in]firstDimensionThe first new dimension
[in]secondDimensionThe second new dimension
[in]fillWhatever fill or not fill the vector with missing values
[in]fillValueThe fill value

◆ ResizeAndFill3DVector()

template<typename T >
void meshkernel::ResizeAndFill3DVector ( std::vector< std::vector< std::vector< T >>> &  v,
UInt const &  firstDimension,
UInt const &  secondDimension,
UInt const &  thirdDim,
bool  fill = false,
const T &  fillValue = {} 
)

Resizes and fills a three dimensional vector.

Template Parameters
TThe type of the vector elements
Parameters
[in]vThe input three dimensional vector
[in]firstDimensionThe first new dimension
[in]secondDimensionThe second new dimension
[in]thirdDimThe third new dimension
[in]fillWhatever fill or not fill the vector with missing values
[in]fillValueThe fill value

◆ sgn()

template<typename T >
int meshkernel::sgn ( val)

Computes the sign of a type.

Template Parameters
TA signed type
Parameters
[in]valthe value to use for computing a sign
Returns
-1 for negatives and +1 for positives

◆ SortedIndices()

template<typename T >
std::vector<UInt> meshkernel::SortedIndices ( const std::vector< T > &  v)

Sort a vector and return the sorted indices.

Parameters
[in]vThe vector to sort
Returns
The indices of elements

◆ SphericalToCartesian3D()

Cartesian3DPoint meshkernel::SphericalToCartesian3D ( const Point sphericalPoint)

Transforms 2D point in spherical coordinates to 3D cartesian coordinates.

Parameters
[in]sphericalPointThe current spherical point (2 coordinates).
Returns
The converted cartesian 3d point.

◆ TranslateSphericalCoordinates()

void meshkernel::TranslateSphericalCoordinates ( std::vector< Point > &  polygon)

For a given polygon the function may shift the input coordinates.

Parameters
[in,out]polygonThe input polygon.
Note
To be called for spherical coordinate system only

◆ Triangulation()

void meshkernel::Triangulation ( int  jatri,
double const *const  xs,
double const *const  ys,
int  ns,
int *const  indx,
int *const  numtri,
int *const  edgeidx,
int *const  numedge,
int *const  triedge,
double *const  xs3,
double *const  ys3,
int *const  ns3,
double  trisize 
)

Function of the Triangle library.

See also
https://www.cs.cmu.edu/~quake/triangle.html

◆ VectorProduct()

Cartesian3DPoint meshkernel::VectorProduct ( const Cartesian3DPoint a,
const Cartesian3DPoint b 
)

Defines vector product for cartesian 3D-space.

Parameters
[in]aThe first cartesian 3D point
[in]bThe second cartesian 3D point
Returns
The vector product

Variable Documentation

◆ HasConversionProjection

template<typename Operation >
concept meshkernel::HasConversionProjection
Initial value:
= requires(Operation op) {{ op.SourceProjection()} -> std::same_as<Projection>;
{ op.TargetProjection()} -> std::same_as<Projection>; }

Ensure any instantiation of the MeshConversion Compute function is able to determine source and target projection types.

◆ InterpolatableType

template<typename T >
concept meshkernel::InterpolatableType
Initial value:
= std::same_as<T, short> ||
std::same_as<T, float> ||
std::same_as<T, double> ||
std::same_as<T, int>

Defines the iterpolatable data types.