HalfspaceIntersection#
- class scipy.spatial.HalfspaceIntersection(halfspaces, interior_point, incremental=False, qhull_options=None)#
N 维半空间交集。
于 0.19.0 版本中添加。
- 参数:
- halfspaces浮点数 ndarray,形状为 (nineq, ndim+1)
如下格式的叠加不等式 Ax + b <= 0 [A; b]
- interior_point浮点数 ndarray,形状为 (ndim,)
半空间定义区域内的明确内点。也称为可行点,可通过线性规划获取。
- incremental布尔值,可选
允许逐步添加新的半空间。这会占用一些额外的资源。
- qhull_options字符串,可选
传给 Qhull 的其他选项。请参阅 Qhull 手册以了解详情。(默认值:ndim > 4 为 “Qx”,否则为 “”)选项“H”总是启用的。
- 引发:
- QhullError
当 Qhull 遇到错误条件(比如未启用解决选项时的几何简并)时引发。
- ValueError
当输入不兼容的数组时引发。
注意
交集使用Qhull 库计算。这再现了 Qhull 的 “qhalf” 功能。
参考
[Qhull][1]S. Boyd,L. Vandenberghe,凸优化,可访问 http://stanford.edu/~boyd/cvxbook/
示例
形成多边形的平面半空间交集
>>> from scipy.spatial import HalfspaceIntersection >>> import numpy as np >>> halfspaces = np.array([[-1, 0., 0.], ... [0., -1., 0.], ... [2., 1., -4.], ... [-0.5, 1., -2.]]) >>> feasible_point = np.array([0.5, 0.5]) >>> hs = HalfspaceIntersection(halfspaces, feasible_point)
将半空间绘制为填充区域和交点
>>> import matplotlib.pyplot as plt >>> fig = plt.figure() >>> ax = fig.add_subplot(1, 1, 1, aspect='equal') >>> xlim, ylim = (-1, 3), (-1, 3) >>> ax.set_xlim(xlim) >>> ax.set_ylim(ylim) >>> x = np.linspace(-1, 3, 100) >>> symbols = ['-', '+', 'x', '*'] >>> signs = [0, 0, -1, -1] >>> fmt = {"color": None, "edgecolor": "b", "alpha": 0.5} >>> for h, sym, sign in zip(halfspaces, symbols, signs): ... hlist = h.tolist() ... fmt["hatch"] = sym ... if h[1]== 0: ... ax.axvline(-h[2]/h[0], label='{}x+{}y+{}=0'.format(*hlist)) ... xi = np.linspace(xlim[sign], -h[2]/h[0], 100) ... ax.fill_between(xi, ylim[0], ylim[1], **fmt) ... else: ... ax.plot(x, (-h[2]-h[0]*x)/h[1], label='{}x+{}y+{}=0'.format(*hlist)) ... ax.fill_between(x, (-h[2]-h[0]*x)/h[1], ylim[sign], **fmt) >>> x, y = zip(*hs.intersections) >>> ax.plot(x, y, 'o', markersize=8)
默认情况下,qhull 并未提供计算内点的方法。可轻松使用线性规划对其进行计算。考虑到形式为 \(Ax + b \leq 0\) 的半空间,编写线性规划
\[ \begin{align}\begin{aligned}max \: y\\s.t. Ax + y ||A_i|| \leq -b\end{aligned}\end{align} \]其中 \(A_i\) 是 A 的行,即每个平面的法线。
将产生一个最深入凸多面体的点 x。准确来说,它是内接在多面体中半径为 y 的最大超球的中心。此点称为多面体的切比雪夫中心(参见 [1] 4.3.1,pp148-149)。Qhull 输出的方程始终经过规范化。
>>> from scipy.optimize import linprog >>> from matplotlib.patches import Circle >>> norm_vector = np.reshape(np.linalg.norm(halfspaces[:, :-1], axis=1), ... (halfspaces.shape[0], 1)) >>> c = np.zeros((halfspaces.shape[1],)) >>> c[-1] = -1 >>> A = np.hstack((halfspaces[:, :-1], norm_vector)) >>> b = - halfspaces[:, -1:] >>> res = linprog(c, A_ub=A, b_ub=b, bounds=(None, None)) >>> x = res.x[:-1] >>> y = res.x[-1] >>> circle = Circle(x, radius=y, alpha=0.3) >>> ax.add_patch(circle) >>> plt.legend(bbox_to_anchor=(1.6, 1.0)) >>> plt.show()
- 属性:
- halfspacesdouble 的 ndarray,形状为 (nineq, ndim+1)
输入半空间。
- interior_point:浮点数的 ndarray,形状为 (ndim,)
输入内点。
- intersectionsdouble 的 ndarray,形状为 (ninter, ndim)
所有半空间的交集。
- dual_pointsdouble 的 ndarray,形状为 (nineq, ndim)
输入半空间的对偶点。
- dual_facetsint 列表列表
形成对偶凸包(不一定简单)阶段的点的索引。
- dual_verticesint 的 ndarray,形状为 (nvertices,)
形成对偶凸包顶点的半空间的索引。对于 2-D 凸包,顶点按逆时针顺序排列。对于其他维度,它们按输入顺序排列。
- dual_equationsdouble 的 ndarray,形状为 (nfacet, ndim+1)
[normal, offset] 形成对偶阶段的超平面方程(有关更多信息,请参阅 Qhull 文档)。
- dual_areafloat
对偶凸包的面积
- dual_volumefloat
对偶凸包的体积
方法
add_halfspaces
(halfspaces[, restart])处理一组额外的全新半空间。
close
()完成增量处理。