pycapacity.algorithms module
Overview
There are three main polytope evaluation algorithms implemented:
hyper_plane_shift_method: Hyper plane shifting method implementation used to find half-space representation of problems of a form
vertex_enumeration_vepoli2: Efficient vertex enumeration algorithm for a problem of a form:
iterative_convex_hull_method: A function calculating the polytopes given by the equations form:
Additionally this module implements different helping functions to half-plane and vertex representation manipulations of the polytopes:
hspace_to_vertex: Function transforming H-representation to V-representation
vertex_to_hspace: Function transforming V-representation to H-representation
vertex_to_faces: Helping function transforming vertices of the polytope to a triangulated faces
- pycapacity.algorithms.chebyshev_ball(A, b)[source]
Calculating chebyshev ball of a polytope given the half-space representation
- Parameters:
A (list) – matrix of half-space representation Ax<b
b (list) – vector of half-space representation Ax<b
- Returns:
returns a chebyshev center of the polytope radius(float): returns a chebyshev radius of the polytope
- Return type:
center(array)
- pycapacity.algorithms.chebyshev_center(A, b)[source]
Calculating chebyshev center of a polytope given the half-space representation
- Parameters:
A (list) – matrix of half-space representation Ax<b
b (list) – vector of half-space representation Ax<b
- Returns:
returns a chebyshev center of the polytope
- Return type:
center(array)
- pycapacity.algorithms.face_index_to_vertex(vertices, indexes)[source]
Helping function for transforming the list of faces with indexes to the vertices
- Parameters:
vertices – list of vertices
indexes – list of vertex indexes forming faces
- Returns:
list of faces composed of vertices
- Return type:
faces
- pycapacity.algorithms.hspace_to_vertex(H, d)[source]
From half-space representation to the vertex representation
- Parameters:
H (list) – matrix of half-space representation Hx<d
d (list) – vector of half-space representation Hx<d
- Returns:
vertices(list) (vertices of the polytope)
face_indexes(list) (indexes of vertices forming triangulated faces of the polytope)
- pycapacity.algorithms.hyper_plane_shift_method(A, x_min, x_max, tol=1e-15)[source]
Hyper plane shifting method implementation used to solve problems of a form:
\[P = \{y~ |~y = Ax, \quad x_{min} \leq x \leq x_{max}\}\]Note
Gouttefarde M., Krut S. (2010) Characterization of Parallel Manipulator Available Wrench Set Facets. In: Lenarcic J., Stanisic M. (eds) Advances in Robot Kinematics: Motion in Man and Machine. Springer, Dordrecht
This algorithm can be used to calculate acceleration polytope, velocity polytoe and even polytope of the joint achievable joint torques based on the muscle forces
- Parameters:
A – projection matrix
x_min – minimal values
x_max – maximal values
- Returns:
H(list) – matrix of half-space representation Hx<d
d(list) – vector of half-space representation Hx<d
- Raises:
ValueError – if the limits are not of the same size as the projection matrix
ValueError – if the limits are not of the same size
- pycapacity.algorithms.iterative_convex_hull_method(A, B, y_min, y_max, tol, P=None, bias=None, G_in=None, h_in=None, G_eq=None, h_eq=None, max_iter=1000, verbose=False)[source]
A function calculating the polytopes of achievable x for equations form:
\[P_x = \{x~ |~ Ax = By,\quad y_{min} \leq y \leq y_{max}\}\](optionally - additional projection matrix)
\[P_x = \{x~ |~Pz = x,~~ Az = By,\quad y_{min} \leq y \leq y_{max}\}\](optionally - additional bias \(b\))
\[P_x = \{x~ |~Ax = By + b, \quad y_{min} \leq y \leq y_{max}\}\](optionally - additional inequality constaints)
\[G_{in} y \leq h_{in}\](optionally - additional equality constaints)
\[G_{eq} y = h_{eq}\]Finally, the most general equation supported is
\[P_x = \{x~ |~Pz = x,~~ Az = By,~~ G_{in} y \leq h_{in} ,~~ G_{eq} y = h_{eq}, ~~ y_{min} \leq y \leq y_{max}\}\]Note
On-line feasible wrench polytope evaluation based on human musculoskeletal models: an iterative convex hull method A.Skuric,V.Padois,N.Rezzoug,D.Daney
- Parameters:
A – matrix
B – matrix
y_min – minimal values
y_max – maximal values
tol – tolerance for the polytope calculation
P – an additional projection matrix
bias – bias in the intermediate space
G_in – matrix - inegality constraint G_in y < h_in
h_in – vector - inegality constraint G_in y < h_in
G_eq – matrix - equality constraint G_eq y = h_eq
h_eq – vector - equality constraint G_eq y = h_eq
max_iter – maximum number of iterations (number of linera programmes solved) - default 1000
verbose – verbose program output - showing execution warnings - default False
- Returns:
x_vert(list) – list of vertices
H(list) – matrix of half-space representation Hx<d
d(list) – vector of half-space representation Hx<d
faces(list) – list of vertex indexes forming polytope faces
y_vert(list) – list of y values producing x_vert
z_vert(list) – list of z values producing x_vert
- Raises:
ValueError – if the dimensions of the input matrices are not valid
ValueError – if the limits dimensions are not valid
ValueError – if the limits min and max are not valid
ValueError – if the inequality/equality constraints dimensions are not valid
ValueError – if the the projection matrix dimensions are not valid
- pycapacity.algorithms.order_index(points)[source]
Order clockwise 2D points
- Parameters:
points – matrix of 2D points
- Returns
indexes(array) : ordered indexes
- pycapacity.algorithms.stack(A, B, dir='v')[source]
Helping function enabling vertical and horizontal stacking of numpy arrays
- pycapacity.algorithms.vertex_enumeration_vepoli2(A, b_max, b_min, b_bias=None)[source]
Efficient vertex enumeration algorithm for a problem of a form:
\[P = \{x~ |~Ax = b, \quad b_{min} \leq b \leq b_{max}\}\]Optional (bias \(b_{bias}\) can be added):
\[P = \{x~ |~Ax = b - b_{bias}, \quad b_{min} \leq b \leq b_{max}\}\]Note
On-line force capability evaluation based on efficient polytope vertex search by A. Skuric, V. Padois, D. Daney
- Parameters:
A – system matrix A
b_max – maximal b
b_min – minimal b
b_bias – b bias vector ( offset from 0 )
- Returns:
x_vertex(list) (vertices of the polytope)
b_vartex(list) (b values producing x_vertex)
- Raises:
ValueError – if the A matrix dimensions are not valid
ValueError – if the b_min and b_max are not of the same size as the A matrix
ValueError – if the b_bias is not of the same size as the A matrix
- pycapacity.algorithms.vertex_to_faces(vertex)[source]
Function grouping the vertices to faces using a ConvexHull algorithm
- Parameters:
vertex (array) – list of vertices
- Returns:
list of triangle faces with vertex indexes which form them
- Return type:
faces(array)
- pycapacity.algorithms.vertex_to_hspace(vertex)[source]
Function transforming vertices to half-space representation using a ConvexHull algorithm
- Parameters:
vertex (array) – list of vertices
- Returns:
H(list) (matrix of half-space representation Hx<d)
d(list) (vector of half-space representation Hx<d)