Пример #1
0
    def load_kernels(self):
        """Load rotated kernels and project to the specific gradient scheme of this subject.
        Dispatch to the proper function, depending on the model.
        """
        if self.model is None:
            ERROR('Model not set; call "set_model()" method first')
        if self.scheme is None:
            ERROR('Scheme not loaded; call "load_data()" first')

        tic = time.time()
        LOG('\n-> Resampling LUT for subject "%s":' %
            self.get_config('subject'))

        # auxiliary data structures
        idx_OUT, Ylm_OUT = amico.lut.aux_structures_resample(
            self.scheme, self.get_config('lmax'))

        # hash table
        self.htable = amico.lut.load_precomputed_hash_table(
            self.get_config('ndirs'))

        # Dispatch to the right handler for each model
        self.KERNELS = self.model.resample(self.get_config('ATOMS_path'),
                                           idx_OUT, Ylm_OUT,
                                           self.get_config('doMergeB0'),
                                           self.get_config('ndirs'))

        LOG('   [ %.1f seconds ]' % (time.time() - tic))
Пример #2
0
    def generate( self, out_path, aux, idx_in, idx_out, ndirs ) :
        if self.scheme.version != 1 :
            ERROR( 'This model requires a "VERSION: STEJSKALTANNER" scheme' )

        scheme_high = amico.lut.create_high_resolution_scheme( self.scheme, b_scale=1E6 )
        filename_scheme = pjoin( out_path, 'scheme.txt' )
        np.savetxt( filename_scheme, scheme_high.raw, fmt='%15.8e', delimiter=' ', header='VERSION: STEJSKALTANNER', comments='' )

        # temporary file where to store "datasynth" output
        filename_signal = pjoin( tempfile._get_default_tempdir(), next(tempfile._get_candidate_names())+'.Bfloat' )

        nATOMS = len(self.Rs) + len(self.d_perps) + len(self.d_isos)
        progress = ProgressBar( n=nATOMS, prefix="   ", erase=False )

        # Cylinder(s)
        for R in self.Rs :
            CMD = 'datasynth -synthmodel compartment 1 CYLINDERGPD %E 0 0 %E -schemefile %s -voxels 1 -outputfile %s 2> /dev/null' % ( self.d_par*1E-6, R, filename_scheme, filename_signal )
            subprocess.call( CMD, shell=True )
            if not exists( filename_signal ) :
                ERROR( 'Problems generating the signal with "datasynth"' )
            signal  = np.fromfile( filename_signal, dtype='>f4' )
            if exists( filename_signal ) :
                remove( filename_signal )

            lm = amico.lut.rotate_kernel( signal, aux, idx_in, idx_out, False, ndirs )
            np.save( pjoin( out_path, 'A_%03d.npy'%progress.i ), lm )
            progress.update()

        # Zeppelin(s)
        for d in self.d_perps :
            CMD = 'datasynth -synthmodel compartment 1 ZEPPELIN %E 0 0 %E -schemefile %s -voxels 1 -outputfile %s 2> /dev/null' % ( self.d_par*1E-6, d*1e-6, filename_scheme, filename_signal )
            subprocess.call( CMD, shell=True )
            if not exists( filename_signal ) :
                ERROR( 'Problems generating the signal with "datasynth"' )
            signal  = np.fromfile( filename_signal, dtype='>f4' )
            if exists( filename_signal ) :
                remove( filename_signal )

            lm = amico.lut.rotate_kernel( signal, aux, idx_in, idx_out, False, ndirs )
            np.save( pjoin( out_path, 'A_%03d.npy'%progress.i ), lm )
            progress.update()

        # Ball(s)
        for d in self.d_isos :
            CMD = 'datasynth -synthmodel compartment 1 BALL %E -schemefile %s -voxels 1 -outputfile %s 2> /dev/null' % ( d*1e-6, filename_scheme, filename_signal )
            subprocess.call( CMD, shell=True )
            if not exists( filename_signal ) :
                ERROR( 'Problems generating the signal with "datasynth"' )
            signal  = np.fromfile( filename_signal, dtype='>f4' )
            if exists( filename_signal ) :
                remove( filename_signal )

            lm = amico.lut.rotate_kernel( signal, aux, idx_in, idx_out, True, ndirs )
            np.save( pjoin( out_path, 'A_%03d.npy'%progress.i ), lm )
            progress.update()
Пример #3
0
    def generate_kernels(self, regenerate=False, lmax=12, ndirs=32761):
        """Generate the high-resolution response functions for each compartment.
        Dispatch to the proper function, depending on the model.

        Parameters
        ----------
        regenerate : boolean
            Regenerate kernels if they already exist (default : False)
        lmax : int
            Maximum SH order to use for the rotation procedure (default : 12)
        ndirs : int
            Number of directions on the half of the sphere representing the possible orientations of the response functions (default : 32761)
         """
        if self.scheme is None:
            ERROR('Scheme not loaded; call "load_data()" first')
        if self.model is None:
            ERROR('Model not set; call "set_model()" method first')
        if not is_valid(ndirs):
            ERROR(
                'Unsupported value for ndirs.\nNote: Supported values for ndirs are [500, 1000, 1500, 2000, 2500, 3000, 3500, 4000, 4500, 5000, 5500, 6000, 6500, 7000, 7500, 8000, 8500, 9000, 9500, 10000, 32761 (default)]'
            )

        # store some values for later use
        self.set_config('lmax', lmax)
        self.set_config('ndirs', ndirs)
        self.model.scheme = self.scheme
        LOG('\n-> Creating LUT for "%s" model:' % self.model.name)

        # check if kernels were already generated
        tmp = glob.glob(pjoin(self.get_config('ATOMS_path'), 'A_*.npy'))
        if len(tmp) > 0 and not regenerate:
            LOG('   [ LUT already computed. USe option "regenerate=True" to force regeneration ]'
                )
            return

        # create folder or delete existing files (if any)
        if not exists(self.get_config('ATOMS_path')):
            makedirs(self.get_config('ATOMS_path'))
        else:
            for f in glob.glob(pjoin(self.get_config('ATOMS_path'), '*')):
                remove(f)

        # auxiliary data structures
        aux = amico.lut.load_precomputed_rotation_matrices(lmax, ndirs)
        idx_IN, idx_OUT = amico.lut.aux_structures_generate(self.scheme, lmax)

        # Dispatch to the right handler for each model
        tic = time.time()
        self.model.generate(self.get_config('ATOMS_path'), aux, idx_IN,
                            idx_OUT, ndirs)
        LOG('   [ %.1f seconds ]' % (time.time() - tic))
Пример #4
0
    def fit( self, y, dirs, KERNELS, params, htable ) :
        nD = dirs.shape[0]
        if nD > 1 : # model works only with one direction
            ERROR( '"%s" model requires exactly 1 orientation' % self.name )

        n1 = len(self.d_perps)
        n2 = len(self.d_isos)
        nATOMS = n1+n2
        if nATOMS == 0 : # empty dictionary
            return [0, 0], None, None, None

        # prepare DICTIONARY from dir and lookup tables
        lut_idx = amico.lut.dir_TO_lut_idx( dirs[0], htable )
        A = np.zeros( (len(y), nATOMS), dtype=np.float64, order='F' )
        A[:,:(nD*n1)] = KERNELS['D'][:,lut_idx,:].T
        A[:,(nD*n1):] = KERNELS['CSF'].T

        # fit
        x = spams.lasso( np.asfortranarray( y.reshape(-1,1) ), D=A, **params ).todense().A1

        # return estimates
        v = x[ :n1 ].sum() / ( x.sum() + 1e-16 )

        # checking that there is more than 1 isotropic compartment
        if self.type == 'Mouse' :
            v_blood = x[ n1 ] / ( x.sum() + 1e-16 )
            v_csf = x[ n1+1 ] / ( x.sum() + 1e-16 )

            return [ v, 1-v, v_blood, v_csf ], dirs, x, A

        else :
            return [ v, 1-v ], dirs, x, A
Пример #5
0
    def resample( self, in_path, idx_out, Ylm_out, doMergeB0, ndirs ) :
        if doMergeB0:
            nS = 1+self.scheme.dwi_count
            merge_idx = np.hstack((self.scheme.b0_idx[0],self.scheme.dwi_idx))
        else:
            nS = self.scheme.nS
            merge_idx = np.arange(nS)
        KERNELS = {}
        KERNELS['model'] = self.id
        KERNELS['D']     = np.zeros( (len(self.d_perps),ndirs,nS), dtype=np.float32 )
        KERNELS['CSF']   = np.zeros( (len(self.d_isos),nS), dtype=np.float32 )

        nATOMS = len(self.d_perps) + len(self.d_isos)
        progress = ProgressBar( n=nATOMS, prefix="   ", erase=False )

        # Tensor compartment(s)
        for i in range(len(self.d_perps)) :
            lm = np.load( pjoin( in_path, 'A_%03d.npy'%progress.i ) )
            if lm.shape[0] != ndirs:
                ERROR( 'Outdated LUT. Call "generate_kernels( regenerate=True )" to update the LUT' )
            KERNELS['D'][i,...] = amico.lut.resample_kernel( lm, self.scheme.nS, idx_out, Ylm_out, False, ndirs )[:,merge_idx]
            progress.update()

        # Isotropic compartment(s)
        for i in range(len(self.d_isos)) :
            lm = np.load( pjoin( in_path, 'A_%03d.npy'%progress.i ) )
            KERNELS['CSF'][i,...] = amico.lut.resample_kernel( lm, self.scheme.nS, idx_out, Ylm_out, True, ndirs )[merge_idx]
            progress.update()

        return KERNELS
