Commit e105516afadc50f51d7f314536b0d2a8b89c38a5

Authored by David Mayerich
1 parent 403139ae

added spherical harmonics file and updated network code to support vectorization

Showing 2 changed files with 154 additions and 0 deletions   Show diff stats
python/network.py
... ... @@ -76,6 +76,7 @@ class Fiber:
76 76  
77 77 #print(volume)
78 78 return volume
  79 +
79 80 '''
80 81 Writes the header given and open file descripion, number of verticies and number of edges.
81 82 '''
... ... @@ -331,3 +332,33 @@ def distancefield(nodeList, edgeList, R=(100, 100, 100)):
331 332 D, I = tree.query(Q)
332 333  
333 334 return D
  335 +
  336 +#returns the number of points in the network
  337 +def npoints(N, F):
  338 + n = 0
  339 + for f in F:
  340 + n = n + len(f.points) - 2
  341 + n = n + len(N)
  342 + return n
  343 +
  344 +#returns the number of linear segments in the network
  345 +def nsegments(N, F):
  346 + n = 0
  347 + for f in F:
  348 + n = n + len(f.points) - 1
  349 + return n
  350 +
  351 +def vectorize(N, F):
  352 + #initialize three coordinate vectors ()
  353 + n = nsegments(N, F)
  354 + X = np.zeros((n))
  355 + Y = np.zeros((n))
  356 + Z = np.zeros((n))
  357 + idx = 0
  358 + for i in range(0, len(F)):
  359 + for j in range(0, len(F[i].points)-1):
  360 + X[idx] = F[i].points[j][0]-F[i].points[j+1][0]
  361 + Y[idx] = F[i].points[j][1]-F[i].points[j+1][1]
  362 + Z[idx] = F[i].points[j][2]-F[i].points[j+1][2]
  363 + idx = idx + 1
  364 + return X, Y, Z
... ...
python/spharmonics.py 0 → 100644
  1 +# -*- coding: utf-8 -*-
  2 +"""
  3 +Created on Mon Dec 18 16:31:36 2017
  4 +
  5 +@author: david
  6 +"""
  7 +
  8 +import numpy
  9 +import scipy
  10 +import matplotlib.pyplot as plt
  11 +from matplotlib import cm, colors
  12 +from mpl_toolkits.mplot3d import Axes3D
  13 +import math
  14 +
  15 +
  16 +def sph2cart(theta, phi, r):
  17 + x = r * numpy.cos(theta) * numpy.sin(phi)
  18 + y = r * numpy.sin(theta) * numpy.sin(phi)
  19 + z = r * numpy.cos(phi)
  20 +
  21 + return x, y, z
  22 +
  23 +def cart2sph(x,y,z):
  24 + r = numpy.sqrt(x**2+y**2+z**2)
  25 + theta = numpy.arctan2(y,x)
  26 + phi = numpy.arccos(z/r)
  27 + #if(x == 0):
  28 + # phi = 0
  29 + #else:
  30 + # phi = numpy.arccos(z/r)
  31 + #phi = numpy.arctan2(numpy.sqrt(x**2 + y**2), z)
  32 + return r, theta, phi
  33 +
  34 +
  35 +def sh(theta, phi, l, m):
  36 +
  37 + if m < 0:
  38 + return numpy.sqrt(2) * (-1)**m * scipy.special.sph_harm(abs(m), l, theta, phi).imag
  39 + elif m > 0:
  40 + return numpy.sqrt(2) * (-1)**m * scipy.special.sph_harm(m, l, theta, phi).real
  41 + else:
  42 + return scipy.special.sph_harm(0, l, theta, phi).real
  43 +
  44 +#calculate a spherical harmonic value from a set of coefficients and a spherical coordinate
  45 +def sh_coeff(theta, phi, C):
  46 +
  47 + s = numpy.zeros(theta.shape)
  48 + c = range(len(C))
  49 +
  50 + for c in range(len(C)):
  51 + l, m = i2lm(c) #get the 2D indices
  52 + s = s + C[c] * sh(theta, phi, l, m)
  53 +
  54 + return s
  55 +
  56 +#plot a spherical harmonic function on a sphere using N points
  57 +def sh_plot(C, N):
  58 + phi = numpy.linspace(0, numpy.pi, N)
  59 + theta = numpy.linspace(0, 2*numpy.pi, N)
  60 + phi, theta = numpy.meshgrid(phi, theta)
  61 +
  62 + # The Cartesian coordinates of the unit sphere
  63 + x = numpy.sin(phi) * numpy.cos(theta)
  64 + y = numpy.sin(phi) * numpy.sin(theta)
  65 + z = numpy.cos(phi)
  66 +
  67 + # Calculate the spherical harmonic Y(l,m) and normalize to [0,1]
  68 + fcolors = sh_coeff(theta, phi, C)
  69 + fmax, fmin = fcolors.max(), fcolors.min()
  70 + fcolors = (fcolors - fmin)/(fmax - fmin)
  71 +
  72 + # Set the aspect ratio to 1 so our sphere looks spherical
  73 + fig = plt.figure(figsize=plt.figaspect(1.))
  74 + ax = fig.add_subplot(111, projection='3d')
  75 + ax.plot_surface(x, y, z, rstride=1, cstride=1, facecolors=cm.seismic(fcolors))
  76 + # Turn off the axis planes
  77 + ax.set_axis_off()
  78 + plt.show()
  79 +
  80 +def i2lm(i):
  81 + l = numpy.floor(numpy.sqrt(i))
  82 + m = i - l *(l + 1)
  83 + return l, m
  84 +
  85 +def lm2i(l, m):
  86 + return l * (l+1) + m
  87 +
  88 +#generates a set of spherical harmonic coefficients from samples using linear least squares
  89 +def linfit(theta, phi, s, nc):
  90 + #allocate space for the matrix
  91 + A = numpy.zeros((nc, nc))
  92 +
  93 + #calculate each of the matrix coefficients
  94 + #(see SH technical report in the vascular_viz repository)
  95 + for i in range(nc):
  96 + li, mi = i2lm(i)
  97 + yi = sh(theta, phi, li, mi)
  98 + for j in range(nc):
  99 + lj, mj = i2lm(j)
  100 + yj = sh(theta, phi, lj, mj)
  101 + A[i, j] = numpy.sum(yi * yj)
  102 +
  103 + #calculate the RHS values
  104 + b = numpy.zeros(nc)
  105 + for j in range(nc):
  106 + lj, mj = i2lm(j)
  107 + yj = sh(theta, phi, lj, mj)
  108 + b[j] = numpy.sum(yj * s)
  109 +
  110 + #solve the system of linear equations
  111 + return numpy.linalg.solve(A, b)
  112 +
  113 +#generate a scatter plot in 3D using spherical coordinates
  114 +def scatterplot3d(theta, phi, r):
  115 + #convert all of the samples to cartesian coordinates
  116 + X, Y, Z = sph2cart(theta, phi, r)
  117 +
  118 + fig = plt.figure()
  119 + ax = fig.add_subplot(111, projection='3d')
  120 + ax.scatter(X, Y, Z)
  121 + plt.show()
  122 +
  123 +
... ...