Crossing invariant fiber response function with FORECAST model#

We show how to obtain a voxel specific response function in the form of axially symmetric tensor and the fODF using the FORECAST model from [Anderson2005] , [Kaden2016] and [Zucchelli2017].

First import the necessary modules:

import numpy as np
import matplotlib.pyplot as plt
from dipy.reconst.forecast import ForecastModel
from dipy.viz import actor, window
from dipy.data import fetch_hbn, get_sphere
import nibabel as nib
import os.path as op
from dipy.core.gradients import gradient_table

Download and read the data for this tutorial. Our implementation of FORECAST requires multi-shell data.fetch_hbn() provides data that was acquired using b-values of 1000 and 2000 as part of the Healthy Brain Network study [Alexander2017] and was preprocessed and quality controlled in the HBN-POD2 dataset [RichieHalford2022].

data_path = fetch_hbn(["NDARAA948VFH"])[1]
dwi_path = op.join(
       data_path, "derivatives", "qsiprep", "sub-NDARAA948VFH",
       "ses-HBNsiteRU", "dwi")

img = nib.load(op.join(
       dwi_path,
       "sub-NDARAA948VFH_ses-HBNsiteRU_acq-64dir_space-T1w_desc-preproc_dwi.nii.gz"))

gtab = gradient_table(
       op.join(dwi_path,
"sub-NDARAA948VFH_ses-HBNsiteRU_acq-64dir_space-T1w_desc-preproc_dwi.bval"),
       op.join(dwi_path,
"sub-NDARAA948VFH_ses-HBNsiteRU_acq-64dir_space-T1w_desc-preproc_dwi.bvec"))

data = np.asarray(img.dataobj)

mask_img = nib.load(
       op.join(dwi_path,
"sub-NDARAA948VFH_ses-HBNsiteRU_acq-64dir_space-T1w_desc-brain_mask.nii.gz"))