Пример #6
0
    def get_table( self) :
        """Return a matrix from the structure.

        The first three columns represent the gradient directions.
        Then, we return two formats to describe each gradient:
            - if the shape of data is Nx4, the 4^ column is the b-value;
            - if the shape of data is Nx7, the last 4 columns are, respectively, the gradient strength, big delta, small delta and TE.
        
        """
        if self.raw is None:
            ERROR( 'The structure has not been created.' )
         
        if self.version == 0 :
            scheme_table = np.zeros([self.b0_count + self.dwi_count, 4])
            for shell in self.shells:
                scheme_table[shell['idx'], 0:3 ]    = shell['grad']
                scheme_table[shell['idx'], 3 ]      = shell['b']

        if self.version == 1 :
            scheme_table = np.zeros([self.b0_count + self.dwi_count, 7])
            for shell in self.shells:
                scheme_table[shell['idx'], 0:3 ]    = shell['grad']
                scheme_table[shell['idx'], 3 ]      = shell['G']
                scheme_table[shell['idx'], 4 ]      = shell['Delta']
                scheme_table[shell['idx'], 5 ]      = shell['delta']
                scheme_table[shell['idx'], 6 ]      = shell['TE']
                    
        return scheme_table
Пример #7
0
def load_precomputed_rotation_matrices(lmax, ndirs):
    """Load precomputed the rotation matrices to rotate the high-resolution kernels.

    Parameters
    ----------
    lmax : int
        Maximum SH order to use for the rotation phase
    ndirs : int
        Number of directions on the half of the sphere representing the possible orientations of the response functions
    """
    filename = pjoin(dipy_home,
                     f'AMICO_aux_matrices_lmax={lmax}_ndirs={ndirs}.pickle')
    if not isfile(filename):
        ERROR('Auxiliary matrices not found; call "amico.setup()" first')
    with open(filename, 'rb') as rotation_matrices_file:
        if sys.version_info.major == 3:
            aux = pickle.load(rotation_matrices_file,
                              fix_imports=True,
                              encoding='bytes')
            # Pickle files written by Python 2 are loaded with byte
            # keys, whereas those written by Python 3 are loaded with
            # str keys, even when both are written using protocol=2 as
            # in precompute_rotation_matrices
            result_aux = {(k.decode() if hasattr(k, "decode") else k): v
                          for k, v in aux.items()}
            return result_aux
        else:
            return pickle.load(rotation_matrices_file)
Пример #8
0
 def set_solver(self, **params):
     """Setup the specific parameters of the solver to fit the model.
     Dispatch to the proper function, depending on the model; a model shoudl provide a "set_solver" function to set these parameters.
     """
     if self.model is None:
         ERROR('Model not set; call "set_model()" method first')
     self.set_config('solver_params', self.model.set_solver(**params))
Пример #9
0
    def __init__( self, data, b0_thr = 0 ) :
        """Initialize the acquisition scheme.

        Parameters
        ----------
        data : string or numpy.ndarray
            The filename of the scheme or a matrix containing the actual values
        b0_thr : float
            The threshold on the b-values to identify the b0 images (default: 0)
        """
        if type(data) is str :
            # try loading from file
            try :
                n = 0 # headers lines to skip to get to the numeric data
                with open(data) as fid :
                    for line in fid :
                        if re_match( r'[+-]?(\d+(\.\d*)?|\.\d+)([eE][+-]?\d+)?', line.strip() ) :
                            break
                        n += 1

                tmp = np.loadtxt( data, skiprows=n )

            except :
                ERROR( 'Unable to open scheme file' )

            self.load_from_table( tmp, b0_thr )
        else :
            # try loading from matrix
            self.load_from_table( data, b0_thr )
Пример #10
0
    def cyl_neuman_le_perp_PGSE( self, d, R, G, delta, smalldel, roots ):

        # When R=0, no need to do any calculation
        if (R == 0.00):
            LE = np.zeros(G.shape) # np.size(R) = 1
            return LE
        else:
            ERROR( '"cyl_neuman_le_perp_PGSE" not yet validated for non-zero values' )
Пример #11
0
    def legendre_gaussian_integral( self, Lpmp, n ):
        if n > 6:
            ERROR( 'The maximum value for n is 6, which corresponds to the 12th order Legendre polynomial' )
        exact = Lpmp>0.05
        approx = Lpmp<=0.05

        mn = n + 1

        I = np.zeros((len(Lpmp),mn))
        sqrtx = np.sqrt(Lpmp[exact])
        I[exact,0] = np.sqrt(np.pi)*scipy.special.erf(sqrtx)/sqrtx
        dx = 1.0/Lpmp[exact]
        emx = -np.exp(-Lpmp[exact])
        for i in range(1,mn):
            I[exact,i] = emx + (i-0.5)*I[exact,i-1]
            I[exact,i] = I[exact,i]*dx

        # Computing the legendre gaussian integrals for large enough Lpmp
        L = np.zeros((len(Lpmp),n+1))
        for i in range(n+1):
            if i == 0:
                L[exact,0] = I[exact,0]
            elif i == 1:
                L[exact,1] = -0.5*I[exact,0] + 1.5*I[exact,1]
            elif i == 2:
                L[exact,2] = 0.375*I[exact,0] - 3.75*I[exact,1] + 4.375*I[exact,2]
            elif i == 3:
                L[exact,3] = -0.3125*I[exact,0] + 6.5625*I[exact,1] - 19.6875*I[exact,2] + 14.4375*I[exact,3]
            elif i == 4:
                L[exact,4] = 0.2734375*I[exact,0] - 9.84375*I[exact,1] + 54.140625*I[exact,2] - 93.84375*I[exact,3] + 50.2734375*I[exact,4]
            elif i == 5:
                L[exact,5] = -(63./256.)*I[exact,0] + (3465./256.)*I[exact,1] - (30030./256.)*I[exact,2] + (90090./256.)*I[exact,3] - (109395./256.)*I[exact,4] + (46189./256.)*I[exact,5]
            elif i == 6:
                L[exact,6] = (231./1024.)*I[exact,0] - (18018./1024.)*I[exact,1] + (225225./1024.)*I[exact,2] - (1021020./1024.)*I[exact,3] + (2078505./1024.)*I[exact,4] - (1939938./1024.)*I[exact,5] + (676039./1024.)*I[exact,6]

        # Computing the legendre gaussian integrals for small Lpmp
        x2=np.power(Lpmp[approx],2)
        x3=x2*Lpmp[approx]
        x4=x3*Lpmp[approx]
        x5=x4*Lpmp[approx]
        x6=x5*Lpmp[approx]
        for i in range(n+1):
            if i == 0:
                L[approx,0] = 2 - 2*Lpmp[approx]/3 + x2/5 - x3/21 + x4/108
            elif i == 1:
                L[approx,1] = -4*Lpmp[approx]/15 + 4*x2/35 - 2*x3/63 + 2*x4/297
            elif i == 2:
                L[approx,2] = 8*x2/315 - 8*x3/693 + 4*x4/1287
            elif i == 3:
                L[approx,3] = -16*x3/9009 + 16*x4/19305
            elif i == 4:
                L[approx,4] = 32*x4/328185
            elif i == 5:
                L[approx,5] = -64*x5/14549535
            elif i == 6:
                L[approx,6] = 128*x6/760543875

        return L
Пример #12
0
    def synth_meas_iso_GPD( self, d, protocol ):
        if protocol['pulseseq'] != 'PGSE' and protocol['pulseseq'] != 'STEAM':
            ERROR( 'synth_meas_iso_GPD() : Protocol %s not translated from NODDI matlab code yet' % protocol['pulseseq'] )

        GAMMA = 2.675987E8
        modQ = GAMMA*protocol['smalldel'].transpose()*protocol['gradient_strength'].transpose()
        modQ_Sq = np.power(modQ,2)
        difftime = protocol['delta'].transpose()-protocol['smalldel']/3.0
        return np.exp(-difftime*modQ_Sq*d)
