def setup_taylor_grad(self): s_manifold = self.mdims['s_manifold'] m_gradients = self.mdims['m_gradients'] r_points = self.mdims['r_points'] l_labels = self.mdims['l_labels'] D = np.eye(r_points) - 1.0 / r_points * np.ones((r_points, r_points)) for j in range(m_gradients): # v : columns correspond to the vertices of triangle j v = self.v[:, self.faces[:, j]] # v0 : centroid of triangle j v0 = self.v_grad[:, j] # E : each row is a tangent vector at v0 pointing to a vertex E = np.zeros((3, 3)) self.log(v0, v.T, E) # C : basis of tangent space at v0, orthonormal wrt. euclidean metric c1 = normalize(E[0, :]).ravel() c2 = normalize(E[1, :] - E[1, :].dot(c1) * c1).ravel() C = np.vstack((c1, c2)).T # M { r_points, s_manifold } M = E.dot(C) self.A[j] = M.T.dot(D).dot(M) self.B[j] = M.T.dot(D)
def __init__(self, refinement=2, vecs=None): """ Setup a triangular grid on the 2-sphere. Args: vecs : array of shape (3, l_labels) refinement : if vecs is None, the triangulation is generated using trisphere(refinement), defined below """ if vecs is None: vecs, tris = trisphere(refinement) vecs = normalize(vecs) else: assert (vecs.shape[0] == 3) vecs = normalize(vecs) sv = SphericalVoronoi(vecs.T) sv.sort_vertices_of_regions() tris = sv._tri.simplices.T assert (tris.shape[0] == 3) s_manifold = 2 m_gradients = tris.shape[1] r_points = 3 l_labels = vecs.shape[1] logging.info("Set up for 2-sphere with {l} labels and {m} gradients." \ .format(s=s_manifold, l=l_labels,m=m_gradients)) self.v = vecs self.faces = tris self.mdims = { 's_manifold': s_manifold, 'm_gradients': m_gradients, 'r_points': r_points, 'l_labels': l_labels } self.v_grad = np.zeros((3, m_gradients), order='C') for j in range(m_gradients): # v : columns correspond to the vertices of triangle j v = self.v[:, self.faces[:, j]] # v0 : centroid of triangle j # taking the mean is a bottle neck here! self.v_grad[:, j] = self.mean(v).ravel() # b { l_labels } # P { m_gradients, r_points } # b : vol. elements, s.t. sum(self.b) converges to 4*PI as n to inf self.b = np.zeros((l_labels, ), order='C') self.b[:] = sphere_areas(vecs) self.b_precond = 1.0 / np.sqrt(np.einsum('i,i->', self.b, self.b)) self.P = np.zeros((m_gradients, r_points), order='C', dtype=np.int64) self.P[:] = tris.T # A { m_gradients, s_manifold, s_manifold } # B { m_gradients, s_manifold, r_points } self.A = np.zeros((m_gradients, s_manifold, s_manifold), order='C') self.B = np.zeros((m_gradients, s_manifold, r_points), order='C') self.setup_taylor_grad()
def embed_barycentric(self, x): """ Embed a 3D vector into the triangulation using barycentric coordinates. Args: x : the vector to be embedded, numpy array of shape (3,) Returns: numpy array of shape (l_labels,) Test code: >>> import numpy as np >>> from qball.sphere import load_sphere >>> n = 2 >>> mf = load_sphere(refinement=n) >>> for theta, phi in zip([23,42,76,155,143,103], [345,220,174,20,18,87]): >>> f_in = np.asarray([ >>> np.sin(theta*np.pi/180.0)*np.cos(phi*np.pi/180.0), >>> np.sin(theta*np.pi/180.0)*np.sin(phi*np.pi/180.0), >>> np.cos(theta*np.pi/180.0) >>> ]) >>> f_out = mf.mean(mf.v, mf.embed_barycentric(f_in)[:,None]).ravel() >>> err = np.sqrt(sum((f_in-f_out)**2)) >>> print("theta: {}, phi: {}, err: {}".format(theta, phi, err)) >>> assert err <= 1e-15 """ xn = normalize(x)[:, 0] labels, tri = self.get_triangle(xn) assert labels is not None _, _, _, alph = barylinear_prep(self, tri, xn) weights = np.zeros((self.mdims['l_labels'], ), order='C') weights[labels] = alph return weights
def barylinear_prep(sph, tri, x): """ Helper function for the calculation of barylinear interpolations in tri at x on the sphere sph Args: sph : instance of qball.sphere.Sphere tri : array of shape (3, 3), columns correspond to vertices v_i x : point at which the interpolation is to be evaluated Returns: N : array of shape (3, 3), each row is a unit tangent vector at x pointing to a vertex of tri C : array of shape (2, 3), orthonormal basis for tangent space at x A : array of shape (2, 2), matrix for computing barycentric coordinates of x; expressed in the basis C alph : array of shape (3,), barycentric coordinates of x """ s_manifold = 2 r_points = 3 # E : each row is a tangent vector at x pointing to a vertex of tri E = np.zeros((3, 3)) sph.log(x, tri.T, E) N = normalize(E.T).T c2 = normalize(N[1] - N[1].dot(N[0]) * N[0]).ravel() C = np.vstack((N[0], c2)) # transform tangent vectors to basis C N = N.dot(C.T) E = E.dot(C.T) A = np.vstack((E[1] - E[0], E[2] - E[0])).T alph = np.hstack((0, np.linalg.solve(A, -E[0]))) alph[0] = 1 - sum(alph[1:]) return N, C, A, alph
def synth_unimodal_odfs(qball_sphere, sphere_vol, directions, cutoff=2 * np.pi, const_width=3, tightness=4.9): verts = qball_sphere.vertices l_labels = verts.shape[0] # `const_width` unimodals for each entry of `directions` val_base = 1e-6 * 290 vals = [tightness * val_base, val_base, val_base] vecs = normalize( np.array( [[np.sin(phi * np.pi / 180.0), np.cos(phi * np.pi / 180.0), 0] for phi in directions]).T).T voxels = [ multi_tensor_odf(verts, np.array((vals, vals)), [v, v], [50, 50]) for v in vecs ] if cutoff < 1.5 * np.pi: # cut off support of distribution functions for v, vox in zip(vecs, voxels): for k in range(l_labels): v_diff1 = verts[k] - v v_diff2 = verts[k] + v if np.einsum('n,n->', v_diff1, v_diff1) > cutoff**2 \ and np.einsum('n,n->', v_diff2, v_diff2) > cutoff**2: vox[k] = 0 fin = np.zeros((l_labels, const_width * len(voxels)), order='C') for i, vox in enumerate(voxels): i1 = i * const_width i2 = (i + 1) * const_width fin[:, i1:i2] = np.tile(vox, (const_width, 1)).T normalize_odf(fin, sphere_vol) return fin
def trisphere(n): """ Calculates a sphere triangulation of 12*(4^n) vertices. The algorithm starts from an icosahedron (12 vertices, 20 triangles) and applies the following procedure `n` times: Each triangular face is split into four triangles using the triangle's edge centers as new vertices. Args: n : grade of refinement; n=0 means 12 vertices (icosahedron). Returns: tuple (verts, tris) verts : numpy array, each column corresponds to a point on the sphere tris : numpy array, each column defines a triangle on the sphere through indices into `verts` """ # Regular icosahedron: # (X, Z) : solution to X/Z = Z/(X + Z) and 1 = X^2 + Z^2 X = ((5 - 5**0.5) / 10)**0.5 Z = (1 - X**2)**0.5 # verts { 3, 12 } verts = np.array([[-X, 0.0, Z], [X, 0.0, Z], [-X, 0.0, -Z], [X, 0.0, -Z], [0.0, Z, X], [0.0, Z, -X], [0.0, -Z, X], [0.0, -Z, -X], [Z, X, 0.0], [-Z, X, 0.0], [Z, -X, 0.0], [-Z, -X, 0.0]]).T # tris { 3, 20 } tris = np.array([[0, 4, 1], [0, 9, 4], [9, 5, 4], [4, 5, 8], [4, 8, 1], [8, 10, 1], [8, 3, 10], [5, 3, 8], [5, 2, 3], [2, 7, 3], [7, 10, 3], [7, 6, 10], [7, 11, 6], [11, 0, 6], [0, 1, 6], [6, 1, 10], [9, 0, 11], [9, 11, 2], [9, 2, 5], [7, 2, 11]]).T # t = np.array([[0, 1, 4]]).T; for i in range(n): vn = verts # numv : current number of verts numv = vn.shape[1] # edgecenters : symmetric matrix of indices into `vn` edgecenters = np.zeros([verts.shape[1]] * 2, dtype=np.int64) # iterate over triangles, adding edge centers to `vn` and # making note of it in `edgecenters` for tri in tris.T: # the if clauses make sure we don't get duplicates if tri[0] < tri[1]: vn = np.hstack( (vn, normalize(0.5 * (verts[:, tri[0]] + verts[:, tri[1]])))) edgecenters[tri[0], tri[1]] = numv edgecenters[tri[1], tri[0]] = numv numv += 1 if tri[1] < tri[2]: vn = np.hstack( (vn, normalize(0.5 * (verts[:, tri[1]] + verts[:, tri[2]])))) edgecenters[tri[1], tri[2]] = numv edgecenters[tri[2], tri[1]] = numv numv += 1 if tri[2] < tri[0]: vn = np.hstack( (vn, normalize(0.5 * (verts[:, tri[2]] + verts[:, tri[0]])))) edgecenters[tri[2], tri[0]] = numv edgecenters[tri[0], tri[2]] = numv numv += 1 # tn { 3, 4*tris.shape[1] } # build up as list then convert to numpy array tn = [] # each triangle is split into four using the edge centers for tri in tris.T: a = edgecenters[tri[1], tri[2]] b = edgecenters[tri[2], tri[0]] c = edgecenters[tri[0], tri[1]] if a == 0 or b == 0 or c == 0: print("Damn") tn.extend([[tri[0], c, b], [tri[1], a, c], [tri[2], b, a], [a, b, c]]) tn = np.asarray(tn).T verts = vn tris = tn return verts, tris
def exp(self, location, vfrom): """ exp_l(v) = cos(|v|) * l + sin(|v|) * v/|v| """ c = np.sqrt(np.einsum('i,i->', vfrom, vfrom)) pto = np.cos(c) * location + np.sin(c) / max(np.spacing(1), c) * vfrom # normalizing prevents errors from accumulating return normalize(pto).reshape(-1)