opencmp.diffuse_interface package
Submodules
opencmp.diffuse_interface.dim module
- class opencmp.diffuse_interface.dim.DIM(DIM_dir, import_dir, t_param)[source]
Bases:
object
Class to hold the diffuse interface methods and parameters.
- _generate_BC_masks()[source]
Function to generate the boundary condition phase field masks.
This function produces numpy arrays, which must later be converted into gridfunctions.
- _generate_DIM_mesh()[source]
Function to get the nonconformal mesh.
This mesh is a simple rectangle or rectangular prism, is constructed of quadrilateral/hexahedral elements and has exterior boundaries “left”, “right”, “top”, “bottom” (and “front” and “back” for 3D). See get_Netgen_nonconformal in mesh_helpers for more details.
- Return type
None
- _generate_phase_field()[source]
Function to generate the phase field.
This function produces a numpy array, which must later be converted into a gridfunction.
- _load_bc_parameters(config, bc_config, quiet=False)[source]
Helper function to load the BC parameters. Should only be called if self.multiple_bcs=True.
- Parameters
config (
ConfigParser
) – The main DIM config file.bc_config (
ConfigParser
) – The DIM BC config file (or the config file containing [VERTICES] and [CENTROIDS]).quiet (
bool
) – If True suppresses warnings about using default parameter values.
- Return type
None
- _load_nonconformal_parameters()[source]
Helper to load parameters needed for generating a nonconformal mesh and phase fields.
- Return type
None
- get_DIM_gridfunctions(mesh, interp_ord)[source]
Function to get all of the phase fields and masks as gridfunctions.
This is either done by loading the phase fields and masks from files or by constructing the numpy arrays and then converting those into gridfunctions.
- Parameters
mesh (
Mesh
) – The mesh for the gridfunctions.interp_ord (
int
) – The interpolant order for the finite element space for the gridfunctions.
opencmp.diffuse_interface.interface module
- opencmp.diffuse_interface.interface.get_binary_2d(boundary_lst, N, scale, offset)[source]
Generate a binary representation of a 2D complex geometry on a numpy array.
This is done by setting the array elements corresponding to points inside the geometry to 1 and all other array elements to 0.
- Parameters
boundary_lst (
List
) – List of coordinates of the geometry’s boundary vertices in counterclockwise order.N (
List
[int
]) – Number of mesh elements in each direction (N+1 nodes).scale (
List
[float
]) – Extent of the meshed domain in each direction ([-2,2] square -> scale=[4,4]).offset (
List
[float
]) – Centers the meshed domain in each direction ([-2,2] square -> offset=[2,2]).
- Return type
ndarray
- Returns
Array containing a binary representation of the complex geometry.
- opencmp.diffuse_interface.interface.get_binary_3d(face_lst, N, scale, offset, mnum=1, close=False)[source]
Generate a binary representation of a 3D complex geometry on a numpy array.
This is done by setting the array elements corresponding to points inside the geometry to 1 and all other array elements to 0.
- Parameters
face_lst (
ndarray
) – List of the vertices and outwards facing normals of the complex geometry’s faces.N (
List
[int
]) – Number of mesh elements in each direction (N+1 nodes).scale (
List
[float
]) – Extent of the meshed domain in each direction ([-2,2] cube -> scale=[4,4,4]).offset (
List
[float
]) – Centers the meshed domain in each direction ([-2,2] cube -> offset=[2,2,2]).mnum (
float
) – Magic number that increases the distance an array element can be from a vertex while still belonging to that vertex. Increase for higher order interpolants if the generated border has gaps.close (
bool
) – If True, wraps the binary hole filling in a binary closing. Use if the generated border has gaps.
- Return type
ndarray
- Returns
Array containing a binary representation of the complex geometry.
- opencmp.diffuse_interface.interface.get_phi(binary, lmbda, N, scale, offset, dim=2)[source]
Generate a phase field from a binary representation of a complex geometry.
The phase field diffuses from 1 inside of the complex geometry to 0 outside of the complex geometry.
- Parameters
binary (
ndarray
) – Array containing binary representation of complex geometry.lmbda (
float
) – Measure of the diffuseness of the phase field boundary.N (
List
[int
]) – Number of mesh elements in each direction (N+1 nodes).scale (
List
[float
]) – Extent of the meshed domain in each direction ([-2,2] square -> scale=[4,4]).offset (
List
[float
]) – Centers the meshed domain in each direction ([-2,2] square -> offset=[2,2]).dim (
int
) – Dimension of the domain (must be 2 or 3).
- Return type
ndarray
- Returns
Array containing the phase field.
- opencmp.diffuse_interface.interface.nonconformal_subdomain_2d(boundary_lst, vertices, N, scale, offset, lmbda_overlap=False, centroid=None)[source]
Function to generate a 2D BC mask.
Generate a numpy array that can be used to mask a NGSolve function based on partitioning a mesh’s interior domain into sections that will have different boundary conditions. Only works in 2D.
- Parameters
boundary_lst (
List
) – List of coordinates of an interior domain’s boundary vertices in counterclockwise order.vertices (
List
) – The coordinates of the two vertices that denote the section of the interior boundary.N (
List
[int
]) – Number of mesh elements in each direction (N+1 nodes).scale (
List
[float
]) – Extent of the meshed domain in each direction ([-2,2] square -> scale=[4,4]).offset (
List
[float
]) – Centers the meshed domain in each direction ([-2,2] square -> offset=[2,2]).lmbda_overlap (
Union
[float
,bool
]) – Measure of the diffuseness of the boundary between sections (sharp boundary if False).centroid (
Optional
[Tuple
[float
,float
],None
]) – Coordinates of the point to use as the centroid of the split (centroid of domain if None).
- Return type
ndarray
- Returns
Mask of domain.
- opencmp.diffuse_interface.interface.nonconformal_subdomain_3d(face_lst, vertices, N, scale, offset, lmbda_overlap=False, centroid=None)[source]
Function to generate a 3D BC mask.
Generate a numpy array that can be used to mask a NGSolve function based on partitioning a mesh’s interior domain into sections that will have different boundary conditions. Only works in 3D.
- Parameters
face_lst (
ndarray
) – List of the vertices and outwards facing normals of the interior domain’s faces.vertices (
str
) – The path to the .msh or .stl file defining the boundary region of the domain.N (
List
[int
]) – Number of mesh elements in each direction (N+1 nodes).scale (
List
[float
]) – Extent of the meshed domain in each direction ([-2,2] cube -> scale=[4,4,4]).offset (
List
[float
]) – Centers the meshed domain in each direction ([-2,2] cube -> offset=[2,2,2]).lmbda_overlap (
Union
[float
,bool
]) – Measure of the diffuseness of the boundary between sections (sharp boundary if False).centroid (
Optional
[Tuple
[float
,float
,float
],None
]) – Coordinates of the point to use as the centroid of the split (centroid of domain if None).
- Return type
ndarray
- Returns
Mask of domain.
- opencmp.diffuse_interface.interface.split_nonconformal_subdomains_2d(boundary_lst, vertices, N, scale, offset, lmbda_overlap=False, remainder=False, centroid={})[source]
Function to generate all BC masks for a 2D diffuse interface simulation.
Generate a dictionary of numpy arrays that can be used to mask a NGSolve function based on partitioning a mesh’s interior domain into sections that will have different boundary conditions. Only works in 2D.
- Parameters
boundary_lst (
List
) – List of coordinates of an interior domain’s boundary vertices in counterclockwise order.vertices (
Dict
) – Dictionary of coordinates of the vertices that denote the different sections of the interior boundary. Vertices must be ordered counterclockwise (unit square with different boundary conditions on each side -> vertices={‘bottom’: [(0,0), (1,0)], ‘right’: [(1,0), (1,1)], ‘top’: [(1,1), (0,1)], ‘left’: [(0,1), (0,0)]}).N (
List
[int
]) – Number of mesh elements in each direction (N+1 nodes).scale (
List
[float
]) – Extent of the meshed domain in each direction ([-2,2] square -> scale=[4,4]).offset (
List
[float
]) – Centers the meshed domain in each direction ([-2,2] square -> offset=[2,2]).lmbda_overlap (
Union
[float
,bool
]) – Measure of the diffuseness of the boundary between sections (sharp boundary if False).remainder (
bool
) – If True a final mask is made from the regions left unmasked by all the previous masks (use if a particular boundary section has poorly defined end vertices).centroid (
Dict
) – Dictionary of the coordinates of the points to use as the centroids of the splits. The dictionary’s keys should correspond to the boundary condition names like those of vertices.
- Return type
Dict
[str
,ndarray
]- Returns
Dictionary of numpy array masks.
- opencmp.diffuse_interface.interface.split_nonconformal_subdomains_3d(face_lst, vertices, N, scale, offset, lmbda_overlap=False, remainder=False, centroid={})[source]
Function to generate all BC masks for a 3D diffuse interface simulation.
Generate a list of numpy arrays that can be used to mask a NGSolve function based on partitioning a mesh’s interior domain into sections that will have different boundary conditions. Only works in 3D.
- Parameters
face_lst (
List
) – List of coordinates of an interior domain’s boundary vertices in counterclockwise order.vertices (
Dict
[str
,str
]) – Dictionary of paths to the .msh or .stl files that denote the different sections of the interior boundary.N (
List
[int
]) – Number of mesh elements in each direction (N+1 nodes).scale (
List
[float
]) – Extent of the meshed domain in each direction ([-2,2] square -> scale=[4,4]).offset (
List
[float
]) – Centers the meshed domain in each direction ([-2,2] square -> offset=[2,2]).lmbda_overlap (
Union
[float
,bool
]) – Measure of the diffuseness of the boundary between sections (sharp boundary if False).remainder (
bool
) – If True a final mask is made from the regions left unmasked by all the previous masks (use if a particular boundary section has poorly defined end vertices).centroid (
Dict
) – Dictionary of the coordinates of the points to use as the centroids of the splits. The dictionary’s keys should correspond to the boundary condition names like those of vertices.
- Return type
Dict
[str
,ndarray
]- Returns
Dictionary of numpy array masks.
opencmp.diffuse_interface.mesh_helpers module
- opencmp.diffuse_interface.mesh_helpers.angle_between(p1_tmp, p2_tmp, p3_tmp)[source]
Calculate the angle formed by p1-p2-p3.
- Parameters
p1_tmp (
List
[float
]) – [x,y] or [x,y,z] coordinate.p2_tmp (
List
[float
]) – [x,y] or [x,y,z] coordinate.p3_tmp (
List
[float
]) – [x,y] or [x,y,z] coordinate.
- Return type
float
- Returns
Calculated angle in radians rescaled to [0,2*pi].
- opencmp.diffuse_interface.mesh_helpers.calc_barycentric(p, v1, v2, v3)[source]
Use barycentric coordinates to determine if point p, which lies on the plane with normal n, also lies within the triangle v1-v2-v3.
- Parameters
p (
ndarray
) – [x,y,z] coordinate.v1 (
ndarray
) – [x,y,z] coordinate of one of the triangle’s vertices.v2 (
ndarray
) – [x,y,z] coordinate of one of the triangle’s vertices.v3 (
ndarray
) – [x,y,z] coordinate of one of the triangle’s vertices.
- Return type
bool
- Returns
True if p lies in the triangle, otherwise False.
- opencmp.diffuse_interface.mesh_helpers.calc_unit_normal(v1, v2, v3)[source]
Calculates the unit normal of the triangle v1-v2-v3.
- Parameters
v1 (
ndarray
) – [x,y,z] coordinate of one of the triangle’s vertices.v2 (
ndarray
) – [x,y,z] coordinate of one of the triangle’s vertices.v3 (
ndarray
) – [x,y,z] coordinate of one of the triangle’s vertices.
- Return type
ndarray
- Returns
Unit normal vector.
- opencmp.diffuse_interface.mesh_helpers.crop_to_mesh_bounds(arr, N, scale, offset, tmp_N, tmp_scale, tmp_offset)[source]
Take an array that exceeds the nonconformal mesh’s boundary and crop it so it fits within an array defined over the nonconformal mesh.
- Parameters
arr (
ndarray
) – The numpy array to crop.N (
List
[int
]) – Number of mesh elements in each direction (N+1 nodes).scale (
List
[float
]) – Extent of the meshed domain in each direction.offset (
List
[float
]) – Centers the meshed domain in each direction.tmp_N (
List
[int
]) – Number of mesh elements used when constructing the phase fields (preserves original dx).tmp_scale (
List
[float
]) – Scale used when constructing the phase fields.tmp_offset (
List
[float
]) – Offset used when constructing the phase fields.
- Return type
ndarray
- Returns
Numpy array covering only the bounds of the mesh and containing portions of arr.
- opencmp.diffuse_interface.mesh_helpers.get_Netgen_nonconformal(N, scale, offset, dim=2, quad=True)[source]
Generate a structured Netgen mesh over a prescribed square/cubic domain.
- Parameters
N (
List
[int
]) – Number of mesh elements in each direction (N+1 nodes).scale (
List
[float
]) – Extent of the meshed domain in each direction ([-2,2] square -> scale=[4,4]).offset (
List
[float
]) – Centers the meshed domain in each direction ([-2,2] square -> offset=[2,2]).dim (
int
) – Dimension of the domain (must be 2 or 3).quad (
bool
) – If True, the mesh consists of quadrilateral/hexahedral elements. Otherwise the mesh consists of triangular/tetrahedral elements.
- Return type
Mesh
- Returns
Structured Netgen mesh.
- opencmp.diffuse_interface.mesh_helpers.get_mesh_boundary_2d(filename)[source]
Get the ordered boundary vertices of a 2D .stl file.
- Parameters
filename (
str
) – Path to .stl file.- Returns
boundary_lst: List of coordinates of .stl file boundary vertices in counterclockwise order.
bounds_lst: List of lists of min and max values of the .stl file boundary vertices in each direction.
- Return type
Tuple[List, List]
- opencmp.diffuse_interface.mesh_helpers.get_mesh_boundary_3d(filename)[source]
Get the edges and boundary edges of a 3D .msh or .stl file.
This function should only be used for surface meshes since the boundary edges are considered to be edges with only one face. Volume meshes do not have such boundary edges. NOTE: Using a .msh file does not always give a comprehensive edge_lst since for some unknown reason ngmesh.Elements1D() does not always give all of the edges in the mesh. This is not necessarily a problem if edge_lst is only being used to produce BC masks (but always check!). However, using a .stl file should always give the full boundary_lst and edge_lst.
- Parameters
filename (
str
) – Path to .msh or .stl file.- Returns
edge_lst: List of coordinates of edge vertices.
boundary_lst: List of coordinates of boundary edge vertices.
- Return type
Tuple[List, List]
- opencmp.diffuse_interface.mesh_helpers.get_mesh_boundary_from_conformal_2d(mesh)[source]
Get the ordered boundary vertices of a 2D Netgen mesh.
- Parameters
mesh (
Mesh
) – 2D Netgen mesh to find the boundary vertices of.- Return type
List
- Returns
List of coordinates of Netgen mesh boundary vertices in counterclockwise order.
- opencmp.diffuse_interface.mesh_helpers.get_mesh_edges_vertices_3d(filename)[source]
Get the edges and vertices of a 3D .stl file.
This is a helper function for calculating mesh quality metrics.
- Parameters
filename (
str
) – Path to the .stl file.- Returns
face_lst: List of the vertices and outwards facing normals of the mesh’s faces.
edge_lst: List of coordinates of edge vertices.
v_lst: List of coordinates of vertices.
v_con_lst: List containing each vertex and an ordered list of its neighbouring vertices (coordinate form).
- Return type
Tuple[List, List, List, List]
- opencmp.diffuse_interface.mesh_helpers.get_new_bounds(bounds_lst, N, scale, offset)[source]
If a .stl file’s boundaries exceed the bounds of the nonconformal mesh, find new bounds for generating the phase fields and masks.
- Parameters
bounds_lst (
List
) – List containing lists of the min/max values in each direction (x,y or x,y,z).N (
List
[int
]) – Number of mesh elements in each direction (N+1 nodes).scale (
List
[float
]) – Extent of the meshed domain in each direction.offset (
List
[float
]) – Centers the meshed domain in each direction.
- Returns
tmp_N: Number of mesh elements to use when constructing the phase fields (preserves original dx).
tmp_scale: Scale to use when constructing the phase fields.
tmp_offset: Offset to use when constructing the phase fields.
- Return type
Tuple[List[int], List[float], List[float]]
- opencmp.diffuse_interface.mesh_helpers.get_stl_faces(filename)[source]
Compile the face vertices and outwards facing normals of a mesh defined in a .stl file into a list.
The mesh will typically be a boundary mesh so face_lst can be used by ray_trace_3d.
- Parameters
filename (
str
) – Path to the .stl file.- Returns
face_lst: Numpy array listing out the vertices and outwards facing normals of the mesh’s faces.
bounds_lst: List of lists of min and max values of the .stl file boundary vertices in each direction.
- Return type
Tuple[ndarray, List]
- opencmp.diffuse_interface.mesh_helpers.index_sublist(lst, val)[source]
Find which sublist of a list contains a given value.
- Parameters
lst (
List
) – List to search through.val (
Any
) – Value to check for.
- Return type
int
- Returns
Index of lst corresponding to the sublist of interest.
- opencmp.diffuse_interface.mesh_helpers.move_vertex(vertex, vertex_lst, vertex_coords)[source]
Find vertex pair [A,B] (or [B,A]) containing vertex A and return vertex B of said pair.
- Parameters
vertex (
List
) – Vertex used to find vertex pair.vertex_lst (
List
) – List of vertex pairs.vertex_coords (
List
) – List of coordinates corresponding to vertex pairs in vertex_lst.
- Returns
v_prime: Vertex B of vertex pair [A,B] (or [B,A]).
p_prime: [x,y] or [x,y,z] coordinates of vertex B.
index: Index of vertex pair [A,B] (or [B,A]) in vertex_lst.
- Return type
Tuple[List, List, int]
- opencmp.diffuse_interface.mesh_helpers.order_ccw(points_lst)[source]
Rearrange an ordered list of polygon vertices to be in counterclockwise order.
- Parameters
points_lst (
List
) – Ordered list of coordinates of polygon vertices (could be in clockwise or counterclockwise order).- Returns
List of coordinates of polygon vertices in counterclockwise order.
- Return type
points_lst
- opencmp.diffuse_interface.mesh_helpers.orient_2d(p1, p2, p3, eps=1e-10)[source]
Check if a set of three 2D points is collinear, ordered counterclockwise, or ordered clockwise.
- Parameters
p1 (
List
[float
]) – [x,y] coordinate.p2 (
List
[float
]) – [x,y] coordinate.p3 (
List
[float
]) – [x,y] coordinate.eps (
float
) – Tolerance on whether or not the points are collinear.
- Return type
str
- Returns
Indicates the ordering of the points.
- opencmp.diffuse_interface.mesh_helpers.orient_3d(p, v1, v2, v3, eps=1e-06)[source]
Determine which side of the plane defined by v1, v2, v3 the point p is on.
This function assumes a point and plane in 3D space.
- Parameters
p (
ndarray
) – [x,y,z] coordinate.v1 (
ndarray
) – [x,y,z] coordinate of one of the plane’s vertices.v2 (
ndarray
) – [x,y,z] coordinate of one of the plane’s vertices.v3 (
ndarray
) – [x,y,z] coordinate of one of the plane’s vertices.eps (
float
) – Tolerance for adding small perturbations to the plane’s vertices if p is degenerate with an edge or vertex of the plane.
- Return type
int
- Returns
Denotes which side of the plane p is on.
- opencmp.diffuse_interface.mesh_helpers.ray_trace_2d(x, y, polygon)[source]
Determine if the given 2D point is located inside the given polygon using the ray-tracing point-in-polygon technique.
Extend a ray from the point along the positive x-axis to infinity and count the number of times the ray intersects the polygon. The point is located inside the polygon if there are an odd number of intersections.
- Parameters
x (
float
) – x-coordinate of point.y (
float
) – y-coordinate of point.polygon (
List
) – List of coordinates of polygon vertices in counterclockwise order.
- Return type
bool
- Returns
Whether the point is inside the polygon.
- opencmp.diffuse_interface.mesh_helpers.reorder_vertices_2d(edge_lst, vertex_lst, vertex_coords)[source]
Rearrange a list of mesh boundary vertices to be in order based on the mesh’s boundary edge connectivities.
- Parameters
edge_lst (
List
) – List of mesh boundary edges.vertex_lst (
List
) – List of pairs of mesh boundary vertices corresponding to the mesh’s boundary edges.vertex_coords (
List
) – List of coordinates corresponding to vertex pairs in vertex_lst.
- Return type
List
- Returns
Ordered list of coordinates of mesh boundary vertices (could be clockwise or counterclockwise order).
opencmp.diffuse_interface.mesh_quality_metrics module
- opencmp.diffuse_interface.mesh_quality_metrics.calc_curvature_3d(v, con_lst)[source]
Estimates the mean curvature at a point by fitting an osculating parabola to the 3D surface surrounding the point.
!!! NOT WORKING PROPERLY !!!
- Parameters
v (
ndarray
) – [x, y, z] coordinates of the vertex of interest.con_lst (
List
) – A list of the coordinates of all the vertices that share an edge with v. The vertices are ordered by their shared edges.
- Return type
float
- Returns
Estimate of the mean curvature at v.
- opencmp.diffuse_interface.mesh_quality_metrics.get_chords_2d(boundary_lst, crossing=True, separation=0)[source]
Get every chord length in a 2D polygon.
A chord length is defined as the Euclidean distance between two different polygon vertices. Chords can be excluded if they cross the polygon’s boundary or if their bounding vertices are too close to each other.
- Parameters
boundary_lst (
List
) – List of coordinates of the polygon boundary vertices in counterclockwise order.crossing (
bool
) – If False, only chords completely contained within the polygon’s boundary are included.separation (
int
) – The minimum number of vertices that must separate a pair of boundary vertices in order to have a chord between those vertices. This can be used to find necks.
- Return type
List
- Returns
List of all of the polygon’s chord lengths.
- opencmp.diffuse_interface.mesh_quality_metrics.get_chords_3d(boundary_lst, face_lst, crossing=True)[source]
Get every chord length in a 3D polygon.
A chord length is defined as the Euclidean distance between two different polygon vertices. Chords can be excluded if they cross the polygon’s boundary.
- Parameters
boundary_lst (
List
) – List of coordinates of the polygon boundary vertices.face_lst (
List
) – List of the vertices and outwards facing normals of the polygon’s faces.crossing (
bool
) – If False, only chords completely contained within the polygon’s boundary are included.
- Return type
List
- Returns
List of all of the polygon’s chord lengths.
- opencmp.diffuse_interface.mesh_quality_metrics.get_radius_curvature_2d(boundary_lst)[source]
Get the radius of curvature at every boundary vertex in a 2D polygon.
Exclude vertices whose connected edges make a 90 degree or acute angle. List these instead as discontinuities with magnitudes equal to the inverse of their angle.
- Parameters
boundary_lst (
List
) – List of coordinates of the polygon boundary vertices in counterclockwise order.- Returns
rc: List of all of the polygon’s radii of curvature.
discont: List of the magnitudes of all of the polygon’s discontinuities.
- Return type
Tuple[List[float], List[float]]
- opencmp.diffuse_interface.mesh_quality_metrics.get_radius_curvature_3d(v_con_lst)[source]
Get the radius of curvature at every boundary vertex in a 3D polygon.
Exclude vertices where the boundary surface has a mean curvature of 90 degrees or less. List these instead asdiscontinuities with magnitudes equal to the inverse of their mean curvature.!!! NOT WORKING PROPERLY !!!- Parameters
v_con_lst (
List
) – List containing each vertex and an ordered list of its neighbouring vertices (coordinate form).- Returns
rc: List of all of the polygon’s radii of curvature.
discont: List of the magnitudes of all of the polygon’s discontinuities.
- Return type
Tuple[List[float], List[float]]
- opencmp.diffuse_interface.mesh_quality_metrics.line_segment_face_intersect_3d(p1_tmp, p2_tmp, v1_tmp, v2_tmp, v3_tmp, n_tmp)[source]
Check if line segment p1-p2 intersects face v1-v2-v3.
- Parameters
p1_tmp (
List
[float
]) – [x,y,z] coordinate of a boundary point of the line segment.p2_tmp (
List
[float
]) – [x,y,z] coordinate of a boundary point of the line segment.v1_tmp (
List
[float
]) – [x,y,z] coordinate of a vertex of the face.v2_tmp (
List
[float
]) – [x,y,z] coordinate of a vertex of the face.v3_tmp (
List
[float
]) – [x,y,z] coordinate of a vertex of the face.
- Return type
bool
- Returns
True if the line segment intersects the face, otherwise False.
- opencmp.diffuse_interface.mesh_quality_metrics.line_segments_intersect_2d(p1, p2, p3, p4)[source]
Check if two 2D line segments intersect.
- Parameters
p1 (
List
[float
]) – [x,y] coordinate of a boundary point of the first line segment.p2 (
List
[float
]) – [x,y] coordinate of a boundary point of the first line segment.p3 (
List
[float
]) – [x,y] coordinate of a boundary point of the second line segment.p4 (
List
[float
]) – [x,y] coordinate of a boundary point of the second line segment.
- Return type
bool
- Returns
True if the line segments intersect, otherwise False.