Пример #13
0
    def fit( self, y, dirs, KERNELS, params, htable ) :
        singleb0 = True if len(y) == (1+self.scheme.dwi_count) else False
        nD = dirs.shape[0]
        if nD != 1 :
            ERROR( '"%s" model requires exactly 1 orientation' % self.name )

        # prepare DICTIONARY from dir and lookup tables
        nWM = len(self.IC_ODs)*len(self.IC_VFs)
        nATOMS = nWM + 1
        if self.isExvivo == True :
            nATOMS += 1
        lut_idx = amico.lut.dir_TO_lut_idx( dirs[0], htable )
        A = np.ones( (len(y), nATOMS), dtype=np.float64, order='F' )
        A[:,:nWM] = KERNELS['wm'][:,lut_idx,:].T
        A[:,-1]  = KERNELS['iso']


        # estimate CSF partial volume (and isotropic restriction, if exvivo) and remove from signal
        x, _ = scipy.optimize.nnls( A, y )
        yy = y - x[-1]*A[:,-1]
        if self.isExvivo == True :
            yy = yy - x[-2]*A[:,-2]
        yy[ yy<0 ] = 0

        # estimate IC and EC compartments and promote sparsity
        if singleb0:
            An = A[1:, :nWM] * KERNELS['norms']
            yy = yy[1:].reshape(-1,1)
        else:
            An = A[ self.scheme.dwi_idx, :nWM ] * KERNELS['norms']
            yy = yy[ self.scheme.dwi_idx ].reshape(-1,1)
        x = spams.lasso( np.asfortranarray(yy), D=np.asfortranarray(An), **params ).todense().A1

        # debias coefficients
        x = np.append( x, 1 )
        if self.isExvivo == True :
            x = np.append( x, 1 )
        idx = x>0
        x[idx], _ = scipy.optimize.nnls( A[:,idx], y )

        # return estimates
        xx = x / ( x.sum() + 1e-16 )
        xWM  = xx[:nWM]
        fISO = xx[-1]
        xWM = xWM / ( xWM.sum() + 1e-16 )
        f1 = np.dot( KERNELS['icvf'], xWM )
        f2 = np.dot( (1.0-KERNELS['icvf']), xWM )
        v = f1 / ( f1 + f2 + 1e-16 )
        k = np.dot( KERNELS['kappa'], xWM )
        od = 2.0/np.pi * np.arctan2(1.0,k)

        if self.isExvivo:
            return [v, od, fISO, xx[-2]], dirs, x, A
        else:
            return [v, od, fISO], dirs, x, A
Пример #14
0
def resample_kernel(KRlm, nS, idx_out, Ylm_out, is_isotropic, ndirs):
    """Project/resample a spherical function to signal space.

    Parameters
    ----------
    KRlm : numpy.array
        Rotated spherical functions (in SH space) to project
    nS : integer
        Number of samples in the subject's acquisition scheme
    idx_out : list of list
        Index of samples in output kernel
    Ylm_out : numpy.array
        Matrix to project back all shells from SH space to signal space (of the subject)
    is_isotropic : boolean
        Indentifies whether Klm is an isotropic function or not
    ndirs : int
        Number of directions on the half of the sphere representing the possible orientations of the response functions

    Returns
    -------
    KR = numpy.array
        Rotated spherical functions projected to signal space of the subject
    """
    if is_isotropic == False:
        KR = np.ones((ndirs, nS), dtype=np.float32)
        try:
            for i in range(ndirs):
                KR[i, idx_out] = np.dot(Ylm_out, KRlm[i, :]).astype(np.float32)
        except:
            ERROR(
                'Outdated LUT. Call "generate_kernels( regenerate=True )" to update the LUT'
            )
    else:
        KR = np.ones(nS, dtype=np.float32)
        try:
            for i in range(ndirs):
                KR[idx_out] = np.dot(Ylm_out, KRlm).astype(np.float32)
        except:
            ERROR(
                'Outdated LUT. Call "generate_kernels( regenerate=True )" to update the LUT'
            )
    return KR
Пример #15
0
def setup( lmax = 12, ndirs = 32761 ) :
    """General setup/initialization of the AMICO framework.
    
    Parameters
    ----------
    lmax : int
        Maximum SH order to use for the rotation phase (default : 12)
    ndirs : int
        Number of directions on the half of the sphere representing the possible orientations of the response functions (default : 32761)
    """
    if not is_valid(ndirs):
        ERROR( 'Unsupported value for ndirs.\nNote: supported values for ndirs are [500, 1000, 1500, 2000, 2500, 3000, 3500, 4000, 4500, 5000, 5500, 6000, 6500, 7000, 7500, 8000, 8500, 9000, 9500, 10000, 32761 (default)]' )
    
    amico.lut.precompute_rotation_matrices( lmax, ndirs )
Пример #16
0
    def set_model( self, model_name ) :
        """Set the model to use to describe the signal contributions in each voxel.

        Parameters
        ----------
        model_name : string
            The name of the model (must match a class name in "amico.models" module)
        """
        # Call the specific model constructor
        if hasattr(amico.models, model_name ) :
            self.model = getattr(amico.models,model_name)()
        else :
            ERROR( 'Model "%s" not recognized' % model_name )

        self.set_config('ATOMS_path', pjoin( self.get_config('study_path'), 'kernels', self.model.id ))

        # setup default parameters for fitting the model (can be changed later on)
        self.set_solver()
Пример #17
0
    def resample( self, in_path, idx_out, Ylm_out, doMergeB0, ndirs ):
        nATOMS = len(self.IC_ODs)*len(self.IC_VFs) + 1
        if doMergeB0:
            nS = 1+self.scheme.dwi_count
            merge_idx = np.hstack((self.scheme.b0_idx[0],self.scheme.dwi_idx))
        else:
            nS = self.scheme.nS
            merge_idx = np.arange(nS)
        KERNELS = {}
        KERNELS['model'] = self.id
        KERNELS['wm']    = np.zeros( (nATOMS-1,ndirs,nS), dtype=np.float32 )
        KERNELS['iso']   = np.zeros( nS, dtype=np.float32 )
        KERNELS['kappa'] = np.zeros( nATOMS-1, dtype=np.float32 )
        KERNELS['icvf']  = np.zeros( nATOMS-1, dtype=np.float32 )
        KERNELS['norms'] = np.zeros( (self.scheme.dwi_count, nATOMS-1) )

        progress = ProgressBar( n=nATOMS, prefix="   ", erase=False )

        # Coupled contributions
        for i in range( len(self.IC_ODs) ):
            for j in range( len(self.IC_VFs) ):
                lm = np.load( pjoin( in_path, 'A_%03d.npy'%progress.i ) )
                if lm.shape[0] != ndirs:
                    ERROR( 'Outdated LUT. Call "generate_kernels( regenerate=True )" to update the LUT' )
                idx = progress.i - 1
                KERNELS['wm'][idx,:,:] = amico.lut.resample_kernel( lm, self.scheme.nS, idx_out, Ylm_out, False, ndirs )[:,merge_idx]
                KERNELS['kappa'][idx] = 1.0 / np.tan( self.IC_ODs[i]*np.pi/2.0 )
                KERNELS['icvf'][idx]  = self.IC_VFs[j]
                if doMergeB0:
                    KERNELS['norms'][:,idx] = 1 / np.linalg.norm( KERNELS['wm'][idx,0,1:] ) # norm of coupled atoms (for l1 minimization)
                else:
                    KERNELS['norms'][:,idx] = 1 / np.linalg.norm( KERNELS['wm'][idx,0,self.scheme.dwi_idx] ) # norm of coupled atoms (for l1 minimization)
                progress.update()

        # Isotropic
        lm = np.load( pjoin( in_path, 'A_%03d.npy'%progress.i ) )
        KERNELS['iso'] = amico.lut.resample_kernel( lm, self.scheme.nS, idx_out, Ylm_out, True, ndirs )[merge_idx]
        progress.update()

        return KERNELS
