Statistical Learning

ants.sparse_decom2(inmatrix, inmask=(None, None), sparseness=(0.01, 0.01), nvecs=3, its=20, cthresh=(0, 0), statdir=None, perms=0, uselong=0, z=0, smooth=0, robust=0, mycoption=0, initialization_list=[], initialization_list2=[], ell1=10, prior_weight=0, verbose=False, rejector=0, max_based=False, version=1)[source]

Decomposes two matrices into paired sparse eigenevectors to maximize canonical correlation - aka Sparse CCA. Note: we do not scale the matrices internally. We leave scaling choices to the user.

ANTsR function: sparseDecom2

Parameters:
  • inmatrix (2-tuple of ndarrays) – input as inmatrix=(mat1,mat2). n by p input matrix and n by q input matrix , spatial variable lies along columns.

  • inmask (2-tuple of ANTsImage python:types (optional - one or both)) – optional pair of image masks

  • sparseness (tuple) – a pair of float values e.g c(0.01,0.1) enforces an unsigned 99 percent and 90 percent sparse solution for each respective view

  • nvecs (python:integer) – number of eigenvector pairs

  • its (python:integer) – number of iterations, 10 or 20 usually sufficient

  • cthresh (2-tuple) – cluster threshold pair

  • statdir (string (optional)) – temporary directory if you want to look at full output

  • perms (python:integer) – number of permutations. settings permutations greater than 0 will estimate significance per vector empirically. For small datasets, these may be conservative. p-values depend on how one scales the input matrices.

  • uselong (boolean) – enforce solutions of both views to be the same - requires matrices to be the same size

  • z (float) – subject space (low-dimensional space) sparseness value

  • smooth (float) – smooth the data (only available when mask is used)

  • robust (boolean) – rank transform input matrices

  • mycoption (python:integer) – enforce 1 - spatial orthogonality, 2 - low-dimensional orthogonality or 0 - both

  • initialization_list (list) – initialization for first view

  • initialization_list2 (list) – initialization for 2nd view

  • ell1 (float) – gradient descent parameter, if negative then l0 otherwise use l1

  • prior_weight (scalar) – Scalar value weight on prior between 0 (prior is weak) and 1 (prior is strong). Only engaged if initialization is used

  • verbose (boolean) – activates verbose output to screen

  • rejector (scalar) – rejects small correlation solutions

  • max_based (boolean) – whether to choose max-based thresholding

Returns:

projectionsndarray

X projections

projections2ndarray

Y projections

eig1ndarray

X components

eig2ndarray

Y components

summarypd.DataFrame

first column is canonical correlations, second column is p-values (these are None if perms > 0)

Return type:

dict w/ following key/value pairs

Example

>>> import numpy as np
>>> import ants
>>> mat = np.random.randn(20, 100)
>>> mat2 = np.random.randn(20, 90)
>>> mydecom = ants.sparse_decom2(inmatrix = (mat,mat2),
                                sparseness=(0.1,0.3), nvecs=3,
                                its=3, perms=0)
ants.eig_seg(mask, img_list, apply_segmentation_to_images=False, cthresh=0, smooth=1)[source]

Segment a mask into regions based on the max value in an image list. At a given voxel the segmentation label will contain the index to the image that has the largest value. If the 3rd image has the greatest value, the segmentation label will be 3 at that voxel.

Parameters:
  • mask (ANTsImage) – D-dimensional mask > 0 defining segmentation region.

  • img_list (collection of ANTsImage or np.ndarray) – images to use

  • apply_segmentation_to_images (boolean) – determines if original image list is modified by the segmentation.

  • cthresh (python:integer) – throw away isolated clusters smaller than this value

  • smooth (float) – smooth the input data first by this value

Return type:

ANTsImage

Example

>>> import ants
>>> mylist = [ants.image_read(ants.get_ants_data('r16')),
              ants.image_read(ants.get_ants_data('r27')),
              ants.image_read(ants.get_ants_data('r85'))]
>>> myseg = ants.eig_seg(ants.get_mask(mylist[0]), mylist)
ants.initialize_eigenanatomy(initmat, mask=None, initlabels=None, nreps=1, smoothing=0)[source]

InitializeEigenanatomy is a helper function to initialize sparseDecom and sparseDecom2. Can be used to estimate sparseness parameters per eigenvector. The user then only chooses nvecs and optional regularization parameters.

Parameters:
  • initmat (np.ndarray or ANTsImage) – input matrix where rows provide initial vector values. alternatively, this can be an antsImage which contains labeled regions.

  • mask (ANTsImage) – mask if available

  • initlabels (list/tuple of python:integers) – which labels in initmat to use as initial components

  • nreps (python:integer) – nrepetitions to use

  • smoothing (float) – if using an initial label image, optionally smooth each roi

Returns:

initlistlist of ANTsImage types

initialization list(s) for sparseDecom(2)

maskANTsImage

mask(s) for sparseDecom(2)

enameslist of strings

string names of components for sparseDecom(2)

Return type:

dict w/ the following key/value pairs

Example

>>> import ants
>>> import numpy as np
>>> mat = np.random.randn(4,100).astype('float32')
>>> init = ants.initialize_eigenanatomy(mat)