brain_mask = mask_img.get_fdata()
  0%|                                                                                                                               | 0/36 [00:00<?, ?it/s]
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/anat/sub-NDARAA948VFH_desc-brain_mask.nii.gz:   0%|    | 0/36 [00:00<?, ?it/s]
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/anat/sub-NDARAA948VFH_desc-brain_mask.nii.gz:   3%| | 1/36 [00:00<00:09,  3.71
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/anat/sub-NDARAA948VFH_desc-preproc_T1w.nii.gz:   3%| | 1/36 [00:00<00:09,  3.7
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/anat/sub-NDARAA948VFH_desc-preproc_T1w.nii.gz:   6%| | 2/36 [00:02<00:40,  1.2
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/anat/sub-NDARAA948VFH_dseg.nii.gz:   6%|▍      | 2/36 [00:02<00:40,  1.20s/it]
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/anat/sub-NDARAA948VFH_dseg.nii.gz:   8%|▌      | 3/36 [00:02<00:24,  1.36it/s]
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/anat/sub-NDARAA948VFH_from-MNI152NLin2009cAsym_to-T1w_mode-image_xfm.h5:   8%|
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/anat/sub-NDARAA948VFH_from-MNI152NLin2009cAsym_to-T1w_mode-image_xfm.h5:  11%|
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/anat/sub-NDARAA948VFH_from-T1w_to-MNI152NLin2009cAsym_mode-image_xfm.h5:  11%|
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/anat/sub-NDARAA948VFH_from-T1w_to-MNI152NLin2009cAsym_mode-image_xfm.h5:  14%|
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/anat/sub-NDARAA948VFH_label-CSF_probseg.nii.gz:  14%|▏| 5/36 [00:08<01:09,  2.
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/anat/sub-NDARAA948VFH_label-CSF_probseg.nii.gz:  17%|▏| 6/36 [00:09<00:48,  1.
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/anat/sub-NDARAA948VFH_label-GM_probseg.nii.gz:  17%|▏| 6/36 [00:09<00:48,  1.6
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/anat/sub-NDARAA948VFH_label-GM_probseg.nii.gz:  19%|▏| 7/36 [00:09<00:35,  1.2
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/anat/sub-NDARAA948VFH_label-WM_probseg.nii.gz:  19%|▏| 7/36 [00:09<00:35,  1.2
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/anat/sub-NDARAA948VFH_label-WM_probseg.nii.gz:  22%|▏| 8/36 [00:10<00:27,  1.0
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/anat/sub-NDARAA948VFH_space-MNI152NLin2009cAsym_desc-brain_mask.nii.gz:  22%|▏
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/anat/sub-NDARAA948VFH_space-MNI152NLin2009cAsym_desc-brain_mask.nii.gz:  25%|▎
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/anat/sub-NDARAA948VFH_space-MNI152NLin2009cAsym_desc-preproc_T1w.nii.gz:  25%|
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/anat/sub-NDARAA948VFH_space-MNI152NLin2009cAsym_desc-preproc_T1w.nii.gz:  28%|
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/anat/sub-NDARAA948VFH_space-MNI152NLin2009cAsym_dseg.nii.gz:  28%|▎| 10/36 [00
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/anat/sub-NDARAA948VFH_space-MNI152NLin2009cAsym_dseg.nii.gz:  31%|▎| 11/36 [00
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/anat/sub-NDARAA948VFH_space-MNI152NLin2009cAsym_label-CSF_probseg.nii.gz:  31%
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/anat/sub-NDARAA948VFH_space-MNI152NLin2009cAsym_label-CSF_probseg.nii.gz:  33%
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/anat/sub-NDARAA948VFH_space-MNI152NLin2009cAsym_label-GM_probseg.nii.gz:  33%|
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/anat/sub-NDARAA948VFH_space-MNI152NLin2009cAsym_label-GM_probseg.nii.gz:  36%|
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/anat/sub-NDARAA948VFH_space-MNI152NLin2009cAsym_label-WM_probseg.nii.gz:  36%|
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/anat/sub-NDARAA948VFH_space-MNI152NLin2009cAsym_label-WM_probseg.nii.gz:  39%|
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/figures/sub-NDARAA948VFH_seg_brainmask.svg:  39%|▍| 14/36 [00:12<00:11,  1.87i
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/figures/sub-NDARAA948VFH_seg_brainmask.svg:  42%|▍| 15/36 [00:13<00:10,  1.92i
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/figures/sub-NDARAA948VFH_ses-HBNsiteRU_acq-64dir_carpetplot.svg:  42%|▍| 15/36
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/figures/sub-NDARAA948VFH_ses-HBNsiteRU_acq-64dir_coreg.svg:  44%|▍| 16/36 [00:
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/figures/sub-NDARAA948VFH_ses-HBNsiteRU_acq-64dir_coreg.svg:  47%|▍| 17/36 [00:
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/figures/sub-NDARAA948VFH_ses-HBNsiteRU_acq-64dir_desc-resampled_b0ref.svg:  47
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/figures/sub-NDARAA948VFH_ses-HBNsiteRU_acq-64dir_desc-resampled_b0ref.svg:  50
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/figures/sub-NDARAA948VFH_ses-HBNsiteRU_acq-64dir_desc-sdc_b0.svg:  50%|▌| 18/3
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/figures/sub-NDARAA948VFH_ses-HBNsiteRU_acq-64dir_desc-sdc_b0.svg:  53%|▌| 19/3
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/figures/sub-NDARAA948VFH_ses-HBNsiteRU_acq-64dir_dwi_denoise_ses_HBNsiteRU_acq
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/figures/sub-NDARAA948VFH_ses-HBNsiteRU_acq-64dir_dwi_denoise_ses_HBNsiteRU_acq
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/figures/sub-NDARAA948VFH_ses-HBNsiteRU_acq-64dir_dwi_denoise_ses_HBNsiteRU_acq
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/figures/sub-NDARAA948VFH_ses-HBNsiteRU_acq-64dir_dwi_denoise_ses_HBNsiteRU_acq
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/figures/sub-NDARAA948VFH_ses-HBNsiteRU_acq-64dir_dwi_denoise_ses_HBNsiteRU_acq
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/figures/sub-NDARAA948VFH_ses-HBNsiteRU_acq-64dir_dwi_denoise_ses_HBNsiteRU_acq
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/figures/sub-NDARAA948VFH_ses-HBNsiteRU_acq-64dir_sampling_scheme.gif:  61%|▌|
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/figures/sub-NDARAA948VFH_ses-HBNsiteRU_acq-64dir_sampling_scheme.gif:  64%|▋|
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/figures/sub-NDARAA948VFH_t1_2_mni.svg:  64%|█▎| 23/36 [00:15<00:05,  2.40it/s]
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/figures/sub-NDARAA948VFH_t1_2_mni.svg:  67%|█▎| 24/36 [00:16<00:04,  2.66it/s]
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/ses-HBNsiteRU/anat/sub-NDARAA948VFH_ses-HBNsiteRU_acq-HCP_from-orig_to-T1w_mod
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/ses-HBNsiteRU/anat/sub-NDARAA948VFH_ses-HBNsiteRU_acq-HCP_from-orig_to-T1w_mod
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/ses-HBNsiteRU/dwi/sub-NDARAA948VFH_ses-HBNsiteRU_acq-64dir_confounds.tsv:  69%
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/ses-HBNsiteRU/dwi/sub-NDARAA948VFH_ses-HBNsiteRU_acq-64dir_confounds.tsv:  72%
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/ses-HBNsiteRU/dwi/sub-NDARAA948VFH_ses-HBNsiteRU_acq-64dir_desc-ImageQC_dwi.cs
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/ses-HBNsiteRU/dwi/sub-NDARAA948VFH_ses-HBNsiteRU_acq-64dir_desc-ImageQC_dwi.cs
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/ses-HBNsiteRU/dwi/sub-NDARAA948VFH_ses-HBNsiteRU_acq-64dir_desc-SliceQC_dwi.js
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/ses-HBNsiteRU/dwi/sub-NDARAA948VFH_ses-HBNsiteRU_acq-64dir_desc-SliceQC_dwi.js
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/ses-HBNsiteRU/dwi/sub-NDARAA948VFH_ses-HBNsiteRU_acq-64dir_dwiqc.json:  78%|▊|
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/ses-HBNsiteRU/dwi/sub-NDARAA948VFH_ses-HBNsiteRU_acq-64dir_dwiqc.json:  81%|▊|
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/ses-HBNsiteRU/dwi/sub-NDARAA948VFH_ses-HBNsiteRU_acq-64dir_space-T1w_desc-brai
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/ses-HBNsiteRU/dwi/sub-NDARAA948VFH_ses-HBNsiteRU_acq-64dir_space-T1w_desc-brai
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/ses-HBNsiteRU/dwi/sub-NDARAA948VFH_ses-HBNsiteRU_acq-64dir_space-T1w_desc-eddy
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/ses-HBNsiteRU/dwi/sub-NDARAA948VFH_ses-HBNsiteRU_acq-64dir_space-T1w_desc-eddy
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/ses-HBNsiteRU/dwi/sub-NDARAA948VFH_ses-HBNsiteRU_acq-64dir_space-T1w_desc-prep
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/ses-HBNsiteRU/dwi/sub-NDARAA948VFH_ses-HBNsiteRU_acq-64dir_space-T1w_desc-prep
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/ses-HBNsiteRU/dwi/sub-NDARAA948VFH_ses-HBNsiteRU_acq-64dir_space-T1w_desc-prep
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/ses-HBNsiteRU/dwi/sub-NDARAA948VFH_ses-HBNsiteRU_acq-64dir_space-T1w_desc-prep
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/ses-HBNsiteRU/dwi/sub-NDARAA948VFH_ses-HBNsiteRU_acq-64dir_space-T1w_desc-prep
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/ses-HBNsiteRU/dwi/sub-NDARAA948VFH_ses-HBNsiteRU_acq-64dir_space-T1w_desc-prep
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/ses-HBNsiteRU/dwi/sub-NDARAA948VFH_ses-HBNsiteRU_acq-64dir_space-T1w_desc-prep
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/ses-HBNsiteRU/dwi/sub-NDARAA948VFH_ses-HBNsiteRU_acq-64dir_space-T1w_desc-prep
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/ses-HBNsiteRU/dwi/sub-NDARAA948VFH_ses-HBNsiteRU_acq-64dir_space-T1w_dwiref.ni
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/ses-HBNsiteRU/dwi/sub-NDARAA948VFH_ses-HBNsiteRU_acq-64dir_space-T1w_dwiref.ni
Downloading C:\Users\skoudoro\.dipy\HBN/derivatives/qsiprep/sub-NDARAA948VFH/ses-HBNsiteRU/dwi/sub-NDARAA948VFH_ses-HBNsiteRU_acq-64dir_space-T1w_dwiref.ni

Let us consider only a single slice for the FORECAST fitting

data_small = data[:, :, 50:51]
mask_small = brain_mask[:, :, 50:51]

Instantiate the FORECAST Model.

“sh_order” is the spherical harmonics order used for the fODF.

dec_alg is the spherical deconvolution algorithm used for the FORECAST basis fitting, in this case we used the Constrained Spherical Deconvolution (CSD) algorithm.

fm = ForecastModel(gtab, sh_order=6, dec_alg='CSD')

Fit the FORECAST to the data

f_fit = fm.fit(data_small, mask_small)
  0%|                                                                                                                           | 0/5148.0 [00:00<?, ?it/s]
  0%|                                                                                                                   | 1/5148.0 [00:00<10:25,  8.22it/s]
  1%|▊                                                                                                                | 36/5148.0 [00:00<00:26, 191.86it/s]
  1%|█▌                                                                                                               | 69/5148.0 [00:00<00:20, 251.03it/s]
  2%|██                                                                                                               | 96/5148.0 [00:00<00:20, 250.08it/s]
  2%|██▋                                                                                                             | 123/5148.0 [00:00<00:19, 254.47it/s]
  3%|███▌                                                                                                            | 166/5148.0 [00:00<00:16, 311.20it/s]
  4%|████▎                                                                                                           | 200/5148.0 [00:00<00:15, 319.82it/s]
  5%|█████▌                                                                                                          | 255/5148.0 [00:00<00:12, 390.79it/s]
  6%|██████▍                                                                                                         | 298/5148.0 [00:00<00:12, 402.10it/s]
  7%|███████▍                                                                                                        | 339/5148.0 [00:01<00:12, 391.85it/s]
  8%|████████▍                                                                                                       | 390/5148.0 [00:01<00:11, 426.28it/s]
  8%|█████████▍                                                                                                      | 433/5148.0 [00:01<00:11, 401.36it/s]
  9%|██████████▎                                                                                                     | 475/5148.0 [00:01<00:11, 401.52it/s]
 10%|███████████▏                                                                                                    | 516/5148.0 [00:01<00:12, 384.63it/s]
 11%|████████████                                                                                                    | 555/5148.0 [00:01<00:13, 352.84it/s]
 11%|████████████▊                                                                                                   | 591/5148.0 [00:01<00:12, 351.47it/s]
 12%|█████████████▋                                                                                                  | 627/5148.0 [00:01<00:13, 334.15it/s]
 13%|██████████████▍                                                                                                 | 661/5148.0 [00:01<00:13, 330.22it/s]
 14%|███████████████                                                                                                 | 695/5148.0 [00:02<00:13, 327.93it/s]
 14%|███████████████▉                                                                                                | 735/5148.0 [00:02<00:12, 342.88it/s]
 15%|████████████████▊                                                                                               | 771/5148.0 [00:02<00:12, 344.09it/s]
 16%|█████████████████▌                                                                                              | 806/5148.0 [00:02<00:12, 343.12it/s]
 16%|██████████████████▎                                                                                             | 841/5148.0 [00:02<00:12, 344.33it/s]
 17%|███████████████████▍                                                                                            | 892/5148.0 [00:02<00:10, 391.67it/s]
 18%|████████████████████▍                                                                                           | 939/5148.0 [00:02<00:10, 413.96it/s]
 19%|█████████████████████▍                                                                                          | 985/5148.0 [00:02<00:09, 426.77it/s]
 20%|██████████████████████▏                                                                                        | 1028/5148.0 [00:02<00:09, 414.62it/s]
 21%|███████████████████████                                                                                        | 1070/5148.0 [00:03<00:10, 382.25it/s]
 22%|███████████████████████▉                                                                                       | 1109/5148.0 [00:03<00:11, 363.40it/s]
 23%|█████████████████████████                                                                                      | 1164/5148.0 [00:03<00:09, 413.28it/s]
 23%|██████████████████████████                                                                                     | 1207/5148.0 [00:03<00:09, 405.59it/s]
 24%|██████████████████████████▉                                                                                    | 1249/5148.0 [00:03<00:10, 377.16it/s]
 25%|███████████████████████████▊                                                                                   | 1288/5148.0 [00:03<00:10, 379.90it/s]
 26%|████████████████████████████▋                                                                                  | 1330/5148.0 [00:03<00:09, 389.88it/s]
 27%|█████████████████████████████▊                                                                                 | 1385/5148.0 [00:03<00:08, 434.38it/s]
 28%|██████████████████████████████▊                                                                                | 1429/5148.0 [00:03<00:09, 400.50it/s]
 29%|███████████████████████████████▊                                                                               | 1474/5148.0 [00:04<00:08, 413.06it/s]
 29%|████████████████████████████████▋                                                                              | 1517/5148.0 [00:04<00:08, 416.91it/s]
 30%|█████████████████████████████████▋                                                                             | 1560/5148.0 [00:04<00:09, 386.15it/s]
 31%|██████████████████████████████████▌                                                                            | 1601/5148.0 [00:04<00:09, 391.68it/s]
 32%|███████████████████████████████████▍                                                                           | 1641/5148.0 [00:04<00:09, 382.34it/s]
 33%|████████████████████████████████████▏                                                                          | 1680/5148.0 [00:04<00:09, 353.83it/s]
 33%|█████████████████████████████████████                                                                          | 1717/5148.0 [00:04<00:09, 357.45it/s]
 34%|█████████████████████████████████████▊                                                                         | 1754/5148.0 [00:04<00:09, 360.34it/s]
 35%|██████████████████████████████████████▉                                                                        | 1804/5148.0 [00:04<00:08, 398.90it/s]
 36%|███████████████████████████████████████▉                                                                       | 1853/5148.0 [00:04<00:07, 424.19it/s]
 37%|████████████████████████████████████████▉                                                                      | 1896/5148.0 [00:05<00:07, 412.72it/s]
 38%|█████████████████████████████████████████▊                                                                     | 1941/5148.0 [00:05<00:07, 422.38it/s]
 39%|██████████████████████████████████████████▊                                                                    | 1984/5148.0 [00:05<00:07, 423.53it/s]
 39%|███████████████████████████████████████████▋                                                                   | 2028/5148.0 [00:05<00:07, 427.25it/s]
 40%|████████████████████████████████████████████▋                                                                  | 2071/5148.0 [00:05<00:07, 414.85it/s]
 41%|█████████████████████████████████████████████▌                                                                 | 2113/5148.0 [00:05<00:07, 415.30it/s]
 42%|██████████████████████████████████████████████▍                                                                | 2155/5148.0 [00:05<00:07, 394.72it/s]
 43%|███████████████████████████████████████████████▌                                                               | 2203/5148.0 [00:05<00:07, 417.94it/s]
 44%|████████████████████████████████████████████████▍                                                              | 2246/5148.0 [00:05<00:06, 418.34it/s]
 45%|█████████████████████████████████████████████████▋                                                             | 2305/5148.0 [00:06<00:06, 467.13it/s]
 46%|██████████████████████████████████████████████████▊                                                            | 2357/5148.0 [00:06<00:05, 481.29it/s]
 47%|███████████████████████████████████████████████████▉                                                           | 2406/5148.0 [00:06<00:06, 420.60it/s]
 48%|████████████████████████████████████████████████████▊                                                          | 2450/5148.0 [00:06<00:06, 413.80it/s]
 48%|█████████████████████████████████████████████████████▊                                                         | 2493/5148.0 [00:06<00:06, 417.16it/s]
 49%|██████████████████████████████████████████████████████▋                                                        | 2539/5148.0 [00:06<00:06, 428.25it/s]
 50%|███████████████████████████████████████████████████████▋                                                       | 2583/5148.0 [00:06<00:05, 430.70it/s]
 52%|█████████████████████████████████████████████████████████▎                                                     | 2660/5148.0 [00:06<00:04, 527.38it/s]
 53%|██████████████████████████████████████████████████████████▊                                                    | 2727/5148.0 [00:06<00:04, 567.37it/s]
 54%|████████████████████████████████████████████████████████████                                                   | 2785/5148.0 [00:07<00:04, 553.33it/s]
 55%|█████████████████████████████████████████████████████████████▎                                                 | 2841/5148.0 [00:07<00:04, 483.93it/s]
 56%|██████████████████████████████████████████████████████████████▎                                                | 2892/5148.0 [00:07<00:04, 453.45it/s]
 57%|███████████████████████████████████████████████████████████████▎                                               | 2939/5148.0 [00:07<00:05, 423.63it/s]
 58%|████████████████████████████████████████████████████████████████▎                                              | 2983/5148.0 [00:07<00:05, 396.54it/s]
 59%|█████████████████████████████████████████████████████████████████▏                                             | 3024/5148.0 [00:07<00:05, 379.53it/s]
 60%|██████████████████████████████████████████████████████████████████▏                                            | 3072/5148.0 [00:07<00:05, 404.31it/s]
 60%|███████████████████████████████████████████████████████████████████▏                                           | 3114/5148.0 [00:07<00:05, 388.84it/s]
 61%|████████████████████████████████████████████████████████████████████                                           | 3157/5148.0 [00:08<00:05, 396.94it/s]
 62%|████████████████████████████████████████████████████████████████████▉                                          | 3198/5148.0 [00:08<00:04, 391.06it/s]
 63%|█████████████████████████████████████████████████████████████████████▊                                         | 3238/5148.0 [00:08<00:04, 390.34it/s]
 64%|██████████████████████████████████████████████████████████████████████▋                                        | 3278/5148.0 [00:08<00:04, 392.40it/s]
 65%|███████████████████████████████████████████████████████████████████████▊                                       | 3328/5148.0 [00:08<00:04, 422.46it/s]
 65%|████████████████████████████████████████████████████████████████████████▋                                      | 3371/5148.0 [00:08<00:04, 424.01it/s]
 66%|█████████████████████████████████████████████████████████████████████████▌                                     | 3414/5148.0 [00:08<00:04, 424.84it/s]
 67%|██████████████████████████████████████████████████████████████████████████▌                                    | 3457/5148.0 [00:08<00:03, 425.38it/s]
 68%|███████████████████████████████████████████████████████████████████████████▍                                   | 3500/5148.0 [00:08<00:04, 401.77it/s]
 69%|████████████████████████████████████████████████████████████████████████████▎                                  | 3541/5148.0 [00:08<00:04, 381.20it/s]
 70%|█████████████████████████████████████████████████████████████████████████████▋                                 | 3601/5148.0 [00:09<00:03, 440.89it/s]
 71%|██████████████████████████████████████████████████████████████████████████████▊                                | 3656/5148.0 [00:09<00:03, 470.77it/s]
 72%|███████████████████████████████████████████████████████████████████████████████▉                               | 3706/5148.0 [00:09<00:03, 478.09it/s]
 73%|████████████████████████████████████████████████████████████████████████████████▉                              | 3755/5148.0 [00:09<00:03, 444.07it/s]
 74%|█████████████████████████████████████████████████████████████████████████████████▉                             | 3801/5148.0 [00:09<00:03, 400.39it/s]
 75%|██████████████████████████████████████████████████████████████████████████████████▊                            | 3843/5148.0 [00:09<00:03, 375.12it/s]
 75%|███████████████████████████████████████████████████████████████████████████████████▋                           | 3882/5148.0 [00:09<00:03, 370.41it/s]
 76%|████████████████████████████████████████████████████████████████████████████████████▌                          | 3924/5148.0 [00:09<00:03, 380.83it/s]
 77%|█████████████████████████████████████████████████████████████████████████████████████▍                         | 3965/5148.0 [00:09<00:03, 387.78it/s]
 78%|███████████████████████████████████████████████████████████████████████████████████████                        | 4036/5148.0 [00:10<00:02, 467.37it/s]
 80%|████████████████████████████████████████████████████████████████████████████████████████▍                      | 4103/5148.0 [00:10<00:01, 522.57it/s]
 81%|█████████████████████████████████████████████████████████████████████████████████████████▌                     | 4156/5148.0 [00:10<00:02, 492.30it/s]
 82%|██████████████████████████████████████████████████████████████████████████████████████████▋                    | 4206/5148.0 [00:10<00:01, 482.59it/s]
 83%|███████████████████████████████████████████████████████████████████████████████████████████▋                   | 4255/5148.0 [00:10<00:02, 432.74it/s]
 84%|████████████████████████████████████████████████████████████████████████████████████████████▋                  | 4300/5148.0 [00:10<00:02, 404.44it/s]
 84%|█████████████████████████████████████████████████████████████████████████████████████████████▌                 | 4342/5148.0 [00:10<00:02, 378.04it/s]
 85%|██████████████████████████████████████████████████████████████████████████████████████████████▋                | 4391/5148.0 [00:10<00:01, 405.52it/s]
 86%|███████████████████████████████████████████████████████████████████████████████████████████████▊               | 4445/5148.0 [00:11<00:01, 440.15it/s]
 87%|████████████████████████████████████████████████████████████████████████████████████████████████▊              | 4491/5148.0 [00:11<00:01, 421.52it/s]
 88%|█████████████████████████████████████████████████████████████████████████████████████████████████▊             | 4535/5148.0 [00:11<00:01, 414.07it/s]
 89%|██████████████████████████████████████████████████████████████████████████████████████████████████▋            | 4577/5148.0 [00:11<00:01, 403.34it/s]
 90%|███████████████████████████████████████████████████████████████████████████████████████████████████▌           | 4618/5148.0 [00:11<00:01, 383.07it/s]
 91%|████████████████████████████████████████████████████████████████████████████████████████████████████▍          | 4659/5148.0 [00:11<00:01, 389.46it/s]
 91%|█████████████████████████████████████████████████████████████████████████████████████████████████████▍         | 4703/5148.0 [00:11<00:01, 402.70it/s]
 92%|██████████████████████████████████████████████████████████████████████████████████████████████████████▍        | 4750/5148.0 [00:11<00:00, 421.15it/s]
 93%|███████████████████████████████████████████████████████████████████████████████████████████████████████▎       | 4793/5148.0 [00:11<00:00, 422.66it/s]
 94%|████████████████████████████████████████████████████████████████████████████████████████████████████████▎      | 4836/5148.0 [00:11<00:00, 423.70it/s]
 95%|█████████████████████████████████████████████████████████████████████████████████████████████████████████▏     | 4879/5148.0 [00:12<00:00, 403.08it/s]
 96%|██████████████████████████████████████████████████████████████████████████████████████████████████████████▏    | 4926/5148.0 [00:12<00:00, 420.92it/s]
 97%|███████████████████████████████████████████████████████████████████████████████████████████████████████████▏   | 4972/5148.0 [00:12<00:00, 431.39it/s]
 97%|████████████████████████████████████████████████████████████████████████████████████████████████████████████▏  | 5018/5148.0 [00:12<00:00, 439.00it/s]
 98%|█████████████████████████████████████████████████████████████████████████████████████████████████████████████▏ | 5066/5148.0 [00:12<00:00, 450.18it/s]
 99%|██████████████████████████████████████████████████████████████████████████████████████████████████████████████▏| 5112/5148.0 [00:12<00:00, 427.12it/s]
100%|███████████████████████████████████████████████████████████████████████████████████████████████████████████████| 5148/5148.0 [00:12<00:00, 405.53it/s]

Calculate the crossing invariant tensor indices [Kaden2016] : the parallel diffusivity, the perpendicular diffusivity, the fractional anisotropy and the mean diffusivity.

d_par = f_fit.dpar
d_perp = f_fit.dperp
fa = f_fit.fractional_anisotropy()
md = f_fit.mean_diffusivity()

Show the indices and save them in FORECAST_indices.png.

fig = plt.figure(figsize=(6, 6))
ax1 = fig.add_subplot(2, 2, 1, title='parallel diffusivity')
ax1.set_axis_off()
ind = ax1.imshow(d_par[:, :, 0].T, interpolation='nearest',
                 origin='lower', cmap=plt.cm.gray)
plt.colorbar(ind, shrink=0.6)
ax2 = fig.add_subplot(2, 2, 2, title='perpendicular diffusivity')
ax2.set_axis_off()
ind = ax2.imshow(d_perp[:, :, 0].T, interpolation='nearest',
                 origin='lower', cmap=plt.cm.gray)
plt.colorbar(ind, shrink=0.6)
ax3 = fig.add_subplot(2, 2, 3, title='fractional anisotropy')
ax3.set_axis_off()
ind = ax3.imshow(fa[:, :, 0].T, interpolation='nearest',
                 origin='lower', cmap=plt.cm.gray)
plt.colorbar(ind, shrink=0.6)
ax4 = fig.add_subplot(2, 2, 4, title='mean diffusivity')
ax4.set_axis_off()
ind = ax4.imshow(md[:, :, 0].T, interpolation='nearest',
                 origin='lower', cmap=plt.cm.gray)
plt.colorbar(ind, shrink=0.6)
plt.savefig('FORECAST_indices.png', dpi=300, bbox_inches='tight')
reconst forecast

FORECAST scalar indices.

Load an ODF reconstruction sphere

sphere = get_sphere('repulsion724')

Compute the fODFs.

odf = f_fit.odf(sphere)
print('fODF.shape (%d, %d, %d, %d)' % odf.shape)
fODF.shape (108, 129, 1, 724)

Display a part of the fODFs

odf_actor = actor.odf_slicer(odf[30:60, 30:60, :], sphere=sphere,
                             colormap='plasma', scale=0.6)
scene = window.Scene()
scene.add(odf_actor)
window.record(scene, out_path='fODFs.png', size=(600, 600), magnification=4)
reconst forecast

Fiber Orientation Distribution Functions, in a small ROI of the brain.

References#

[Anderson2005]

Anderson A. W., “Measurement of Fiber Orientation Distributions Using High Angular Resolution Diffusion Imaging”, Magnetic Resonance in Medicine, 2005.

[Kaden2016] (1,2)

Kaden E. et al., “Quantitative Mapping of the Per-Axon Diffusion Coefficients in Brain White Matter”, Magnetic Resonance in Medicine, 2016.

[Zucchelli2017]

Zucchelli E. et al., “A generalized SMT-based framework for Diffusion MRI microstructural model estimation”, MICCAI Workshop on Computational DIFFUSION MRI (CDMRI), 2017.

[Alexander2017]

Alexander LM, Escalera J, Ai L, et al. An open resource for transdiagnostic research in pediatric mental health and learning disorders. Sci Data. 2017;4:170181.

[RichieHalford2022]

Richie-Halford A, Cieslak M, Ai L, et al. An analysis-ready and quality controlled resource for pediatric brain white-matter research. Scientific Data. 2022;9(1):1-27.

Total running time of the script: (1 minutes 16.956 seconds)

Gallery generated by Sphinx-Gallery