Пример #18
0
    def save_results(self, path_suffix=None):
        """Save the output (directions, maps etc).

        Parameters
        ----------
        path_suffix : string
            Text to be appended to the output path (default : None)
        """
        if self.RESULTS is None:
            ERROR('Model not fitted to the data; call "fit()" first')
        if self.get_config('OUTPUT_path') is None:
            RESULTS_path = pjoin('AMICO', self.model.id)
            if path_suffix:
                RESULTS_path = RESULTS_path + '_' + path_suffix
            self.RESULTS['RESULTS_path'] = RESULTS_path
            LOG('\n-> Saving output to "%s/*":' % RESULTS_path)

            # delete previous output
            RESULTS_path = pjoin(self.get_config('DATA_path'), RESULTS_path)
        else:
            RESULTS_path = self.get_config('OUTPUT_path')
            if path_suffix:
                RESULTS_path = RESULTS_path + '_' + path_suffix
            self.RESULTS['RESULTS_path'] = RESULTS_path
            LOG('\n-> Saving output to "%s/*":' % RESULTS_path)

        if not exists(RESULTS_path):
            makedirs(RESULTS_path)
        else:
            for f in glob.glob(pjoin(RESULTS_path, '*')):
                remove(f)

        # configuration
        print('\t- configuration', end=' ')
        with open(pjoin(RESULTS_path, 'config.pickle'), 'wb+') as fid:
            pickle.dump(self.CONFIG, fid, protocol=2)
        print(' [OK]')

        # estimated orientations
        print('\t- FIT_dir.nii.gz', end=' ')
        niiMAP_img = self.RESULTS['DIRs']
        affine = self.niiDWI.affine if nibabel.__version__ >= '2.0.0' else self.niiDWI.get_affine(
        )
        niiMAP = nibabel.Nifti1Image(niiMAP_img, affine)
        niiMAP_hdr = niiMAP.header if nibabel.__version__ >= '2.0.0' else niiMAP.get_header(
        )
        niiMAP_hdr['cal_min'] = -1
        niiMAP_hdr['cal_max'] = 1
        niiMAP_hdr['scl_slope'] = 1
        niiMAP_hdr['scl_inter'] = 0
        nibabel.save(niiMAP, pjoin(RESULTS_path, 'FIT_dir.nii.gz'))
        print(' [OK]')

        # fitting error
        if self.get_config('doComputeNRMSE'):
            print('\t- FIT_nrmse.nii.gz', end=' ')
            niiMAP_img = self.RESULTS['NRMSE']
            niiMAP = nibabel.Nifti1Image(niiMAP_img, affine)
            niiMAP_hdr = niiMAP.header if nibabel.__version__ >= '2.0.0' else niiMAP.get_header(
            )
            niiMAP_hdr['cal_min'] = 0
            niiMAP_hdr['cal_max'] = 1
            niiMAP_hdr['scl_slope'] = 1
            niiMAP_hdr['scl_inter'] = 0
            nibabel.save(niiMAP, pjoin(RESULTS_path, 'FIT_nrmse.nii.gz'))
            print(' [OK]')

        if self.get_config('doSaveCorrectedDWI'):
            if self.model.name == 'Free-Water':
                print('\t- dwi_fw_corrected.nii.gz', end=' ')
                niiMAP_img = self.RESULTS['DWI_corrected']
                niiMAP = nibabel.Nifti1Image(niiMAP_img, affine)
                niiMAP_hdr = niiMAP.header if nibabel.__version__ >= '2.0.0' else niiMAP.get_header(
                )
                niiMAP_hdr['cal_min'] = 0
                niiMAP_hdr['cal_max'] = 1
                nibabel.save(niiMAP,
                             pjoin(RESULTS_path, 'dwi_fw_corrected.nii.gz'))
                print(' [OK]')
            else:
                WARNING(
                    '"doSaveCorrectedDWI" option not supported for "%s" model'
                    % self.model.name)

        # voxelwise maps
        for i in range(len(self.model.maps_name)):
            print('\t- FIT_%s.nii.gz' % self.model.maps_name[i], end=' ')
            niiMAP_img = self.RESULTS['MAPs'][:, :, :, i]
            niiMAP = nibabel.Nifti1Image(niiMAP_img, affine)
            niiMAP_hdr = niiMAP.header if nibabel.__version__ >= '2.0.0' else niiMAP.get_header(
            )
            niiMAP_hdr['descrip'] = self.model.maps_descr[i]
            niiMAP_hdr['cal_min'] = niiMAP_img.min()
            niiMAP_hdr['cal_max'] = niiMAP_img.max()
            niiMAP_hdr['scl_slope'] = 1
            niiMAP_hdr['scl_inter'] = 0
            nibabel.save(
                niiMAP,
                pjoin(RESULTS_path, 'FIT_%s.nii.gz' % self.model.maps_name[i]))
            print(' [OK]')

        LOG('   [ DONE ]')
Пример #19
0
    def fit(self):
        """Fit the model to the data iterating over all voxels (in the mask) one after the other.
        Call the appropriate fit() method of the actual model used.
        """
        if self.niiDWI is None:
            ERROR('Data not loaded; call "load_data()" first')
        if self.model is None:
            ERROR('Model not set; call "set_model()" first')
        if self.KERNELS is None:
            ERROR(
                'Response functions not generated; call "generate_kernels()" and "load_kernels()" first'
            )
        if self.KERNELS['model'] != self.model.id:
            ERROR('Response functions were not created with the same model')

        self.set_config('fit_time', None)
        totVoxels = np.count_nonzero(self.niiMASK_img)
        LOG('\n-> Fitting "%s" model to %d voxels:' %
            (self.model.name, totVoxels))

        # setup fitting directions
        peaks_filename = self.get_config('peaks_filename')
        if peaks_filename is None:
            DIRs = np.zeros([
                self.get_config('dim')[0],
                self.get_config('dim')[1],
                self.get_config('dim')[2], 3
            ],
                            dtype=np.float32)
            nDIR = 1
            if self.get_config('doMergeB0'):
                gtab = gradient_table(
                    np.hstack((0, self.scheme.b[self.scheme.dwi_idx])),
                    np.vstack((np.zeros(
                        (1, 3)), self.scheme.raw[self.scheme.dwi_idx, :3])))
            else:
                gtab = gradient_table(self.scheme.b, self.scheme.raw[:, :3])
            DTI = dti.TensorModel(gtab)
        else:
            niiPEAKS = nibabel.load(
                pjoin(self.get_config('DATA_path'), peaks_filename))
            DIRs = niiPEAKS.get_data().astype(np.float32)
            nDIR = np.floor(DIRs.shape[3] / 3)
            print('\t* peaks dim = %d x %d x %d x %d' % DIRs.shape[:4])
            if DIRs.shape[:3] != self.niiMASK_img.shape[:3]:
                ERROR('PEAKS geometry does not match with DWI data')

        # setup other output files
        MAPs = np.zeros([
            self.get_config('dim')[0],
            self.get_config('dim')[1],
            self.get_config('dim')[2],
            len(self.model.maps_name)
        ],
                        dtype=np.float32)

        if self.get_config('doComputeNRMSE'):
            NRMSE = np.zeros([
                self.get_config('dim')[0],
                self.get_config('dim')[1],
                self.get_config('dim')[2]
            ],
                             dtype=np.float32)

        if self.get_config('doSaveCorrectedDWI'):
            DWI_corrected = np.zeros(self.niiDWI.shape, dtype=np.float32)

        # fit the model to the data
        # =========================
        t = time.time()
        progress = ProgressBar(n=totVoxels, prefix="   ", erase=False)
        for iz in range(self.niiMASK_img.shape[2]):
            for iy in range(self.niiMASK_img.shape[1]):
                for ix in range(self.niiMASK_img.shape[0]):
                    if self.niiMASK_img[ix, iy, iz] == 0:
                        continue

                    # prepare the signal
                    y = self.niiDWI_img[ix, iy, iz, :].astype(np.float64)
                    y[y < 0] = 0  # [NOTE] this should not happen!

                    # fitting directions
                    if peaks_filename is None:
                        dirs = DTI.fit(y).directions[0]
                    else:
                        dirs = DIRs[ix, iy, iz, :]

                    # dispatch to the right handler for each model
                    MAPs[ix, iy,
                         iz, :], DIRs[ix, iy, iz, :], x, A = self.model.fit(
                             y, dirs.reshape(-1, 3), self.KERNELS,
                             self.get_config('solver_params'), self.htable)

                    # compute fitting error
                    if self.get_config('doComputeNRMSE'):
                        y_est = np.dot(A, x)
                        den = np.sum(y**2)
                        NRMSE[ix, iy,
                              iz] = np.sqrt(np.sum(
                                  (y - y_est)**2) / den) if den > 1e-16 else 0

                    if self.get_config('doSaveCorrectedDWI'):

                        if self.model.name == 'Free-Water':
                            n_iso = len(self.model.d_isos)

                            # keep only FW components of the estimate
                            x[0:x.shape[0] - n_iso] = 0

                            # y_fw_corrected below is the predicted signal by the anisotropic part (no iso part)
                            y_fw_part = np.dot(A, x)

                            # y is the original signal
                            y_fw_corrected = y - y_fw_part
                            y_fw_corrected[
                                y_fw_corrected <
                                0] = 0  # [NOTE] this should not happen!

                            if self.get_config('doNormalizeSignal'
                                               ) and self.scheme.b0_count > 0:
                                y_fw_corrected = y_fw_corrected * self.mean_b0s[
                                    ix, iy, iz]

                            if self.get_config('doKeepb0Intact'
                                               ) and self.scheme.b0_count > 0:
                                # put original b0 data back in.
                                y_fw_corrected[self.scheme.b0_idx] = y[
                                    self.scheme.b0_idx] * self.mean_b0s[ix, iy,
                                                                        iz]

                            DWI_corrected[ix, iy, iz, :] = y_fw_corrected

                    progress.update()

        self.set_config('fit_time', time.time() - t)
        LOG('   [ %s ]' % (time.strftime(
            "%Hh %Mm %Ss", time.gmtime(self.get_config('fit_time')))))

        # store results
        self.RESULTS = {}
        self.RESULTS['DIRs'] = DIRs
        self.RESULTS['MAPs'] = MAPs
        if self.get_config('doComputeNRMSE'):
            self.RESULTS['NRMSE'] = NRMSE
        if self.get_config('doSaveCorrectedDWI'):
            self.RESULTS['DWI_corrected'] = DWI_corrected
Пример #20
0
 def set_solver( self ) :
     ERROR( 'Not implemented' )
