pyoints package

Submodules

pyoints.about module

pyoints.assertion module

Functions to ensure the properties of frequently used data structures.

pyoints.assertion.ensure_coords(coords, by_col=False, dim=None, min_dim=2, max_dim=inf)[source]

Ensures required properties of an array associated with coordinates.

Parameters:
  • coords (array_like(Number, shape=(n, k))) – Represents n data points of k dimensions in a Cartesian coordinate system.
  • by_col (optional, bool) – Indicates whether or not the coordinates are provided column by column instead of row by row.
  • dim,min_dim,max_dim (optional, positive int) – See ensure_dim.
Returns:

coords – Coordinates with guaranteed properties.

Return type:

np.ndarray(Number, shape=(n, k))

Raises:

TypeError, ValueError

Examples

Coordinates provided row by row.

>>> coords = ensure_coords([(3, 2), (2, 4), (-1, 2), (9, 3)])
>>> print(isinstance(coords, np.ndarray))
True
>>> print_rounded(coords)
[[ 3  2]
 [ 2  4]
 [-1  2]
 [ 9  3]]

Coordinates provided column by column.

>>> coords = ensure_coords([(3, 2, -1, 9), (2, 4, 2, 3)], by_col=True)
>>> print_rounded(coords)
[[ 3  2]
 [ 2  4]
 [-1  2]
 [ 9  3]]

See also

ensure_polar()

pyoints.assertion.ensure_dim(value, dim=None, min_dim=2, max_dim=inf)[source]

Ensure a dimension value to be in a specific range.

Parameters:
  • value (int) – Value representing a dimension.
  • min_dim, max_dim (dim,) – Minimum and maximum allowed dimensions. If dim is provided, the check_dim has to be exactly dim. If not, the check_dim must be in range [min_dim, max_dim].
Returns:

Dimension value with ensured properties.

Return type:

int

Raises:

ValueError

pyoints.assertion.ensure_indices(v, min_value=0, max_value=inf)[source]

Ensures an index array to be in a specific range.

Parameters:
  • v (array_like(int, shape=(n))) – Array of indices to check.
  • max_value (min_value,) – Minimum and maximum allowed value of v.
Returns:

Array of indices.

Return type:

np.ndarray(int, shape=(n))

Raises:

TypeError, ValueError

pyoints.assertion.ensure_json(js)[source]

Ensures the properties of a serializable json object.

Parameters:js (dict) – Dictionary to convert to a serializable json object.
Returns:Serializable json object.
Return type:dict
pyoints.assertion.ensure_length(value, length=None, min_length=0, max_length=inf)[source]

Ensure a length value to be in a specific range.

Parameters:
  • value (int) – Length value to check.
  • length,min_length,max_length (optional, positive int) – Minimum and maximum allowed length. If length is provided, check_length has to be exactly length. If not, the check_length must be in range [min_length, max_length].
Returns:

Length value with ensured properties.

Return type:

int

Raises:

ValueError

pyoints.assertion.ensure_numarray(arr, shape=None)[source]

Ensures the properties of an numeric numpy ndarray.

Parameters:arr (array_like(Number)) – Array like numeric object.
Returns:Array with guaranteed properties.
Return type:np.ndarray(Number)
Raises:TypeError, ValueError

Examples

>>> print_rounded(ensure_numarray([0,1,2]))
[0 1 2]
>>> print_rounded(ensure_numarray((-4,-5)))
[-4 -5]
pyoints.assertion.ensure_numvector(v, length=None, min_length=1, max_length=inf)[source]

Ensures the properties of a numeric vector.

Parameters:
  • v (array_like(Number, shape=(k))) – Vector of length n.
  • length,min_length,max_length (optional, positive int) – See ensure_length
Returns:

v – Vector with guaranteed properties.

Return type:

np.ndarray(Number, shape=(n))

Examples

Check a valid vector.

>>> v = (3, 2, 4, 4)
>>> v = ensure_numvector(v)
>>> print_rounded(v)
[3 2 4 4]

Vector of insufficient length.

>>> try:
...     ensure_numvector(v, length=5)
... except ValueError as e:
...     print(e)
length 5 required
Raises:TypeError, ValueError
pyoints.assertion.ensure_polar(pcoords, by_col=False, dim=None, min_dim=2, max_dim=inf)[source]

Ensures the properties of polar coordinates.

Parameters:
  • pcoords (array_like(Number, shape=(n,k))) – Represents n data points of k dimensions in a polar coordinate system.
  • by_col (optional, bool) – Defines whether or not the coordinates are provided column by column instead of row by row.
  • dim,min_dim,max_dim (optional, positive int) – See ensure_dim.
Raises:

TypeError, ValueError

Returns:

pcoords – Polar coordinates with guaranteed properties.

Return type:

np.ndarray(Number, shape=(n,k))

See also

ensure_coords()

pyoints.assertion.ensure_shape(shape, dim=None, min_dim=1, max_dim=inf)[source]

Ensures the properties of an array shape.

Parameters:shape (array_like(int, shape=(k))) – Shape of k dimensions to validate.
Returns:Shape with ensured properties.
Return type:np.ndarray(int, shape=(k))
Raises:ValueError, TypeError
pyoints.assertion.ensure_tmatrix(T, dim=None, min_dim=2, max_dim=inf)[source]

Ensures the properties of transformation matrix.

Parameters:
  • T (array_like(Number, shape=(k+1,k+1))) – Transformation matrix.
  • dim,min_dim,max_dim (optional, positive int) – See ensure_dim.
Returns:

T – Transformation matrix with guaranteed properties.

Return type:

np.matrix(Number, shape=(k+1,k+1))

Raises:

TypeError, ValueError

See also

transformation.matrix()

pyoints.assertion.iscoord(coord)[source]

Checks if a value can be associated with a coordinate.

Parameters:coord (array_like) – Value associated with a coordinate.
Returns:Indicates whether or not the value is a coordinate.
Return type:bool
pyoints.assertion.isnumeric(value, min_th=-inf, max_th=inf)[source]

Checks if a value is numeric.

Parameters:
  • value (Number) – Value to validate.
  • min_th,max_th (optional, Number) – Minimum and maximum value allowed range.
Returns:

Indicates whether or not the value is numeric.

Return type:

bool

pyoints.assign module

Module to find pairs of points

class pyoints.assign.KnnMatcher(coords, radii)[source]

Bases: pyoints.assign.Matcher

Find pairs of points. Each point is assigned to k closest points within a predefined sphere. Duplicate assignents are possible.

See also

Matcher

Examples

Create coordinates and matcher.

>>> A = np.array([(0, 0), (0, 0.1), (1, 1), (1, 0), (0.5, 0.5), (-1, -2)])
>>> B = np.array([(0.4, 0.4), (0.2, 0), (0.1, 1.2), (2, 1), (-1.1, -1.2)])
>>> matcher = KnnMatcher(A, [0.5, 0.5])

Try to assign one neighbour.

>>> pairs = matcher(B)
>>> print_rounded(pairs)
[[4 0]
 [0 1]]

Try to assign two neighbours.

>>> pairs = matcher(B, k=2)
>>> print_rounded(pairs)
[[4 0]
 [0 1]
 [1 1]]
class pyoints.assign.Matcher(coords, radii)[source]

Bases: object

Base class to simplify point matching. Points of a reference point set A are assigned to points of a point set B.

Parameters:
  • A (array_like(Number, shape=(n, k))) – Represents n points of k dimensions. These points are used as a reference point set.
  • radii (array_like(Number, shape=(k))) – Defines the sphere within the points can get assigned.
class pyoints.assign.PairMatcher(coords, radii)[source]

Bases: pyoints.assign.Matcher

Find unique pairs of points. A point a of point set A is assigned to its closest point b of point set B if a is also the nearest neighbour to b. Thus, duplicate assignment is not possible.

See also

Matcher

Examples

Create point sets.

>>> A = np.array([(0, 0), (0, 0.1), (1, 1), (1, 0), (0.5, 0.5), (-1, -2)])
>>> B = np.array([(0.4, 0.4), (0.2, 0), (0.1, 1.2), (2, 1), (-1.1, -1.2)])

Match points.

>>> matcher = PairMatcher(A, [0.3, 0.2])
>>> pairs = matcher(B)
>>> print_rounded(pairs)
[[4 0]
 [0 1]]
>>> print_rounded(A[pairs[:, 0], :] - B[pairs[:, 1], :])
[[ 0.1  0.1]
 [-0.2  0. ]]
class pyoints.assign.SphereMatcher(coords, radii)[source]

Bases: pyoints.assign.Matcher

Find pairs of points. Each point is assigned to all the points within a previously defined sphere. Duplicate assignments are possible.

See also

Matcher

Examples

Create point sets.

>>> A = np.array([(0, 0), (0, 0.1), (1, 1), (1, 0), (0.5, 0.5), (-1, -2)])
>>> B = np.array([(0.4, 0.4), (0.2, 0), (0.1, 1.2), (2, 1), (-1.1, -1.2)])

Match points.

>>> matcher = SphereMatcher(A, [0.3, 0.2])
>>> pairs = matcher(B)
>>> print_rounded(pairs)
[[4 0]
 [0 1]
 [1 1]]
>>> print_rounded(A[pairs[:, 0], :] - B[pairs[:, 1], :])
[[ 0.1  0.1]
 [-0.2  0. ]
 [-0.2  0.1]]

pyoints.classification module

Collection of functions to classify or reclassify values or cluster values.

pyoints.classification.classes_to_dict(classification, ids=None, min_size=1, max_size=inf, missing_value=-1)[source]

Converts a list of class indices to a dictionary of grouped classes.

Parameters:
  • classification (array_like(shape=(n))) – Array of class indices.
  • ids (optional, array_like(int, shape=(n))) – Indices to specify a subset of classification. If none, the indices are numbered consecutively.
  • min_size,max_size (optional, positive int) – Minimum and maximum desired size of a class to be kept in the result.
  • missing_value (optional, object) – Default value for unclassified values.
Returns:

Dictionary of class indices. The dictionary keys represent the class ids, while the values represent the indices in the original array.

Return type:

dict

See also

dict_to_classes()
Dictionary representation of classification.

Examples

>>> classes = ['cat', 'cat', 'dog', 'bird', 'dog', 'bird', 'cat', 'dog']
>>> class_dict = classes_to_dict(classes)
>>> print(sorted(class_dict))
['bird', 'cat', 'dog']
>>> print_rounded(class_dict['cat'])
[0 1 6]
>>> classes = [0, 0, 1, 2, 1, 0, 3, 3, 5, 3, 2, 1, 0]
>>> print(classes_to_dict(classes))
{0: [0, 1, 5, 12], 1: [2, 4, 11], 2: [3, 10], 3: [6, 7, 9], 5: [8]}
pyoints.classification.dict_to_classes(classes_dict, n, min_size=1, max_size=inf, missing_value=-1)[source]

Converts a dictionary of classes into a list of classes.

Parameters:
  • classes_dict (dict) – Dictionary of class indices.
  • n (positive int) – Desired size of the output array. It must be at least the size of the maximum class index.
  • min_size,max_size (optional, positive int) – Minimum and maximum desired size of a class to be kept in the result.
  • missing_value (optional, object) – Default value for unclassified values.
Returns:

Array representation of classes_dict.

Return type:

np.ndarray(int, shape=(n))

Notes

Only a minimal input validation is provided.

Examples

Alphanumeric classes.

>>> classes_dict = {'bird': [0, 1, 5, 4], 'dog': [3, 6, 8], 'cat': [7]}
>>> print(dict_to_classes(classes_dict, 10, missing_value=''))
['bird' 'bird' '' 'dog' 'bird' 'bird' 'dog' 'cat' 'dog' '']

Omit small classes.

>>> print(dict_to_classes(classes_dict, 10, min_size=2))
['bird' 'bird' -1 'dog' 'bird' 'bird' 'dog' -1 'dog' -1]

Numeric classes.

>>> classes_dict = {0: [0, 1, 5], 1: [3, 6], 2: [7, 2]}
>>> print(dict_to_classes(classes_dict, 9))
[0 0 2 1 -1 0 1 2 -1]
pyoints.classification.majority(classes, empty_value=-1)[source]

Finds most frequent class or value in an array.

Parameters:
  • classes (array_like(object, shape=(n))) – Classes or values to check.
  • empty_value (optional, object) – Class value in case that no decision can be made.
Returns:

Most frequent class.

Return type:

object

Notes

Only a limited input validation is provided.

Examples

Find majority class.

>>> classes =['cat', 'dog', 'dog', 'bird', 'cat', 'dog']
>>> print(majority(classes))
dog
>>> classes =[1, 8, 9, 0, 0, 2, 4, 2, 4, 3, 2, 3, 5, 6]
>>> print_rounded(majority(classes))
2

No decision possible.

>>> classes =[1, 2, 3, 4, 4, 3]
>>> print_rounded(majority(classes))
-1
pyoints.classification.rename_dict(d, ids=None)[source]

Assigns new key names to a dictionary.

Parameters:
  • d (dict) – Dictionary to rename.
  • ids (optional, array_like(shape=(len(d)))) – Desired key names. If none, the keys are numbered consecutively.
Returns:

Dictionary with new names.

Return type:

dict

Examples

>>> d = {1: [0, 1], 2: None, 3: 'text'}
>>> renamed_dict = rename_dict(d, ['first', 'second', 'last'])
>>> print(sorted(renamed_dict))
['first', 'last', 'second']
pyoints.classification.split_by_breaks(values, breaks)[source]

Classifies values by ranges.

Parameters:
  • values (array_like(Number, shape=(n))) – Values to classify.
  • breaks (array_like(Number, shape=(m))) – Series of value ranges.
Returns:

