def test_equivalence_csd_and_parametric_fod( odi=0.15, mu=[0., 0.], lambda_par=1.7e-9): stick = cylinder_models.C1Stick() watsonstick = distribute_models.SD1WatsonDistributed( [stick]) params = {'SD1Watson_1_odi': odi, 'SD1Watson_1_mu': mu, 'C1Stick_1_lambda_par': lambda_par} data = watsonstick(scheme, **params) sh_mod = modeling_framework.MultiCompartmentSphericalHarmonicsModel( [stick]) sh_mod.set_fixed_parameter('C1Stick_1_lambda_par', lambda_par) sh_fit = sh_mod.fit(scheme, data) fod = sh_fit.fod(sphere.vertices) watson = distributions.SD1Watson(mu=[0., 0.], odi=0.15) sf = watson(sphere.vertices) assert_array_almost_equal(fod[0], sf, 2) fitted_signal = sh_fit.predict() assert_array_almost_equal(data, fitted_signal[0], 4)
def test_watson_kappa(): # test for Wn(k2) > Wn(k1) when k2>k1 along mu random_mu = np.random.rand(2) random_mu_cart = utils.sphere2cart(np.r_[1., random_mu]) odi1 = .8 odi2 = .6 watson = distributions.SD1Watson(mu=random_mu) Wn1 = watson(n=random_mu_cart, odi=odi1) Wn2 = watson(n=random_mu_cart, odi=odi2) assert_equal(Wn2 > Wn1, True)
def test_bingham_equal_to_watson(beta_fraction=0): # test if bingham with beta=0 equals watson distribution mu_ = np.random.rand(2) n_cart = utils.sphere2cart(np.r_[1., mu_]) psi_ = np.random.rand() * np.pi odi_ = np.random.rand() bingham = distributions.SD2Bingham(mu=mu_, psi=psi_, odi=odi_, beta_fraction=beta_fraction) watson = distributions.SD1Watson(mu=mu_, odi=odi_) Bn = bingham(n=n_cart) Wn = watson(n=n_cart) assert_almost_equal(Bn, Wn, 3)
def test_construction_observation_matrix( odi=0.15, mu=[0., 0.], lambda_par=1.7e-9, lmax=8): stick = cylinder_models.C1Stick(lambda_par=lambda_par) watsonstick = distribute_models.SD1WatsonDistributed( [stick]) params = {'SD1Watson_1_odi': odi, 'SD1Watson_1_mu': mu, 'C1Stick_1_lambda_par': lambda_par} data = watsonstick(scheme, **params) watson = distributions.SD1Watson(mu=mu, odi=odi) sh_watson = watson.spherical_harmonics_representation(lmax) stick_rh = stick.rotational_harmonics_representation(scheme) A = construct_model_based_A_matrix(scheme, stick_rh, lmax) data_approximated = A.dot(sh_watson) np.testing.assert_array_almost_equal(data_approximated, data, 4)
def test_watson_integral_unity(): # test for integral to unity for isotropic concentration (k=0) odi = 1. # isotropic distribution # first test for one orientation n random_n = np.random.rand(3) random_n /= np.linalg.norm(random_n) random_mu = np.random.rand(2) watson = distributions.SD1Watson(mu=random_mu, odi=odi) Wn = watson(n=random_n) # in random direction spherical_integral = Wn * 4 * np.pi # for isotropic distribution assert_equal(spherical_integral, 1.) # third test for unity when k>0 sphere = get_sphere('repulsion724') n_sphere = sphere.vertices odi = np.random.rand() Wn_sphere = watson(n=n_sphere, mu=random_mu, odi=odi) spherical_integral = sum(Wn_sphere) / n_sphere.shape[0] * 4 * np.pi assert_almost_equal(spherical_integral, 1., 3)
def test_watson_orienting(): # test for orienting the axis of the Watson distribution mu for k>0 # first test to see if Wn is highest along mu sphere = get_sphere('repulsion724') n = sphere.vertices indices = np.array(range(n.shape[0])) np.random.shuffle(indices) mu_index = indices[0] mu_cart = n[mu_index] mu_sphere = utils.cart2sphere(mu_cart)[1:] odi = np.random.rand() watson = distributions.SD1Watson(mu=mu_sphere, odi=odi) Wn_vector = watson(n=n) assert_almost_equal(Wn_vector[mu_index], max(Wn_vector)) # second test to see if Wn is lowest prependicular to mu mu_perp = utils.perpendicular_vector(mu_cart) Wn_perp = watson(n=mu_perp) assert_equal(np.all(Wn_perp < Wn_vector), True)
def test_spherical_convolution_watson_sh(sh_order=4): sphere = get_sphere('symmetric724') n = sphere.vertices bval = np.tile(1e9, len(n)) scheme = acquisition_scheme_from_bvalues(bval, n, delta, Delta) indices_sphere_orientations = np.arange(sphere.vertices.shape[0]) np.random.shuffle(indices_sphere_orientations) mu_index = indices_sphere_orientations[0] mu_watson = sphere.vertices[mu_index] mu_watson_sphere = utils.cart2sphere(mu_watson)[1:] watson = distributions.SD1Watson(mu=mu_watson_sphere, odi=.3) f_sf = watson(n=sphere.vertices) f_sh = sf_to_sh(f_sf, sphere, sh_order) lambda_par = 2e-9 stick = cylinder_models.C1Stick(mu=[0, 0], lambda_par=lambda_par) k_sf = stick(scheme) sh_matrix, m, n = real_sym_sh_mrtrix(sh_order, sphere.theta, sphere.phi) sh_matrix_inv = np.linalg.pinv(sh_matrix) k_sh = np.dot(sh_matrix_inv, k_sf) k_rh = k_sh[m == 0] fk_convolved_sh = sh_convolution(f_sh, k_rh) fk_convolved_sf = sh_to_sf(fk_convolved_sh, sphere, sh_order) # assert if spherical mean is the same between kernel and convolved kernel assert_almost_equal(abs(np.mean(k_sf) - np.mean(fk_convolved_sf)), 0., 2) # assert if the lowest signal attenuation (E(b,n)) is orientation along # the orientation of the watson distribution. min_position = np.argmin(fk_convolved_sf) if min_position == mu_index: assert_equal(min_position, mu_index) else: # then it's the opposite direction sphere_positions = np.arange(sphere.vertices.shape[0]) opposite_index = np.all( np.round(sphere.vertices - mu_watson, 2) == 0, axis=1 ) min_position_opposite = sphere_positions[opposite_index] assert_equal(min_position_opposite, mu_index)
def test_multi_voxel_parametric_to_sm_to_sh_fod_watson(): stick = cylinder_models.C1Stick() zeppelin = gaussian_models.G2Zeppelin() watsonstick = distribute_models.SD1WatsonDistributed([stick, zeppelin]) watsonstick.set_equal_parameter('G2Zeppelin_1_lambda_par', 'C1Stick_1_lambda_par') watsonstick.set_tortuous_parameter('G2Zeppelin_1_lambda_perp', 'G2Zeppelin_1_lambda_par', 'partial_volume_0') mc_mod = modeling_framework.MultiCompartmentModel([watsonstick]) parameter_dict = { 'SD1WatsonDistributed_1_SD1Watson_1_mu': np.random.rand(10, 2), 'SD1WatsonDistributed_1_partial_volume_0': np.linspace(0.1, 0.9, 10), 'SD1WatsonDistributed_1_G2Zeppelin_1_lambda_par': np.linspace(1.5, 2.5, 10) * 1e-9, 'SD1WatsonDistributed_1_SD1Watson_1_odi': np.linspace(0.3, 0.7, 10) } data = mc_mod.simulate_signal(scheme, parameter_dict) sm_mod = modeling_framework.MultiCompartmentSphericalMeanModel( [stick, zeppelin]) sm_mod.set_equal_parameter('G2Zeppelin_1_lambda_par', 'C1Stick_1_lambda_par') sm_mod.set_tortuous_parameter('G2Zeppelin_1_lambda_perp', 'G2Zeppelin_1_lambda_par', 'partial_volume_0', 'partial_volume_1') sf_watson = [] for mu, odi in zip( parameter_dict['SD1WatsonDistributed_1_SD1Watson_1_mu'], parameter_dict['SD1WatsonDistributed_1_SD1Watson_1_odi']): watson = distributions.SD1Watson(mu=mu, odi=odi) sf_watson.append(watson(sphere.vertices)) sf_watson = np.array(sf_watson) sm_fit = sm_mod.fit(scheme, data) sh_mod = sm_fit.return_spherical_harmonics_fod_model() sh_fit_auto = sh_mod.fit(scheme, data) # will pick tournier fod_tournier = sh_fit_auto.fod(sphere.vertices) assert_array_almost_equal(fod_tournier, sf_watson, 1) sh_fit_tournier = sh_mod.fit(scheme, data, solver='csd_tournier07', unity_constraint=False) fod_tournier = sh_fit_tournier.fod(sphere.vertices) assert_array_almost_equal(fod_tournier, sf_watson, 1) sh_fit_cvxpy = sh_mod.fit(scheme, data, solver='csd_cvxpy', unity_constraint=True, lambda_lb=0.) fod_cvxpy = sh_fit_cvxpy.fod(sphere.vertices) assert_array_almost_equal(fod_cvxpy, sf_watson, 2) sh_fit_cvxpy = sh_mod.fit(scheme, data, solver='csd_cvxpy', unity_constraint=False, lambda_lb=0.) fod_cvxpy = sh_fit_cvxpy.fod(sphere.vertices) assert_array_almost_equal(fod_cvxpy, sf_watson, 2)