MeshKernel
|
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< UInt > | SortedIndices (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 > | |
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< Point > | 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). More... | |
std::vector< Point > | ComputeEdgeCenters (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< Point > | ComputePolyLineDiscretization (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. | |
Contains the logic of the C++ static library.
|
strong |
enum meshkernel::ExitCode |
Enumeration of exit codes.
|
strong |
|
strong |
Mesh locations enumeration.
Enumerator | |
---|---|
Faces | Faces. |
Nodes | Nodes. |
Edges | Edges. |
Unknown | Unknown. |
|
strong |
Calculates the absolute difference between to Index
numbers.
[in] | number_1 | The first number |
[in] | number_2 | The second number |
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.
[in] | normal | The increment direction. |
[in] | increment | The increment to use. |
[in] | referencePoint | The reference point containing the reference latitude to use. |
[in] | projection | The coordinate system projection. |
[in,out] | point | The point to be incremented. |
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)
[in] | firstSegmentFirstPoint | The first point of the first segment |
[in] | firstSegmentSecondPoint | The second point of the first segment |
[in] | secondSegmentFirstPoint | The first point of the second segment |
[in] | secondSegmentSecondPoint | The second point of the second segment |
[in] | adimensionalCrossProduct | Whether to compute the dimensionless cross product |
[in] | projection | The coordinate system projection |
Point meshkernel::Cartesian3DToSpherical | ( | const Cartesian3DPoint & | cartesianPoint, |
double | referenceLongitude | ||
) |
Transforms 3D cartesian coordinates to 2D point in spherical coordinates.
[in] | cartesianPoint | The 3d cartesian point |
[in] | referenceLongitude | The reference longitude |
Point meshkernel::CircumcenterOfTriangle | ( | const Point & | firstNode, |
const Point & | secondNode, | ||
const Point & | thirdNode, | ||
const Projection & | projection | ||
) |
Computes the circumcenter of a triangle.
[in] | firstNode | The first triangle node |
[in] | secondNode | The second triangle node |
[in] | thirdNode | The third triangle node |
[in] | projection | The coordinate system projection |
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.
[in] | v | The vector of points |
[in] | projection | The projection to use. |
Point meshkernel::ComputeAverageCoordinate | ( | const std::vector< Point > & | points, |
const Projection & | projection | ||
) |
Given a series of point computes the average coordinate.
[in] | points | The point series. |
[in] | projection | The projection to use. |
double meshkernel::ComputeDistance | ( | const Point & | firstPoint, |
const Point & | secondPoint, | ||
const Projection & | projection | ||
) |
Computes the distance between two points (dbdistance)
[in] | firstPoint | The first point. |
[in] | secondPoint | The second point. |
[in] | projection | The coordinate system projection. |
std::vector<Point> meshkernel::ComputeEdgeCenters | ( | const std::vector< Point > & | nodes, |
const std::vector< Edge > & | edges | ||
) |
Computes the edge centers.
[in] | nodes | The vector of edge nodes. |
[in] | edges | The vector of edge indices. |
Point meshkernel::ComputeMiddlePoint | ( | const Point & | firstPoint, |
const Point & | secondPoint, | ||
const Projection & | projection | ||
) |
Computes the middle point.
[in] | firstPoint | The first point of the segment. |
[in] | secondPoint | The second point of the segment. |
[in] | projection | The coordinate system projection. |
Point meshkernel::ComputeMiddlePointAccountingForPoles | ( | const Point & | firstPoint, |
const Point & | secondPoint, | ||
const Projection & | projection | ||
) |
Computes the middle point (account for poles, latitudes close to 90 degrees)
[in] | firstPoint | The first point of the segment. |
[in] | secondPoint | The second point of the segment. |
[in] | projection | The coordinate system projection. |
Vector meshkernel::ComputeNormalToline | ( | const Point & | start, |
const Point & | end, | ||
const Projection | projection | ||
) |
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)
[in] | coordinates | The spline node coordinates |
[in] | coordinatesDerivatives | The spline nodal derivatives |
[in] | pointAdimensionalCoordinate | The adimensinal coordinate where to perform the interpolation |
std::vector<Point> meshkernel::ComputePolyLineDiscretization | ( | std::vector< Point > const & | polyline, |
std::vector< double > & | chainages, | ||
Projection | projection | ||
) |
Computes the discretization points along a polyline.
polyline | A polyline described by its nodes |
chainages | The chainages used for dicretizing the current polyline |
projection | The projection to use |
std::vector<double> meshkernel::ComputePolyLineEdgesLengths | ( | std::vector< Point > const & | polyline, |
Projection | projection | ||
) |
Computes the lengths of each polyline segment.
polyline | A polyline described by its nodes |
projection | The projection to use |
std::vector<double> meshkernel::ComputePolyLineNodalChainages | ( | std::vector< Point > const & | polyline, |
Projection | projection | ||
) |
Computes the chainages of each polyline node.
polyline | A polyline described by its nodes |
projection | The projection to use |
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.
[in] | firstPoint | The first point. |
[in] | secondPoint | The second point. |
[in] | projection | The coordinate system projection. |
Compute the cross product of two points (as vectors).
This is mainly a convenience function
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)
[in] | firstSegmentFirstPoint | The first point of the first segment |
[in] | firstSegmentSecondPoint | The second point of the first segment |
[in] | secondSegmentFirstPoint | The second point of the second segment |
[in] | secondSegmentSecondPoint | The second point of the second segment |
[in] | projection | The coordinate system projection |
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).
[in] | leftDiscretization | The discretization of the left side. |
[in] | rightDiscretization | The discretization of the right side. |
[in] | bottomDiscretization | The discretization of the bottom. |
[in] | upperDiscretization | The discretization of the upper side. |
[in] | projection | The projection to use. |
[in] | numM | The number of columns to generate. |
[in] | numN | The number of rows to generate. |
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)
[in] | point | The point to consider in the distance calculation. |
[in] | firstNode | The first point of the segment. |
[in] | secondNode | The second point of the segment. |
[in] | projection | The coordinate system projection. |
Compute the dot product of a point and a vector.
This is mainly a convenience function
Compute the dot product of two points.
This is mainly a convenience function
Compute the dot product of a point and a vector.
This is mainly a convenience function
double meshkernel::FindFunctionRootWithGoldenSectionSearch | ( | F | func, |
double | min, | ||
double | max | ||
) |
Algorithm performing the zero's search using the golden section algorithm's.
[in] | func | Function to search for a root |
[in] | min | The minimum value of the interval |
[in] | max | The maximum value of the interval |
UInt meshkernel::FindIndex | ( | const std::vector< T > & | vec, |
T | el | ||
) |
Find index of a certain element.
[in] | vec | The vector to search in |
[in] | el | The element to search for |
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.
[in] | vec | The vector with separator |
[in] | start | The start of the range to search for |
[in] | end | The end of the range to search for |
[in] | separator | The value of the separator |
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
Vector meshkernel::GetDelta | ( | const Point & | firstPoint, |
const Point & | secondPoint, | ||
const Projection & | projection | ||
) |
Get the delta (dx, dy) for the given projection.
[in] | firstPoint | |
[in] | secondPoint | |
[in] | projection | The coordinate system projection. |
double meshkernel::GetDx | ( | const Point & | firstPoint, |
const Point & | secondPoint, | ||
const Projection & | projection | ||
) |
Gets dx for the given projection.
[in] | firstPoint | |
[in] | secondPoint | |
[in] | projection | The coordinate system projection. |
double meshkernel::GetDy | ( | const Point & | firstPoint, |
const Point & | secondPoint, | ||
const Projection & | projection | ||
) |
Gets dy for the given projection.
[in] | firstPoint | |
[in] | secondPoint | |
[in] | projection | The coordinate system projection. |
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
double meshkernel::InnerProduct | ( | const Cartesian3DPoint & | a, |
const Cartesian3DPoint & | b | ||
) |
Defines inner product in cartesian 3D-space.
[in] | a | The first cartesian 3D point |
[in] | b | The second cartesian 3D point |
double meshkernel::InnerProductTwoSegments | ( | const Point & | firstPointFirstSegment, |
const Point & | secondPointFirstSegment, | ||
const Point & | firstPointSecondSegment, | ||
const Point & | secondPointSecondSegment, | ||
const Projection & | projection | ||
) |
Inner product of two segments (dprodin)
[in] | firstPointFirstSegment | The first point of the first segment |
[in] | secondPointFirstSegment | The second point of the first segment |
[in] | firstPointSecondSegment | The first point of the second segment |
[in] | secondPointSecondSegment | The second point of the second segment |
[in] | projection | The coordinate system projection |
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.
[in] | points | The array of points |
[in] | start | The start of the range in which to search |
[in] | end | The end of the range in which to search |
Test points for equality upto a tolerance.
[in] | p1 | First point to compare |
[in] | p2 | Second point to compare |
[in] | epsilon | Relative tolerance to which the values are compared. |
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.
[in] | point | The point to check |
[in] | polygonNodes | A series of closed polygons |
[in] | startNode | The start index in polygonNodes |
[in] | endNode | The end index in polygonNodes |
[in] | projection | The coordinate system projection. |
[in] | polygonCenter | A coordinate needed in case of sphericalAccurate projection |
bool meshkernel::IsPointOnPole | ( | const Point & | point | ) |
Determines if a point is close to the poles (latitude close to 90 degrees).
[in] | point | The current point. |
|
inline |
Compute the dot product of a point with itself.
This is mainly a convenience function
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.
[in] | interpolationPoint | The point where to interpolate. |
[in] | polygon | The polygon containing the triangle nodes. |
[in] | values | The values at each node. |
[in] | projection | The projection to use. |
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
Get the next backward index, for a range in 0 .. size-1.
[in] | currentIndex | The current index. |
[in] | size | The size of the vector. |
Get the next forward index, for a range in 0 .. size-1.
[in] | currentIndex | The current index. |
[in] | size | The size of the vector. |
|
inline |
Return the normalised vector.
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)
[in] | firstPointFirstSegment | The first point of the first segment |
[in] | secondPointFirstSegment | The second point of the first segment |
[in] | firstPointSecondSegment | The first point of the second segment |
[in] | secondPointSecondSegment | The second point of the second segment |
[in] | projection | The coordinate system projection |
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.
[in] | firstPoint | The first point of the segment. |
[in] | secondPoint | The second point of the segment. |
[in] | insidePoint | The inside point of the segment |
[in] | projection | The coordinate system projection. |
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.
Point meshkernel::NormalVectorOutside | ( | const Point & | firstPoint, |
const Point & | secondPoint, | ||
const Projection & | projection | ||
) |
Computes the normal vector outside (normalout)
Test points for equality using a default tolerance.
|
inline |
Multiply point by a scalar.
|
inline |
Multiply point by a scalar.
|
inline |
Multiply point by a point.
Computes a point-wise product.
|
inline |
Add points and scalar.
|
inline |
Add points and scalar.
|
inline |
Add two points.
|
inline |
Add vector to a point.
|
inline |
Subtract point from scalar.
|
inline |
Subtract scalar from point.
|
inline |
Subtract two points.
|
inline |
Unary minus.
|
inline |
Subtract vector from a point.
|
inline |
Unary minus.
|
inline |
Divide point by a scalar.
|
inline |
Divide point by a point.
Computes a point-wise division.
Test points for equality using a default tolerance.
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.
firstNode | The first node of the segment |
secondNode | The second node of the segment |
pointToProject | The point to project |
|
inline |
Compute the point at some position along the line connecting start- and end-point.
\( r = (1-\lambda) s + \lambda e \)
Point meshkernel::ReferencePoint | ( | const std::vector< Point > & | nodes, |
const std::vector< UInt > & | polygonIndices, | ||
const Projection & | projection | ||
) |
For a given polygon compute a reference point.
[in] | nodes | The input nodes |
[in] | polygonIndices | The polygon node indices. |
[in] | projection | The coordinate system projection. |
Point meshkernel::ReferencePoint | ( | const std::vector< Point > & | polygon, |
const Projection & | projection | ||
) |
For a given polygon compute a reference point.
[in] | polygon | The input polygon. |
[in] | projection | The coordinate system projection. |
auto meshkernel::ReorderVector | ( | const std::vector< T > & | v, |
const std::vector< UInt > & | order | ||
) |
Reorder vector accordingly to a specific order.
[in] | v | The vector to reorder |
[in] | order | The order to use |
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.
T | The type of the vector elements |
[in] | v | The input two dimensional vector |
[in] | firstDimension | The first new dimension |
[in] | secondDimension | The second new dimension |
[in] | fill | Whatever fill or not fill the vector with missing values |
[in] | fillValue | The fill value |
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.
T | The type of the vector elements |
[in] | v | The input three dimensional vector |
[in] | firstDimension | The first new dimension |
[in] | secondDimension | The second new dimension |
[in] | thirdDim | The third new dimension |
[in] | fill | Whatever fill or not fill the vector with missing values |
[in] | fillValue | The fill value |
int meshkernel::sgn | ( | T | val | ) |
Computes the sign of a type.
T | A signed type |
[in] | val | the value to use for computing a sign |
std::vector<UInt> meshkernel::SortedIndices | ( | const std::vector< T > & | v | ) |
Sort a vector and return the sorted indices.
[in] | v | The vector to sort |
Cartesian3DPoint meshkernel::SphericalToCartesian3D | ( | const Point & | sphericalPoint | ) |
Transforms 2D point in spherical coordinates to 3D cartesian coordinates.
[in] | sphericalPoint | The current spherical point (2 coordinates). |
void meshkernel::TranslateSphericalCoordinates | ( | std::vector< Point > & | polygon | ) |
For a given polygon the function may shift the input coordinates.
[in,out] | polygon | The input polygon. |
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.
Cartesian3DPoint meshkernel::VectorProduct | ( | const Cartesian3DPoint & | a, |
const Cartesian3DPoint & | b | ||
) |
Defines vector product for cartesian 3D-space.
[in] | a | The first cartesian 3D point |
[in] | b | The second cartesian 3D point |
concept meshkernel::HasConversionProjection |
Ensure any instantiation of the MeshConversion Compute function is able to determine source and target projection types.
concept meshkernel::InterpolatableType |
Defines the iterpolatable data types.