30#include <MeshKernel/Parameters.hpp>
32#include <MeshKernelApi/Contacts.hpp>
33#include <MeshKernelApi/CurvilinearGrid.hpp>
34#include <MeshKernelApi/GeometryList.hpp>
35#include <MeshKernelApi/GriddedSamples.hpp>
36#include <MeshKernelApi/Mesh1D.hpp>
37#include <MeshKernelApi/Mesh2D.hpp>
40#if !defined(MKERNEL_API)
41#define MKERNEL_API __declspec(dllexport)
44#define MKERNEL_API __attribute__((visibility("default")))
131 const int* oneDNodeMask,
133 double searchRadius);
142 const int* oneDNodeMask);
153 const int* oneDNodeMask,
155 double projectionFactor);
165 const int* oneDNodeMask,
176 const int* oneDNodeMask,
291 double xPointCoordinate,
292 double yPointCoordinate);
404 double repulsionParameter,
405 double xFirstNodeOnTheLine,
406 double yFirstNodeOnTheLine,
407 double xSecondNodeOnTheLine,
408 double ySecondNodeOnTheLine,
409 double xLowerLeftCorner,
410 double yLowerLeftCorner,
411 double xUpperRightCorner,
412 double yUpperRightCorner);
424 double mirroringFactor,
426 double xFirstGridLineNode,
427 double yFirstGridLineNode,
428 double xSecondGridLineNode,
429 double ySecondGridLineNode);
487 double xFromCoordinate,
488 double yFromCoordinate,
489 double xToCoordinate,
490 double yToCoordinate);
503 double xLowerLeftCorner,
504 double yLowerLeftCorner,
505 double xUpperRightCorner,
506 double yUpperRightCorner,
534 double xLowerLeftCorner,
535 double yLowerLeftCorner,
536 double xUpperRightCorner,
537 double yUpperRightCorner);
553 double xLowerLeftCorner,
554 double yLowerLeftCorner,
555 double xUpperRightCorner,
556 double yUpperRightCorner);
581 double xFirstGridLineNode,
582 double yFirstGridLineNode,
583 double xSecondGridLineNode,
584 double ySecondGridLineNode,
597 double& xFirstFrozenLineCoordinate,
598 double& yFirstFrozenLineCoordinate,
599 double& xSecondFrozenLineCoordinate,
600 double& ySecondFrozenLineCoordinate);
607 int& numFrozenLines);
614 int* frozenLinesIds);
624 double xFirstGridLineNode,
625 double yFirstGridLineNode,
626 double xSecondGridLineNode,
627 double ySecondGridLineNode);
638 int smoothingIterations,
639 double xLowerLeftCorner,
640 double yLowerLeftCorner,
641 double xUpperRightCorner,
642 double yUpperRightCorner);
657 int smoothingIterations,
658 double xFirstGridlineNode,
659 double yFirstGridlineNode,
660 double xSecondGridLineNode,
661 double ySecondGridLineNode,
662 double xLowerLeftCornerSmoothingArea,
663 double yLowerLeftCornerSmoothingArea,
664 double xUpperRightCornerSmootingArea,
665 double yUpperRightCornerSmootingArea);
679 double sectionControlPoint1x,
680 double sectionControlPoint1y,
681 double sectionControlPoint2x,
682 double sectionControlPoint2y,
698 double sectionControlPoint1x,
699 double sectionControlPoint1y,
700 double sectionControlPoint2x,
701 double sectionControlPoint2y,
875 int numberOfPointsBetweenNodes);
927 int averagingMethodType,
928 double relativeSearchSize,
929 size_t minNumSamples,
970 int startSplineIndex,
1027 const double minimumRefinementDepth);
1039 int projectToLandBoundaryOption,
1094 int& numberOfMeshNodes);
1108 double smallFlowEdgesLengthThreshold,
1109 int& numSmallFlowEdges);
1120 int invertDeletion);
1135 double xLowerLeftBoundingBox,
1136 double yLowerLeftBoundingBox,
1137 double xUpperRightBoundingBox,
1138 double yUpperRightBoundingBox);
1170 double smallFlowEdgesThreshold,
1171 double minFractionalAreaTriangles);
1188 int isTriangulationRequired,
1189 int projectToLandBoundaryRequired,
1206 double xCoordinateIn,
1207 double yCoordinateIn,
1208 double searchRadius,
1209 double xLowerLeftBoundingBox,
1210 double yLowerLeftBoundingBox,
1211 double xUpperRightBoundingBox,
1212 double yUpperRightBoundingBox,
1213 double& xCoordinateOut,
1214 double& yCoordinateOut);
1263 double xLowerLeftBoundingBox,
1264 double yLowerLeftBoundingBox,
1265 double xUpperRightBoundingBox,
1266 double yUpperRightBoundingBox,
1309 int& geometryListDimension);
1339 int& locationIndex);
1364 double searchRadius,
1365 double xLowerLeftBoundingBox,
1366 double yLowerLeftBoundingBox,
1367 double xUpperRightBoundingBox,
1368 double yUpperRightBoundingBox,
1380 int* selectedNodes);
1417 double smallFlowEdgesThreshold,
1437 int projectToLandBoundaryOption,
1465 int& firstNodeIndex,
1466 int& secondNodeIndex,
1493 double* edgeDistances,
1494 double* segmentDistances,
1495 int* segmentIndexes,
1498 int* faceEdgeIndex);
1598 bool useNodalRefinement);
1622 double relativeSearchRadius,
1623 int minimumNumSamples,
1640 double relativeSearchRadius,
1641 int minimumNumSamples,
1642 int numberOfSmoothingIterations,
1704 double* fixedChainages,
1705 int sizeFixedChainages,
1707 double fixedChainagesOffset);
1740 int& numberOfPolygonNodes);
1758 int& numberOfPolygonNodes);
1774 int& numberOfPolygonNodes);
1817 int secondNodeIndex,
1818 double targetEdgeLength,
1833 int secondNodeIndex,
Contains all structs and functions exposed at the API level.
Definition BoundingBox.hpp:33
MKERNEL_API int mkernel_mesh2d_count_hanging_edges(int meshKernelId, int &numEdges)
Count the number of hanging edges in a mesh2d. An hanging edge is an edge where one of the two nodes ...
MKERNEL_API int mkernel_expunge_state(int meshKernelId)
Deallocate mesh state and remove it completely, no undo for this meshKernelId will be possible after ...
MKERNEL_API int mkernel_mesh2d_merge_nodes(int meshKernelId, const GeometryList &geometryListIn)
Merges the mesh2d nodes within a distance of 0.001 m, effectively removing all small edges.
MKERNEL_API int mkernel_mesh2d_delete_edge_by_index(int meshKernelId, int edgeIndex)
Deletes a mesh2d edge given the index of the edge. The coordinates of the edge middle points are used...
MKERNEL_API int mkernel_mesh2d_get_filtered_face_polygons(int meshKernelId, int propertyValue, double minValue, double maxValue, const GeometryList &facePolygons)
Gets the geometry list containing the face polygons within the filtering range.
MKERNEL_API int mkernel_get_interpolation_type_short(int &type)
Get the integer indicating the interpolation type short.
MKERNEL_API int mkernel_get_averaging_method_min_absolute_value(int &method)
Gets an int indicating the minimum absolute value averaging method type.
MKERNEL_API int mkernel_mesh2d_delete_edge(int meshKernelId, double xCoordinate, double yCoordinate, double xLowerLeftBoundingBox, double yLowerLeftBoundingBox, double xUpperRightBoundingBox, double yUpperRightBoundingBox)
Deletes the closest mesh2d edge to a point. The coordinates of the edge middle points are used for ca...
MKERNEL_API int mkernel_curvilinear_get_location_index(int meshKernelId, double xCoordinate, double yCoordinate, int locationType, const BoundingBox &boundingBox, int &locationIndex)
Gets the grid location closet to a specific coordinate.
MKERNEL_API double mkernel_get_separator()
Gets the double value used in the back-end library as separator and missing value.
MKERNEL_API int mkernel_get_splines(const GeometryList &geometryListIn, GeometryList &geometryListOut, int numberOfPointsBetweenNodes)
Get the computed spline points between two corner nodes.
MKERNEL_API int mkernel_network1d_to_mesh1d(int meshKernelId, double minFaceSize)
Convert network chainages to mesh1d nodes and edges.
MKERNEL_API int mkernel_curvilinear_compute_circular_grid(int meshKernelId, const meshkernel::MakeGridParameters ¶meters)
Compute a rectangular or circular curvilinear grid.
MKERNEL_API int mkernel_mesh2d_compute_orthogonalization(int meshKernelId, int projectToLandBoundaryOption, const meshkernel::OrthogonalizationParameters &orthogonalizationParameters, const GeometryList &selectingPolygon, const GeometryList &landBoundaries)
MKERNEL_API int mkernel_mesh2d_get_data(int meshKernelId, Mesh2D &mesh2d)
Gets the Mesh2D dimensions data.
MKERNEL_API int mkernel_undo_state_count_for_id(int meshKernelId, int &committedCount, int &restoredCount)
Count the number of undo actions for a particular meshKernelId.
MKERNEL_API int mkernel_curvilinear_initialize_orthogonal_grid_from_splines(int meshKernelId, const GeometryList &geometryList, const meshkernel::CurvilinearParameters &curvilinearParameters, const meshkernel::SplinesToCurvilinearParameters &splinesToCurvilinearParameters)
Generates a curvilinear grid from splines with the advancing front method. Initialization step (inter...
MKERNEL_API int mkernel_curvilinear_orthogonalize(int meshKernelId, const meshkernel::OrthogonalizationParameters &orthogonalizationParameters, double xLowerLeftCorner, double yLowerLeftCorner, double xUpperRightCorner, double yUpperRightCorner)
Define a block on the curvilinear grid where to perform orthogonalization.
MKERNEL_API int mkernel_mesh1d_get_data(int meshKernelId, Mesh1D &mesh1d)
Gets the Mesh1D data.
MKERNEL_API int mkernel_curvilinear_insert_face(int meshKernelId, double xCoordinate, double yCoordinate)
Inserts a new face on a curvilinear grid. The new face will be inserted on top of the closest edge by...
MKERNEL_API int mkernel_polygon_count_linear_refine(int meshKernelId, const GeometryList &polygonToRefine, int firstIndex, int secondIndex, int &numberOfPolygonNodes)
Counts the number of polygon nodes resulting from polygon refinement with mkernel_polygon_linear_refi...
MKERNEL_API int mkernel_mesh2d_averaging_interpolation(int meshKernelId, const GeometryList &samples, int locationType, int averagingMethodType, double relativeSearchSize, size_t minNumSamples, GeometryList &results)
AveragingInterpolation interpolation (ec_module)
MKERNEL_API int mkernel_mesh2d_count_nodes_in_polygons(int meshKernelId, const GeometryList &geometryListIn, int inside, int &numberOfMeshNodes)
Counts the number of selected mesh node indices.
MKERNEL_API int mkernel_mesh2d_get_mesh_boundaries_as_polygons(int meshKernelId, GeometryList &selectionPolygon, GeometryList &boundaryPolygons)
Retrieves the boundaries of a mesh as a series of separated polygons.
MKERNEL_API int mkernel_contacts_compute_boundary(int meshKernelId, const int *oneDNodeMask, const GeometryList &polygons, double searchRadius)
Computes 1d-2d contacts, where 1d nodes are connected to the closest 2d faces at the boundary (ggeo_m...
MKERNEL_API int mkernel_mesh2d_delete_orthogonalization(int meshKernelId)
Cleans the orthogonalization algorithm state, allocated in mkernel_mesh2d_initialize_orthogonalizatio...
MKERNEL_API int mkernel_mesh2d_get_hanging_edges(int meshKernelId, int *edges)
Gets the indices of hanging edges. An hanging edge is an edge where one of the two nodes is not conne...
MKERNEL_API int mkernel_get_projection_cartesian(int &projection)
Gets an int indicating the cartesian projection.
MKERNEL_API int mkernel_polygon_count_offset(int meshKernelId, const GeometryList &geometryListIn, int innerPolygon, double distance, int &numberOfPolygonNodes)
Counts the number of polygon nodes resulting from polygon offset.
MKERNEL_API int mkernel_curvilinear_compute_rectangular_grid_on_extension(int meshKernelId, const meshkernel::MakeGridParameters &makeGridParameters)
Computes a rectangular curvilinear grid on a defined extension.
MKERNEL_API int mkernel_curvilinear_convert_to_mesh2d(int meshKernelId)
Converts a curvilinear grid to an unstructured mesh.
MKERNEL_API int mkernel_mesh2d_connect_meshes(int meshKernelId, const Mesh2D &mesh2d, double searchFraction, bool connect)
Connect two or more disconnected regions along boundary.
MKERNEL_API int mkernel_curvilinear_frozen_line_add(int meshKernelId, double xFirstGridLineNode, double yFirstGridLineNode, double xSecondGridLineNode, double ySecondGridLineNode, int &frozenLineId)
Sets a new frozen line in the meshkernel state.
MKERNEL_API int mkernel_mesh2d_delete_node(int meshKernelId, int nodeIndex)
Deletes a mesh2d node.
MKERNEL_API int mkernel_curvilinear_compute_orthogonal_grid_from_splines(int meshKernelId, const GeometryList &geometryList, const meshkernel::CurvilinearParameters &curvilinearParameters, const meshkernel::SplinesToCurvilinearParameters &splinesToCurvilinearParameters)
Generates curvilinear grid from splines with the advancing front method.
MKERNEL_API int mkernel_clear_undo_state()
Clear the undo state for all mesh kernel ids, no undo is possible after this.
MKERNEL_API double mkernel_get_inner_outer_separator()
Gets the double value used to separate the inner part of a polygon from its outer part.
MKERNEL_API int mkernel_mesh2d_make_triangular_mesh_from_samples(int meshKernelId, const GeometryList &samples)
Makes a triangular mesh from a set of samples, triangulating the sample points.
MKERNEL_API int mkernel_get_interpolation_type_float(int &type)
Get the integer indicating the interpolation type float.
MKERNEL_API int mkernel_curvilinear_get_boundaries_as_polygons(int meshKernelId, int lowerLeftN, int lowerLeftM, int upperRightN, int upperRightM, GeometryList &boundaryPolygons)
Gets the boundary polygon of a curvilinear grid, nodes with invalid coordinates are excluded.
MKERNEL_API int mkernel_network1d_compute_fixed_chainages(int meshKernelId, double *fixedChainages, int sizeFixedChainages, double minFaceSize, double fixedChainagesOffset)
Compute the network chainages from fixed point locations.
MKERNEL_API int mkernel_polygon_get_offset(int meshKernelId, const GeometryList &geometryListIn, int inWard, double distance, GeometryList &geometryListOut)
Generate a new polygon from an existing one by offsetting the perimeter by a given distance.
MKERNEL_API int mkernel_mesh2d_casulli_derefinement_on_polygon(int meshKernelId, const GeometryList &polygons)
De-refine mesh using the Casulli de-refinement algorithm.
MKERNEL_API int mkernel_polygon_count_refine(int meshKernelId, const GeometryList &polygonToRefine, int firstIndex, int secondIndex, double distance, int &numberOfPolygonNodes)
Counts the number of polygon nodes resulting from polygon refinement with mkernel_polygon_refine.
MKERNEL_API int mkernel_get_exit_code_mesh_geometry_error(int &exitCode)
Gets the exit code of an exception of type MeshGeometryexitCode.
MKERNEL_API int mkernel_mesh2d_count_mesh_boundaries_as_polygons(int meshKernelId, GeometryList &selectionPolygon, int &numberOfPolygonNodes)
Counts the number of polygon nodes contained in the mesh boundary polygons computed in function mkern...
MKERNEL_API int mkernel_mesh2d_prepare_outer_iteration_orthogonalization(int meshKernelId)
Prepares an outer orthogonalization iteration, computing the new orthogonalization and smoothing weig...
MKERNEL_API int mkernel_mesh2d_get_edge(int meshKernelId, double xCoordinate, double yCoordinate, double xLowerLeftBoundingBox, double yLowerLeftBoundingBox, double xUpperRightBoundingBox, double yUpperRightBoundingBox, int &edgeIndex)
Gets the closest mesh2d edge to a point in a bounding box.
MKERNEL_API int mkernel_redo_state(bool &redone, int &meshKernelId)
Attempt to redo by one undo-action.
MKERNEL_API int mkernel_mesh2d_delete_small_flow_edges_and_small_triangles(int meshKernelId, double smallFlowEdgesThreshold, double minFractionalAreaTriangles)
Deletes all small mesh2d flow edges and small triangles. The flow edges are the edges connecting face...
MKERNEL_API int mkernel_curvilinear_snap_to_spline(int meshKernelId, const GeometryList &spline, double sectionControlPoint1x, double sectionControlPoint1y, double sectionControlPoint2x, double sectionControlPoint2y, double regionControlPointX=mkernel_get_separator(), double regionControlPointY=mkernel_get_separator())
Sets the curvilinear grid.
MKERNEL_API int mkernel_mesh2d_intersections_from_polygon(int meshKernelId, const GeometryList &boundaryPolygon, int *edgeNodes, int *edgeIndex, double *edgeDistances, double *segmentDistances, int *segmentIndexes, int *faceIndexes, int *faceNumEdges, int *faceEdgeIndex)
Gets the edges intersected by a polygon, with additional information on the intersections.
MKERNEL_API int mkernel_mesh2d_count_small_flow_edge_centers(int meshKernelId, double smallFlowEdgesLengthThreshold, int &numSmallFlowEdges)
Counts the number of small mesh2d flow edges. The flow edges are the edges connecting faces circumcen...
MKERNEL_API int mkernel_get_projection_spherical_accurate(int &projection)
Gets an int indicating the spherical accurate projection.
MKERNEL_API int mkernel_mesh2d_get_obtuse_triangles_mass_centers(int meshKernelId, GeometryList &result)
Gets the mass centers of obtuse mesh2d triangles. Obtuse triangles are those having one edge longer t...
MKERNEL_API int mkernel_get_version(char *version)
Gets pointer to version string.
MKERNEL_API int mkernel_curvilinear_set_line_line_shift(int meshKernelId, double xFirstGridLineNode, double yFirstGridLineNode, double xSecondGridLineNode, double ySecondGridLineNode)
Sets the start and end nodes of the line to shift.
MKERNEL_API int mkernel_mesh2d_get_closest_node(int meshKernelId, double xCoordinateIn, double yCoordinateIn, double searchRadius, double xLowerLeftBoundingBox, double yLowerLeftBoundingBox, double xUpperRightBoundingBox, double yUpperRightBoundingBox, double &xCoordinateOut, double &yCoordinateOut)
Gets the closest mesh2d node coordinates to a point, searching within a radius.
MKERNEL_API int mkernel_mesh2d_merge_two_nodes(int meshKernelId, int firstNode, int secondNode)
Merges two mesh2d nodes into one.
MKERNEL_API int mkernel_get_exit_code_stdlib_exception(int &exitCode)
Gets the exit code of an exception of type std::exception.
MKERNEL_API int mkernel_mesh2d_get_small_flow_edge_centers(int meshKernelId, double smallFlowEdgesThreshold, GeometryList &result)
Gets the small mesh2d flow edges. The flow edges are the edges connecting faces circumcenters.
MKERNEL_API int mkernel_mesh2d_remove_disconnected_regions(int meshKernelId)
Remove any disconnected regions from a mesh2d.
MKERNEL_API int mkernel_curvilinear_compute_grid_from_splines(int meshKernelId, const GeometryList &geometryList, const meshkernel::CurvilinearParameters &curvilinearParameters)
Generates curvilinear grid from splines.
MKERNEL_API int mkernel_curvilinear_compute_smoothness(int meshKernelId, int direction, double *smoothness)
Computes the smoothness of a curvilinear grid.
MKERNEL_API int mkernel_curvilinear_line_attraction_repulsion(int meshKernelId, double repulsionParameter, double xFirstNodeOnTheLine, double yFirstNodeOnTheLine, double xSecondNodeOnTheLine, double ySecondNodeOnTheLine, double xLowerLeftCorner, double yLowerLeftCorner, double xUpperRightCorner, double yUpperRightCorner)
Attracts/repulses grid lines in a block towards another set grid line.
MKERNEL_API int mkernel_deallocate_state(int meshKernelId)
Deallocate mesh state.
MKERNEL_API int mkernel_mesh2d_get_property_dimension(int meshKernelId, int propertyValue, int &dimension)
The dimension of a specified property of a 2D mesh.
MKERNEL_API int mkernel_get_averaging_method_inverse_distance_weighting(int &method)
Gets an int indicating the inverse distance weights averaging method type.
MKERNEL_API int mkernel_get_geometry_error(int &invalidIndex, int &type)
Gets the index of the erroneous entity.
MKERNEL_API int mkernel_curvilinear_initialize_line_shift(int meshKernelId)
Initializes the curvilinear line shift algorithm.
MKERNEL_API int mkernel_get_averaging_method_closest_point(int &method)
Gets an int indicating the closest point averaging method type.
MKERNEL_API int mkernel_clear_state()
Clear all internal mesh kernel state and undo actions, no undo will be possible after this.
MKERNEL_API int mkernel_curvilinear_compute_rectangular_grid_from_polygon(int meshKernelId, const meshkernel::MakeGridParameters &makeGridParameters, const GeometryList &geometryList)
Computes a rectangular curvilinear grid from polygon.
MKERNEL_API int mkernel_mesh2d_make_rectangular_mesh_from_polygon(int meshKernelId, const meshkernel::MakeGridParameters &makeGridParameters, const GeometryList &geometryList)
Makes a rectangular mesh from a polygon.
MKERNEL_API int mkernel_curvilinear_get_data(int meshKernelId, CurvilinearGrid &curvilinearGrid)
Gets the curvilinear grid data as a CurvilinearGrid struct (converted as set of edges and nodes)
MKERNEL_API int mkernel_mesh2d_get_dimensions(int meshKernelId, Mesh2D &mesh2d)
Gets the Mesh2D dimensions.
MKERNEL_API int mkernel_mesh2d_make_rectangular_mesh(int meshKernelId, const meshkernel::MakeGridParameters &makeGridParameters)
Makes a rectangular mesh.
MKERNEL_API int mkernel_get_exit_code_algorithm_error(int &exitCode)
Gets the exit code of an exception of type AlgorithmexitCode.
MKERNEL_API int mkernel_mesh2d_make_rectangular_mesh_on_extension(int meshKernelId, const meshkernel::MakeGridParameters &makeGridParameters)
Makes a rectangular mesh on a defined extension.
MKERNEL_API int mkernel_curvilinear_move_node_line_shift(int meshKernelId, double xFromCoordinate, double yFromCoordinate, double xToCoordinate, double yToCoordinate)
Moves a node of the line to shift, the operation can be performed multiple times.
MKERNEL_API int mkernel_mesh2d_flip_edges(int meshKernelId, int isTriangulationRequired, int projectToLandBoundaryRequired, const GeometryList &selectingPolygon, const GeometryList &landBoundaries)
Flips mesh2d edges, to optimize the mesh smoothness. This operation is usually performed after mkerne...
MKERNEL_API int mkernel_mesh2d_get_face_polygons(int meshKernelId, int numEdges, const GeometryList &facePolygons)
Gets the faces polygons with a number of edges larger or equal to numEdges.
MKERNEL_API int mkernel_mesh2d_make_triangular_mesh_from_polygon(int meshKernelId, const GeometryList &polygonPoints)
Generates a triangular mesh2d grid within a polygon. The size of the triangles is determined from the...
MKERNEL_API int mkernel_mesh2d_casulli_refinement_on_polygon(int meshKernelId, const GeometryList &polygons)
Refine mesh using the Casulli refinement algorithm.
MKERNEL_API int mkernel_mesh2d_get_node_index(int meshKernelId, double xCoordinate, double yCoordinate, double searchRadius, double xLowerLeftBoundingBox, double yLowerLeftBoundingBox, double xUpperRightBoundingBox, double yUpperRightBoundingBox, int &nodeIndex)
Finds the mesh2d node closest to a point, within a search radius.
MKERNEL_API int mkernel_mesh2d_move_node(int meshKernelId, double xCoordinate, double yCoordinate, int nodeIndex)
Moves a mesh2d node to a new position.
MKERNEL_API int mkernel_get_interpolation_type_int(int &type)
Get the integer indicating the interpolation type int.
MKERNEL_API int mkernel_mesh2d_convert_to_curvilinear(int meshKernelId, double xPointCoordinate, double yPointCoordinate)
Converts a mesh to a curvilinear mesh.
MKERNEL_API int mkernel_contacts_get_dimensions(int meshKernelId, Contacts &contacts)
Gets the number of 1d-2d contacts.
MKERNEL_API int mkernel_get_nodes_location_type(int &type)
Gets an int indicating the node location type.
MKERNEL_API int mkernel_mesh1d_get_dimensions(int meshKernelId, Mesh1D &mesh1d)
Gets the Mesh1D data dimensions.
MKERNEL_API int mkernel_get_exit_code_constraint_error(int &exitCode)
Gets the exit code of an exception of type ConstraintexitCode.
MKERNEL_API int mkernel_mesh2d_refine_based_on_samples(int meshKernelId, const GeometryList &samples, double relativeSearchRadius, int minimumNumSamples, const meshkernel::MeshRefinementParameters &meshRefinementParameters)
Refines a mesh2d based on samples. Refinement is achieved by successive splits of the face edges.
MKERNEL_API int mkernel_get_projection_spherical(int &projection)
Gets an int indicating the spherical projection.
MKERNEL_API int mkernel_mesh2d_get_property(int meshKernelId, int propertyValue, int locationId, const GeometryList &geometrylist)
Retrieves a specified property of a 2D mesh.
MKERNEL_API int mkernel_curvilinear_get_dimensions(int meshKernelId, CurvilinearGrid &curvilinearGrid)
Gets the curvilinear grid dimensions as a CurvilinearGrid struct (converted as set of edges and nodes...
MKERNEL_API int mkernel_mesh2d_casulli_derefinement(int meshKernelId)
De-refine mesh using the Casulli de-refinement algorithm.
MKERNEL_API int mkernel_mesh2d_delete_hanging_edges(int meshKernelId)
Deletes all hanging edges. An hanging edge is an edge where one of the two nodes is not connected.
MKERNEL_API int mkernel_mesh2d_translate(int meshKernelId, double translationX, double translationY)
Translate a mesh2d.
MKERNEL_API int mkernel_mesh2d_casulli_refinement(int meshKernelId)
Refine mesh using the Casulli refinement algorithm.
MKERNEL_API int mkernel_mesh2d_refine_ridges_based_on_gridded_samples(int meshKernelId, const GriddedSamples &griddedSamples, double relativeSearchRadius, int minimumNumSamples, int numberOfSmoothingIterations, const meshkernel::MeshRefinementParameters &meshRefinementParameters)
Refines a mesh2d based on samples with ridge refinement. This method automatically detects the ridges...
MKERNEL_API int mkernel_curvilinear_line_mirror(int meshKernelId, double mirroringFactor, int numRowsToMirror, double xFirstGridLineNode, double yFirstGridLineNode, double xSecondGridLineNode, double ySecondGridLineNode)
Mirrors a boundary gridline outwards. The boundary grid line is defined by its starting and ending po...
MKERNEL_API int mkernel_get_exit_code_range_error(int &exitCode)
Gets the exit code of an exception of type RangeexitCode.
MKERNEL_API int mkernel_curvilinear_delete_node(int meshKernelId, double xPointCoordinate, double yPointCoordinate)
Delete the node closest to a point.
MKERNEL_API int mkernel_mesh2d_convert_projection(int meshKernelId, int projectionType, const char *const zoneString)
Converts the projection of a mesh2d.
MKERNEL_API int mkernel_mesh2d_insert_edge_from_coordinates(int meshKernelId, double firstNodeX, double firstNodeY, double secondNodeX, double secondNodeY, int &firstNodeIndex, int &secondNodeIndex, int &edgeIndex)
Insert a new mesh2d edge from 2 coordinates.
MKERNEL_API int mkernel_curvilinear_iterate_orthogonal_grid_from_splines(int meshKernelId, int layer)
One advancement of the front in curvilinear grid from splines (interactive)
MKERNEL_API int mkernel_allocate_state(int projectionType, int &meshKernelId)
Creates a new mesh state and returns the generated meshKernelId.
MKERNEL_API int mkernel_mesh2d_insert_node(int meshKernelId, double xCoordinate, double yCoordinate, int &nodeIndex)
Insert a new mesh2d node at a specific coordinate.
MKERNEL_API int mkernel_curvilinear_snap_to_landboundary(int meshKernelId, const GeometryList &land, double sectionControlPoint1x, double sectionControlPoint1y, double sectionControlPoint2x, double sectionControlPoint2y, double regionControlPointX=mkernel_get_separator(), double regionControlPointY=mkernel_get_separator())
Sets the curvilinear grid.
MKERNEL_API int mkernel_mesh2d_rotate(int meshKernelId, double centreX, double centreY, double theta)
Rotate a mesh2d about a point.
MKERNEL_API int mkernel_mesh2d_triangulation_interpolation(int meshKernelId, const GeometryList &samples, int locationType, GeometryList &results)
Triangle interpolation (ec_module)
MKERNEL_API int mkernel_get_interpolation_type_double(int &type)
Get the integer indicating the interpolation type double.
MKERNEL_API int mkernel_curvilinear_refine(int meshKernelId, double xLowerLeftCorner, double yLowerLeftCorner, double xUpperRightCorner, double yUpperRightCorner, int refinement)
Directional curvilinear grid refinement or de-refinement. Additional gridlines are added or removed p...
MKERNEL_API int mkernel_curvilinear_compute_transfinite_from_triangle(int meshKernelId, const GeometryList &polygon, int firstNode, int secondNode, int thirdNode)
Computes a curvilinear grid in a triangle. 3 separate polygon nodes need to be selected.
MKERNEL_API int mkernel_curvilinear_delete_interior(int meshKernelId, const BoundingBox &boundingBox)
Delete the interior part of a curvilinear grid.
MKERNEL_API int mkernel_contacts_compute_with_points(int meshKernelId, const int *oneDNodeMask, const GeometryList &points)
Computes 1d-2d contacts, where 1d nodes are connected to the 2d faces mass centers containing the inp...
MKERNEL_API int mkernel_curvilinear_set_block_line_shift(int meshKernelId, double xLowerLeftCorner, double yLowerLeftCorner, double xUpperRightCorner, double yUpperRightCorner)
Defines a block on the curvilinear where the shifting is distributed.
MKERNEL_API int mkernel_contacts_compute_with_polygons(int meshKernelId, const int *oneDNodeMask, const GeometryList &polygons)
Computes 1d-2d contacts, where a 2d face per polygon is connected to the closest 1d node (ggeo_make1D...
MKERNEL_API int mkernel_get_averaging_method_min(int &method)
Gets an int indicating the minimum averaging method type.
MKERNEL_API int mkernel_curvilinear_finalize_line_shift(int meshKernelId)
Resets the instance of the line shift algorithm in MeshKernelState.
MKERNEL_API int mkernel_curvilinear_delete_exterior(int meshKernelId, const BoundingBox &boundingBox)
Delete the exterior part of a curvilinear grid.
MKERNEL_API int mkernel_deallocate_property(int meshKernelId, int propertyId)
Deallocate property calculator.
MKERNEL_API int mkernel_mesh1d_set(int meshKernelId, const Mesh1D &mesh1d)
Sets the meshkernel::Mesh1D state.
MKERNEL_API int mkernel_mesh2d_insert_edge(int meshKernelId, int startNode, int endNode, int &newEdgeIndex)
Insert a new mesh2d edge connecting two nodes.
MKERNEL_API int mkernel_mesh2d_compute_inner_ortogonalization_iteration(int meshKernelId)
Performs inner orthogonalization iteration, by slowly moving the mesh nodes to new optimal positions ...
MKERNEL_API int mkernel_mesh2d_refine_based_on_gridded_samples(int meshKernelId, const GriddedSamples &griddedSamples, const meshkernel::MeshRefinementParameters &meshRefinementParameters, bool useNodalRefinement)
Refine based on gridded samples.
MKERNEL_API int mkernel_mesh2d_refine_based_on_polygon(int meshKernelId, const GeometryList &geometryList, const meshkernel::MeshRefinementParameters &meshRefinementParameters)
Refines a mesh2d within a polygon. Refinement is achieved by splitting the edges contained in the pol...
MKERNEL_API int mkernel_curvilinear_compute_rectangular_grid(int meshKernelId, const meshkernel::MakeGridParameters &makeGridParameters)
Computes a rectangular curvilinear grid.
MKERNEL_API int mkernel_mesh2d_count_obtuse_triangles(int meshKernelId, int &numObtuseTriangles)
Gets the number of obtuse mesh2d triangles. Obtuse triangles are those having one edge longer than th...
MKERNEL_API int mkernel_curvilinear_frozen_line_is_valid(int meshKernelId, int frozenLineId, bool &isValid)
Checks if a curvilinear frozen line is valid.
MKERNEL_API int mkernel_mesh2d_set(int meshKernelId, const Mesh2D &mesh2d)
Sets the meshkernel::Mesh2D state.
MKERNEL_API int mkernel_mesh2d_get_location_index(int meshKernelId, double xCoordinate, double yCoordinate, int locationType, const BoundingBox &boundingBox, int &locationIndex)
Gets the mesh location closet to a specific coordinate.
MKERNEL_API int mkernel_mesh2d_merge_nodes_with_merging_distance(int meshKernelId, const GeometryList &geometryListIn, double mergingDistance)
Merges the mesh2d nodes within a distance of 0.001 m, effectively removing all small edges.
MKERNEL_API int mkernel_polygon_snap_to_landboundary(int meshKernelId, const GeometryList &land, GeometryList &polygon, int startIndex, int endIndex)
Snaps the polygon to the land boundary.
MKERNEL_API int mkernel_curvilinear_move_node(int meshKernelId, double xFromPoint, double yFromPoint, double xToPoint, double yToPoint)
Moves a point of a curvilinear grid from one location to another.
MKERNEL_API int mkernel_get_faces_location_type(int &type)
Gets an int indicating the faces location type.
MKERNEL_API int mkernel_mesh2d_delete(int meshKernelId, const GeometryList &polygon, int deletionOption, int invertDeletion)
Deletes a mesh in a polygon using several options.
MKERNEL_API int mkernel_contacts_compute_single(int meshKernelId, const int *oneDNodeMask, const GeometryList &polygons, double projectionFactor)
Computes 1d-2d contacts, where each single 1d node is connected to one mesh2d face circumcenter (ggeo...
MKERNEL_API int mkernel_curvilinear_compute_transfinite_from_splines(int meshKernelId, const GeometryList &splines, const meshkernel::CurvilinearParameters &curvilinearParameters)
Generates curvilinear grid from splines with transfinite interpolation.
MKERNEL_API int mkernel_mesh2d_set_property(int meshKernelId, const meshkernel::InterpolationParameters &interpolationParameters, const GeometryList &sampleData, int &propertyId)
Sets the property data for the mesh, the sample data points do not have to match the mesh2d nodes.
MKERNEL_API int mkernel_curvilinear_frozen_lines_get_ids(int meshKernelId, int *frozenLinesIds)
Gets the ids of the frozen lines.
MKERNEL_API int mkernel_contacts_set(int meshKernelId, const Contacts &contacts)
Sets the 1d-2d contacts.
MKERNEL_API int mkernel_mesh2d_get_filtered_face_polygons_dimension(int meshKernelId, int propertyValue, double minValue, double maxValue, int &geometryListDimension)
Retrieves the dimension of the geometry list containing the face polygons within the filtering range.
MKERNEL_API int mkernel_mesh2d_casulli_derefinement_elements(int meshKernelId, GeometryList &elements)
Get list of elements that will be removed after the Casulli de-refinement algorithm.
MKERNEL_API int mkernel_mesh2d_snap_to_landboundary(int meshKernelId, const GeometryList &selectingPolygon, const GeometryList &landBoundaries)
Snaps a mesh to a land boundary.
MKERNEL_API int mkernel_contacts_get_data(int meshKernelId, Contacts &contacts)
Gets the 1d-2d contacts indices (from index / to indices)
MKERNEL_API int mkernel_network1d_set(int meshKernelId, const GeometryList &polylines)
Sets the meshkernel::Network1D state.
MKERNEL_API int mkernel_get_averaging_method_simple_averaging(int &method)
Gets an int indicating the simple averaging method type.
MKERNEL_API int mkernel_curvilinear_refresh_orthogonal_grid_from_splines(int meshKernelId)
Converts curvilinear grid to mesh and refreshes the state (interactive)
MKERNEL_API int mkernel_contacts_compute_multiple(int meshKernelId, const int *oneDNodeMask)
Computes 1d-2d contacts, where a single 1d node is connected to multiple 2d face circumcenters (ggeo_...
MKERNEL_API int mkernel_mesh2d_get_node_edge_data(int meshKernelId, Mesh2D &mesh2d)
Gets only the node and edge Mesh2D data.
MKERNEL_API int mkernel_mesh2d_initialize_orthogonalization(int meshKernelId, int projectToLandBoundaryOption, meshkernel::OrthogonalizationParameters &orthogonalizationParameters, const GeometryList &selectingPolygon, const GeometryList &landBoundaries)
Initialization of the meshkernel::OrthogonalizationAndSmoothing algorithm.
MKERNEL_API int mkernel_mesh2d_is_valid_property(int meshKernelId, const int propertyId, const int locationId, bool &propertyIsAvailable)
Determine if the property data for the mesh can be computed.
MKERNEL_API int mkernel_mesh2d_casulli_derefinement_elements_on_polygon(int meshKernelId, const GeometryList &polygonGeometry, GeometryList &elements)
Get list of elements that will be removed after the Casulli de-refinement algorithm.
MKERNEL_API int mkernel_mesh2d_get_orthogonality_property_type(int &type)
Gets an int indicating the orthogonality property type for mesh2d.
MKERNEL_API int mkernel_set_undo_size(int undoStackSize)
Set the maximum size of the undo stack.
MKERNEL_API int mkernel_mesh1d_add(int meshKernelId, const Mesh1D &mesh1d)
Adds a mesh to the meshkernel::Mesh1D state.
MKERNEL_API int mkernel_curvilinear_set(int meshKernelId, const CurvilinearGrid &grid)
Sets the curvilinear grid.
MKERNEL_API int mkernel_curvilinear_compute_transfinite_from_polygon(int meshKernelId, const GeometryList &polygons, int firstNode, int secondNode, int thirdNode, int useFourthSide)
Computes a curvilinear grid in a polygon. 3 separate polygon nodes need to be selected.
MKERNEL_API int mkernel_curvilinear_smoothing_directional(int meshKernelId, int smoothingIterations, double xFirstGridlineNode, double yFirstGridlineNode, double xSecondGridLineNode, double ySecondGridLineNode, double xLowerLeftCornerSmoothingArea, double yLowerLeftCornerSmoothingArea, double xUpperRightCornerSmootingArea, double yUpperRightCornerSmootingArea)
Smooths a curvilinear grid along the direction specified by a segment.
MKERNEL_API int mkernel_network1d_compute_offsetted_chainages(int meshKernelId, double offset)
Compute the network chainages at a regular offset.
MKERNEL_API int mkernel_get_exit_code_not_implemented_error(int &exitCode)
Gets the exit code of an exception of type NotImplementedCode.
MKERNEL_API int mkernel_mesh2d_finalize_inner_ortogonalization_iteration(int meshKernelId)
Finalizes the orthogonalization outer iteration, computing the new coefficients for grid adaption and...
MKERNEL_API int mkernel_is_valid_state(int meshKernelId, bool &isValid)
Determine if the meshKernelId is valid.
MKERNEL_API int mkernel_curvilinear_line_shift(int meshKernelId)
Computes the new grid, shifting the line towards the moved nodes and distributing the shifting in blo...
MKERNEL_API int mkernel_curvilinear_frozen_line_delete(int meshKernelId, int frozenLineId)
Deletes an existing frozen line in the meshkernel state.
MKERNEL_API int mkernel_mesh2d_casulli_refinement_based_on_depths(int meshKernelId, const GeometryList &polygons, int propertyId, const meshkernel::MeshRefinementParameters &meshRefinementParameters, const double minimumRefinementDepth)
Refine mesh using the Casulli refinement algorithm based on the depth values.
MKERNEL_API int mkernel_polygon_linear_refine(int meshKernelId, const GeometryList &polygonToRefine, int firstNodeIndex, int secondNodeIndex, GeometryList &refinedPolygon)
Linear refines the polygon perimeter between two nodes.
MKERNEL_API int mkernel_curvilinear_frozen_lines_get_count(int meshKernelId, int &numFrozenLines)
Gets the number of stored frozen lines in the meshkernel state.
MKERNEL_API int mkernel_mesh2d_get_smoothness(int meshKernelId, GeometryList &geometryList)
Gets the smoothness, expressed as the ratio between the values of two neighboring faces areas.
MKERNEL_API int mkernel_get_exit_code_linear_algebra_error(int &exitCode)
Gets the exit code of an exception of type LinearAlgebraexitCode.
MKERNEL_API int mkernel_mesh2d_get_nodes_in_polygons(int meshKernelId, const GeometryList &geometryListIn, int inside, int *selectedNodes)
Gets the indices of the mesh2d nodes selected with a polygon.
MKERNEL_API int mkernel_get_exit_code_success(int &exitCode)
Gets the success exit code.
MKERNEL_API int mkernel_curvilinear_compute_curvature(int meshKernelId, int direction, double *curvature)
Computes the curvature of a curvilinear grid.
MKERNEL_API int mkernel_get_exit_code_meshkernel_error(int &exitCode)
Gets the exit code of an exception of type MeshKernelError.
MKERNEL_API int mkernel_curvilinear_frozen_line_get(int meshKernelId, int frozenLineId, double &xFirstFrozenLineCoordinate, double &yFirstFrozenLineCoordinate, double &xSecondFrozenLineCoordinate, double &ySecondFrozenLineCoordinate)
Sets a new frozen line in the meshkernel state.
MKERNEL_API int mkernel_curvilinear_delete_orthogonal_grid_from_splines(int meshKernelId)
Finalizes curvilinear grid from splines algorithm.
MKERNEL_API int mkernel_get_projection(int meshKernelId, int &projection)
Gets the coordinate projection of the meshkernel state.
MKERNEL_API int mkernel_get_edges_location_type(int &type)
Gets an int indicating the edge location type.
MKERNEL_API int mkernel_curvilinear_count_boundaries_as_polygons(int meshKernelId, int lowerLeftN, int lowerLeftM, int upperRightN, int upperRightM, int &numberOfPolygonNodes)
Count the number of nodes in curvilinear grid boundary polygons.
MKERNEL_API int mkernel_polygon_refine(int meshKernelId, const GeometryList &polygonToRefine, int firstNodeIndex, int secondNodeIndex, double targetEdgeLength, GeometryList &refinedPolygon)
Refines the polygon perimeter between two nodes. This interval is refined to achieve a target edge le...
MKERNEL_API int mkernel_clear_undo_state_for_id(int meshKernelId)
Clear the undo state for particular mesh kernel id, no undo for the id is possible after this.
MKERNEL_API int mkernel_get_averaging_method_max(int &method)
Gets an int indicating the max value averaging method type.
MKERNEL_API int mkernel_mesh2d_add(int meshKernelId, const Mesh2D &mesh2d)
Adds a mesh to the meshkernel::Mesh2D state.
MKERNEL_API int mkernel_mesh2d_get_face_polygons_dimension(int meshKernelId, int numEdges, int &geometryListDimension)
Gets the dimension of faces polygons with a number of edges larger or equal to numNodes.
MKERNEL_API int mkernel_mesh2d_get_orthogonality(int meshKernelId, GeometryList &geometryList)
Gets the mesh orthogonality, expressed as the ratio between the edges and the segments connecting the...
MKERNEL_API int mkernel_undo_state_count(int &committedCount, int &restoredCount)
Count the number of undo actions.
MKERNEL_API int mkernel_splines_snap_to_landboundary(int meshKernelId, const GeometryList &land, GeometryList &splines, int startSplineIndex, int endSplineIndex)
Snaps the spline (or splines) to the land boundary.
MKERNEL_API int mkernel_undo_state(bool &undone, int &meshKernelId)
Attempt to undo by one undo-action.
MKERNEL_API int mkernel_mesh2d_split_row(int meshKernelId, int firstNode, int secondNode)
An-isotropically refines the elements along a row or column, given a starting edge.
MKERNEL_API int mkernel_get_exit_code_unknown_exception(int &exitCode)
Gets the exit code of an exception of unknown type.
MKERNEL_API int mkernel_curvilinear_smoothing(int meshKernelId, int smoothingIterations, double xLowerLeftCorner, double yLowerLeftCorner, double xUpperRightCorner, double yUpperRightCorner)
Smooths a curvilinear grid.
MKERNEL_API int mkernel_mesh2d_make_global(int meshKernelId, int numLongitudeNodes, int numLatitudeNodes)
Compute the global mesh with a given number of points along the longitude and latitude directions.
MKERNEL_API int mkernel_get_error(char *errorMessage)
Gets pointer to error message.
MKERNEL_API int mkernel_curvilinear_full_refine(int meshKernelId, int mRefinement, int nRefinement)
Curvilinear grid refinement. Additional gridlines are added in both directions, over the entire grid.
MKERNEL_API int mkernel_polygon_get_included_points(int meshKernelId, const GeometryList &selectingPolygon, const GeometryList &polygonToSelect, GeometryList &selectionResults)
Selects the polygon nodes within another polygon.
A struct used to describe parameters for generating a curvilinear grid in a C-compatible manner.
Definition Parameters.hpp:111
Parameters used by the sample interpolation.
Definition Parameters.hpp:270
This struct describes the necessary parameters to create a new curvilinear grid in a C-compatible man...
Definition Parameters.hpp:38
A struct used to describe the mesh refinement parameters in a C-compatible manner.
Definition Parameters.hpp:187
A struct used to describe the orthogonalization parameters in a C-compatible manner.
Definition Parameters.hpp:238
A struct used to describe the spline to curvilinear grid parameters in a C-compatible manner.
Definition Parameters.hpp:139
A struct describing a bounding box.
Definition BoundingBox.hpp:36
A struct used to describe the values of a curvilinear grid in a C-compatible manner.
Definition CurvilinearGrid.hpp:34
A struct used to describe a list of geometries in a C-compatible manner.
Definition GeometryList.hpp:34
A struct describing gridded samples.
Definition GriddedSamples.hpp:34
A struct used to describe the values of a mesh 1d in a C-compatible manner.
Definition Mesh1D.hpp:34
A struct used to describe the values of an unstructured, two-dimensional mesh in a C-compatible manne...
Definition Mesh2D.hpp:34