fractopo.tval package

Submodules

fractopo.tval.proximal_traces module

Determine traces that could be integrated together.

determine_proximal_traces takes an input of GeoSeries or GeoDataFrame of LineString geometries and returns a GeoDataFrame with a new column Merge which has values of True or False depending on if nearby proximal traces were found.

fractopo.tval.proximal_traces.determine_proximal_traces(traces: GeoSeries | GeoDataFrame, buffer_value: float | int | Real, azimuth_tolerance: float | int | Real) GeoDataFrame

Determine proximal traces.

Takes an input of GeoSeries or GeoDataFrame of LineString geometries and returns a GeoDataFrame with a new column Merge which has values of True or False depending on if nearby proximal traces were found.

E.g.

>>> lines = [
...     LineString([(0, 0), (0, 3)]),
...     LineString([(1, 0), (1, 3)]),
...     LineString([(5, 0), (5, 3)]),
...     LineString([(0, 0), (-3, -3)]),
... ]
>>> traces = gpd.GeoDataFrame({"geometry": lines})
>>> buffer_value = 1.1
>>> azimuth_tolerance = 10
>>> determine_proximal_traces(traces, buffer_value, azimuth_tolerance)
                   geometry  Merge
0    LINESTRING (0 0, 0 3)   True
1    LINESTRING (1 0, 1 3)   True
2    LINESTRING (5 0, 5 3)  False
3  LINESTRING (0 0, -3 -3)  False
fractopo.tval.proximal_traces.is_similar_azimuth(trace: LineString, other: LineString, tolerance: float | int | Real) bool

Determine if azimuths of trace and other are close.

Checks both the start – end -azimuth and regression-based azimuth.

E.g.

>>> trace = LineString([(0, 0), (0, 3)])
>>> other = LineString([(0, 0), (0, 4)])
>>> is_similar_azimuth(trace, other, 1)
True
>>> trace = LineString([(0, 0), (1, 1)])
>>> other = LineString([(0, 0), (0, 4)])
>>> is_similar_azimuth(trace, other, 40)
False
>>> is_similar_azimuth(trace, other, 50)
True
>>> is_similar_azimuth(trace, other, 45)
False
fractopo.tval.proximal_traces.is_within_buffer_distance(trace: LineString, other: LineString, buffer_value: float | int | Real) bool

Determine if trace and other are within buffer distance.

Threshold distance is buffer_value (both are buffered with half of buffer_value) of each other.

E.g.

>>> trace = LineString([(0, 0), (0, 3)])
>>> other = LineString([(0, 0), (0, 3)])
>>> is_within_buffer_distance(trace, other, 0.1)
True
>>> line = LineString([(0, 0), (0, 3)])
>>> other = LineString([(3, 0), (3, 3)])
>>> is_within_buffer_distance(trace, other, 2)
False
>>> is_within_buffer_distance(trace, other, 4)
True

fractopo.tval.trace_validation module

Contains main entrypoint class for validating trace data, Validation.

Create Validation objects from traces and their target areas to validate the traces for further analysis (fractopo.analysis.network.Network).

class fractopo.tval.trace_validation.Validation(traces: GeoDataFrame, area: GeoDataFrame, name: str, allow_fix: bool, SNAP_THRESHOLD: float = 0.01, SNAP_THRESHOLD_ERROR_MULTIPLIER: float = 1.1, AREA_EDGE_SNAP_MULTIPLIER: float = 1.5, TRIANGLE_ERROR_SNAP_MULTIPLIER: float = 10.0, OVERLAP_DETECTION_MULTIPLIER: float = 50.0, STACKED_DETECTOR_BUFFER_MULTIPLIER: float = 5.0, SHARP_AVG_THRESHOLD: float = 135.0, SHARP_PREV_SEG_THRESHOLD: float = 100.0, ERROR_COLUMN: str = 'VALIDATION_ERRORS', GEOMETRY_COLUMN: str = 'geometry', determine_validation_nodes: bool = True)

Bases: object

Validate traces data delineated by target area(s).

If allow_fix is True, some automatic fixing will be done to e.g. convert MultiLineStrings to LineStrings.