classification – Desired class affiliation of values. A value of classification[i] equal to k means that ‘values[i]’ is in range [breaks[k], breaks[k][

Return type:

np.ndarray(int, shape=(n))

Examples

>>> values = np.arange(10)
>>> breaks = [0.5, 5, 7.5]
>>> classes = split_by_breaks(values, breaks)
>>> print_rounded(classes)
[0 1 1 1 1 2 2 2 3 3]

pyoints.clustering module

Clustering algorithms to assign group points.

pyoints.clustering.clustering(indexKD, r, get_class, order=None, clusters=None, auto_set=True)[source]

Generic clustering based on spatial neighbourhood.

Parameters:
  • indexKD (IndexKD) – Spatial index with n points.
  • r (positive float) – Radius to identify the cluster affiliation of neighbored points.
  • get_class (callable) – Function to define the cluster id (affiliation) of a point. It receives a list of cluster ids of neigboured points to define the cluster id of selected point. It returns -1 if the point is not associated with any cluster.
  • order (optional, array_like(int)) – Defines the order to apply the clustering algorithm. It can also be used to subsample points for clustering. If None, the order is defined by decreasing point density.
  • clusters (optional, array_like(int, shape=(n))) – List of n integers. Each element represents the preliminary cluster id of a point in indexKD. A cluster id of -1 represents no class.
  • auto_set (optional, bool) – Defines whether or not a cluster id is set automatically if -1 (no class) was returned by get_class. If True, a new cluster id is set to max(clusters) + 1.
Returns:

Dictionary of clusters. The keys correspond to the class ids. The values correspond to the point indices associated with the cluster.

Return type:

dict

pyoints.clustering.dbscan(indexKD, min_pts, epsilon=None, quantile=0.8, factor=3)[source]

Variant of the DBSCAN algorithm [1] with automatic estimation of the epsilon parameter using point density. Useful for automatic outlier identification.

Parameters:
  • indexKD (IndexKD) – Spatial index with n points to cluster.
  • min_pts (int) – Corresponds to the min_pts parameter of the DBSCAN algorithm.
  • epsilon (optional, positive float) – Corresponds to the epsilon parameter of DBSCAN algorithm. If None, a suitable value is estimated by investigating the nearest neighbour distances dists of all points in indexKD with `epsilon = np.percentile(dists, quantile * 100) * factor`.
  • quantile (optional, positive float) – Used to calculate epsilon.
  • factor (optional, positive float) – Used to calculate epsilon.

References

[1] M. Ester, et al. (1996): “A Density-Based Algorithm for Discovering Clusters in Large Spatial Databases with Noise”, KDD-96 Proceedings, pp. 226-231.

Examples

>>> coords = [(0, 0), (0, 1), (1, 1), (0, 0.5), (2, 2), (2, 2.5), (19, 29)]
>>> indexKD = IndexKD(coords)

User defined epsilon.

>>> clusters = dbscan(indexKD, 1, epsilon=1)
>>> print_rounded(clusters)
[0 0 0 0 1 1 2]

Automatic epsilon estimation for outlier removal.

>>> clusters = dbscan(indexKD, 2)
>>> print_rounded(clusters)
[ 0  0  0  0  0  0 -1]

Adjust automatic epsilon estimation to achieve small clusters.

>>> clusters = dbscan(indexKD, 1, quantile=0.7, factor=1)
>>> print_rounded(clusters)
[0 0 1 0 2 2 3]
pyoints.clustering.majority_clusters(indexKD, r, **kwargs)[source]

Clustering by majority voting. The algorithm assigns points iteratively to the most dominant class within a given radius.

Parameters:
  • indexKD (IndexKD) – Spatial index with n points.
  • r (positive float) – Radius to identify the cluster affiliation of neighbored points.
  • **kwargs (optional) – Optional arguments of the clustering function.

See also

clustering()

Examples

>>> coords = [(0, 0), (1, 1), (2, 1), (3, 3), (0, 1), (2, 3), (3, 4)]
>>> clusters = majority_clusters(IndexKD(coords), 2)
>>> print_rounded(clusters)
[ 1  1 -1  0  1  0  0]
pyoints.clustering.weight_clusters(indexKD, r, weights=None, **kwargs)[source]

Clustering by class weight.

Parameters:
  • indexKD (IndexKD) – Spatial index with n points.
  • r (positive float) – Radius to identify the cluster affiliation of neighbored points.
  • weights (optional, array_like(Number, shape=(len(indexKD)))) – Weighting of each point. The class with highest weight wins. If None, all weights are set to 1, which results in similar results than majority_clusters.
  • **kwargs (optional) – Optional arguments passed to clustering.

Examples

Clustering with equal weights.

>>> coords = [(0, 0), (0, 1), (1, 1), (0, 0.5), (2, 2), (2, 2.5), (2.5, 2)]
>>> indexKD = IndexKD(coords)
>>> initial_clusters = np.arange(len(coords), dtype=int)
>>> clusters = weight_clusters(indexKD, 1.5, clusters=initial_clusters)
>>> print_rounded(clusters)
[0 0 4 3 6 5 5]

Clustering with differing weights.

>>> weights = np.arange(len(coords))
>>> clusters = weight_clusters(
...     indexKD,
...     1.5,
...     weights=weights,
...     clusters=initial_clusters
... )
>>> print_rounded(clusters)
[3 6 6 3 5 5 5]

pyoints.coords module

Data structures to handle multi-dimensional point data.

class pyoints.coords.Coords[source]

Bases: numpy.ndarray, object

Class to represent coordinates. It provides a spatial index to conveniently optimize spatial neighborhood queries. The spatial index is calculated on demand only and deleted, if coordinates have been changed.

Parameters:coords (array_like(Number)) – Represents n data points of k dimensions in a Cartesian coordinate system. Any desired shape of at least length two is allowed to enable the representation point, voxel or raster data. The last shape element always represents the coordinate dimension.
dim

Number of coordinate dimensions.

Type:positive int
flattened

List representation of the coordinates.

Type:array_like(Number, shape=(n, k))

Examples

Create some 2D points.

>>> coords = Coords([(0, 1), (2, 1), (2, 3), (0, -1)])
>>> print_rounded(coords)
[[ 0  1]
 [ 2  1]
 [ 2  3]
 [ 0 -1]]

Get Extent.

>>> print_rounded(coords.extent())
[ 0 -1  2  3]
>>> print_rounded(coords.extent(dim=1))
[0 2]

Use IndexKD, which updates automatically, if a coordinate has changed.

>>> coords.indexKD().ball([0, 0], 2)
[0, 3]
>>> coords[1, 0] = -1
>>> coords.indexKD().ball([0, 0], 2)
[0, 1, 3]
dim
extent(dim=None)[source]

Provides the spatial extent of the coordinates.

Parameters:dim (optional, positive int) – Define how many coordinate dimensions to use. If None, all dimensions are used.
Returns:extent – Spatial extent of the coordinates.
Return type:Extent

Notes

The extents are calculated on demand and are cached automatically. Setting new coordinates clears the cache.

See also

Extent()

Examples

>>> coords = Coords([(2, 3), (3, 2), (0, 1), (-1, 2.2), (9, 5)])
>>> print_rounded(coords.extent())
[-1.  1.  9.  5.]
>>> print_rounded(coords.extent(dim=1))
[-1.  9.]
>>> print_rounded(coords.indexKD().ball([0, 0], 2))
[2]
flattened
indexKD(dim=None)[source]

Get a spatial index of the coordinates.

Parameters:dim (optional, positive int) – Desired dimension of the spatial index. If None, the all coordinate dimensions are used.
Returns:Spatial index of the coordinates of desired dimension.
Return type:IndexKD

Notes

The spatial indices are generated on demand and are cached automatically. Setting new coordinates clears the cache.

See also

poynts.IndexKD()

Examples

>>> coords = Coords([(2, 3, 1), (3, 2, 3), (0, 1, 0), (9, 5, 4)])
>>> print_rounded(coords.indexKD().dim)
3
>>> print_rounded(coords.indexKD(dim=2).dim)
2
transform(T)[source]

Transform coordinates.

Parameters:T (array_like(Number, shape=(self.dim+1, self.dim+1))) – Transformation matrix to apply.
Returns:transformed coords.
Return type:Coords(shape=self.shape)

Examples

Transform structured coordinates.

>>> coords = Coords(
...             [[(2, 3), (2, 4), (3, 2)], [(0, 0), (3, 5), (9, 4)]])
>>> print_rounded(coords)
[[[2 3]
  [2 4]
  [3 2]]
<BLANKLINE>
 [[0 0]
  [3 5]
  [9 4]]]
>>> T = transformation.matrix(t=[10, 20], s=[0.5, 1])
>>> tcoords = coords.transform(T)
>>> print_rounded(tcoords)
[[[ 11.   23. ]
  [ 11.   24. ]
  [ 11.5  22. ]]
<BLANKLINE>
 [[ 10.   20. ]
  [ 11.5  25. ]
  [ 14.5  24. ]]]

pyoints.distance module

Various distance metrics.

pyoints.distance.dist(p, coords)[source]

Calculates the distances between points.

Parameters:
  • p (array_like(Number, shape=(n, k)) or array_like(Number, shape=(k))) – Represents n points or a single point of k dimensions.
  • coords (array_like(Number, shape=(n, k))) – Represents n points of k dimensions.
Returns:

Normed values.

Return type:

Number or array_like(Number, shape=(n))

See also

sdist()

Examples

Point to points distance.

>>> p = (1, 2)
>>> coords = [(2, 2), (1, 1), (1, 2), (9, 8)]
>>> print_rounded(dist(p, coords))
[  1.   1.   0.  10.]

Points to points distance.

>>> A = [(2, 2), (1, 1), (1, 2)]
>>> B = [(4, 2), (2, 1), (9, 8)]
>>> print_rounded(dist(A, B))
[  2.   1.  10.]
pyoints.distance.idw(dists, p=2)[source]

Calculates the weights for Inverse Distance Weighting method.

Parameters:
  • dists (Number or array_like(Number, shape=(n))) – Represent n distance values.
  • p (optional, Number) – Weighting power.
Returns:

Weights according to Inverse Distance Weighting.

Return type:

Number or array_like(Number, shape=(n))

Examples

>>> dists = [0, 1, 4]
>>> print_rounded(idw(dists))
[ 1.    0.25  0.04]
>>> print_rounded(idw(dists, p=1))
[ 1.   0.5  0.2]
pyoints.distance.norm(coords)[source]

Normalization of coordinates.

Parameters:coords (array_like(Number, shape=(k)) or array_like(Number, shape=(n, k))) – Represents n points or a single point of k dimensions.
Returns:Normed values.
Return type:array_like(shape=(n, ))

See also

snorm()

Examples

>>> coords = [(3, 4), (0, 1), (4, 3), (0, 0), (8, 6)]
>>> print_rounded(norm(coords))
[  5.   1.   5.   0.  10.]
pyoints.distance.rmse(A, B=None)[source]

Calculates the Root Mean Squared Error of corresponding data sets.

Parameters:B (A,) – Represent n points or a single point of k dimensions.
Returns:Root Mean Squared Error.
Return type:Number

Examples

>>> A = [(2, 2), (1, 1), (1, 2)]
>>> B = [(2.2, 2), (0.9, 1.1), (1, 2.1)]
>>> print_rounded(rmse(A, B))
0.15
pyoints.distance.sdist(p, coords)[source]

Calculates the squared distances between points.

Parameters:
  • p (array_like(Number, shape=(n, k)) or array_like(Number, shape=(k))) – Represents n points or a single point of k dimensions.
  • coords (array_like(Number, shape=(n, k))) – Represents n points of k dimensions.
Returns:

Squared distances between the points.

Return type:

Number or array_like(Number, shape=(n))

See also

dist()

Examples

Squared point to points distance.

>>> p = (1, 2)
>>> coords = [(2, 4), (1, 1), (1, 2), (9, 8)]
>>> print_rounded(sdist(p, coords))
[  5   1   0 100]

Squared points to points distance.

>>> A = [(2, 2), (1, 1), (1, 2)]
>>> B = [(4, 2), (2, 1), (9, 8)]
>>> print_rounded(sdist(A, B))
[  4   1 100]
pyoints.distance.snorm(coords)[source]

Squared normalization of coordinates.

Parameters:coords (array_like(Number, shape=(k)) or array_like(Number, shape=(n, k))) – Represents n points or a single point of k dimensions.
Returns:Squared normed values.
Return type:Number or array_like(Number, shape=(n))

See also

norm()

Examples

>>> coords = [(3, 4), (0, 1), (4, 3), (0, 0), (8, 6)]
>>> print_rounded(snorm(coords))
[ 25   1  25   0 100]

pyoints.extent module

Handles spatial extents.

class pyoints.extent.Extent[source]

Bases: numpy.recarray, object

Specifies spatial extent (or bounding box) of coordinates in k dimensions.

Parameters:ext (array_like(Number, shape=(2 * k)) or array_like(Number, shape=(n, k))) – Defines spatial extent of k dimensions as either minimum corner and maximum corner or as a set of n points. If a set of points is given, the bounding box of these coordinates is calculated.
dim

Number of coordinate dimensions.

Type:positive int
ranges

Ranges between each coordinate dimension.

Type:np.ndarray(Number, shape=(self.dim))
min_corner,max_corner

Minimum and maximum values in each coordinate dimension.

Type:array_like(Number, shape=(self.dim))
center

Focal point of the extent.

Type:array_like(Number, shape=(self.dim))

Examples

Derive the extent of a list of points.

>>> points = [(0, 0), (1, 4), (0, 1), (1, 0.5), (0.5, 0.7)]
>>> ext = Extent(points)
>>> print(ext)
[ 0.  0.  1.  4.]

Create a extent based on minimum and maximum values.

>>> ext = Extent([-1, 0, 1, 4, ])
>>> print(ext)
[-1  0  1  4]

Get some properties.

>>> print_rounded(ext.dim)
2
>>> print_rounded(ext.min_corner)
[-1  0]
>>> print_rounded(ext.max_corner)
[1 4]
>>> print_rounded(ext.ranges)
[2 4]
>>> print(ext.center)
[ 0.  2.]
>>> print_rounded(ext.corners)
[[-1  0]
 [ 1  0]
 [ 1  4]
 [-1  4]]
center
corners

Provides each corner of the extent box.

Returns:corners – Corners of the extent.
Return type:np.ndarray(Number, shape=(2**self.dim, self.dim))

Examples

Two dimensional case.

>>> ext = Extent([-1, -2, 1, 2])
>>> print_rounded(ext.corners)
[[-1 -2]
 [ 1 -2]
 [ 1  2]
 [-1  2]]

Three dimensional case.

>>> ext = Extent([-1, -2, -3, 1, 2, 3])
>>> print_rounded(ext.corners)
[[-1 -2 -3]
 [ 1 -2 -3]
 [ 1  2 -3]
 ...,
 [ 1  2  3]
 [ 1 -2  3]
 [-1 -2  3]]
dim
intersection(coords, dim=None)[source]

Tests if coordinates are located within the extent.

Parameters:
  • coords (array_like(Number, shape=(n, k)) or) –
  • shape=(k)) (array_like(Number,) – Represents n data points of k dimensions.
  • dim (positive int) – Desired number of dimensions to consider.
Returns:

indices – Indices of coordinates which are within the extent. If just a single point is given, a boolean value is returned.

Return type:

np.ndarray(int, shape=(n)) or np.ndarray(bool, shape=(n))

Examples

Point within extent?

>>> ext = Extent([0, 0.5, 1, 4])
>>> print(ext.intersection([(0.5, 1)]))
True

Points within extent?

>>> print_rounded(ext.intersection([(1, 2), (-1, 1), (0.5, 1)]))
[0 2]

Corners are located within the extent.

>>> print_rounded(ext.intersection(ext.corners))
[0 1 2 3]
max_corner
min_corner
ranges
split()[source]

Splits the extent into the minimum and maximum corners.

Returns:min_corner,max_corner – Minimum and maximum values in each coordinate dimension.
Return type:np.ndarray(Number, shape=(self.dim))

pyoints.filters module

Point filters.

pyoints.filters.ball(indexKD, r, order=None, inverse=False, axis=-1, min_pts=1)[source]

Filters coordinates by radius. This algorithm is suitable to remove duplicate points or to get an almost uniform point density.

Parameters:
  • indexKD (IndexKD) – IndexKD containing n points to filter.
  • r (positive float or array_like(float, shape=(n))) – Ball radius or radii to select neighbouring points.
  • order (optional, array_like(int, shape=(m))) – Order to proceed. If m < n, only a subset of points is investigated. If none, ordered by axis.
  • axis (optional, int) – Coordinate axis to use to generate the order.
  • inverse (bool) – Indicates whether or not the order is inversed.
  • min_pts (optional, int) – Specifies how many neighbouring points within radius r shall be required to yield a filtered point. This parameter can be used to filter noisy point sets.
Yields:

positive int – Indices of the filtered points.

Notes

Within a dense point cloud, the filter guarantees the distance of neighboured points in a range of ]r, 2*r[.

pyoints.filters.dem_filter(coords, r, max_angle=70)[source]

Selects points suitable for generating a digital elevation model.

Parameters:
  • coords (array_like(Number, shape=(n, k))) – Represents n points of k dimensions to filter.
  • r (Number or array_like(Number, shape=(n))) – Radius or radii to apply.
  • max_angle (Number) – Maximum allowed slope of a simplex.
Returns:

Desired indices of points suitable for generating a surface model.

Return type:

array_like(int, shape=(m))

Notes

Applies a k-1 dimensional ball filter to identify a digital elevation model. The point order is defined by the last coordinate dimension. The filtering radius r is defined by the user. To optimize the quality of the dem, all points with less than 5 neighbours within radius 2 * r are removed from further analysis. In addition a TIN (Triangulated Irregular Network) is generated to model the surface. Each facet with a slope larger than max_angle is removed iteratively.

pyoints.filters.extrema(indexKD, attributes, r, inverse=False)[source]

Finds local maxima or minima of given point set.

Parameters:
  • indexKD (IndexKD) – Spatial index of n points to filter.
  • attributes (array_like(Number, shape=(n))) – Attributes to search for extrema. If none, the last coordinate dimension is used as the attribute.
  • r (positive Number) – Maximum distance between two points to be considered as neighbours.
  • inverse (optional, bool) – Indicates if local maxima (False) or local minima (True) shall be yielded.
Yields:

positive int – Indices of local maxima or minima.

Examples

Find local maxima.

>>> indexKD = IndexKD([(0, 0), (0, 1), (1, 1), (1, 0), (0.5, 0.5) ])
>>> attributes = [2, 0.1, 1, 0, 0.5]
>>> fIds = list(extrema(indexKD, attributes, 1.1))
>>> print_rounded(fIds)
[0 2]

Find local minima.

>>> fIds = list(extrema(indexKD, attributes, 1.1, inverse=True))
>>> print_rounded(fIds)
[3 1]
>>> fIds = list(extrema(indexKD, attributes, 1.5, inverse=True))
>>> print_rounded(fIds)
[3]
pyoints.filters.has_neighbour(indexKD, r)[source]

Filters points which have neighbours within a specific radius.

Parameters:
  • indexKD (IndexKD) – Spatial index of n points to filter.
  • r (positive Number) – Maximum distance of a point to a neighbouring point to be still considered as isolated.
Yields:

positive int – Index of a point with at least one neighbouring point with radius r.

See also

is_isolated()

Examples

>>> coords = [(0, 0), (0.5, 0.5), (0, 1), (0.7, 0.5), (-1, -1)]
>>> indexKD = IndexKD(coords)
>>> print_rounded(list(has_neighbour(indexKD, 0.7)))
[1 3]
pyoints.filters.in_convex_hull(hull_coords, coords)[source]

Decides whether or not points are within a convex hull.

Parameters:
  • hull_coords (array_like(Number, shape=(m, k))) – Represents m hull points of k dimensions.
  • coords (array_like(Number, shape=(n, k))) – Represents n data points to check whether or not they are located within the convex hull.
Returns:

Indicates whether or not the points are located within the convex hull.

Return type:

array_like(bool, shape=(n))

Examples

>>> hull = [(0, 0), (1, 0), (1, 2)]
>>> coords = [(0, 0), (0.5, 0.5), (0, 1), (0.7, 0.5), (-1, -1)]
>>> print(in_convex_hull(hull, coords))
[ True  True False  True False]
pyoints.filters.is_isolated(indexKD, r)[source]

Filters points which have no neighbours within a specific radius.

Parameters:
  • indexKD (IndexKD) – Spatial index of n points to filter.
  • r (positive Number) – Maximum distance of a point to a neighbouring point to be still considered as isolated.
Yields:

positive int – Index of an isolated point with no neighbouring point within radius r.

See also

has_neighbour()

Examples

>>> coords = [(0, 0), (0.5, 0.5), (0, 1), (0.7, 0.5), (-1, -1)]
>>> indexKD = IndexKD(coords)
>>> print_rounded(list(is_isolated(indexKD, 0.7)))
[0 2 4]
pyoints.filters.min_filter(indexKD, attributes, r, inverse=False)[source]

Finds minima or maxima within a specified radius for all points. For each point the neighbouring point with extreme attribute is marked.

Parameters:
  • indexKD (IndexKD) – Spatial index of n points to filter.
  • attributes (array_like(Number, shape=(n))) – Attributes to search for extrema. If none, the last coordinate dimension is used as the attribute.
  • r (positive Number) – Local radius to search for a minimum.
  • inverse (optional, bool) – Indicates if minima (False) or maxima (True) shall be identfied.
Returns:

Filtered point indices.

Return type:

np.ndarray(int, shape=(n))

pyoints.filters.radial_dem_filter(coords, angle_res, center=None, max_angle=70)[source]

Filters surface points based on distance to the center. The algorithm is designed to create digital elevation models of terrestrial laser scans. Terrestrial laser scans are characterized by decreasing point densities with increasing distance to the scanner position. Thus, the radius to identify neighbouring points is adjusted to the distance of the scanner.

Parameters:
  • coords (array_like(Number, shape=(n, k))) – Points to filter.
  • angle_res (positive Number) – Filter resolution expressed as an angle.
  • center (optional, array_like(Number, shape=(k))) – Desired center.
  • max_angle (Number) – Maximum allowed angle of a simplex.
Returns:

Indices of filtered points.

Return type:

array_like(int, shape=(m))

See also

dem_filter()

pyoints.filters.surface(indexKD, r, order=None, inverse=False, axis=-1)[source]

Filters points associated with a surface.

Parameters:
  • indexKD (IndexKD) – IndexKD containing n points to filter.
  • r (positive float or array_like(float, shape=(n))) – Ball radius or radii to apply.
  • order (optional, array_like(int, shape=(m))) – Order to proceed. If m < n, only a subset of points is investigated.
  • inverse (optional, bool) – Indicates whether or not to inverse the order.
  • axis (optional, int) – Axis to use for generating the order.
Yields:

positive int – Indices of the filtered points.

pyoints.fit module

Fits shapes or functions to points.

pyoints.fit.fit_cylinder(coords, vec=None)[source]

Fits a cylinder to points.

Parameters:
  • coords (array_like(Number, shape=(n, k))) – Represents n data points of k dimensions.
  • vec (optional, array_like(Number, shape(k))) – Estimated orientation of the cylinder axis.
Returns:

  • vec (vector.Vector) – Orientaton vector.
  • r (Number) – Radius of the cylinder.
  • resid (Number) – Remaining residuals.

Examples

Prepare roto-translated half cylinder.

>>> r = 2.5
>>> x = np.arange(-1, 0, 0.01) * r
>>> y = np.sqrt(r**2 - x**2)
>>> y[::2] = - y[::2]
>>> z = np.repeat(5, len(x))
>>> z[::2] = -5
>>> T = transformation.matrix(t=[10, 20, 30], r=[0.3, 0.2, 0.0])
>>> coords = transformation.transform(np.array([x, y, z]).T, T)

Get cylinder.

>>> vec, r, residuals = fit_cylinder(coords, vec=[0, 0, 1])
>>> print_rounded(r)
2.5
>>> print_rounded(vec.origin)
[ 10.  20.  30.]

Check distances to vector.

>>> dists = vec.distance(coords)
>>> print_rounded([np.min(dists), np.max(dists)])
[ 2.5  2.5]
pyoints.fit.fit_sphere(coords, weights=1.0)[source]

Least square fitting of a sphere to a set of points.

Parameters:
  • coords (array_like(Number, shape=(n, k))) – Represents n data points of k dimensions.
  • weights ((k+1,k+1), array_like) – Transformation matrix.
Returns:

  • center (np.ndarray(Number, shape=(k))) – Center of the sphere.
  • r (Number) – Radius of the sphere.

Notes

Idea taken from [1].

References

[1] A. Bruenner (2001): URL http://www.arndt-bruenner.de/mathe/scripts/kreis3p.htm

Examples

Draw points on a half circle with radius 5 and center (2, 4) and try to determine the circle parameters.

>>> x = np.arange(-1, 1, 0.1)
>>> y = np.sqrt(5**2 - x**2)
>>> coords = np.array([x,y]).T + [2,4]
>>> center, r, residuals = fit_sphere(coords)
>>> print_rounded(center)
[ 2.  4.]
>>> print_rounded(r)
5.0

pyoints.georecords module

Generic data structure to handle point data.

class pyoints.georecords.GeoRecords[source]

Bases: numpy.recarray, object

Abstraction class to ease handling point sets as well as structured matrices of point like objects. This gives the opportunity to handle unstructured point sets the same way as rasters or voxels. The class also provides an IndexKD object on demand to speed up neigbourhood analyses.

Parameters:
  • proj (projection.Proj) – Projection object to provide the coordinate reference system.
  • rec (np.recarray) – A numpy record array provides coordinates and attributes of n points. The field coords is required and represents coordinates of k dimensions.
  • T (optional, array_like(Number, shape=(k+1, k+1))) – Linear transformation matrix to transform the coordinates. Set to a identity matrix if not given. This eases handling structured point sets like rasters or voxels.
proj

Projection of the coordinates.

Type:projection.Proj
t

Linear transformation matrix to transform the coordinates.

Type:np.ndarray(Number, shape=(k+1, k+1))
dim

Number of coordinate dimensions of the coords field.

Type:positive int
count

Number of objects within the data structure. E.g. number of points of a point cloud or number of cells of a raster.

Type:positive int
keys

Keys or indices of the data structure.

Type:np.ndarray(int, shape=self.shape)
date

Date of capture.

Type:datetime

Examples

>>> data = {
...    'coords': [(2, 3), (3, 2), (0, 1), (-1, 2.2), (9, 5)],
...    'values': [1, 3, 4, 0, 6]
... }
>>> geo = GeoRecords(projection.Proj(), data)
>>> print_rounded(geo.values)
[1 3 4 0 6]
>>> print_rounded(geo.coords)
[[ 2.   3. ]
 [ 3.   2. ]
 [ 0.   1. ]
 [-1.   2.2]
 [ 9.   5. ]]

Set new coordinates.

>>> geo['coords'] = [(1, 2), (9, 2), (8, 2), (-7, 3), (7, 8)]
>>> print_rounded(geo.coords)
[[ 1.  2.]
 [ 9.  2.]
 [ 8.  2.]
 [-7.  3.]
 [ 7.  8.]]

Use structured data (two dimensional matrix).

>>> data = {
...    'coords': [
...                 [(2, 3.2), (-3, 2.2)],
...                 [(0, 1.1), (-1, 2.2)],
...                 [(-7, -1), (9.2, -5)]
...             ],
...    'values': [[1, 3], [4, 0], [-4, 2]]
... }
>>> data = nptools.recarray(data,dim=2)
>>> geo = GeoRecords(None, data)
>>> print_rounded(geo.shape)
(3, 2)
>>> print_rounded(geo.coords)
[[[ 2.   3.2]
  [-3.   2.2]]
<BLANKLINE>
 [[ 0.   1.1]
  [-1.   2.2]]
<BLANKLINE>
 [[-7.  -1. ]
  [ 9.2 -5. ]]]
add_fields(dtypes, data=None)[source]

Adds data fields to the georecords.

Parameters:
  • dtypes (np.dtype) – Data types of the new fields.
  • data (optional, list of arrays) – Data values of the new fields.
Returns:

geo – New GeoRecords with additional fields.

Return type:

self.__class__

See also

nptools.add_fields()

Examples

Unstructured data.

>>> data = {
...    'coords': [(2, 3, 1), (3, 2, 3), (0, 1, 0), (9, 5, 4)],
... }
>>> geo = GeoRecords(None, data)
>>> geo2 = geo.add_fields(
...     [('A', int), ('B', object)],
...     data=[[1, 2, 3, 4], ['a', 'b', 'c', 'd']]
... )
>>> print(geo2.dtype.names)
('coords', 'A', 'B')
>>> print(geo2.B)
['a' 'b' 'c' 'd']
apply_function(func, dtypes=[<class 'object'>])[source]

Applies or maps a function to each element of the data array.

Parameters:
  • func (function) – This function is applied to record of the array.
  • dtypes (optional, np.dtype) – Desired data type of the output array.
Returns:

geo – New GeoRecords with shape self.shape.

Return type:

self.__class__

See also

nptools.apply_function()

coords
count
date
dim
extent(*args, **kwargs)[source]
indexKD(*args, **kwargs)[source]
keys

Keys or indices of the data structure.

Returns:Array of indices. Each entry provides an index tuple e.g. to receive data elements with __getitem__.
Return type:np.ndarray(int, shape=(self.count))

Examples

Unstructured data.

>>> data = {
...    'coords': [(2, 3, 1), (3, 2, 3), (0, 1, 0), (9, 5, 4)],
...    'values': [1, 3, 4, 0]
... }
>>> geo = GeoRecords(None, data)
>>> print_rounded(geo.keys)
[[0]
 [1]
 [2]
 [3]]

Structured data (two dimensional matrix).

>>> data = np.ones(
...         (4,3), dtype=[('coords', float, 2)]).view(np.recarray)
>>> geo = GeoRecords(None, data)
>>> print_rounded(geo.keys)
[[[0 0]
  [0 1]
  [0 2]]
<BLANKLINE>
 [[1 0]
  [1 1]
  [1 2]]
<BLANKLINE>
 [[2 0]
  [2 1]
  [2 2]]
<BLANKLINE>
 [[3 0]
  [3 1]
  [3 2]]]
merge(rec)[source]

Merges a record array with the georecords.

Parameters:rec (array_like) – Record array with same fields as self.
Returns:geo – New GeoRecords.
Return type:self.__class__

See also

nptools.merge()

proj
project(proj)[source]

Projects the coordinates to a different coordinate system.

Parameters:proj (Proj) – Desired output projection system.
Returns:geo – New GeoRecords with updated spatial reference.
Return type:self.__class__

See also

Proj()

Examples

TODO

records()[source]

Provides the flattened data records. Useful if structured data (like matrices) are used.

Returns:Flattened
Return type:array_like(shape=(self.count), dtype=self.dtype)

Examples

>>> data = {
...    'coords': [
...                 [(2, 3.2), (-3, 2.2)],
...                 [(0, 1.1), (-1, 2.2)],
...                 [(-7, -1), (9.2, -5)]
...             ],
... }
>>> data = nptools.recarray(data, dim=2)
>>> geo = GeoRecords(None, data)
>>> print_rounded(geo.shape)
(3, 2)
>>> print_rounded(geo.records().coords)
[[ 2.   3.2]
 [-3.   2.2]
 [ 0.   1.1]
 [-1.   2.2]
 [-7.  -1. ]
 [ 9.2 -5. ]]
>>> geo.coords[0, 0] = (1, 1)
t
transform(T)[source]

Transforms coordinates.

Parameters:T (array_like(Number, shape=(self.dim+1, self.dim+1))) – Transformation matrix to apply.
Returns:
Return type:self

See also

Coords.transform()

Examples

>>> data = {
...    'coords': [(2, 3), (3, 2), (0, 1), (9, 5)],
...    'values': [1, 3, 4, 0]
... }
>>> geo = GeoRecords(None, data)
>>> T = transformation.matrix(t=[10, 20], s=[0.5, 1])
>>> _ = geo.transform(T)
>>> print_rounded(geo.coords)
[[11 23]
 [11 22]
 [10 21]
 [14 25]]
class pyoints.georecords.LasRecords[source]

Bases: pyoints.georecords.GeoRecords

Data structure extending GeoRecords to provide an optimized API for LAS data.

last_return

Array indicating if a point is a last return point.

Type:np.ndarray(bool)
first_return

Array indicating if a point is a first return point.

Type:np.ndarray(bool)
only_return

Array indicating if a point is the only returned point.

Type:np.ndarray(bool)
total_gps_seconds

Total gps time in seconds. Requires the fields ‘gps_time’ and ‘date’ to be set.

Type:np.ndarray(np.float32) or None

See also

GeoRecords

CUSTOM_FIELDS = [('coords', <class 'float'>, 3), ('classification', <class 'numpy.uint8'>), ('num_returns', <class 'numpy.uint8'>), ('return_num', <class 'numpy.uint8'>), ('synthetic', <class 'bool'>), ('keypoint', <class 'bool'>), ('withheld', <class 'bool'>), ('scan_direction_flag', <class 'bool'>), ('edge_of_flight_line', <class 'bool'>)]
EXTRA_FIELDS = [('normals', <class 'float'>, 3)]
STANDARD_FIELDS = [('user_data', <class 'numpy.uint8'>), ('intensity', <class 'numpy.uint8'>), ('pt_src_id', <class 'numpy.uint16'>), ('gps_time', <class 'numpy.float64'>), ('scan_angle_rank', <class 'numpy.int8'>), ('red', <class 'numpy.uint8'>), ('green', <class 'numpy.uint8'>), ('blue', <class 'numpy.uint8'>), ('nir', <class 'numpy.uint8'>)]
activate(field_name)[source]

Activates a desired field on demand.

Parameters:field_name (String) – Name of the field to activate.
static available_fields()[source]
class_indices(*classes)[source]

Filters by classes.

Parameters:*classes (int) – Classes to filter by.
Returns:Filtered record indices.
Return type:np.ndarray(int)
first_return
grd()[source]

Filters by points classified as ground.

Returns:Filtered records.
Return type:LasRecords
last_return
only_return
total_gps_seconds
veg()[source]

Filters by points classified as vegetation.

Returns:Filtered records.
Return type:LasRecords

pyoints.indexkd module

A generic spatial index.

class pyoints.indexkd.IndexKD(coords, T=None, leafsize=16, quickbuild=True, copy=True)[source]

Bases: object

Wrapper class for several spatial indices to speed up spatial queries and ease the usage.

Parameters:
  • coords (array_like(Number, shape=(n, k))) – Represents n data points of k dimensions in a Cartesian coordinate system.
  • T (optional, array_like(Number, shape=(k+1, k+1))) – Represents any kind of transformation matrix applied to the coordinates before index computation.
  • leafsize (optional, positive int) – Leaf size of KD-Tree.
  • quickbuild (optional, bool) – Indicates whether or not the spatial index shall be optimized for quick building (True) or quick spatial queries (False).
  • copy (optional, bool) – Indicates whether or not to copy the coordinate array.
dim

Number of coordinate dimensions k.

Type:positive int
t

Transformation matrix.

Type:np.matrix(Number, shape=(k+1, k+1))
coords

Coordinates of the spatial index.

Type:np.ndarray(Number, shape=(n, k))
kd_tree

KD-tree for rapid neighbourhood queries. Generated on first demand.

Type:scipy.spatial.cKDTree
r_tree

R-tree for rapid box queries. Generated on first demand.

Type:rtree.Rtree

Notes

Most spatial index operations are time critical. So it is usually avoided to check each input parameter in detail.

Examples

Create a simple spatial index.

>>> coords = np.indices((5, 10)).reshape((2, 50)).T
>>> indexKD = IndexKD(coords)
>>> len(indexKD)
50
>>> indexKD.dim
2

Query points within a sphere.

>>> nids = indexKD.ball([0, 0], 1.1)
>>> print_rounded(nids)
[ 0  1 10]
>>> print_rounded(indexKD.coords[nids, :])
[[0 0]
 [0 1]
 [1 0]]

Scale the coordinates using a transformation matrix to enable queries in the shape of an ellipsoid.

>>> T = [(0.5, 0, 0), (0, 0.8, 0), (0, 0, 1)]
>>> indexKD = IndexKD(coords, T)
>>> nids = indexKD.ball([0, 0], 1.1)
>>> print_rounded(nids)
[ 0  1 20 10 11]
>>> print_rounded(indexKD.coords[nids, :])
[[ 0.   0. ]
 [ 0.   0.8]
 [ 1.   0. ]
 [ 0.5  0. ]
 [ 0.5  0.8]]
ball(coords, r, bulk=100000, **kwargs)[source]

Finds all points within distance r of point or points coords.

Parameters:
  • coords (array_like(Number, shape=(n, k))) – Represents n data points of k dimensions.
  • r (positive float) – Radius of ball.
  • bulk (optional, positive int) – Reduces required memory by performing bulk queries.
  • **kwargs (optional) – Additional parameters passed to scipy.spatial.cKDTree.query_ball_point
Returns:

nIds – If coords is a single point, this returns a list of neighbours. If coords is a list of points, this returns a list containing the lists of neighbours.

Return type:

list or array of lists

Examples

>>> coords = np.indices((5, 10)).reshape((2, 50)).T
>>> indexKD = IndexKD(coords)
>>> indexKD.ball((0, 0), 1)
[0, 1, 10]
>>> indexKD.ball(np.array([(0, 0), (1, 1)]), 1)
[[0, 1, 10], [1, 10, 11, 12, 21]]
ball_count(r, coords=None, **kwargs)[source]

Counts numbers of neighbours within radius.

Parameters:
  • r (float or iterable of float) – Iterable radii to query.
  • coords (optional, array_like(Number, shape=(n, k)) or iterable) – Represents n points of k dimensions. If none, it is set to self.coords.
  • **kwargs (optional) – Additional parameters passed to scipy.spatial.cKDTree.query_ball_point
Returns:

Number of neigbours for each point.

Return type:

numpy.ndarray(int, shape=(n))

Examples

>>> coords = [(0, 0), (0, 1), (1, 1), (2, 1), (1, 0.5), (0.5, 1)]
>>> indexKD = IndexKD(coords)
>>> counts = indexKD.ball_count(1)
>>> print_rounded(counts)
[2 4 5 2 3 4]
>>> counts = indexKD.ball_count(1, coords=[0.5, 0.5])
>>> print_rounded(counts)
5
ball_count_iter(r, coords=None, **kwargs)[source]

Counts numbers of neighbours within radius.

Parameters:
  • r (float or iterable of float) – Iterable radii to query.
  • coords (optional, array_like(Number, shape=(n, k)) or iterable) – Represents n points of k dimensions. If none, it is set to self.coords.
  • **kwargs (optional) – Additional parameters passed to scipy.spatial.cKDTree.query_ball_point
Yields:

int – Number of neigbours for each point.

ball_iter(coords, r, bulk=10000, **kwargs)[source]

Similar to ball, but yields lists of neighbours.

Yields:nIds (list of int) – Lists of indices of neighbouring points.

See also

ball(), balls_iter()

balls_iter(coords, radii, **kwargs)[source]

Similar to ball_iter, but with differing radii.

Parameters:radii (iterable of float) – Radii to query.
Yields:nIds (list) – Lists of indices of neighbouring points.

See also

ball(), ball_iter()

box(extent)[source]

Selects points within a given extent.

Parameters:extent (array_like(Number, shape=(2 * self.dim))) – Specifies the points to return. A point p is returned, if np.all(p <= extent[0: dim]) and np.all(p >= extent[dim+1: 2*dim])
Returns:Indices of points within the extent.
Return type:list of int

See also

ball(), cube(), slice()

box_count(extent)[source]

Counts all points within a given extent.

Parameters:extent (array_like(Number, shape=(2*self.dim))) – Specifies the points to return. A point p is returned, if p <= extent[0:dim] and p >= extent[dim+1:2*dim]
Returns:Number of points within the extent.
Return type:list of int

See also

box()

closest(ids, **kwargs)[source]

Provides nearest neighbour of a point.

Parameters:
  • ids (int or array_like(int, shape=(m))) – Index of point or indices of points in self.coords.
  • **kwargs (optional) – Additional parameters similar to scipy.spatial.cKDTree.query.
  • Returns
  • --------
  • distances (np.ndarray(Number, shape=(n))) – Distance to closest point.
  • indices (np.ndarray(int, shape=(n))) – Index of closest point.

See also

knn()

Examples

>>> coords = np.indices((5, 10)).reshape((2, 50)).T
>>> indexKD = IndexKD(coords)
>>> indexKD.closest(3)
(1.0, 4)
>>> print_rounded(indexKD.closest([0, 2, 5, 3])[1])
[1 1 6 4]
coords
cube(coords, r, **kwargs)[source]

Provides points within a cube.

Notes

Wrapper for self.ball with p=np.inf.

See also

ball()

dim
kd_tree
knn(coords, k=1, bulk=100000, **kwargs)[source]

Query for k nearest neighbours.

Parameters:
  • coords (array_like(Number, shape=(n, k))) – Represents n points of k dimensions.
  • k (optional, positive int) – Number of nearest neighbours to return.
  • **kwargs (optional) – Additional parameters passed to scipy.spatial.cKDTree.query_ball_point
Returns:

  • dists (np.ndarray(Number, shape=(n, k))) – Distances to nearest neighbours.
  • indices (np.ndarray(int, shape=(n, k))) – Indices of nearest neighbours.

Examples

>>> coords = [(0, 0), (0, 1), (1, 1), (2, 1), (1, 0.5), (0.5, 1)]
>>> indexKD = IndexKD(coords)
>>> dists, nids = indexKD.knn((0.5, 1), 2)
>>> print_rounded(dists)
[ 0.   0.5]
>>> print_rounded(nids)
[5 2]
>>> dists, nids = indexKD.knn([(0.5, 1), (1.5, 1)], 2)
>>> print_rounded(dists)
[[ 0.   0.5]
 [ 0.5  0.5]]
>>> print_rounded(nids)
[[5 2]
 [3 2]]
>>> dists, nids = indexKD.knn([(0.5, 1), (1.5, 1), (1, 1)], [3, 1, 2])
>>> print(dists)
(array([ 0. ,  0.5,  0.5]), array([ 0.5]), array([ 0. ,  0.5]))
>>> print(nids)
(array([5, 1, 2]), array([2]), array([2, 4]))

See also

knn(), knns_iter(), scipy.spatial.cKDTree.query()

knn_iter(coords, k=1, bulk=100000, **kwargs)[source]

Similar to knn, but yields lists of neighbours.

Yields:
  • dists (np.ndarray(Number, shape=(k))) – List of distances to nearest neighbours.
  • indices (np.ndarray(int, shape=(k))) – List of ??? and corresponding point indices.

See also

knn(), knns_iter()

knns_iter(coords, ks, **kwargs)[source]

Similar to knn, but yields lists of neighbours.

Parameters:

ks (iterable of int) – Iterable numbers of neighbours to query.

Yields:
  • dists (np.ndarray(Number, shape=(k))) – List of distances to nearest neighbours.
  • indices (np.ndarray(int, shape=(k))) – List of ??? and corresponding point indices.

See also

knn(), knn_iter()

nn

Provides the nearest neighbours for each point.

Returns:
  • distances (np.ndarray(Number, shape=(n))) – Distances to each nearest neighbour.
  • indices (np.ndarray(int, shape=(n))) – Indices of nearest neighbours.
r_tree
slice(min_th=-inf, max_th=inf, axis=-1)[source]

Selects points with coordinate value of axis axis within the range of [min_th, max_th].

Parameters:
  • min_th,max_th (Number) – A point p is returned if min_th <= p[axis] <= max_th.
  • axis (int) – Axis to evaluate.
Returns:

Indices of points within the slice.

Return type:

list of int

See also

box(), cube(), slice()

sphere(coord, r_min, r_max, **kwargs)[source]

Counts numbers of neighbours within radius.

Parameters:
  • r_min (float) – Inner radius of the sphere.
  • r_max (float) – Outer radius of the sphere.
  • coord (array_like(Number, shape=(k))) – Center of sphere.
  • **kwargs (optional) – Additional parameters passed to scipy.spatial.cKDTree.query_ball_point
Returns:

Indices of points within sphere.

Return type:

list of int

Examples

>>> coords = np.indices((5, 10)).reshape((2, 50)).T
>>> indexKD = IndexKD(coords)
>>> print_rounded(indexKD.ball((3, 3), 1))
[32 34 33 43 23]
>>> print_rounded(indexKD.ball((3, 3), 1.5))
[22 42 32 34 33 43 44 23 24]
>>> print_rounded(indexKD.sphere((3, 3), 1, 1.5))
[23 32 33 34 43]
t

pyoints.interpolate module

Spatial interpolation methods.

class pyoints.interpolate.Interpolator(coords, values)[source]

Bases: object

Abstract generic multidimensional interpolation class.

Parameters:
  • coords (array_like(Number, shape=(n, k))) – Represents n data points of k dimensions.
  • values (array_like(Number, shape=(n)) or array_like(Number, shape=(n, m))) – One dimensional or m dimensional values to interpolate.
coords

Provided coordinates.

Type:np.ndarray(Number, shape=(n, k))
dim

Number of coordinate dimensions.

Type:positive int
coords
dim
class pyoints.interpolate.KnnInterpolator(coords, values, k=None, max_dist=None)[source]

Bases: pyoints.interpolate.Interpolator

Nearest neighbor interpolation.

Parameters:
  • coords (array_like(Number, shape=(n, k))) – Represents n data points of k dimensions.
  • values (array_like(Number, shape=(n)) or array_like(Number, shape=(n, m))) – One dimensional or m dimensional values to interpolate.
  • k (optional, positive int) – Number of neighbors used for interpolation.
  • max_dist (optional, positive Number) – Maximum distance of a neigbouring point to be used for interpolation.

See also

Interpolator, sklearn.neighbors.KNeighborsRegressor

Examples

Interpolation of one-dimensional values.

>>> coords = [(0, 0), (0, 2), (2, 1)]
>>> values = [0, 3, 6]
>>> interpolator = KnnInterpolator(coords, values, k=2, max_dist=1)
>>> print(interpolator([0, 1]))
1.5
>>> print(interpolator([(0, 1), (0, 0), (0, -1)]))
[ 1.5  0.   0. ]

Interpolation of multi-dimensional values.

>>> coords = [(0, 0), (0, 2), (2, 3)]
>>> values = [(0, 1, 3), (2, 3, 8), (6, 4, 4)]
>>> interpolator = KnnInterpolator(coords, values, k=2)
>>> print(interpolator([0, 1]))
[ 1.   2.   5.5]
>>> print(interpolator([(0, 1), (0, 0), (0, -1)]))
[[ 1.    2.    5.5 ]
 [ 0.    1.    3.  ]
 [ 0.5   1.5   4.25]]
class pyoints.interpolate.LinearInterpolator(coords, values)[source]

Bases: pyoints.interpolate.Interpolator

Linear interpolation using Delaunay triangilation.

Parameters:
  • coords (array_like(Number, shape=(n, k))) – Represents n data points of k dimensions.
  • values (array_like(Number, shape=(n)) or array_like(Number, shape=(n, m))) – One dimensional or m dimensional values to interpolate.

See also

Interpolator, scipy.interpolate.LinearNDInterpolator

Examples

Interpolation of one-dimensional values.

>>> coords = [(0, 0), (0, 2), (2, 1)]
>>> values = [0, 3, 6]
>>> interpolator = LinearInterpolator(coords, values)
>>> print_rounded(interpolator([0, 1]))
1.5
>>> print_rounded(interpolator([(0, 1), (0, 0), (0, -1)]))
[ 1.5  0.   nan]

Interpolation of multi-dimensional values.

>>> coords = [(0, 0), (0, 2), (2, 3)]
>>> values = [(0, 1, 3), (2, 3, 8), (6, 4, 4)]
>>> interpolator = LinearInterpolator(coords, values)
>>> print_rounded(interpolator([0, 1]))
[ 1.   2.   5.5]
>>> print_rounded(interpolator([(0, 1), (0, 0), (0, -1)]))
[[ 1.   2.   5.5]
 [ 0.   1.   3. ]
 [ nan  nan  nan]]
delaunay
simplex(coords)[source]

Finds the corresponding simplex to a coordinate.

Parameters:coords (array_like(Number, shape=(n, k))) – Represents n points of k dimensions.
Returns:Simplex indices. If no simplex is found, an empty list is returned.
Return type:list(int) or list of np.ndarray(int shape=(k+1))

Examples

>>> coords = [(0, 0), (0, 2), (3, 0), (3, 2), (1, 1)]
>>> values = [0, 3, 6, 5, 9]
>>> interpolator = LinearInterpolator(coords, values)
>>> simplex_ids = interpolator.simplex([0, 1])
>>> print(simplex_ids)
[4 1 0]
>>> print(interpolator.coords[simplex_ids, :])
[[1 1]
 [0 2]
 [0 0]]
>>> simplex_ids = interpolator.simplex([(0.5, 0.5), (-1, 1), (3, 2)])
>>> print(interpolator.coords[simplex_ids[0], :])
[[1 1]
 [0 2]
 [0 0]]
>>> print(interpolator.coords[simplex_ids[1], :])
[]
>>> print(interpolator.coords[simplex_ids[2], :])
[[3 2]
 [1 1]
 [3 0]]
class pyoints.interpolate.PolynomInterpolator(coords, values, deg=2, weights=None, interaction_only=False)[source]

Bases: pyoints.interpolate.Interpolator

Polynomial interpolation.

Parameters:
  • coords (array_like(Number, shape=(n, k))) – Represents n data points of k dimensions.
  • values (array_like(Number, shape=(n)) or array_like(Number, shape=(n, m))) – One dimensional or m dimensional values to interpolate.
  • deg (optional, positive int) – The degree of the polynomial features.
  • weights (optional, array_like(shape=(n))) – Weights for each sample.
  • interaction_only (optional, bool) – Indicates whether or not to calculate interaction only.

See also

Interpolator, sklearn.preprocessing.PolynomialFeatures, sklearn.linear_model.LinearRegression

Examples

Interpolation of one-dimensional values.

>>> coords = [(0, 0), (0, 2), (2, 1)]
>>> values = [0, 3, 6]
>>> interpolator = PolynomInterpolator(coords, values, deg=1)
>>> print_rounded(interpolator([0, 1]))
1.5
>>> print_rounded(interpolator([(0, 1), (0, 0), (0, -1)]))
[ 1.5  0.  -1.5]

Interpolation of multi-dimensional values.

>>> coords = [(0, 0), (0, 2), (2, 3)]
>>> values = [(0, 1, 3), (2, 3, 8), (6, 4, 4)]
>>> interpolator = PolynomInterpolator(coords, values, deg=1)
>>> print_rounded(interpolator([0, 1]))
[ 1.   2.   5.5]
>>> print_rounded(interpolator([(0, 1), (0, 0), (0, -1)]))
[[ 1.   2.   5.5]
 [ 0.   1.   3. ]
 [-1.   0.   0.5]]
coef
intercept

pyoints.misc module

Some random functions, which ease development.

pyoints.misc.get_size(obj, seen=None)[source]

Recursively finds size of objects

pyoints.misc.list_licences(requirements_file)[source]
pyoints.misc.print_object_size(obj)[source]

Get the size of cached objects.

Parameters:obj (object) – Object to determine size from.
pyoints.misc.print_rounded(values, decimals=2)[source]

Prints rounded values.

Parameters:
  • values (Number or array_like(Number)) – Values to display.
  • decimals (optional, int) – Number of decimals passed to ‘numpy.round’.

Notes

Sets negative values close to zero to zero.

pyoints.misc.sizeof_fmt(num, suffix='B')[source]

Notes

Taken form [1]. Originally posted by [2].

References

[1] jan-glx (2018), https://stackoverflow.com/questions/24455615/python-how-to-display-size-of-all-variables, (acessed: 2018-08-16) [2] Fred Cirera, https://stackoverflow.com/a/1094933/1870254 (acessed: 2018-08-16)

pyoints.misc.tic()[source]
pyoints.misc.toc()[source]

pyoints.normals module

Derive normals of point clouds.

pyoints.normals.approximate_normals(coords, r=inf, k=None, preferred=None)[source]

Approximates normals of points by selecting nearest neighbours and assigning the derived normal to all neighbours.

Parameters:
  • coords (array_like(Number, shape=(n, dim))) – Represents n points of dim dimensions.
  • r (optional, positive Number) – Maximum radius to select neighbouring points.
  • k (optional, positive int) – Specifies the number of neighbours to select. If specified, at least dim neighbours are required.
  • preferred (optional, array_like(Number, shape=(k)) or array_like(Number, shape=(n, k))) – Preferred orientation of the normals.
Returns:

Approximated normals for the coordinates coords.

Return type:

array_like(Number, shape=(n, dim))

Examples

>>> coords = [(0, 0), (1, 1), (2, 3), (3, 3), (4, 2), (5, 1), (5, 0)]

Approximate two normals using all neighbours within radius n.

>>> normals = approximate_normals(coords, 2.5, preferred=(1, 0))
>>> print_rounded(normals)
[[ 0.71 -0.71]
 [ 0.45 -0.89]
 [ 0.45 -0.89]
 [ 0.45 -0.89]
 [ 0.88  0.47]
 [ 0.88  0.47]
 [ 0.88  0.47]]

Approximate two normals using k nearest neighbours.

>>> normals = approximate_normals(coords, k=4, preferred=(1, 0))
>>> print_rounded(normals)
[[ 0.76 -0.65]
 [ 0.76 -0.65]
 [ 0.59  0.81]
 [ 0.81  0.59]
 [ 0.81  0.59]
 [ 0.81  0.59]
 [ 0.81  0.59]]

Approximate two normals using k nearest neighbours within radius r.

>>> normals = approximate_normals(coords, k=4, r=2.5, preferred=(1, 0))
>>> print_rounded(normals)
[[ 0.71 -0.71]
 [ 0.45 -0.89]
 [ 0.45 -0.89]
 [ 0.45 -0.89]
 [ 0.88  0.47]
 [ 0.88  0.47]
 [ 0.88  0.47]]
pyoints.normals.fit_normals(coords, r=inf, k=None, indices=None, preferred=None)[source]

Fits normals to points by selecting nearest neighbours.

Parameters:
  • coords (array_like(Number, shape=(n, dim))) – Represents n points of dim dimensions.
  • r (optional, positive Number) – Maximum radius to select neighbouring points.
  • k (optional, positive int) – Specifies the number of neighbours to select. If specified, at least dim neighbours are required.
  • indices (optional, array_like(int, shape=(m))) – Vector of point indices to subsample the point cloud (m <= n). If None, indices is set to range(n).
  • preferred (optional, array_like(Number, shape=(k)) or array_like(Number, shape=(n, k))) – Preferred orientation of the normals.
Returns:

Desired normals of coordinates coords.

Return type:

array_like(Number, shape=(m, k))

Examples

>>> coords = [(0, 0), (1, 1), (2, 3), (3, 3), (4, 2), (5, 1), (5, 0)]

Fit normals using k nearest neighbours.

>>> normals = fit_normals(coords, k=2, preferred=[1, 0])
>>> print_rounded(normals, 2)
[[ 0.71 -0.71]
 [ 0.71 -0.71]
 [ 0.    1.  ]
 [ 0.    1.  ]
 [ 0.71  0.71]
 [ 1.    0.  ]
 [ 1.    0.  ]]

Fit normals a using all nearest neighbours within radius r.

>>> normals = fit_normals(coords, r=2.5, preferred=[1, 0])
>>> print_rounded(normals, 2)
[[ 0.71 -0.71]
 [ 0.84 -0.54]
 [ 0.45 -0.89]
 [ 0.47  0.88]
 [ 0.71  0.71]
 [ 0.88  0.47]
 [ 0.88  0.47]]

Fit normals using k nearest neighbours within radius r.

>>> normals = fit_normals(coords, r=2.5, k=3, preferred=[1, 0])
>>> print_rounded(normals)
[[ 0.71 -0.71]
 [ 0.84 -0.54]
 [ 0.76 -0.65]
 [ 0.47  0.88]
 [ 0.71  0.71]
 [ 0.88  0.47]
 [ 0.88  0.47]]
pyoints.normals.normalize(vectors, dim=None, n=None)[source]

Make a vector or vectors a set of normals with a given shape.

Parameters:
  • vectors (array_like(Number, shape=(k)) or array_like(Number, shape=(n, k))) – Orientation of the normals.
  • shape ((n, k)) – Desired shape of the output normals.
Returns:

Normals with ensured properties.

Return type:

np.ndarray(Number, shape=(n, k))

Examples

Normalize two dimensional vectors.

>>> vectors = [(3, 4), (8, 6), (2, 0), (0, 0)]
>>> normals = normalize(vectors)
>>> print_rounded(normals)
[[ 0.6  0.8]
 [ 0.8  0.6]
 [ 1.   0. ]
 [ 0.   0. ]]
>>> print_rounded(distance.norm(normals))
[ 1.  1.  1.  0.]

Normalize three dimensional vectors.

>>> vectors = [(3, 0, 4), (2, 0, 0), (0, 0, 0)]
>>> normals = normalize(vectors)
>>> print_rounded(normals)
[[ 0.6  0.   0.8]
 [ 1.   0.   0. ]
 [ 0.   0.   0. ]]
>>> print_rounded(distance.norm(normals))
[ 1.  1.  0.]

Normalize individual vectors.

>>> print_rounded(normalize([3, 4]))
[ 0.6  0.8]
>>> print_rounded(normalize((3, 0, 4)))
[ 0.6  0.   0.8]
pyoints.normals.prefer_orientation(normals, preferred)[source]

Orients normals using a prefered normals.

Parameters:
  • normals (array_like(Number, shape=(n, k))) – Normals of n points with k dimensions.
  • preferred (array_like(Number, shape=(k)) or array_like(Number, shape=(n, k))) – Preferred normal orientation for each normal in normals.
Returns:

Oriented normals. If the angle between a normal in normals and the corresponding normal in preferred is greater than 90 degree, the it is flipped.

Return type:

np.ndarray(Number, shape=(n, k))

pyoints.nptools module

Functions for convenient handling of numpy arrays.

pyoints.nptools.add_fields(arr, dtypes, data=None)[source]

Adds additional fields to a numpy record array.

Parameters:
  • arr (np.recarray) – Numpy record array to add fields to.
  • dtypes (np.dtype) – Data types of the new fields.
  • data (optional, list of array_like) – Data values of the new fields. The shape of each array has to be compatible to arr.
Returns:

Record array similar to A, but with additional fields of type dtypes and values of data.

Return type:

np.recarray

Examples

>>> A = recarray({'a': [0, 1, 2, 3]})
>>> C = add_fields(A, [('b', float, 2), ('c', int)])
>>> print(sorted(C.dtype.names))
['a', 'b', 'c']
>>> D = add_fields(A, [('d', int), ('e', str)], data=[[1, 2, 3, 4], None])
>>> print(D)
[(0, 1, '') (1, 2, '') (2, 3, '') (3, 4, '')]
pyoints.nptools.apply_function(arr, func, dtype=None)[source]

Applies a function to each record of a numpy array.

Parameters:
  • arr (np.ndarray or np.recarray) – Numpy array to apply function to.
  • func (function) – Function to apply to each record.
  • dtypes (optional, np.dtype) – Desired data type of the output array.
Returns:

Record array similar to input array, but with function applied to.

Return type:

np.recarray

Examples

Apply a function to a numpy ndarray.

>>> arr = np.ones((2, 3), dtype=[('a', int), ('b', int)])
>>> func = lambda item: item[0] + item[1]
>>> print(apply_function(arr, func))
[[2 2 2]
 [2 2 2]]

Aggregate a one dimensional numpy record array.

>>> data = { 'a': [0, 1, 2, 3], 'b': [1, 2, 3, 4] }
>>> arr = recarray(data)
>>> func = lambda record: record.a + record.b
>>> print(apply_function(arr, func))
[1 3 5 7]

Two dimensional case.

>>> data = { 'a': [[0, 1], [2, 3]], 'b': [[1, 2], [3, 4]] }
>>> arr = recarray(data, dim=2)
>>> func = lambda record: record.a ** record.b
>>> print(apply_function(arr, func))
[[ 0  1]
 [ 8 81]]

Specify the output data type.

>>> func = lambda record: (record.a + record.b, record.a ** record.b)
>>> print(apply_function(arr, func, dtype=[('c', float), ('d', int)]))
[[( 1.,  0) ( 3.,  1)]
 [( 5.,  8) ( 7., 81)]]

Specify a multi-dimensional output data type.

>>> func = lambda record: (record.a + 2, [record.a ** 2, record.b * 3])
>>> print(apply_function(arr, func, dtype=[('c', float), ('d', int, 2)]))
[[( 2., [ 0,  3]) ( 3., [ 1,  6])]
 [( 4., [ 4,  9]) ( 5., [ 9, 12])]]
>>> func = lambda record: ([record.a ** 2, record.b * 3],)
>>> print(apply_function(arr, func, dtype=[('d', int, 2)]))
[[([ 0,  3],) ([ 1,  6],)]
 [([ 4,  9],) ([ 9, 12],)]]
pyoints.nptools.colzip(arr)[source]

Splits a two dimensional numpy array into a list of columns.

Parameters:arr (np.ndarray(shape=(n, k)) or np.recarray(shape=(n, ))) – Numpy array with n rows and k columns.
Returns:columns – List of k numpy arrays.
Return type:list of np.ndarray
Raises:TypeError, ValueError

Examples

>>> arr = np.eye(3, dtype=int)
>>> cols = colzip(arr)
>>> len(cols)
3
>>> print(cols[0])
[1 0 0]
pyoints.nptools.dtype_subset(dtype, names)[source]

Creates a subset of a numpy type object.

Parameters:
  • dtype (list or np.dtype) – Numpy data type.
  • names (list of str) – Fields to select.
Raises:

TypeError

Returns:

Desired subset of numpy data type descriptions.

Return type:

list

Examples

>>> dtypes = [('coords', float, 3), ('values', int), ('text', '<U0')]
>>> print(dtype_subset(dtypes, ['text', 'coords']))
[('text', '<U0'), ('coords', '<f8', (3,))]
pyoints.nptools.flatten_dtypes(np_dtypes)[source]

Extract name, datatype and shape information from a numpy data type.

Parameters:np_dtypes (np.dtype) – Numpy data types to flatten.
Returns:
  • names (list of str) – Names of fields.
  • dtypes (list of dtypes) – Data types of fields.
  • shapes (list of tuples) – Shapes of fields.

Examples

>>> dtype = np.dtype([
...     ('simple', np.int8),
...     ('multidimensional', np.float32, 3),
... ])
>>> names, dtypes, shapes = flatten_dtypes(dtype)
>>> names
['simple', 'multidimensional']
>>> dtypes
[dtype('int8'), dtype('float32')]
>>> shapes
[0, 3]
pyoints.nptools.fuse(*recarrays)[source]

Fuses multiple numpy record arrays of identical shape to one array.

Parameters:*recarrays (np.recarray) – Numpy record arrays to fuse.
Returns:Record array with all fields of recarrays.
Return type:np.recarray

Examples

Fuse one dimensional arrays.

>>> A = recarray({'a': [0, 1, 2, 3]})
>>> B = recarray({'b': [4, 5, 6, 7]})
>>> C = fuse(A, B)
>>> print(C.shape)
(4,)
>>> print(C.dtype.names)
('a', 'b')

Fuse multiple two dimensional arrays.

>>> A = recarray({'a': [[0, 1], [2, 3]]}, dim = 2)
>>> B = recarray({'b': [[4, 5], [6, 7]]}, dim = 2)
>>> C = recarray({
...         'c1': [['c1', 'c2'], ['c3', 'c4']],
...         'c2': [[0.1, 0.2], [0.3, 0.3]],
...     }, dim = 2)
>>> D = fuse(A, B, C)
>>> print(sorted(D.dtype.names))
['a', 'b', 'c1', 'c2']
>>> print(D.shape)
(2, 2)
>>> print_rounded(D.a)
[[0 1]
 [2 3]]
>>> print(D.c1)
[['c1' 'c2']
 ['c3' 'c4']]
>>> print_rounded(D.c2)
[[ 0.1  0.2]
 [ 0.3  0.3]]
pyoints.nptools.haskeys(d)[source]

Checks if an object has keys and can be treated like a dictionary.

Parameters:d (object) – Object to be checked.
Returns:Indicates whether or not the object has accessable keys.
Return type:bool

Examples

>>> haskeys({'a': 5, 'b': 3})
True
>>> haskeys([5, 6])
False
>>> haskeys(np.recarray(3, dtype=[('a', int)]))
False
pyoints.nptools.indices(shape, flatten=False)[source]

Create keys or indices of a numpy ndarray.

Parameters:shape (array_like(int)) – Shape of desired output array.
Returns:Array of indices with desired shape. Each entry provides an index tuple to access the array entries.
Return type:np.ndarray(int, shape=(*shape, len(shape)))

Examples

One dimensional case.

>>> keys = indices(9)
>>> print(keys.shape)
(9,)
>>> print(keys)
[0 1 2 3 4 5 6 7 8]

Two dimensional case.

>>> keys = indices((3, 4))
>>> keys.shape
(3, 4, 2)
>>> print(keys)
[[[0 0]
  [0 1]
  [0 2]
  [0 3]]
<BLANKLINE>
 [[1 0]
  [1 1]
  [1 2]
  [1 3]]
<BLANKLINE>
 [[2 0]
  [2 1]
  [2 2]
  [2 3]]]

Get iterable of indices.

>>> keys = indices((3, 4), flatten=True)
>>> print(keys)
[[0 0]
 [0 1]
 [0 2]
 ...,
 [2 1]
 [2 2]
 [2 3]]
pyoints.nptools.isarray(o)[source]

Checks whether or nor an object is an array.

Parameters:o (object) – Some object.
Returns:Indicates whether or not the object is an array.
Return type:bool

Examples

>>> isarray([1, 2, 3])
True
>>> isarray('text')
False
pyoints.nptools.isnumeric(arr, dtypes=[<class 'numpy.uint8'>, <class 'numpy.uint16'>, <class 'numpy.uint32'>, <class 'numpy.uint64'>, <class 'numpy.int8'>, <class 'numpy.int16'>, <class 'numpy.int32'>, <class 'numpy.int64'>, <class 'numpy.float16'>, <class 'numpy.float32'>, <class 'numpy.float64'>])[source]

Checks if the data type of an array is numeric.

Parameters:
  • arr (array_like) – Numpy array to check.
  • dtypes (optional, tuple) – Tuple of allowed numeric data types.
Returns:

Indicates whether or not the array is numeric.

Return type:

bool

Raises:

TypeError

Examples

>>> isnumeric([1, 2, 3])
True
>>> isnumeric(['1', '2', '3'])
False
>>> isnumeric([1, 2, None])
False
pyoints.nptools.max_value_range(dtype)[source]

Returns the maximum value range of a numeric numpy data type.

Parameters:dtype (np.dtype) – Numeric data type to check
Returns:min_value,max_value – Minimum and maximum value
Return type:int

Examples

>>> value_range = max_value_range(np.dtype(np.uint8))
>>> print(value_range)
(0, 255)
>>> value_range = max_value_range(np.dtype(np.uint16))
>>> print(value_range)
(0, 65535)
>>> value_range = max_value_range(np.dtype(np.int8))
>>> print(value_range)
(-128, 127)
>>> value_range = max_value_range(np.dtype(np.int16))
>>> print(value_range)
(-32768, 32767)
>>> value_range = max_value_range(np.dtype(np.float16))
>>> print(value_range)
(-65504.0, 65504.0)
pyoints.nptools.merge(arrays, strategy=<built-in function concatenate>)[source]

Merges multiple arrays with similar fields.

Parameters:
  • arrays (list of np.recarray) – Numpy arrays to merge.
  • strategy (optional, function) – Aggregate function to apply during merging. Suggested values: np.conatenate, np.hstack, np.vstack, np.dstack.
Returns:

Merged numpy record array of same data type as the first input array.

Return type:

np.recarray

Raises:

TypeError

Examples

One dimensional arrays.

>>> A = recarray({'a': [(0, 1), (2, 3), (4, 5)], 'b': ['e', 'f', 'g']})
>>> B = recarray({'a': [(6, 7), (8, 9), (0, 1)], 'b': ['h', 'i', 'j']})
>>> C = recarray({'a': [(2, 3), (4, 5), (6, 7)], 'b': ['k', 'l', 'm']})
>>> D = merge((A, B, C))
>>> print(sorted(D.dtype.names))
['a', 'b']
>>> print(D.b)
['e' 'f' 'g' 'h' 'i' 'j' 'k' 'l' 'm']
>>> print(D.shape)
(9,)
>>> D = merge((A, B, C), strategy=np.hstack)
>>> print(D.shape)
(9,)
>>> D = merge((A, B, C), strategy=np.vstack)
>>> print(D.shape)
(3, 3)
>>> D = merge((A, B, C), strategy=np.dstack)
>>> print(D.shape)
(1, 3, 3)

Merge two dimensional arrays.

>>> A = recarray({
...     'a': [(0, 1), (2, 3)], 'b': [('e', 'f'), ('g', 'h')]
... }, dim=2)
>>> B = recarray({
...     'a': [(4, 5), (6, 7)], 'b': [('i', 'j'), ('k', 'l')]
... }, dim=2)
>>> C = recarray({
...     'a': [(1, 3), (7, 2)], 'b': [('m', 'n'), ('o', 'p')]
... }, dim=2)
>>> D = merge((A, B, C))
>>> print(sorted(D.dtype.names))
['a', 'b']
>>> print(D.b)
[['e' 'f']
 ['g' 'h']
 ['i' 'j']
 ['k' 'l']
 ['m' 'n']
 ['o' 'p']]
>>> print(D.shape)
(6, 2)
>>> D = merge((A, B, C), strategy=np.hstack)
>>> print(D.shape)
(2, 6)
>>> D = merge((A, B, C), strategy=np.vstack)
>>> print(D.shape)
(6, 2)
>>> D = merge((A, B, C), strategy=np.dstack)
>>> print(D.shape)
(2, 2, 3)
>>> D = merge((A, B, C), strategy=np.concatenate)
>>> print(D.shape)
(6, 2)
>>> A = np.recarray(1, dtype=[('a', object, 2), ('b', str)])
>>> B = np.recarray(2, dtype=[('a', object, 2), ('b', str)])
>>> D = merge((A, B), strategy=np.concatenate)
>>> print(D)
[([None, None], '') ([None, None], '') ([None, None], '')]
pyoints.nptools.minimum_numeric_dtype(arr)[source]

Determines the minimum required data type of a numpy without loosing accuracy.

Parameters:arr (np.ndarray(Number)) – Numeric array to find minimum data type for.
Returns:Minimum required data type.
Return type:np.dtype

Examples

Find minimum data type for integer arrays.

>>> arr = np.array([0, 255], dtype=np.int32)
>>> print(arr.dtype)
int32
>>> print(minimum_numeric_dtype(arr))
uint8
>>> arr = np.array([0, 256])
>>> print(minimum_numeric_dtype(arr))
uint16
>>> arr = np.array([-5, 127])
>>> print(minimum_numeric_dtype(arr))
int8
>>> arr = np.array([-5, 128])
>>> print(minimum_numeric_dtype(arr))
int16
>>> arr = np.array([-5, 214748364])
>>> print(minimum_numeric_dtype(arr))
int32

Find minimum data type for floating point arrays.

>>> arr = np.array([-5.2, 100.3])
>>> print(arr.dtype)
float64
>>> print(minimum_numeric_dtype(arr))
float16
pyoints.nptools.missing(data)[source]

Find missing values in an array.

Parameters:data (array_like) – A array like object which might contain missing values. Missing values are assumed to be either None or NaN.
Returns:Boolean values indicate missing values.
Return type:array_like(bool, shape=data.shape)
Raises:ValueError

Examples

Finding missing values in a list.

>>> arr = ['str', 1, None, np.nan, np.NaN]
>>> print(missing(arr))
[False False  True  True  True]

Finding missing values in a multi-dimensional array.

>>> arr = np.array([(0, np.nan), (None, 1), (2, 3)], dtype=float)
>>> print(missing(arr))
[[False  True]
 [ True False]
 [False False]]
pyoints.nptools.range_filter(arr, min_value=-inf, max_value=inf)[source]

Filter values by range.

Parameters:
  • arr (array_like(Number)) – Numeric array to filter.
  • min_value,max_value (Number) – Minimum and maximum values to define the desired value range [min_value, max_value] of arr.
Returns:

Indices of all values of arr in the desired range.

Return type:

np.ndarray(int)

Examples

Filter a one dimensional array.

>>> a = [0, 2, 1, -1, 5, 7, 9, 4, 3, 2, -2, -11]
>>> indices = range_filter(a, min_value=0)
>>> print(indices)
[0 1 2 4 5 6 7 8 9]
>>> indices = range_filter(a, max_value=5)
>>> print(indices)
[ 0  1  2  3  4  7  8  9 10 11]
>>> idx = range_filter(a, min_value=0, max_value=5)
>>> print(idx)
[0 1 2 4 7 8 9]
>>> print(np.array(a)[idx])
[0 2 1 5 4 3 2]

Filter a multi-dimensional array.

>>> a = [(1, 0), (-2, -1), (3, -5), (4, 2), (-7, 9), (0.5, 2)]
>>> idx = range_filter(a, min_value=2)
>>> print(idx)
((2, 3, 3, 4, 5), (0, 0, 1, 1, 1))
>>> print(np.array(a)[idx])
[ 3.  4.  2.  9.  2.]
>>> idx = range_filter(a, min_value=2, max_value=5)
>>> print(idx)
((2, 3, 3, 5), (0, 0, 1, 1))
>>> print(np.array(a)[idx])
[ 3.  4.  2.  2.]
pyoints.nptools.recarray(data_dict, dtype=[], dim=1)[source]

Converts a dictionary of array like objects to a numpy record array. This function is mostly used for convenience.

Parameters:
  • data_dict (dict) – Dictionary of array like objects to convert to a numpy record array. Each key in data_dict represents a field name of the output record array. Each item in data_dict represents the corresponding values. Thus, at least shape[0:dim] of all input arrays in data_dict have to match.
  • dtype (optional, numpy.dtype) – Describes the desired data types of specific fields. If the data type of a field is not given, the data type is estimated by numpy.
  • dim (positive int) – Desired dimension of the resulting numpy record array.
Returns:

Numpy record array build from input dictionary.

Return type:

np.recarray

Raises:

TypeError, ValueError

Examples

Creation of a one dimensional numpy record array using a dictionary.

>>> rec = recarray({
...    'coords': [ (3, 4), (3, 2), (0, 2), (5, 2)],
...    'text': ['text1', 'text2', 'text3', 'text4'],
...    'n':  [1, 3, 1, 2],
...    'missing':  [None, None, 'str', None],
... })
>>> print(sorted(rec.dtype.names))
['coords', 'missing', 'n', 'text']
>>> print(rec.coords)
[[3 4]
 [3 2]
 [0 2]
 [5 2]]

Create a two dimensional record array.

>>> data = {
...    'coords': [
...                 [(2, 3.2, 1), (-3, 2.2, 4)],
...                 [(0, 1.1, 2), (-1, 2.2, 5)],
...                 [(-7, -1, 3), (9.2, -5, 6)]
...             ],
...    'values': [[1, 3], [4, 0], [-4, 2]]
... }
>>> rec = recarray(data, dim=2)
>>> print(rec.shape)
(3, 2)
>>> print(rec.coords)
[[[ 2.   3.2  1. ]
  [-3.   2.2  4. ]]
<BLANKLINE>
 [[ 0.   1.1  2. ]
  [-1.   2.2  5. ]]
<BLANKLINE>
 [[-7.  -1.   3. ]
  [ 9.2 -5.   6. ]]]
>>> print(rec.values)
[[ 1  3]
 [ 4  0]
 [-4  2]]
>>> print(sorted(rec.dtype.names))
['coords', 'values']
pyoints.nptools.unnest(arr, deep=False)[source]

Unnest a numpy record array. This function recursively splits a nested numpy array to a list of arrays.

Parameters:
  • rec (np.recarray or np.ndarray) – Numpy array to unnest.
  • deep (bool) – Indicates whether or not numpy ndarrays shall be splitted into individual colums or not.
Raises:

TypeError

Returns:

List of unnested fields.

Return type:

list

Examples

>>> dtype = [
...    ('regular', np.int, 1),
...    ('nested', [
...         ('child1', np.str),
...         ('child2', np.float, 2)
...    ])
... ]
>>> rec = np.ones(2, dtype=dtype).view(np.recarray)
>>> print(rec.nested.child2)
[[ 1.  1.]
 [ 1.  1.]]
>>> unnested = unnest(rec)
>>> print(unnested[0])
[1 1]
>>> print(unnested[1])
['' '']
>>> print(unnested[2])
[[ 1.  1.]
 [ 1.  1.]]

pyoints.polar module

Handling of polar coordinates.

pyoints.polar.coords_to_polar(coords)[source]

Converts Cartesian coordinates to polar coordinates.

Parameters:coords (array_like(Number, shape=(n, k))) – Represents n data points of k dimensions in a Cartesian coordinate system.
Returns:pcoords – Represents n data points of k dimensions in a polar coordinate system. First column represents the distance to the origin of the coordinate system. All other columns represent the corresponding axes angles.
Return type:array_like(Number, shape=(n, k))

Examples

2D coordinates.

>>> coords = [(0, 0), (0, 1), (1, 0), (1, 1), (-1, 1), (2, -5)]
>>> pcoords = coords_to_polar(coords)
>>> print_rounded(pcoords, 3)
[[ 0.     0.   ]
 [ 1.     1.571]
 [ 1.     0.   ]
 [ 1.414  0.785]
 [ 1.414  2.356]
 [ 5.385 -1.19 ]]

3D coordinates.

>>> coords = [(0, 0, 0), (1, 1, 0), (-1, -1, -1), (2, -5, 9)]
>>> pcoords = coords_to_polar(coords)
>>> print_rounded(pcoords, 3)
[[  0.      0.      0.   ]
 [  1.414   0.785   1.571]
 [  1.732  -2.356   2.186]
 [ 10.488  -1.19    0.539]]
pyoints.polar.polar_to_coords(pcoords)[source]

Converts polar coordinates to Cartesian coordinates.

Parameters:pcoords (array_like(Number, shape=(n, k))) – Represents n data points of k dimensions in a polar coordinate system. First column represents the distance to the origin of the coordinate system. All other columns represent the corresponding axes angles.
Returns:coords – Represents n data points of k dimensions in a Cartesian coordinate system.
Return type:array_like(Number, shape=(n, k))

Examples

2D coordinates.

>>> pcoords = [(0, 0), (3, 0), (3, np.pi), (4, -0.5*np.pi), (1, 0.5)]
>>> coords = polar_to_coords(pcoords)
>>> print_rounded(coords, 3)
[[ 0.     0.   ]
 [ 3.     0.   ]
 [-3.     0.   ]
 [ 0.    -4.   ]
 [ 0.878  0.479]]

3D coordinates.

>>> pcoords = [(0, 0, 0), (2, 0, 0),(4, 0, np.pi), (4, 0.5*np.pi, 0.5)]
>>> coords = polar_to_coords(pcoords)
>>> print_rounded(coords, 3)
[[ 0.     0.     0.   ]
 [ 0.     0.     2.   ]
 [ 0.     0.    -4.   ]
 [ 0.     1.918  3.51 ]]

pyoints.projection module

Coordinate Reference Systems and two dimensional geograpic coordinate transformations.

class pyoints.projection.GeoTransform(from_proj, to_proj)[source]

Bases: object

Provides a coordinate transformation between different spatial reference systems.

Parameters:from_proj,to_proj (Proj) – Define the coordinate transformation from the origin projection system from_proj to the target projection system to_proj.

Examples

Transform coordinates.

>>> wgs84 = Proj.from_epsg(4326)
>>> gk2 = Proj.from_epsg(31466)
>>> coords = [
...     (6.842, 49.971),
...     (6.847, 49.969),
...     (6.902, 49.991),
...     (6.922, 50.101)
... ]
>>> geoTransfrom = GeoTransform(wgs84, gk2)
>>> tCoords = geoTransfrom(coords)
>>> print_rounded(tCoords, 3)
[[ 2560446.801  5537522.386]
 [ 2560808.009  5537303.984]
 [ 2564724.211  5539797.116]
 [ 2566007.32   5552049.646]]

Reverse transformation.

>>> print_rounded(geoTransfrom(tCoords, reverse=True), 3)
[[  6.842  49.971]
 [  6.847  49.969]
 [  6.902  49.991]
 [  6.922  50.101]]
from_proj
to_proj
class pyoints.projection.Proj(proj4='+proj=latlong +datum=WGS84 +to +proj=latlong +datum=WGS84 +units=m +no_defs')[source]

Bases: object

Wrapper class for different commonly coordinate reference system formats.

Parameters:proj4 (optional, str) – Coordinate reference system definition in Proj4 format. WGS84, if None or empty string.
proj4

Projection in Proj4 format.

Type:str
wkt

Projection in Well Known Text format.

Type:str
pyproj

Projection as pyproj.Proj instance.

Type:pyproj.Proj

Examples

Create from EPSG code.

>>> proj = Proj.from_epsg(4326)
>>> print('AUTHORITY["EPSG","4326"]' in proj.wkt)
True
>>> print('WGS84' in proj.proj4)
True

Create from Proj4 string.

>>> proj = Proj.from_proj4('+proj=longlat +datum=WGS84 +no_defs')
>>> print(proj.proj4)
+proj=longlat +datum=WGS84 +no_defs
>>> print('AUTHORITY["EPSG","4326"]' in proj.wkt)
True

Create from Well Known Text.

>>> proj = Proj.from_wkt(proj.wkt)
>>> print('WGS84' in proj.proj4)
True
classmethod from_epsg(epsg)[source]

Create Proj object from EPSG code.

Parameters:epsg (int) – Coordinate projection definition in EPSG format.
classmethod from_proj4(proj4)[source]

Create Proj object from Proj4 format.

Parameters:proj4 (str) – Coordinate projection definition in Proj4 format.
classmethod from_wkt(wkt)[source]

Create Proj object from Well Known Text.

Parameters:wkt (str) – Coordinate projection definition in Well Known Text format.
osr
proj4
pyproj
wkt
pyoints.projection.project(coords, from_proj, to_proj)[source]

Applies the coordinate transformation.

Parameters:
  • coords (array_like(Number, shape=(n, k))) – Represents n points of k dimensions to transform.
  • from_proj,to_proj (Proj) – Define the coordinate transformation from the origin projection system from_proj to the target projection system to_proj.

See also

GeoTransform()

Returns:Transformed coordinates.
Return type:np.ndarray(Number, shape=(n, k))

pyoints.smoothing module

Collection of algorithms to smooth point clouds.

pyoints.smoothing.mean_ball(coords, r, num_iter=1, update_pairs=False, f=<function <lambda>>)[source]

Smoothing of spatial structures by iterative averaging the coordinates of neighboured points.

Parameters:
  • coords (array_like(Number, shape=(n, k))) – Array representing n points of k dimensions.
  • r (Number) – Maximum distance to nearby points used to average the coordinates.
  • num_iter (optional, positive int) – Number of iterations.
  • update_pairs (optional, bool) – Specifies weather or not point pairs are updated on each iteration.
  • f (callable) – Aggregate function used for smoothing. It receives the original point coordinate and the coordinates of neighboured points as an argument and returns a smoothed coordinate.

See also

mean_knn()

Examples

Create a three dimensional irregular surface of points.

>>> coords = np.ones((100, 3), dtype=float)
>>> coords[:, 0:2] = np.vstack(np.mgrid[0:10, 0:10].T)
>>> coords[:, 2] = np.tile([1.05, 0.95], 50)

Get value range in each coordinate dimension.

>>> print_rounded(np.ptp(coords, axis=0))
[ 9.   9.   0.1]

Smooth coordinates to get a more regular surface. But the first two coordinate dimensions are affected, too.

>>> scoords = mean_ball(coords, 1.5)
>>> print_rounded(np.ptp(scoords, axis=0), 3)
[ 8.     8.     0.033]

Modify the aggregation function to smooth the third coordinate axis only.

>>> def aggregate_function(coord, ncoords):
...     coord[2] = ncoords[:, 2].mean(0)
...     return coord
>>> scoords = mean_ball(coords, 1.5, f=aggregate_function)
>>> print_rounded(np.ptp(scoords, axis=0), 3)
[ 9.     9.     0.026]

Increase number of iterations to get a smoother result.

>>> scoords = mean_ball(coords, 1.5, num_iter=3, f=aggregate_function)
>>> print_rounded(np.ptp(scoords, axis=0), 3)
[ 9.    9.    0.01]
pyoints.smoothing.mean_knn(coords, k, num_iter=1, update_pairs=False, f=<function <lambda>>)[source]

Smoothing of spatial structures by averaging neighboured point coordinates.

Parameters:
  • coords (array_like(Number, shape=(n, l))) – Array representing n points with l dimensions.
  • k (float) – Number of nearest points used to average the coordinates.
  • num_iter (optional, int) – Number of iterations.
  • update_pairs (optional, bool) – Specifies weather or not point pairs are updated on each iteration.
  • f (callable) – Aggregate function used for smoothing. It receives the original point coordinate and the coordinates of neighboured points as an argument and returns a smoothed coordinate.

See also

mean_ball()

Examples

Create a three dimensional irregular surface of points.

>>> coords = np.ones((100, 3), dtype=float)
>>> coords[:, 0:2] = np.vstack(np.mgrid[0:10, 0:10].T)
>>> coords[:, 2] = np.tile([1.05, 0.95], 50)

Get value range in each coordinate dimension.

>>> print_rounded(np.ptp(coords, axis=0))
[ 9.   9.   0.1]

Smooth coordinates to get a more regular surface. But the first two coordinate dimensions are affected, too.

>>> scoords = mean_knn(coords, 5)
>>> print_rounded(np.ptp(scoords, axis=0), 3)
[ 8.2   8.2   0.02]

Modify the aggregation function to smooth the third coordinate axis only.

>>> def aggregate_function(coord, ncoords):
...     coord[2] = ncoords[:, 2].mean(0)
...     return coord
>>> scoords = mean_knn(coords, 5, f=aggregate_function)
>>> print_rounded(np.ptp(scoords, axis=0), 3)
[ 9.     9.     0.033]

pyoints.surface module

Create surface models.

class pyoints.surface.Surface(coords, method=<class 'pyoints.interpolate.KnnInterpolator'>, **kwargs)[source]

Bases: object

Creates a surface model based on representative points.

Parameters:
  • coords (array_like(Number, shape=(n, k))) – Represents n data points of k dimensions representing a surface.
  • method (optional, Interpolator) – Interpolation method to use.
  • **kwargs (optional) – Arguments passed to the interpolation method.

Examples

>>> method = interpolate.LinearInterpolator
>>> surface = Surface([(0, 0, 0), (0, 2, 0), (2, 1, 4)], method=method)
>>> print(surface([(1, 1)]))
2.0
>>> print(surface([(1, 1), (0.5, 1)]))
[ 2.  1.]
>>> print(surface([(1, 1, 2), (0.5, 1, 9)]))
[ 2.  1.]

pyoints.transformation module

Multidimensional transformation matrices and coordinate transformations.

class pyoints.transformation.LocalSystem[source]

Bases: numpy.ndarray, object

Defines a local coordinate system based on a transformation matrix.

Parameters:T (array_like(Number, shape=(k+1, k+1))) – Transformation matrix of k coordinate dimensions.
dim

Number of coordinate dimensions.

Type:positive int
origin

Global origin of the local coordinate system.

Type:np.ndarray(Number, shape=(k))

Examples

>>> L = LocalSystem([(2, -0.5, 2), (0.5, 1, 3), (0, 0, 1)])
>>> print_rounded(L)
[[ 2.  -0.5  2. ]
 [ 0.5  1.   3. ]
 [ 0.   0.   1. ]]
>>> print_rounded(L.dim)
2
>>> print_rounded(L.origin)
[ 2.  3.]

See also

matrix

decomposition()[source]

Determinates some characteristic parameters of the transformation matrix.

Returns:Decomposition values.
Return type:tuple

See also

decomposition()

det
dim
explained_variance(global_coords)[source]

Calculate the variance of global coordinates explained by the coordinate axes.

Parameters:global_coords (array_like(Number, shape=(n, k))) – Represents ‘n’ points of k dimensions in the global coordinate system.
Returns:Total amount of explained variance by a specific coordinate axis.
Return type:np.ndarray(Number, shape=(self.dim))

Examples

>>> T = s_matrix([1, 2])
>>> print_rounded(T)
[[ 1.  0.  0.]
 [ 0.  2.  0.]
 [ 0.  0.  1.]]
>>> e = T.explained_variance([(2, 1), (0, 0), (1, 1), (2, 3)])
>>> print_rounded(e, 3)
[ 0.688  4.75 ]
explained_variance_ratio(global_coords)[source]

Calculate the relative variance of global coordinates explained by the coordinate axes.

Parameters:global_coords (array_like(Number, shape=(n, k))) – Represents ‘n’ points of k dimensions in the global coordinate system.
Returns:Relative amount of variance explained by each coordinate axis.
Return type:np.ndarray(Number, shape=(self.dim))

Examples

>>> T = s_matrix([1, 2])
>>> print_rounded(T)
[[ 1.  0.  0.]
 [ 0.  2.  0.]
 [ 0.  0.  1.]]
>>> e = T.explained_variance_ratio([(2, 1), (0, 0), (1, 1), (2, 3)])
>>> print_rounded(e, 3)
[ 0.126  0.874]
inv
origin
reflect(axis=0)[source]

Reflects a specific coordinate axis.

Parameters:axis (positive int) – Coordinate axis to reflect.

Examples

Reflect a two dimensional transformation matrix.

>>> L = LocalSystem(matrix(t=[1, 2], s=[0.5, 2]))
>>> print_rounded(L)
[[ 0.5  0.   1. ]
 [ 0.   2.   2. ]
 [ 0.   0.   1. ]]
>>> print_rounded(L.to_local([2, 3]))
[ 2.  8.]
>>> L.reflect()
>>> print_rounded(L)
[[-0.5  0.  -1. ]
 [ 0.   2.   2. ]
 [ 0.   0.   1. ]]
>>> print_rounded(L.to_local([2, 3]))
[-2.  8.]

Reflect a three dimensional transformation matrix.

>>> L = LocalSystem(matrix(t=[1, 2, 3], s=[0.5, -2, 1]))
>>> print_rounded(L)
[[ 0.5  0.   0.   1. ]
 [ 0.  -2.   0.   2. ]
 [ 0.   0.   1.   3. ]
 [ 0.   0.   0.   1. ]]
>>> print_rounded(L.to_local([1, 2, 3]))
[ 1.5 -2.   6. ]
>>> L.reflect(axis=1)
>>> print_rounded(L)
[[ 0.5  0.   0.   1. ]
 [ 0.   2.   0.  -2. ]
 [ 0.   0.   1.   3. ]
 [ 0.   0.   0.   1. ]]
>>> print_rounded(L.to_local([1, 2, 3]))
[ 1.5  2.   6. ]
to_global(local_coords)[source]

Transforms local coordinates into global coordinates.

Parameters:local_coords (array_like(Number, shape=(n, k))) – Represents n points of k dimensions in the local coordinate system.
Returns:Global coordinates.
Return type:np.ndarray(Number, shape=(n, k))

Examples

>>> T = matrix(t=[2, 3], s=[0.5, 10])
>>> gcoords = T.to_global([(2, 3), (2, 13), (2.5, 3), (1.5, -7)])
>>> print_rounded(gcoords)
[[ 0.  0.]
 [ 0.  1.]
 [ 1.  0.]
 [-1. -1.]]
to_local(global_coords)[source]

Transforms global coordinates into local coordinates.

Parameters:global_coords (array_like(Number, shape=(n, k))) – Represents n points of k dimensions in the global coordinate system.
Returns:Local coordinates.
Return type:np.ndarray(Number, shape=(n, k))

Examples

>>> T = matrix(t=[2, 3], s=[0.5, 10])
>>> lcoords = T.to_local([(0, 0), (0, 1), (1, 0), (-1, -1)])
>>> print_rounded(lcoords)
[[  2.    3. ]
 [  2.   13. ]
 [  2.5   3. ]
 [  1.5  -7. ]]
class pyoints.transformation.PCA[source]

Bases: pyoints.transformation.LocalSystem

Principal Component Analysis (PCA).

Parameters:coords (array_like(Number, shape=(n, k))) – Represents n data points of k dimensions. These coordinates are used to fit a PCA.
eigen_values

Characteristic Eigenvalues of the PCA.

Type:np.ndarray(Number, shape=(k))

Notes

Implementation idea taken from [1].

References

[1] https://stackoverflow.com/questions/13224362/principal-component-analysis-pca-in-python

Examples

>>> coords = [(0, 0), (3, 4)]
>>> T = PCA(coords)
>>> print_rounded(T)
[[ 0.6  0.8 -2.5]
 [-0.8  0.6  0. ]
 [ 0.   0.   1. ]]
>>> print_rounded(np.linalg.inv(T))
[[ 0.6 -0.8  1.5]
 [ 0.8  0.6  2. ]
 [ 0.   0.   1. ]]
>>> print_rounded(transform(coords, T))
[[-2.5  0. ]
 [ 2.5  0. ]]
eigen_values
pc(k)[source]

Select a specific principal component.

Parameters:k (positive int) – k-th principal component to return.
Returns:k-th principal component.
Return type:np.ndarray(Number, shape=(self.dim))

Examples

>>> coords = [(-3, -4), (0, 0), (12, 16)]
>>> T = PCA(coords)
>>> print_rounded(T)
[[ 0.6  0.8 -5. ]
 [-0.8  0.6  0. ]
 [ 0.   0.   1. ]]
>>> print_rounded(T.pc(1))
[ 0.6  0.8]
>>> print_rounded(T.pc(2))
[-0.8  0.6]
pyoints.transformation.add_dim(T)[source]

Adds a dimension to a transformation matrix.

Parameters:T (LocalSystem(Number, shape=(k+1, k+1))) – Transformation matrix of k coordinate dimensions.
Returns:Transformation matrix with an additional dimension.
Return type:np.ndarray(float, shape=(k+2, k+2))

Examples

Two dimensional case.

>>> T = matrix(t=[2, 3 ], s=[0.5, 2])
>>> print_rounded(T, 3)
[[ 0.5  0.   2. ]
 [ 0.   2.   3. ]
 [ 0.   0.   1. ]]
>>> T = add_dim(T)
>>> print_rounded(T, 3)
[[ 0.5  0.   0.   2. ]
 [ 0.   2.   0.   3. ]
 [ 0.   0.   1.   0. ]
 [ 0.   0.   0.   1. ]]
pyoints.transformation.decomposition(T, ignore_warnings=False)[source]

Determines some characteristic parameters of a transformation matrix.

Parameters:T (array_like(Number, shape=(k+1, k+1))) – Transformation matrix of k coordinate dimensions.
Returns:
  • t,r,s (optional, np.ndarray(Number, shape=(k))) – Transformation coefficients for translation t, rotation r and scale s.
  • det (float) – Determinant det indicates a distorsion of T.

Examples

>>> T = matrix(t=[2, 3], r=0.2, s=[0.5, 2])
>>> t, r, s, det = decomposition(T)
>>> print_rounded(t)
[ 2.  3.]
>>> print_rounded(r)
0.2
>>> print_rounded(s)
[ 0.5  2. ]
>>> print_rounded(det)
1.0

See also

matrix()

Notes

Idea taken from [1], [2] and [3].

References

[1] https://math.stackexchange.com/questions/237369/given-this-transformation-matrix-how-do-i-decompose-it-into-translation-rotati [2] https://math.stackexchange.com/questions/13150/extracting-rotation-scale-values-from-2d-transformation-matrix/13165#13165 [3] https://www.learnopencv.com/rotation-matrix-to-euler-angles/

pyoints.transformation.eigen(coords)[source]

Fit eigenvectors to coordinates.

Parameters:coords (array_like(Number, shape=(n, k))) – Represents n data points of k dimensions to fit the eigenvectors to.
Returns:
  • eigen_vectors (np.ndarray(Number, shape=(k, k))) – Columnwise normalized eigenvectors in descending order of eigenvalue.
  • eigen_values (np.ndarray(Number, shape=(k))) – Eigenvalues in descending order of magnitude.

Notes

Implementation idea taken from [1].

References

[1] https://stackoverflow.com/questions/13224362/principal-component-analysis-pca-in-python

See also

PCA()

Examples

>>> coords = [(0, 0), (3, 4)]
>>> eigen_vectors, eigen_values = eigen(coords)
>>> print_rounded(eigen_vectors)
[[ 0.6 -0.8]
 [ 0.8  0.6]]
>>> print_rounded(eigen_values)
[ 12.5   0. ]
pyoints.transformation.homogenious(coords, value=1)[source]

Create homogeneous coordinates by adding an additional dimension.

Parameters:
  • coords (array_like(Number, shape=(n, k))) – Represents n data points of k dimensions.
  • value (optional, Number) – Desired values of additional dimension.
Returns:

hcoords – Homogeneous coordinates.

Return type:

np.ndarray(Number, shape=(n, k+1))

pyoints.transformation.i_matrix(dim)[source]

Creates an identity transformation matrix.

Parameters:dim (positive int) – Desired dimension of the transformation matrix. At least a value of two is reqired.
Returns:Identity transformation matrix.
Return type:LocalSystem(int, shape=(dim+1, dim+1))

See also

LocalSystem()

Examples

>>> print_rounded(i_matrix(3))
[[ 1.  0.  0.  0.]
 [ 0.  1.  0.  0.]
 [ 0.  0.  1.  0.]
 [ 0.  0.  0.  1.]]
pyoints.transformation.matrix(t=None, r=None, s=None, order='srt')[source]

Creates a transformation matrix based on translation, rotation and scale coefficients.

Parameters:
  • t,r,s (optional, np.ndarray(Number, shape=(k))) – Transformation coefficients for translation t, rotation r and scale s.
  • order (optional, string) – Order to compute the matrix multiplications. For example, an order of trs means to first translate t, then rotate r, and finally scale s.
Returns:

Transformation matrix according to arguments t, r, s and order order.

Return type:

LocalSystem(shape=(k+1, k+1))

pyoints.transformation.matrix_from_gdal(t)[source]

Creates a transformation matrix using a gdal geotransfom array.

Parameters:t (array_like(Number, shape=(6))) – Gdal geotransform array.
Returns:Matrix representation of the gdal geotransform array.
Return type:LocalSystem(Number, shape=(3, 3))

Examples

>>> T = matrix_from_gdal([-28493, 9, 2, 4255884, -2, -9.0])
>>> print_rounded(T.astype(int))
[[      9       2  -28493]
 [     -2      -9 4255884]
 [      0       0       1]]
pyoints.transformation.matrix_to_gdal(T)[source]

Creates gdal geotransfom array based on a transformation matrix.

Parameters:T (array_like(Number, shape=(3, 3))) – Matrix representation of the gdal geotransform array.
Returns:t – Gdal geotransform array.
Return type:array_like(Number, shape=(6))

Examples

>>> T = np.array([(9, -2, -28493), (2, -9, 4255884), (0, 0, 1)])
>>> t = matrix_to_gdal(T)
>>> print_rounded(t)
(-28493.0, 9.0, -2.0, 4255884.0, 2.0, -9.0)
pyoints.transformation.q_matrix(q)[source]

Creates an rotation matrix based on quaternion values.

Parameters:q (array_like(Number, shape=(4))) – Quaternion parameters (w, x, y, z).
Returns:Rotation matrix derived from quaternions.
Return type:LocalSystem(int, shape=(4, 4))

Examples

>>> T = q_matrix([0.7071, 0.7071, 0, 0])
>>> print_rounded(T, 2)
[[ 1.  0.  0.  0.]
 [ 0.  0. -1.  0.]
 [ 0.  1.  0.  0.]
 [ 0.  0.  0.  1.]]
>>> t, r, s, det = decomposition(T)
>>> print_rounded(r)
[ 1.57  0.    0.  ]
>>> print_rounded(t)
[ 0.  0.  0.]
pyoints.transformation.r_matrix(a, order='xyz')[source]

Creates a rotation matrix.

Parameters:
  • r (Number or array_like(Number, shape=(k))) – Rotation coefficients for each coordinate dimension. At least two coefficients are required.
  • order (optional, String) – For at least three axes, order specifies the order of rotations. For example, an order of zxy means first to translate according to z axis, then rotate according to x-axis, and finally rotate according to y-axis.
Returns:

Rotation matrix.

Return type:

LocalSystem(shape=(k+1, k+1))

See also

LocalSystem()

Notes

Currently only two and tree dimensions are supported yet.

Examples

Two dimensional case.

>>> R = r_matrix(np.pi/4)
>>> print_rounded(R, 3)
[[ 0.707 -0.707  0.   ]
 [ 0.707  0.707  0.   ]
 [ 0.     0.     1.   ]]

Three dimensional case.

>>> R = r_matrix([np.pi/2, 0, np.pi/4])
>>> print_rounded(R, 3)
[[ 0.707  0.     0.707  0.   ]
 [ 0.707  0.    -0.707  0.   ]
 [ 0.     1.     0.     0.   ]
 [ 0.     0.     0.     1.   ]]
pyoints.transformation.s_matrix(s)[source]

Creates a scaling matrix.

Parameters:s (array_like(Number, shape=(k))) – Scaling coefficients for each coordinate dimension. At least two coefficients are required.
Returns:Scaling matrix.
Return type:LocalSystem(shape=(k+1, k+1))

See also

LocalSystem()

Examples

>>> print_rounded(s_matrix([0.5, 1, -3]))
[[ 0.5  0.   0.   0. ]
 [ 0.   1.   0.   0. ]
 [ 0.   0.  -3.   0. ]
 [ 0.   0.   0.   1. ]]
pyoints.transformation.t_matrix(t)[source]

Creates a translation matrix.

Parameters:t (array_like(Number, shape=(k))) – Translation coefficients for each coordinate dimension. At least two coefficients are required.
Returns:Translation matrix.
Return type:LocalSystem(shape=(k+1, k+1))

See also

LocalSystem()

Examples

>>> print_rounded(t_matrix([1, 2, 3]))
[[ 1.  0.  0.  1.]
 [ 0.  1.  0.  2.]
 [ 0.  0.  1.  3.]
 [ 0.  0.  0.  1.]]
pyoints.transformation.transform(coords, T, inverse=False, precise=False)[source]

Performs a linear transformation to coordinates using a transformation matrix.

Parameters:
  • coords (array_like(Number, shape=(n, k)) or array_like(Number, shape=(k))) – Represents n data points of k dimensions.
  • T (array_like(Number, shape=(k+1, k+1))) – Transformation matrix.
  • precise (optional, bool) – Indicates whether or not to calculate the transformaton (at the expend of increased computation time) extra precise.
Returns:

coords – Transformed coordinates.

Return type:

np.ndarray(Number, shape=(n, k)) or np.ndarray(Number, shape=(k))

Examples

Transform coordinates forth and back again.

>>> coords = [(0, 0), (1, 2), (2, 3), (0, 3)]
>>> T = [(1, 0, 5), (0, 1, 3), (0, 0, 1)]
>>> tcoords = transform(coords, T)
>>> print_rounded(tcoords)
[[ 5.  3.]
 [ 6.  5.]
 [ 7.  6.]
 [ 5.  6.]]
>>> print_rounded(transform(tcoords, T, inverse=True))
[[ 0.  0.]
 [ 1.  2.]
 [ 2.  3.]
 [ 0.  3.]]

pyoints.vector module

Various vector operations.

class pyoints.vector.Plane(origin, *vec)[source]

Bases: object

Class to represent two hyperplanes and handle them conveniently.

Parameters:
  • origin (np.ndarray(Number, shape=(k))) – Defines the origin of the planes k dimensional local coordinate system.
  • *vec (np.ndarray(Number, shape=(k))) – The k-1 vectors of vec define the orientation of the plane. The missing axis (perpenticular to the vectors) are calculated automatically using principal component analysis. So, parallel vectors in vec are valid, but might result in unexpected results.
origin

Origin of the plane.

Type:np.ndarray(Number, shape=(k))
vec

Vectors defining the orientation of the plane.

Type:np.ndarray(Number, shape=(k-1, k))
target

Points the vectors are targeting at.

Type:np.ndarray(Number, shape=(k-1, k))
dim

Number of coordinate dimensions of the vector.

Type:positive int
t

Roto-translation matrix representation of the local coordinate system defined by the plane.

Type:PCA(Number, shape=(k+1, k+1))

Examples

Creation of a plane object using an origin and two vectors defining the plane.

>>> origin = [1, 2, 3]
>>> v = [0, 0, 2]
>>> w = [1, 0, 0]
>>> plane = Plane(origin, v, w)
>>> print(plane)
origin: [1 2 3]; vec: [0 0 2], [1 0 0]

Get some properties.

>>> print_rounded(plane.dim)
3
>>> print_rounded(plane.t.inv)
[[ 0.  1.  0.  1.]
 [ 0.  0.  1.  2.]
 [ 1.  0.  0.  3.]
 [ 0.  0.  0.  1.]]
>>> print_rounded(plane.t.eigen_values)
[ 8.  2.  0.]

Transformation of global coordinates to the plane coordinate system.

>>> global_coords = [
...     plane.origin,
...     plane.origin + plane.vec[0, :],
...     -plane.vec[1, :],
...     plane.vec[0, :] - 3,
...     plane.origin + 2 * plane.vec[0, :] + 3 * plane.vec[1, :]
... ]
>>> print_rounded(np.array(global_coords))
[[ 1  2  3]
 [ 1  2  5]
 [-1  0  0]
 [-3 -3 -1]
 [ 4  2  7]]
>>> local_coords = plane.t.to_local(global_coords)
>>> print_rounded(local_coords)
[[ 0.  0.  0.]
 [ 2.  0.  0.]
 [-3. -2. -2.]
 [-4. -4. -5.]
 [ 4.  3.  0.]]
>>> print_rounded(plane.t.to_global(local_coords))
[[ 1.  2.  3.]
 [ 1.  2.  5.]
 [-1.  0.  0.]
 [-3. -3. -1.]
 [ 4.  2.  7.]]

Calculation of the distance of the global points to the plane.

>>> print_rounded(plane.distance(global_coords))
[ 0.  0.  2.  5.  0.]

Creation of the special case of a line in a two dimensional space.

>>> plane = Plane([1, 2], [3, 4])
>>> print(plane)
origin: [1 2]; vec: [3 4]
>>> print_rounded(plane.t.inv)
[[ 0.6 -0.8  1. ]
 [ 0.8  0.6  2. ]
 [ 0.   0.   1. ]]
>>> print_rounded(plane.t.eigen_values)
[ 50.   0.]
>>> print_rounded(plane.dim)
2

Transformation of global coordinates to the plane coordinate system.

>>> global_coords = [
...     plane.origin,
...     plane.origin + plane.vec[0, :],
...     plane.vec[0, :] - 3,
...     plane.origin + 2 * plane.vec[0, :]
... ]
>>> print_rounded(np.array(global_coords))
[[ 1  2]
 [ 4  6]
 [ 0  1]
 [ 7 10]]
>>> local_coords = plane.t.to_local(global_coords)
>>> print_rounded(local_coords)
[[  0.    0. ]
 [  5.    0. ]
 [ -1.4   0.2]
 [ 10.    0. ]]
>>> print_rounded(plane.t.to_global(local_coords))
[[  1.   2.]
 [  4.   6.]
 [  0.   1.]
 [  7.  10.]]

Calculation of the distance of the global points to the plane.

>>> print_rounded(plane.distance(global_coords))
[ 0.   0.   0.2  0. ]
dim
distance(global_coords)[source]

Calculate the distance between points and the plane.

Parameters:global_coords (array_like(Number, shape=(n, self))) – Represents n data points.
Returns:Distances of the points to the plane.
Return type:np.ndarray(Number, shape=(n))
origin
t
target
transform(T)[source]

Transforms the vector using a transformation matrix.

Parameters:T (array_like(Number, shape=(self.dim+1, self.dim+1))) – Transformation matrix to apply.
Returns:
Return type:self

Examples

Create a plane and a tranformation matrix.

>>> plane = Plane([1, 1], [1, 0])
>>> print_rounded(plane.origin)
[1 1]
>>> print_rounded(plane.vec)
[[1 0]]
>>> r = 45 * np.pi / 180.0
>>> t = [1, 2]
>>> T = transformation.matrix(t=t, r=r, order='rts')
>>> print_rounded(T)
[[ 0.71 -0.71  1.  ]
 [ 0.71  0.71  2.  ]
 [ 0.    0.    1.  ]]
>>> plane = plane.transform(T)
>>> print_rounded(plane.origin)
[ 1.    3.41]
>>> print_rounded(plane.vec)
[[ 0.71  0.71]]
vec
class pyoints.vector.Vector(origin, vec)[source]

Bases: object

Class to represent vectors and handle them conveniently.

Parameters:origin,vec (array_like(Number, shape=(k))) – The arrays origin and vec define the location and orientation of the vector in a k dimensional vector space.
origin,target,vec

The vector vec starts at point origin and points to target.

Type:np.ndarray(Number, shape=(k))
length

Length of the vector vec.

Type:positive float
dim

Number of coordinate dimensions of the vector.

Type:positive int
t

Roto-translation matrix representation of the local coordinate system defined by the vector.

Type:PCA(Number, shape=(k+1, k+1))

Examples

Two dimensional case.

>>> v = Vector((5, 7), (3, 4))
>>> print(v)
origin: [5 7]; vec: [3 4]
>>> print_rounded(v.target)
[ 8 11]
>>> print_rounded(v.length)
5.0

Three dimensional case.

>>> v = Vector((1, 1, 1), (2, 3, 4))
>>> print(v)
origin: [1 1 1]; vec: [2 3 4]
>>> print_rounded(v.target)
[3 4 5]

Edit the vector.

>>> v = Vector((1, 1), (3, 4))
>>> print(v)
origin: [1 1]; vec: [3 4]
>>> v.vec = (5, 2)
>>> print(v)
origin: [1 1]; vec: [5 2]
>>> v.origin = (-1, -2)
>>> print(v)
origin: [-1 -2]; vec: [5 2]
>>> v.target = (5, 4)
>>> print(v)
origin: [-1 -2]; vec: [6 6]
angles(deg=False)[source]

Calculates the angles of the vector in relation to the coordinate axes.

Parameters:deg (bool) – Indicates whether or not to provide the angles in degree.
Returns:Angles according to the coordinate axes.
Return type:np.ndarray(Number, shape=(self.dim))

Examples

>>> v = Vector((1, 1, 1), (2, 3, 4))
>>> angles = v.angles(deg=True)
>>> print_rounded(angles, 3)
[ 68.199  56.145  42.031]
dim
distance(global_coords)[source]

Calculate the distance between points and the vector.

Parameters:global_coords (array_like(Number, shape=(n, self))) – Represents n data points.
Returns:Distances of the points to the vector.
Return type:np.ndarray(Number, shape=(n))

Examples

>>> v = Vector((1, -1), (1, 2))
>>> dist = v.distance([(1, -1), (0, -3), (2, 1), (2, -2), (0, 0)])
>>> print_rounded(dist, 3)
[ 0.     0.     0.     1.342  1.342]
>>> print_rounded(np.linalg.inv(v.t).origin)
[ 1. -1.]
>>> print_rounded(v.t.pc(1) * v.length)
[ 1.  2.]
>>> print_rounded(v.t.pc(2) * v.length)
[-2.  1.]
classmethod from_anges(origin, rot)[source]
classmethod from_coords(coords)[source]

Vector creation from coordinates using Principal Component Analysis.

Parameters:coords (array_like(Number, shape=(n, k))) – Represents n data points of k dimensions. These coordinates are used to fit a PCA.

See also

transformation.PCA()

Examples

>>> v = Vector.from_coords([(1, 1), (1, 2), (1, 6)])
>>> print(v)
origin: [ 1.  3.]; vec: [ 0.  3.]
>>> print_rounded(v.t)
[[ 0.  1. -3.]
 [ 1.  0. -1.]
 [ 0.  0.  1.]]
k(global_coords)[source]

Calculates the relative position of points in vector direction.

Parameters:global_coords (array_like(Number, shape=(n, k))) – Represents n points of k dimensions in a global coordinate system.
Returns:Relative relative position of points in vector direction. The origin of the vector defines zero and the target defines one.
Return type:np.ndarray(Number, shape=(k))

Examples

>>> v = Vector((1, 1, 1), (2, 3, 4))
>>> ks = v.k([v.origin - v.vec, v.target, v.origin - 2 * v.vec])
>>> print_rounded(ks)
[-1.  1. -2.]
length
origin
t
target
transform(T)[source]

Transforms the vector using a transformation matrix.

Parameters:T (array_like(Number, shape=(self.dim+1, self.dim+1))) – Transformation matrix to apply.
Returns:
Return type:self

Examples

Create a vector and a transformation matix.

>>> vec = Vector([1, 1], [1, 0])
>>> r = 45 * np.pi / 180.0
>>> t = [1, 2]
>>> T = transformation.matrix(t=t, r=r, order='rts')
>>> print_rounded(T)
[[ 0.71 -0.71  1.  ]
 [ 0.71  0.71  2.  ]
 [ 0.    0.    1.  ]]
>>> vec = vec.transform(T)
>>> print_rounded(vec.origin)
[ 1.    3.41]
>>> print_rounded(vec.vec)
[ 0.71  0.71]
vec
pyoints.vector.angle(v, w, deg=False)[source]

Calculate angle between two vectors.

Parameters:
  • w (v,) – Vector or n vectors of k dimensions.
  • deg (optional, bool) – Indicates whether or not the angle is returned in degree.
Returns:

Angle between vectors v and w.

Return type:

Number

Examples

Angle between two dimensional vectors.

>>> angle([0, 1], [1, 0], deg=True)
90.0
>>> round(angle([0, 1], [1, 1], deg=True), 1)
45.0
>>> angle([2, 2], [1, 1], deg=True)
0.0
>>> angle([0, 0], [1, 1], deg=True)
inf

Angle between three dimensional vectors.

>>> angle([1, 0, 1], [0, 1, 0], deg=True)
90.0

4D >>> angle([1, 0, 0, 0], [0, 1, 1, 0], deg=True) 90.0

Multiple vectors at once.

>>> print_rounded(angle([[0, 1], [1, 1]], [[1, 0], [2, 0]], deg=True))
[ 90.  45.]
pyoints.vector.axes_angles(v, deg=False)[source]

Calculates the angles of a vetor to all coordinate axes.

Parameters:
  • v (array_like(Number, shape=(k))) – Vector of k dimensions.
  • deg (optional, bool) – Indicates whether or not the angles are returned in degree.
Returns:

Rotation angles.

Return type:

np.ndarray(Number, shape=(k))

Examples

>>> v = [1, 1]
>>> print_rounded(axes_angles(v, deg=True))
[ 45.  45.]
>>> v = [0, 1, 0]
>>> print_rounded(axes_angles(v, deg=True))
[ 90.   0.  90.]
pyoints.vector.azimuth(v, deg=False)[source]

Calculates the azimuth angle of a vector or vectors.

Parameters:v (array_like(Number, shape=(k)) or array_like(Number, shape=(n, k))) – Vector or vectors of k dimensions.
Returns:Azimuth angle for each vector.
Return type:Number or np.ndarray(Number, shape=(n))

Examples

>>> v = [1, 1]
>>> print_rounded(azimuth(v, deg=True))
45.0
>>> v = [1, 1, 5]
>>> print_rounded(azimuth(v, deg=True))
45.0
>>> v = [(0, 0), (0, 1), (1, 1), (1, 0), (2, -2), (0, -1), (-1, 1)]
>>> print_rounded(azimuth(v, deg=True))
[  nan    0.   45.   90.  135.  180.  315.]
pyoints.vector.basis(vec, origin=None)[source]

Generates a local coordinate system based on a vector. The local coordinate system is represented by a rotation matrix.

Parameters:
  • vec (array_like(Number, shape=(k))) – Vector of k dimensions to defines the direction of the first coordinate axis of the new coordinate system. The other k-1 axes are build perpendicular to it and each other.
  • origin (optional, array_like(Number, shape=(k))) – Defines the origin of the local coordinate system. If None, no shift is assumed.
Returns:

Transformation matrix representing the desired coordinate system. Since the norm of vector vec has no influence on the scale, the norm of all basic vectors is one.

Return type:

np.ndarray(Number, shape=(k+1, k+1))

Examples

Create some two dimensional coordinates.

>>> coords = [(0, 0), (0, 1), (1, 0), (1, 1), (0.5, 0.5), (-1, -1)]

Flip the basic axes.

>>> b = basis([0, 1])
>>> print_rounded(b)
[[ 0.  1.  0.]
 [ 1.  0.  0.]
 [ 0.  0.  1.]]
>>> local_coords = transformation.transform(coords, b)
>>> print_rounded(local_coords)
[[ 0.   0. ]
 [ 1.   0. ]
 [ 0.   1. ]
 [ 1.   1. ]
 [ 0.5  0.5]
 [-1.  -1. ]]

Keep the original orientation, but set a new origin.

>>> b = basis([2, 0], origin=[2, 3])
>>> print_rounded(b)
[[ 1.  0. -2.]
 [ 0.  1. -3.]
 [ 0.  0.  1.]]
>>> local_coords = transformation.transform(coords, b)
>>> print_rounded(local_coords)
[[-2.  -3. ]
 [-2.  -2. ]
 [-1.  -3. ]
 [-1.  -2. ]
 [-1.5 -2.5]
 [-3.  -4. ]]

Use a diagonal basis.

>>> b = basis([3, 4])
>>> print_rounded(b)
[[ 0.6  0.8  0. ]
 [-0.8  0.6  0. ]
 [ 0.   0.   1. ]]
>>> local_coords = transformation.transform(coords, b)
>>> print_rounded(local_coords)
[[ 0.   0. ]
 [ 0.8  0.6]
 [ 0.6 -0.8]
 [ 1.4 -0.2]
 [ 0.7 -0.1]
 [-1.4  0.2]]

Three dimensional case.

>>> b = basis([3, -4, 0], origin=[1, 2, 3])
>>> print_rounded(b)
[[-0.6  0.8  0.  -1. ]
 [ 0.   0.   1.  -3. ]
 [-0.8 -0.6  0.   2. ]
 [ 0.   0.   0.   1. ]]
pyoints.vector.deg_to_rad(deg)[source]

Converts angles from degree to radiant.

Parameters:deg (Number or array_like(Number, shape=(k))) – Angle or angles in degree.
Returns:Angle or angles in radiant.
Return type:Number or np.ndarray(Number, shape=(k))

See also

rad_to_deg()

Examples

>>> print_rounded(deg_to_rad(90), 3)
1.571
>>> rad = deg_to_rad([0, 45, 180, 360])
>>> print_rounded(rad, 3)
[ 0.     0.785  3.142  6.283]
pyoints.vector.direction(v, deg=False)[source]

Calculate the direction angles of a vector. This direction can be used to create a rotation matrix.

Parameters:
  • v (array_like(Number, shape=(k))) – Vector of k dimensions.
  • deg (optional, bool) – Indicates whether or not the direction is returned in degree.
Returns:

Rotation angles.

Return type:

Number or np.ndarray(Number, shape=(k))

Examples

>>> v = [0, 1]
>>> print_rounded(direction(v, deg=True))
90.0
>>> v = [0, 1, 1]
>>> print_rounded(direction(v, deg=True))
[ 90.  45.]
pyoints.vector.orthogonal(v, w)[source]

Check whether or not two vectors are orthogonal.

Parameters:v,w (array_like(Number, shape=(k))) – Vector of k dimensions.
Returns:True, if v is orthogonal to w.
Return type:boolean

Examples

>>> orthogonal([1, 1], [1, -1])
True
>>> orthogonal([1, 1], [1, 0])
False
pyoints.vector.rad_to_deg(rad)[source]

Convert angles from radiant to degree.

Parameters:rad (Number or array_like(Number, shape=(k))) – Angle or angles in radiant.
Returns:Angle or angles in degree.
Return type:Number or np.ndarray(Number, shape=(k))

See also

deg_to_rad()

Examples

>>> rad_to_deg(0.5*np.pi)
90.0
>>> print_rounded(rad_to_deg([0, np.pi/4, np.pi, 2*np.pi]))
[   0.   45.  180.  360.]
pyoints.vector.scalarproduct(v, w)[source]

Calculates the scalar product or dot product of two vectors.

Parameters:v,w (array_like(Number, shape=(k))) – Vector of k dimensions.
Returns:Scalar product or dot product of the vectors.
Return type:float

Examples

>>> scalarproduct([1, 2], [3, 4])
11
>>> scalarproduct([1, 2, 3], [4, 5, 6])
32
>>> scalarproduct([1, 1], [1, -1])
0
pyoints.vector.vector_plane_intersection(vec, plane)[source]

Calculates the intersection point of a k dimensional vector and a k dimensional plane.

Parameters:
  • vec (Vector) – Vector to calculate the intersection point for.
  • plane (Plane) – Plane to calculate the intersection point for.
Returns:

coord – Intersection point of the vector and the plane.

Return type:

np.array_like(Number, shape=(k))

Notes

The algorithm solves the linear equation system: plane.origin - vec.origin | k * vec.vec - s * plane.vec

Examples

Create a plane and a vector in three dimensional space and calculate the intersection point.

>>> vec = Vector((1, 1, -1), (0, 0, 3))
>>> plane = Plane((-1, 2, 5), (0, 2, 0), (1, 2, 2))
>>> intersection = vector_plane_intersection(vec, plane)
>>> print(intersection)
[ 1.  1.  9.]
>>> print_rounded(vec.distance([intersection]))
[ 0.]

Create a plane and a vector in two dimensional space and calculate the intersection point.

>>> vec = Vector((1, 1), (1, 2))
>>> plane = Plane((-1, 5), (0, 1))
>>> intersection = vector_plane_intersection(vec, plane)
>>> print_rounded(intersection)
[-1. -3.]
>>> print_rounded(vec.distance([intersection]))
[ 0.]
pyoints.vector.vector_surface_intersection(vec, surface, eps=0.001, max_iter=20)[source]

Approximates the intersection point between a k dimensional vector and a k dimensional surface iteratively.

Parameters:
  • vec (Vector) – Vector to calculate the intersection point for.
  • surface (callable) – Surface model. The model needs to receive coordinates as an argument and needs to return the distance to the surface.
Returns:

coord – Approximate intersection point between the vector and the surface.

Return type:

np.array_like(Number, shape=(k))

Examples

>>> from pyoints import surface, interpolate

Create a callable surface object and a vector.

>>> surface = surface.Surface(
...         [(0, 0, 0), (0, 2, 0), (2, 1, 4)],
...         method=interpolate.LinearInterpolator
...     )
>>> vec = Vector((1, 1, -1), (0, 0, 1))

Calculate the intersection point.

>>> print_rounded(vector_surface_intersection(vec, surface))
[ 1.  1.  2.]
pyoints.vector.zenith(v, axis=-1, deg=False)[source]

Angle between a vector or vectors in relation to a specific coordinate axes.

Parameters:
  • v (array_like(Number, shape=(k)) or array_like(Number, shape=(n, k))) – Vector or vectors of k dimensions.
  • axis (optional, int) – Defines which axis to compare the vector with. If not provided, the last dimension is used.
  • deg (optional, bool) – Provide the angle in degree.
Returns:

Angle between the vector and the selected coordinate axis.

Return type:

Number or np.ndarray(Number, shape=(n))

Examples

>>> v = [1, 0]
>>> print_rounded(zenith(v, deg=True))
90.0
>>> v = [0, 0, 1]
>>> print_rounded(zenith(v, deg=True))
0.0
>>> v = [(0, 0), (1, 0), (0, 1), (1, 1)]
>>> print_rounded(zenith(v, deg=True))
[ nan  90.   0.  45.]
>>> v = [1, 0, 1]
>>> print_rounded(zenith([1, 0, 1], axis=2, deg=True))
45.0

Module contents

Pyoints: A Python package for point cloud, voxel and raster processing.