Пример #21
0
    def fit(self):
        """Fit the model to the data iterating over all voxels (in the mask) one after the other.
        Call the appropriate fit() method of the actual model used.
        """
        @delayed
        @wrap_non_picklable_objects
        def fit_voxel(self, ix, iy, iz, dirs, DTI):
            """Perform the fit in a single voxel.
            """
            # prepare the signal
            y = self.niiDWI_img[ix, iy, iz, :].astype(np.float64)
            y[y < 0] = 0  # [NOTE] this should not happen!

            # fitting directions if not
            if not self.get_config('doDirectionalAverage') and DTI is not None:
                dirs = DTI.fit(y).directions[0].reshape(-1, 3)

            # dispatch to the right handler for each model
            results, dirs, x, A = self.model.fit(
                y, dirs, self.KERNELS, self.get_config('solver_params'),
                self.htable)

            # compute fitting error
            if self.get_config('doComputeNRMSE'):
                y_est = np.dot(A, x)
                den = np.sum(y**2)
                NRMSE = np.sqrt(np.sum(
                    (y - y_est)**2) / den) if den > 1e-16 else 0
            else:
                NRMSE = 0.0

            y_fw_corrected = None
            if self.get_config('doSaveCorrectedDWI'):

                if self.model.name == 'Free-Water':
                    n_iso = len(self.model.d_isos)

                    # keep only FW components of the estimate
                    x[0:x.shape[0] - n_iso] = 0

                    # y_fw_corrected below is the predicted signal by the anisotropic part (no iso part)
                    y_fw_part = np.dot(A, x)

                    # y is the original signal
                    y_fw_corrected = y - y_fw_part
                    y_fw_corrected[y_fw_corrected <
                                   0] = 0  # [NOTE] this should not happen!

                    if self.get_config(
                            'doNormalizeSignal') and self.scheme.b0_count > 0:
                        y_fw_corrected = y_fw_corrected * self.mean_b0s[ix, iy,
                                                                        iz]

                    if self.get_config(
                            'doKeepb0Intact') and self.scheme.b0_count > 0:
                        # put original b0 data back in.
                        y_fw_corrected[self.scheme.b0_idx] = y[
                            self.scheme.b0_idx] * self.mean_b0s[ix, iy, iz]

            return results, dirs, NRMSE, y_fw_corrected

        if self.niiDWI is None:
            ERROR('Data not loaded; call "load_data()" first')
        if self.model is None:
            ERROR('Model not set; call "set_model()" first')
        if self.KERNELS is None:
            ERROR(
                'Response functions not generated; call "generate_kernels()" and "load_kernels()" first'
            )
        if self.KERNELS['model'] != self.model.id:
            ERROR('Response functions were not created with the same model')
        n_jobs = self.get_config('parallel_jobs')
        if n_jobs == -1:
            n_jobs = cpu_count()
        elif n_jobs == 0 or n_jobs < -1:
            ERROR('Number of parallel jobs must be positive or -1')

        self.set_config('fit_time', None)
        totVoxels = np.count_nonzero(self.niiMASK_img)
        LOG(f'\n-> Fitting "{self.model.name}" model to {totVoxels} voxels (using {n_jobs} job{"s" if n_jobs>1 else ""}):'
            )

        # setup fitting directions
        peaks_filename = self.get_config('peaks_filename')
        if peaks_filename is None:
            DIRs = np.zeros([
                self.get_config('dim')[0],
                self.get_config('dim')[1],
                self.get_config('dim')[2], 3
            ],
                            dtype=np.float32)
            nDIR = 1
            if self.get_config('doMergeB0'):
                gtab = gradient_table(
                    np.hstack((0, self.scheme.b[self.scheme.dwi_idx])),
                    np.vstack((np.zeros(
                        (1, 3)), self.scheme.raw[self.scheme.dwi_idx, :3])))
            else:
                gtab = gradient_table(self.scheme.b, self.scheme.raw[:, :3])
            DTI = dti.TensorModel(gtab)
        else:
            if not isfile(pjoin(self.get_config('DATA_path'), peaks_filename)):
                ERROR('PEAKS file not found')
            niiPEAKS = nibabel.load(
                pjoin(self.get_config('DATA_path'), peaks_filename))
            DIRs = niiPEAKS.get_data().astype(np.float32)
            nDIR = np.floor(DIRs.shape[3] / 3)
            PRINT('\t* peaks dim = %d x %d x %d x %d' % DIRs.shape[:4])
            if DIRs.shape[:3] != self.niiMASK_img.shape[:3]:
                ERROR('PEAKS geometry does not match with DWI data')
            DTI = None

        # setup other output files
        MAPs = np.zeros([
            self.get_config('dim')[0],
            self.get_config('dim')[1],
            self.get_config('dim')[2],
            len(self.model.maps_name)
        ],
                        dtype=np.float32)

        # fit the model to the data
        # =========================
        t = time.time()

        ix, iy, iz = np.nonzero(self.niiMASK_img)
        n_per_thread = np.floor(totVoxels / n_jobs)
        idx = np.arange(0, totVoxels + 1, n_per_thread, dtype=np.int32)
        idx[-1] = totVoxels

        estimates = Parallel(n_jobs=n_jobs)(
            fit_voxel(self, ix[i], iy[i], iz[i], DIRs[ix[i], iy[i],
                                                      iz[i], :], DTI)
            for i in tqdm(range(totVoxels),
                          ncols=70,
                          bar_format='   |{bar}| {percentage:4.1f}%',
                          disable=(get_verbose() < 3)))

        self.set_config('fit_time', time.time() - t)
        LOG('   [ %s ]' % (time.strftime(
            "%Hh %Mm %Ss", time.gmtime(self.get_config('fit_time')))))

        # store results
        self.RESULTS = {}

        for i in range(totVoxels):
            MAPs[ix[i], iy[i], iz[i], :] = estimates[i][0]
            DIRs[ix[i], iy[i], iz[i], :] = estimates[i][1]
        self.RESULTS['DIRs'] = DIRs
        self.RESULTS['MAPs'] = MAPs

        if self.get_config('doComputeNRMSE'):
            NRMSE = np.zeros([
                self.get_config('dim')[0],
                self.get_config('dim')[1],
                self.get_config('dim')[2]
            ],
                             dtype=np.float32)
            for i in range(totVoxels):
                NRMSE[ix[i], iy[i], iz[i]] = estimates[i][2]
            self.RESULTS['NRMSE'] = NRMSE

        if self.get_config('doSaveCorrectedDWI'):
            DWI_corrected = np.zeros(self.niiDWI.shape, dtype=np.float32)
            for i in range(totVoxels):
                DWI_corrected[ix[i], iy[i], iz[i], :] = estimates[i][3]
            self.RESULTS['DWI_corrected'] = DWI_corrected