AREA_EDGE_SNAP_MULTIPLIER: float = 1.5
ERROR_COLUMN: str = 'VALIDATION_ERRORS'
GEOMETRY_COLUMN: str = 'geometry'
OVERLAP_DETECTION_MULTIPLIER: float = 50.0
SHARP_AVG_THRESHOLD: float = 135.0
SHARP_PREV_SEG_THRESHOLD: float = 100.0
SNAP_THRESHOLD: float = 0.01
SNAP_THRESHOLD_ERROR_MULTIPLIER: float = 1.1
STACKED_DETECTOR_BUFFER_MULTIPLIER: float = 5.0
TRIANGLE_ERROR_SNAP_MULTIPLIER: float = 10.0
allow_fix: bool
area: GeoDataFrame
determine_validation_nodes: bool = True
property endpoint_nodes: list[tuple[Point, ...]]

Get endpoints of all traces.

Returned as a list of tuples wherein each tuple represents the nodes of a trace in traces i.e. endpoint_nodes[index] are the nodes for traces[index].

property faulty_junctions: set[int] | None

Determine indexes with Multi Junctions.

property intersect_nodes: list[tuple[Point, ...]]

Get intersection nodes of all traces.

Returned as a list of tuples wherein each tuple represents the nodes of a trace in traces i.e. intersect_nodes[index] are the nodes for traces[index].

name: str
run_validation(first_pass: bool = True, choose_validators: tuple[type[BaseValidator]] | None = None, allow_empty_area: bool = True) GeoDataFrame

Run validation.

Main entrypoint for validation. Returns validated traces GeoDataFrame.

set_general_nodes()

Set _intersect_nodes and _endpoint_nodes attributes.

property spatial_index: Any | None

Get geopandas spatial index of traces.

traces: GeoDataFrame
property vnodes: set[int] | None

Determine indexes with V-Nodes.

fractopo.tval.trace_validation_utils module

Direct utilities of trace validation.

fractopo.tval.trace_validation_utils.determine_middle_in_triangle(segments: list[LineString], snap_threshold: float | int | Real, snap_threshold_error_multiplier: float | int | Real) list[LineString]

Determine the middle segment within a triangle error.

The middle segment always intersects the other two.

fractopo.tval.trace_validation_utils.determine_trace_candidates(geom: LineString, idx: int, traces: GeoDataFrame, spatial_index: SpatialIndex | None) GeoSeries

Determine potentially intersecting traces with spatial index.

fractopo.tval.trace_validation_utils.is_underlapping(geom: LineString, trace: LineString, endpoint: Point, snap_threshold: float | int | Real, snap_threshold_error_multiplier: float | int | Real) bool | None

Determine if a geom is underlapping.

fractopo.tval.trace_validation_utils.linestring_segment(linestring: LineString, dist: float | int | Real, threshold_length: float | int | Real) tuple[tuple[float | int | Real, float | int | Real], tuple[float | int | Real, float | int | Real]]

Get LineString segment from dist to dist + threshold_length.

fractopo.tval.trace_validation_utils.segment_within_buffer(linestring: LineString, multilinestring: MultiLineString, snap_threshold: float | int | Real, snap_threshold_error_multiplier: float | int | Real, overlap_detection_multiplier: float | int | Real, stacked_detector_buffer_multiplier: float | int | Real) bool

Check if segment is within buffer of multilinestring.

First check if given linestring completely overlaps any part of multilinestring and if it does, returns True.

Next it starts to segmentize the multilinestring to smaller linestrings and consequently checks if these segments are completely within a buffer made of the given linestring. It also checks that the segment size is reasonable.

TODO: segmentize_linestring is very inefficient.

fractopo.tval.trace_validation_utils.segmentize_linestring(linestring: LineString, threshold_length: float | int | Real) list[tuple[tuple[float | int | Real, float | int | Real], tuple[float | int | Real, float | int | Real]]]

Segmentize LineString to smaller parts.

Resulting parts are not guaranteed to be mergeable back to the original.

fractopo.tval.trace_validation_utils.split_to_determine_triangle_errors(trace: LineString, splitter_trace: LineString, snap_threshold: float | int | Real, triangle_error_snap_multiplier: float | int | Real) bool

Split trace with splitter_trace to determine triangle intersections.

fractopo.tval.trace_validators module

Contains Validator classes which each have their own error to handle and mark.

class fractopo.tval.trace_validators.BaseValidator

Bases: object

Base validator that all classes inherit.

TODO: Validation requires overhaul at some point.

ERROR = 'BASE ERROR'
INTERACTION_NODES_COLUMN = 'IP'
LINESTRING_ONLY = True
static fix_method(**_) LineString | None

Abstract fix method.

static validation_method(**_) bool

Abstract validation method.

class fractopo.tval.trace_validators.EmptyTargetAreaValidator

Bases: BaseValidator

Stub validator for empty target area.

Validation for this error occurs at the start of run_validation.

ERROR = 'EMPTY TARGET AREA'
class fractopo.tval.trace_validators.GeomNullValidator

Bases: BaseValidator

Validate the geometry for NULL GEOMETRY errors.

ERROR = 'NULL GEOMETRY'
LINESTRING_ONLY = False
static validation_method(geom: Any, **_) bool

Validate for empty and null geometries.

E.g. some validations are handled by GeomTypeValidator

>>> GeomNullValidator.validation_method(Point(1, 1))
True

Empty geometries are not valid.

>>> GeomNullValidator.validation_method(LineString())
False
>>> GeomNullValidator.validation_method(LineString([(-1, 1), (1, 1)]))
True
class fractopo.tval.trace_validators.GeomTypeValidator

Bases: BaseValidator

Validates the geometry type.

Validates that all traces are LineStrings. Tries to use shapely.ops.linemerge to merge MultiLineStrings into LineStrings.

ERROR = 'GEOM TYPE MULTILINESTRING'
LINESTRING_ONLY = False
static fix_method(geom: Any, **_) LineString | None

Fix mergeable MultiLineStrings to LineStrings.

E.g. mergeable MultiLineString

>>> mls = MultiLineString([((0, 0), (1, 1)), ((1, 1), (2, 2))])
>>> GeomTypeValidator.fix_method(geom=mls).wkt
'LINESTRING (0 0, 1 1, 2 2)'

Unhandled types will just return as None.

>>> GeomTypeValidator.fix_method(geom=Point(1, 1)) is None
True
static validation_method(geom: Any, **_) bool

Validate geometries.

E.g. Anything but LineString

>>> GeomTypeValidator.validation_method(Point(1, 1))
False

With LineString:

>>> GeomTypeValidator.validation_method(LineString([(0, 0), (1, 1)]))
True
class fractopo.tval.trace_validators.MultiJunctionValidator

Bases: BaseValidator

Validates that junctions consists of a maximum of two lines crossing.

ERROR = 'MULTI JUNCTION'
static determine_faulty_junctions(all_nodes: list[tuple[Point, ...]], snap_threshold: float | int | Real, snap_threshold_error_multiplier: float | int | Real) set[int]

Determine when a point of interest represents a multi junction.

I.e. when there are more than 2 points within the buffer distance of each other.

Two points is the limit because an intersection point between two traces is added to the list of interest points twice, once for each trace.

Faulty junction can also represent a slightly overlapping trace i.e. a snapping error.

>>> all_nodes = [
...     (Point(0, 0), Point(1, 1)),
...     (Point(1, 1),),
...     (Point(5, 5),),
...     (Point(0, 0), Point(1, 1)),
... ]
>>> snap_threshold = 0.01
>>> snap_threshold_error_multiplier = 1.1
>>> MultiJunctionValidator.determine_faulty_junctions(
...     all_nodes, snap_threshold, snap_threshold_error_multiplier
... )
{0, 1, 3}
static validation_method(idx: int, faulty_junctions: set[int], **_) bool

Validate for multi junctions between traces.

>>> MultiJunctionValidator.validation_method(
...     idx=1, faulty_junctions=set([1, 2])
... )
False
class fractopo.tval.trace_validators.MultipleCrosscutValidator