Пример #22
0
    def load_data(self,
                  dwi_filename='DWI.nii',
                  scheme_filename='DWI.scheme',
                  mask_filename=None,
                  b0_thr=0):
        """Load the diffusion signal and its corresponding acquisition scheme.

        Parameters
        ----------
        dwi_filename : string
            The file name of the DWI data, relative to the subject folder (default : 'DWI.nii')
        scheme_filename : string
            The file name of the corresponding acquisition scheme (default : 'DWI.scheme')
        mask_filename : string
            The file name of the (optional) binary mask (default : None)
        b0_thr : float
            The threshold below which a b-value is considered a b0 (default : 0)
        """

        # Loading data, acquisition scheme and mask (optional)
        LOG('\n-> Loading data:')
        tic = time.time()

        PRINT('\t* DWI signal')
        if not isfile(pjoin(self.get_config('DATA_path'), dwi_filename)):
            ERROR('DWI file not found')
        self.set_config('dwi_filename', dwi_filename)
        self.niiDWI = nibabel.load(
            pjoin(self.get_config('DATA_path'), dwi_filename))
        self.niiDWI_img = self.niiDWI.get_data().astype(np.float32)
        hdr = self.niiDWI.header if nibabel.__version__ >= '2.0.0' else self.niiDWI.get_header(
        )
        self.set_config('dim', self.niiDWI_img.shape[:3])
        self.set_config('pixdim', tuple(hdr.get_zooms()[:3]))
        PRINT('\t\t- dim    = %d x %d x %d x %d' % self.niiDWI_img.shape)
        PRINT('\t\t- pixdim = %.3f x %.3f x %.3f' % self.get_config('pixdim'))
        # Scale signal intensities (if necessary)
        if (np.isfinite(hdr['scl_slope']) and np.isfinite(hdr['scl_inter'])
                and hdr['scl_slope'] != 0
                and (hdr['scl_slope'] != 1 or hdr['scl_inter'] != 0)):
            PRINT('\t\t- rescaling data ', end='')
            self.niiDWI_img = self.niiDWI_img * hdr['scl_slope'] + hdr[
                'scl_inter']
            PRINT('[OK]')

        PRINT('\t* Acquisition scheme')
        if not isfile(pjoin(self.get_config('DATA_path'), scheme_filename)):
            ERROR('SCHEME file not found')
        self.set_config('scheme_filename', scheme_filename)
        self.set_config('b0_thr', b0_thr)
        self.scheme = amico.scheme.Scheme(
            pjoin(self.get_config('DATA_path'), scheme_filename), b0_thr)
        PRINT(
            f'\t\t- {self.scheme.nS} samples, {len(self.scheme.shells)} shells'
        )
        PRINT(f'\t\t- {self.scheme.b0_count} @ b=0', end=' ')
        for i in range(len(self.scheme.shells)):
            PRINT(
                f', {len(self.scheme.shells[i]["idx"])} @ b={self.scheme.shells[i]["b"]:.1f}',
                end=' ')
        PRINT()

        if self.scheme.nS != self.niiDWI_img.shape[3]:
            ERROR('Scheme does not match with DWI data')

        PRINT('\t* Binary mask')
        if mask_filename is not None:
            if not isfile(pjoin(self.get_config('DATA_path'), mask_filename)):
                ERROR('MASK file not found')
            self.niiMASK = nibabel.load(
                pjoin(self.get_config('DATA_path'), mask_filename))
            self.niiMASK_img = self.niiMASK.get_data().astype(np.uint8)
            niiMASK_hdr = self.niiMASK.header if nibabel.__version__ >= '2.0.0' else self.niiMASK.get_header(
            )
            PRINT('\t\t- dim    = %d x %d x %d' % self.niiMASK_img.shape[:3])
            PRINT('\t\t- pixdim = %.3f x %.3f x %.3f' %
                  niiMASK_hdr.get_zooms()[:3])
            if self.niiMASK.ndim != 3:
                ERROR('The provided MASK if 4D, but a 3D dataset is expected')
            if self.get_config('dim') != self.niiMASK_img.shape[:3]:
                ERROR('MASK geometry does not match with DWI data')
        else:
            self.niiMASK = None
            self.niiMASK_img = np.ones(self.get_config('dim'))
            PRINT('\t\t- not specified')
        PRINT(f'\t\t- voxels = {np.count_nonzero(self.niiMASK_img)}')

        LOG(f'   [ {time.time() - tic:.1f} seconds ]')

        # Preprocessing
        LOG('\n-> Preprocessing:')
        tic = time.time()

        if self.get_config('doDebiasSignal'):
            PRINT('\t* Debiasing signal... ', end='')
            sys.stdout.flush()
            if self.get_config('DWI-SNR') == None:
                ERROR(
                    "Set noise variance for debiasing (eg. ae.set_config('RicianNoiseSigma', sigma))"
                )
            self.niiDWI_img = debiasRician(self.niiDWI_img,
                                           self.get_config('DWI-SNR'),
                                           self.niiMASK_img, self.scheme)
            PRINT(' [OK]')

        if self.get_config('doNormalizeSignal'):
            PRINT('\t* Normalizing to b0... ', end='')
            sys.stdout.flush()
            if self.scheme.b0_count > 0:
                self.mean_b0s = np.mean(self.niiDWI_img[:, :, :,
                                                        self.scheme.b0_idx],
                                        axis=3)
            else:
                ERROR('No b0 volume to normalize signal with')
            norm_factor = self.mean_b0s.copy()
            idx = self.mean_b0s <= 0
            norm_factor[idx] = 1
            norm_factor = 1 / norm_factor
            norm_factor[idx] = 0
            for i in range(self.scheme.nS):
                self.niiDWI_img[:, :, :, i] *= norm_factor
            PRINT(
                f'[ min={self.niiDWI_img.min():.2f},  mean={self.niiDWI_img.mean():.2f}, max={self.niiDWI_img.max():.2f} ]'
            )

        if self.get_config('doMergeB0'):
            PRINT('\t* Merging multiple b0 volume(s)')
            mean = np.expand_dims(np.mean(self.niiDWI_img[:, :, :,
                                                          self.scheme.b0_idx],
                                          axis=3),
                                  axis=3)
            self.niiDWI_img = np.concatenate(
                (mean, self.niiDWI_img[:, :, :, self.scheme.dwi_idx]), axis=3)
        else:
            PRINT('\t* Keeping all b0 volume(s)')

        if self.get_config('doDirectionalAverage'):
            PRINT(
                '\t* Performing the directional average on the signal of each shell... '
            )
            numShells = len(self.scheme.shells)
            dir_avg_img = self.niiDWI_img[:, :, :, :(numShells + 1)]
            scheme_table = np.zeros([numShells + 1, 7])

            id_bval = 0
            dir_avg_img[:, :, :,
                        id_bval] = np.mean(self.niiDWI_img[:, :, :,
                                                           self.scheme.b0_idx],
                                           axis=3)
            scheme_table[id_bval, :] = np.array([1, 0, 0, 0, 0, 0, 0])

            bvals = []
            for shell in self.scheme.shells:
                bvals.append(shell['b'])

            sort_idx = np.argsort(bvals)

            for shell_idx in sort_idx:
                shell = self.scheme.shells[shell_idx]
                id_bval = id_bval + 1
                dir_avg_img[:, :, :,
                            id_bval] = np.mean(self.niiDWI_img[:, :, :,
                                                               shell['idx']],
                                               axis=3)
                scheme_table[id_bval, :] = np.array([
                    1, 0, 0, shell['G'], shell['Delta'], shell['delta'],
                    shell['TE']
                ])

            self.niiDWI_img = dir_avg_img.astype(np.float32)
            self.set_config('dim', self.niiDWI_img.shape[:3])
            PRINT('\t\t- dim    = %d x %d x %d x %d' % self.niiDWI_img.shape)
            PRINT('\t\t- pixdim = %.3f x %.3f x %.3f' %
                  self.get_config('pixdim'))

            PRINT('\t* Acquisition scheme')
            self.scheme = amico.scheme.Scheme(scheme_table, b0_thr)
            PRINT(
                f'\t\t- {self.scheme.nS} samples, {len(self.scheme.shells)} shells'
            )
            PRINT(f'\t\t- {self.scheme.b0_count} @ b=0', end=' ')
            for i in range(len(self.scheme.shells)):
                PRINT(
                    f', {len(self.scheme.shells[i]["idx"])} @ b={self.scheme.shells[i]["b"]:.1f}',
                    end=' ')
            PRINT()

            if self.scheme.nS != self.niiDWI_img.shape[3]:
                ERROR('Scheme does not match with DWI data')

        LOG(f'   [ {time.time() - tic:.1f} seconds ]')
Пример #23
0
    def load_data(self,
                  dwi_filename='DWI.nii',
                  scheme_filename='DWI.scheme',
                  mask_filename=None,
                  b0_thr=0):
        """Load the diffusion signal and its corresponding acquisition scheme.

        Parameters
        ----------
        dwi_filename : string
            The file name of the DWI data, relative to the subject folder (default : 'DWI.nii')
        scheme_filename : string
            The file name of the corresponding acquisition scheme (default : 'DWI.scheme')
        mask_filename : string
            The file name of the (optional) binary mask (default : None)
        b0_thr : float
            The threshold below which a b-value is considered a b0 (default : 0)
        """

        # Loading data, acquisition scheme and mask (optional)
        LOG('\n-> Loading data:')
        tic = time.time()

        print('\t* DWI signal')
        self.set_config('dwi_filename', dwi_filename)
        self.niiDWI = nibabel.load(
            pjoin(self.get_config('DATA_path'), dwi_filename))
        self.niiDWI_img = self.niiDWI.get_data().astype(np.float32)
        hdr = self.niiDWI.header if nibabel.__version__ >= '2.0.0' else self.niiDWI.get_header(
        )
        self.set_config('dim', self.niiDWI_img.shape[:3])
        self.set_config('pixdim', tuple(hdr.get_zooms()[:3]))
        print('\t\t- dim    = %d x %d x %d x %d' % self.niiDWI_img.shape)
        print('\t\t- pixdim = %.3f x %.3f x %.3f' % self.get_config('pixdim'))
        # Scale signal intensities (if necessary)
        if (np.isfinite(hdr['scl_slope']) and np.isfinite(hdr['scl_inter'])
                and hdr['scl_slope'] != 0
                and (hdr['scl_slope'] != 1 or hdr['scl_inter'] != 0)):
            print('\t\t- rescaling data ', end='')
            self.niiDWI_img = self.niiDWI_img * hdr['scl_slope'] + hdr[
                'scl_inter']
            print('[OK]')

        print('\t* Acquisition scheme')
        self.set_config('scheme_filename', scheme_filename)
        self.set_config('b0_thr', b0_thr)
        self.scheme = amico.scheme.Scheme(
            pjoin(self.get_config('DATA_path'), scheme_filename), b0_thr)
        print('\t\t- %d samples, %d shells' %
              (self.scheme.nS, len(self.scheme.shells)))
        print('\t\t- %d @ b=0' % (self.scheme.b0_count), end=' ')
        for i in range(len(self.scheme.shells)):
            print(', %d @ b=%.1f' % (len(
                self.scheme.shells[i]['idx']), self.scheme.shells[i]['b']),
                  end=' ')
        print()

        if self.scheme.nS != self.niiDWI_img.shape[3]:
            ERROR('Scheme does not match with DWI data')

        print('\t* Binary mask')
        if mask_filename is not None:
            self.niiMASK = nibabel.load(
                pjoin(self.get_config('DATA_path'), mask_filename))
            self.niiMASK_img = self.niiMASK.get_data().astype(np.uint8)
            niiMASK_hdr = self.niiMASK.header if nibabel.__version__ >= '2.0.0' else self.niiMASK.get_header(
            )
            print('\t\t- dim    = %d x %d x %d' % self.niiMASK_img.shape[:3])
            print('\t\t- pixdim = %.3f x %.3f x %.3f' %
                  niiMASK_hdr.get_zooms()[:3])
            if self.niiMASK.ndim != 3:
                ERROR('The provided MASK if 4D, but a 3D dataset is expected')
            if self.get_config('dim') != self.niiMASK_img.shape[:3]:
                ERROR('MASK geometry does not match with DWI data')
        else:
            self.niiMASK = None
            self.niiMASK_img = np.ones(self.get_config('dim'))
            print('\t\t- not specified')
        print('\t\t- voxels = %d' % np.count_nonzero(self.niiMASK_img))

        LOG('   [ %.1f seconds ]' % (time.time() - tic))

        # Preprocessing
        LOG('\n-> Preprocessing:')
        tic = time.time()

        if self.get_config('doDebiasSignal'):
            print('\t* Debiasing signal... ', end='')
            sys.stdout.flush()
            if self.get_config('DWI-SNR') == None:
                ERROR(
                    "Set noise variance for debiasing (eg. ae.set_config('RicianNoiseSigma', sigma))"
                )
            self.niiDWI_img = debiasRician(self.niiDWI_img,
                                           self.get_config('DWI-SNR'),
                                           self.niiMASK_img, self.scheme)
            print(' [OK]')

        if self.get_config('doNormalizeSignal'):
            print('\t* Normalizing to b0... ', end='')
            sys.stdout.flush()
            if self.scheme.b0_count > 0:
                self.mean_b0s = np.mean(self.niiDWI_img[:, :, :,
                                                        self.scheme.b0_idx],
                                        axis=3)
            else:
                ERROR('No b0 volume to normalize signal with')
            norm_factor = self.mean_b0s.copy()
            idx = self.mean_b0s <= 0
            norm_factor[idx] = 1
            norm_factor = 1 / norm_factor
            norm_factor[idx] = 0
            for i in range(self.scheme.nS):
                self.niiDWI_img[:, :, :, i] *= norm_factor
            print('[ min=%.2f,  mean=%.2f, max=%.2f ]' %
                  (self.niiDWI_img.min(), self.niiDWI_img.mean(),
                   self.niiDWI_img.max()))

        if self.get_config('doMergeB0'):
            print('\t* Merging multiple b0 volume(s)')
            mean = np.expand_dims(np.mean(self.niiDWI_img[:, :, :,
                                                          self.scheme.b0_idx],
                                          axis=3),
                                  axis=3)
            self.niiDWI_img = np.concatenate(
                (mean, self.niiDWI_img[:, :, :, self.scheme.dwi_idx]), axis=3)
        else:
            print('\t* Keeping all b0 volume(s)')

        LOG('   [ %.1f seconds ]' % (time.time() - tic))
Пример #24
0
    def save_results(self, path_suffix=None, save_dir_avg=False):
        """Save the output (directions, maps etc).

        Parameters
        ----------
        path_suffix : string
            Text to be appended to the output path (default : None)
        save_dir_avg : boolean
            If true and the option doDirectionalAverage is true 
            the directional average signal and the scheme 
            will be saved in files (default : False)
        """
        if self.RESULTS is None:
            ERROR('Model not fitted to the data; call "fit()" first')
        if self.get_config('OUTPUT_path') is None:
            RESULTS_path = pjoin('AMICO', self.model.id)
            if path_suffix:
                RESULTS_path = RESULTS_path + '_' + path_suffix
            self.RESULTS['RESULTS_path'] = RESULTS_path
            LOG('\n-> Saving output to "%s/*":' % RESULTS_path)

            # delete previous output
            RESULTS_path = pjoin(self.get_config('DATA_path'), RESULTS_path)
        else:
            RESULTS_path = self.get_config('OUTPUT_path')
            if path_suffix:
                RESULTS_path = RESULTS_path + '_' + path_suffix
            self.RESULTS['RESULTS_path'] = RESULTS_path
            LOG('\n-> Saving output to "%s/*":' % RESULTS_path)

        if not exists(RESULTS_path):
            makedirs(RESULTS_path)
        else:
            for f in glob.glob(pjoin(RESULTS_path, '*')):
                remove(f)

        # configuration
        print('\t- configuration', end=' ')
        with open(pjoin(RESULTS_path, 'config.pickle'), 'wb+') as fid:
            pickle.dump(self.CONFIG, fid, protocol=2)
        print(' [OK]')

        affine = self.niiDWI.affine if nibabel.__version__ >= '2.0.0' else self.niiDWI.get_affine(
        )
        hdr = self.niiDWI.header if nibabel.__version__ >= '2.0.0' else self.niiDWI.get_header(
        )
        hdr['datatype'] = 16
        hdr['bitpix'] = 32

        # estimated orientations
        if not self.get_config('doDirectionalAverage'):
            print('\t- FIT_dir.nii.gz', end=' ')
            niiMAP_img = self.RESULTS['DIRs']
            niiMAP = nibabel.Nifti1Image(niiMAP_img, affine, hdr)
            niiMAP_hdr = niiMAP.header if nibabel.__version__ >= '2.0.0' else niiMAP.get_header(
            )
            niiMAP_hdr['cal_min'] = -1
            niiMAP_hdr['cal_max'] = 1
            niiMAP_hdr['scl_slope'] = 1
            niiMAP_hdr['scl_inter'] = 0
            nibabel.save(niiMAP, pjoin(RESULTS_path, 'FIT_dir.nii.gz'))
            print(' [OK]')

        # fitting error
        if self.get_config('doComputeNRMSE'):
            print('\t- FIT_nrmse.nii.gz', end=' ')
            niiMAP_img = self.RESULTS['NRMSE']
            niiMAP = nibabel.Nifti1Image(niiMAP_img, affine, hdr)
            niiMAP_hdr = niiMAP.header if nibabel.__version__ >= '2.0.0' else niiMAP.get_header(
            )
            niiMAP_hdr['cal_min'] = 0
            niiMAP_hdr['cal_max'] = 1
            niiMAP_hdr['scl_slope'] = 1
            niiMAP_hdr['scl_inter'] = 0
            nibabel.save(niiMAP, pjoin(RESULTS_path, 'FIT_nrmse.nii.gz'))
            print(' [OK]')

        if self.get_config('doSaveCorrectedDWI'):
            if self.model.name == 'Free-Water':
                print('\t- dwi_fw_corrected.nii.gz', end=' ')
                niiMAP_img = self.RESULTS['DWI_corrected']
                niiMAP = nibabel.Nifti1Image(niiMAP_img, affine, hdr)
                niiMAP_hdr = niiMAP.header if nibabel.__version__ >= '2.0.0' else niiMAP.get_header(
                )
                niiMAP_hdr['cal_min'] = 0
                niiMAP_hdr['cal_max'] = 1
                nibabel.save(niiMAP,
                             pjoin(RESULTS_path, 'dwi_fw_corrected.nii.gz'))
                print(' [OK]')
            else:
                WARNING(
                    '"doSaveCorrectedDWI" option not supported for "%s" model'
                    % self.model.name)

        # voxelwise maps
        for i in range(len(self.model.maps_name)):
            print('\t- FIT_%s.nii.gz' % self.model.maps_name[i], end=' ')
            niiMAP_img = self.RESULTS['MAPs'][:, :, :, i]
            niiMAP = nibabel.Nifti1Image(niiMAP_img, affine, hdr)
            niiMAP_hdr = niiMAP.header if nibabel.__version__ >= '2.0.0' else niiMAP.get_header(
            )
            niiMAP_hdr['descrip'] = self.model.maps_descr[
                i] + ' (AMICO v%s)' % self.get_config('version')
            niiMAP_hdr['cal_min'] = niiMAP_img.min()
            niiMAP_hdr['cal_max'] = niiMAP_img.max()
            niiMAP_hdr['scl_slope'] = 1
            niiMAP_hdr['scl_inter'] = 0
            nibabel.save(
                niiMAP,
                pjoin(RESULTS_path, 'FIT_%s.nii.gz' % self.model.maps_name[i]))
            print(' [OK]')

        # Directional average signal
        if save_dir_avg:
            if self.get_config('doDirectionalAverage'):
                print('\t- dir_avg_signal.nii.gz', end=' ')
                niiMAP = nibabel.Nifti1Image(self.niiDWI_img, affine, hdr)
                niiMAP_hdr = niiMAP.header if nibabel.__version__ >= '2.0.0' else niiMAP.get_header(
                )
                niiMAP_hdr[
                    'descrip'] = 'Directional average signal of each shell' + ' (AMICO v%s)' % self.get_config(
                        'version')
                nibabel.save(niiMAP,
                             pjoin(RESULTS_path, 'dir_avg_signal.nii.gz'))
                print(' [OK]')

                print('\t- dir_avg.scheme', end=' ')
                np.savetxt(pjoin(RESULTS_path, 'dir_avg.scheme'),
                           self.scheme.get_table(),
                           fmt="%.06f",
                           delimiter="\t",
                           header="VERSION: {}".format(self.scheme.version),
                           comments='')
                print(' [OK]')
            else:
                WARNING(
                    'The directional average signal was not created (The option doDirectionalAverage is False).'
                )

        LOG('   [ DONE ]')