Bases: BaseValidator

Find traces that cross-cut each other multiple times.

This also indicates the possibility of duplicate traces.

ERROR = 'MULTIPLE CROSSCUTS'
static validation_method(geom: LineString, trace_candidates: GeoSeries, **_) bool

Validate for multiple crosscuts.

>>> geom = LineString([Point(-3, -4), Point(-3, -1)])
>>> trace_candidates = gpd.GeoSeries(
...     [
...         LineString(
...             [Point(-4, -3), Point(-2, -3), Point(-4, -2), Point(-2, -1)]
...         )
...     ]
... )
>>> MultipleCrosscutValidator.validation_method(geom, trace_candidates)
False
>>> geom = LineString([Point(-3, -4), Point(-3, -4.5)])
>>> trace_candidates = gpd.GeoSeries(
...     [
...         LineString(
...             [Point(-4, -3), Point(-2, -3), Point(-4, -2), Point(-2, -1)]
...         )
...     ]
... )
>>> MultipleCrosscutValidator.validation_method(geom, trace_candidates)
True
class fractopo.tval.trace_validators.SharpCornerValidator

Bases: BaseValidator

Find sharp cornered traces.

ERROR = 'SHARP TURNS'
static validation_method(geom: LineString, sharp_avg_threshold: float | int | Real, sharp_prev_seg_threshold: float | int | Real, **_) bool

Validate for sharp cornered traces.

class fractopo.tval.trace_validators.SimpleGeometryValidator

Bases: BaseValidator

Use shapely is_simple and is_ring attributes to validate LineString.

Checks that trace does not cut itself.

ERROR = 'CUTS ITSELF'
static validation_method(geom: LineString, **_) bool

Validate for self-intersections.

class fractopo.tval.trace_validators.StackedTracesValidator

Bases: BaseValidator

Find stacked traces and small triangle intersections.

ERROR = 'STACKED TRACES'
static validation_method(geom: LineString, trace_candidates: GeoSeries, snap_threshold: float | int | Real, snap_threshold_error_multiplier: float | int | Real, overlap_detection_multiplier: float | int | Real, triangle_error_snap_multiplier: float | int | Real, stacked_detector_buffer_multiplier: float | int | Real, **_) bool

Validate for stacked traces and small triangle intersections.

>>> geom = LineString([(0, 0), (0, 1)])
>>> trace_candidates = gpd.GeoSeries([LineString([(0, -1), (0, 2)])])
>>> snap_threshold = 0.01
>>> snap_threshold_error_multiplier = 1.1
>>> overlap_detection_multiplier = 50
>>> triangle_error_snap_multiplier = 10
>>> stacked_detector_buffer_multiplier = 5
>>> StackedTracesValidator.validation_method(
...     geom,
...     trace_candidates,
...     snap_threshold,
...     snap_threshold_error_multiplier,
...     overlap_detection_multiplier,
...     triangle_error_snap_multiplier,
...     stacked_detector_buffer_multiplier,
... )
False
>>> geom = LineString([(10, 0), (10, 1)])
>>> trace_candidates = gpd.GeoSeries([LineString([(0, -1), (0, 2)])])
>>> snap_threshold = 0.01
>>> snap_threshold_error_multiplier = 1.1
>>> overlap_detection_multiplier = 50
>>> triangle_error_snap_multiplier = 10
>>> stacked_detector_buffer_multiplier = 10
>>> StackedTracesValidator.validation_method(
...     geom,
...     trace_candidates,
...     snap_threshold,
...     snap_threshold_error_multiplier,
...     overlap_detection_multiplier,
...     triangle_error_snap_multiplier,
...     stacked_detector_buffer_multiplier,
... )
True
class fractopo.tval.trace_validators.TargetAreaSnapValidator

Bases: BaseValidator

Validator for traces that underlap the target area.

ERROR = 'TRACE UNDERLAPS TARGET AREA'
static is_candidate_underlapping(endpoint: Point, geom: LineString, area_polygon: Polygon | MultiPolygon, snap_threshold: float | int | Real) bool

Determine if endpoint is candidate for Underlapping error.

static simple_underlapping_checks(endpoint: Point, geom: LineString, area_polygon: Polygon | MultiPolygon, snap_threshold: float | int | Real) bool | None

Perform simple underlapping checks.

static validation_method(geom: LineString, area: GeoDataFrame | GeoSeries, snap_threshold: float | int | Real, snap_threshold_error_multiplier: float | int | Real, area_edge_snap_multiplier: float | int | Real, **_) bool

Validate for trace underlaps.

>>> geom = LineString([(0, 0), (0, 1)])
>>> area = gpd.GeoDataFrame(
...     geometry=[Polygon([(0, 0), (0, 1), (1, 1), (1, 0)])]
... )
>>> snap_threshold = 0.01
>>> snap_threshold_error_multiplier = 1.1
>>> area_edge_snap_multiplier = 1.0
>>> TargetAreaSnapValidator.validation_method(
...     geom,
...     area,
...     snap_threshold,
...     snap_threshold_error_multiplier,
...     area_edge_snap_multiplier,
... )
True
>>> geom = LineString([(0.5, 0.5), (0.5, 0.98)])
>>> area = gpd.GeoDataFrame(
...     geometry=[Polygon([(0, 0), (0, 1), (1, 1), (1, 0)])]
... )
>>> snap_threshold = 0.01
>>> snap_threshold_error_multiplier = 1.1
>>> area_edge_snap_multiplier = 10
>>> TargetAreaSnapValidator.validation_method(
...     geom,
...     area,
...     snap_threshold,
...     snap_threshold_error_multiplier,
...     area_edge_snap_multiplier,
... )
False
class fractopo.tval.trace_validators.UnderlappingSnapValidator

Bases: BaseValidator

Find snapping errors of underlapping traces.

Uses a multiple of the given snap_threshold.

ERROR = 'UNDERLAPPING SNAP'
classmethod validation_method(geom: LineString, trace_candidates: GeoSeries, snap_threshold: float | int | Real, snap_threshold_error_multiplier: float | int | Real, **_) bool

Validate for UnderlappingSnapValidator errors.

>>> snap_threshold = 0.01
>>> snap_threshold_error_multiplier = 1.1
>>> geom = LineString(
...     [
...         (0, 0),
...         (0, 1 + snap_threshold * snap_threshold_error_multiplier * 0.99),
...     ]
... )
>>> trace_candidates = gpd.GeoSeries([LineString([(-1, 1), (1, 1)])])
>>> UnderlappingSnapValidator.validation_method(
...     geom, trace_candidates, snap_threshold, snap_threshold_error_multiplier
... )
False
>>> snap_threshold = 0.01
>>> snap_threshold_error_multiplier = 1.1
>>> geom = LineString([(0, 0), (0, 1)])
>>> trace_candidates = gpd.GeoSeries([LineString([(-1, 1), (1, 1)])])
>>> UnderlappingSnapValidator.validation_method(
...     geom, trace_candidates, snap_threshold, snap_threshold_error_multiplier
... )
True
class fractopo.tval.trace_validators.VNodeValidator

Bases: BaseValidator

Finds V-nodes within trace data.

ERROR = 'V NODE'
static determine_v_nodes(endpoint_nodes: list[tuple[Point, ...]], snap_threshold: float | int | Real, snap_threshold_error_multiplier: float | int | Real) set[int]

Determine V-node errors.

>>> endpoint_nodes = [
...     (Point(0, 0), Point(1, 1)),
...     (Point(1, 1),),
... ]
>>> snap_threshold = 0.01
>>> snap_threshold_error_multiplier = 1.1
>>> VNodeValidator.determine_v_nodes(
...     endpoint_nodes, snap_threshold, snap_threshold_error_multiplier
... )
{0, 1}
static validation_method(idx: int, vnodes: set[int], **_) bool

Validate for V-nodes.

>>> VNodeValidator.validation_method(idx=1, vnodes=set([1, 2]))
False
>>> VNodeValidator.validation_method(idx=5, vnodes=set([1, 2]))
True

Module contents

Package with trace validation utility.