Пример #25
0
 def fit( self, y, dirs, KERNELS, params ) :
     ERROR( 'Not implemented' )
Пример #26
0
    def load_from_table( self, data, b0_thr = 0 ) :
        """Build the structure from an input matrix.

        The first three columns represent the gradient directions.
        Then, we accept two formats to describe each gradient:
            - if the shape of data is Nx4, the 4^ column is the b-value;
            - if the shape of data is Nx7, the last 4 columns are, respectively, the gradient strength, big delta, small delta and TE.

        Parameters
        ----------
        data : numpy.ndarray
            Matrix containing tall the values.
        b0_thr : float
            The threshold on the b-values to identify the b0 images (default: 0)
        """
        if data.ndim == 1 :
            data = np.expand_dims( data, axis=0 )
        self.raw = data


        # number of samples
        # self.nS = self.raw.shape[0] JL: incomplete getter/setter incompatible with 3.6; this is never used any as getter always returns derived value

        # set/calculate the b-values
        if self.raw.shape[1] == 4 :
            self.version = 0
            self.b = self.raw[:,3]
        elif self.raw.shape[1] == 7 :
            self.version = 1
            self.b = ( 267.513e6 * self.raw[:,3] * self.raw[:,5] )**2 * (self.raw[:,4] - self.raw[:,5]/3.0) * 1e-6 # in mm^2/s
        else :
            ERROR( 'Unrecognized scheme format' )

        # store information about the volumes
        self.b0_thr    = b0_thr
        self.b0_idx    = np.where( self.b <= b0_thr )[0]
        self.b0_count  = len( self.b0_idx )
        self.dwi_idx   = np.where( self.b > b0_thr )[0]
        self.dwi_count = len( self.dwi_idx )

        # ensure the directions are in the spherical range [0,180]x[0,180]
        idx = np.where( self.raw[:,1] < 0 )[0]
        self.raw[idx,0:3] = -self.raw[idx,0:3]

        # store information about each shell in a dictionary
        self.shells = []

        tmp = np.ascontiguousarray( self.raw[:,3:] )
        schemeUnique, schemeUniqueInd = np.unique( tmp.view([('', tmp.dtype)]*tmp.shape[1]), return_index=True )
        schemeUnique = schemeUnique.view(tmp.dtype).reshape((schemeUnique.shape[0], tmp.shape[1]))
        schemeUnique = [tmp[index] for index in sorted(schemeUniqueInd)]
        bUnique = [self.b[index] for index in sorted(schemeUniqueInd)]
        for i in range(len(schemeUnique)) :
            if bUnique[i] <= b0_thr :
                continue
            shell = {}
            shell['b'] = bUnique[i]
            if self.version == 0 :
                shell['G']     = None
                shell['Delta'] = None
                shell['delta'] = None
                shell['TE']    = None
            else :
                shell['G']     = schemeUnique[i][0]
                shell['Delta'] = schemeUnique[i][1]
                shell['delta'] = schemeUnique[i][2]
                shell['TE']    = schemeUnique[i][3]

            shell['idx']  = np.where((tmp == schemeUnique[i]).all(axis=1))[0]
            shell['grad'] = self.raw[shell['idx'],0:3]
            self.shells.append( shell )
Пример #27
0
    def watson_SH_coeff( self, kappa ):

        if isinstance(kappa,np.ndarray):
            ERROR( '"watson_SH_coeff()" not implemented for multiple kappa input yet' )

        # In the scope of AMICO only a single value is used for kappa
        n = 6

        C = np.zeros((n+1))
        # 0th order is a constant
        C[0] = 2*np.sqrt(np.pi)

        # Precompute the special function values
        sk = np.sqrt(kappa)
        sk2 = sk*kappa
        sk3 = sk2*kappa
        sk4 = sk3*kappa
        sk5 = sk4*kappa
        sk6 = sk5*kappa
        sk7 = sk6*kappa
        k2 = np.power(kappa,2)
        k3 = k2*kappa
        k4 = k3*kappa
        k5 = k4*kappa
        k6 = k5*kappa
        k7 = k6*kappa

        erfik = scipy.special.erfi(sk)
        ierfik = 1/erfik
        ek = np.exp(kappa)
        dawsonk = 0.5*np.sqrt(np.pi)*erfik/ek

        if kappa > 0.1:

            # for large enough kappa
            C[1] = 3*sk - (3 + 2*kappa)*dawsonk
            C[1] = np.sqrt(5)*C[1]*ek
            C[1] = C[1]*ierfik/kappa

            C[2] = (105 + 60*kappa + 12*k2)*dawsonk
            C[2] = C[2] -105*sk + 10*sk2
            C[2] = .375*C[2]*ek/k2
            C[2] = C[2]*ierfik

            C[3] = -3465 - 1890*kappa - 420*k2 - 40*k3
            C[3] = C[3]*dawsonk
            C[3] = C[3] + 3465*sk - 420*sk2 + 84*sk3
            C[3] = C[3]*np.sqrt(13*np.pi)/64/k3
            C[3] = C[3]/dawsonk

            C[4] = 675675 + 360360*kappa + 83160*k2 + 10080*k3 + 560*k4
            C[4] = C[4]*dawsonk
            C[4] = C[4] - 675675*sk + 90090*sk2 - 23100*sk3 + 744*sk4
            C[4] = np.sqrt(17)*C[4]*ek
            C[4] = C[4]/512/k4
            C[4] = C[4]*ierfik

            C[5] = -43648605 - 22972950*kappa - 5405400*k2 - 720720*k3 - 55440*k4 - 2016*k5
            C[5] = C[5]*dawsonk
            C[5] = C[5] + 43648605*sk - 6126120*sk2 + 1729728*sk3 - 82368*sk4 + 5104*sk5
            C[5] = np.sqrt(21*np.pi)*C[5]/4096/k5
            C[5] = C[5]/dawsonk

            C[6] = 7027425405 + 3666482820*kappa + 872972100*k2 + 122522400*k3  + 10810800*k4 + 576576*k5 + 14784*k6
            C[6] = C[6]*dawsonk
            C[6] = C[6] - 7027425405*sk + 1018467450*sk2 - 302630328*sk3 + 17153136*sk4 - 1553552*sk5 + 25376*sk6
            C[6] = 5*C[6]*ek
            C[6] = C[6]/16384/k6
            C[6] = C[6]*ierfik

        # for very large kappa
        if kappa>30:
            lnkd = np.log(kappa) - np.log(30)
            lnkd2 = lnkd*lnkd
            lnkd3 = lnkd2*lnkd
            lnkd4 = lnkd3*lnkd
            lnkd5 = lnkd4*lnkd
            lnkd6 = lnkd5*lnkd
            C[1] = 7.52308 + 0.411538*lnkd - 0.214588*lnkd2 + 0.0784091*lnkd3 - 0.023981*lnkd4 + 0.00731537*lnkd5 - 0.0026467*lnkd6
            C[2] = 8.93718 + 1.62147*lnkd - 0.733421*lnkd2 + 0.191568*lnkd3 - 0.0202906*lnkd4 - 0.00779095*lnkd5 + 0.00574847*lnkd6
            C[3] = 8.87905 + 3.35689*lnkd - 1.15935*lnkd2 + 0.0673053*lnkd3 + 0.121857*lnkd4 - 0.066642*lnkd5 + 0.0180215*lnkd6
            C[4] = 7.84352 + 5.03178*lnkd - 1.0193*lnkd2 - 0.426362*lnkd3 + 0.328816*lnkd4 - 0.0688176*lnkd5 - 0.0229398*lnkd6
            C[5] = 6.30113 + 6.09914*lnkd - 0.16088*lnkd2 - 1.05578*lnkd3 + 0.338069*lnkd4 + 0.0937157*lnkd5 - 0.106935*lnkd6
            C[6] = 4.65678 + 6.30069*lnkd + 1.13754*lnkd2 - 1.38393*lnkd3 - 0.0134758*lnkd4 + 0.331686*lnkd5 - 0.105954*lnkd6

        if kappa <= 0.1:
            # for small kappa
            C[1] = 4/3*kappa + 8/63*k2
            C[1] = C[1]*np.sqrt(np.pi/5)

            C[2] = 8/21*k2 + 32/693*k3
            C[2] = C[2]*(np.sqrt(np.pi)*0.2)

            C[3] = 16/693*k3 + 32/10395*k4
            C[3] = C[3]*np.sqrt(np.pi/13)

            C[4] = 32/19305*k4
            C[4] = C[4]*np.sqrt(np.pi/17)

            C[5] = 64*np.sqrt(np.pi/21)*k5/692835

            C[6] = 128*np.sqrt(np.pi)*k6/152108775

        return C