API

Preprocessing Workflows

qsiprep base processing workflows

qsiprep.workflows.base.init_qsiprep_wf(subject_list, run_uuid, work_dir, output_dir, bids_dir, ignore, debug, low_mem, anat_only, dwi_only, longitudinal, b0_threshold, hires, anatomical_contrast, denoise_before_combining, dwi_denoise_window, denoise_method, unringing_method, b1_biascorrect_stage, no_b0_harmonization, output_resolution, infant_mode, combine_all_dwis, distortion_group_merge, pepolar_method, omp_nthreads, bids_filters, force_spatial_normalization, skull_strip_template, skull_strip_fixed_seed, freesurfer, hmc_model, impute_slice_threshold, hmc_transform, shoreline_iters, eddy_config, write_local_bvecs, template, motion_corr_to, b0_to_t1w_transform, intramodal_template_iters, intramodal_template_transform, prefer_dedicated_fmaps, fmap_bspline, fmap_demean, use_syn, force_syn, raw_image_sdc)[source]

This workflow organizes the execution of qsiprep, with a sub-workflow for each subject.

(Source code)

Parameters:

subject_listlist

List of subject labels

run_uuidstr

Unique identifier for execution instance

work_dirstr

Directory in which to store workflow execution state and temporary files

output_dirstr

Directory in which to save derivatives

bids_dirstr

Root directory of BIDS dataset

anatomical_contrststr

Which contrast to use for the anatomical reference

ignorelist

Preprocessing steps to skip (may include “slicetiming”, “fieldmaps”)

low_membool

Write uncompressed .nii files in some cases to reduce memory usage

anat_onlybool

Disable diffusion workflows

dwi_onlybool

Disable the anatomical (T1w/T2w) workflows

longitudinalbool

Treat multiple sessions as longitudinal (may increase runtime) See sub-workflows for specific differences

infant_modebool

Run the pipeline for infant DWIs

b0_thresholdint

Images with b-values less than this value will be treated as a b=0 image.

dwi_denoise_windowint

window size in voxels for image-based denoising. Must be odd. If 0, ‘ ‘denoising will not be run’

denoise_methodstr

Either ‘dwidenoise’, ‘patch2self’ or ‘none’

unringing_methodstr

algorithm to use for removing Gibbs ringing. Options: none, mrdegibbs

b1_biascorrect_stagestr

‘final’, ‘none’ or ‘legacy’

no_b0_harmonizationbool

skip rescaling dwi scans to have matching b=0 intensities across scans

denoise_before_combiningbool

‘run dwidenoise before combining dwis. Requires combine_all_dwis

combine_all_dwisbool

Combine all dwi sequences within a session into a single data set

distortion_group_mergestr

How to combine multiple distortion groups after correction. ‘concat’, ‘average’ or ‘none’

pepolar_methodstr

Either ‘DRBUDDI’, ‘TOPUP’ or ‘TOPUP+DRBUDDI’. The method for SDC when EPI fieldmaps are used.

omp_nthreadsint

Maximum number of threads an individual process may use

skull_strip_templatestr

Name of ANTs skull-stripping template (‘OASIS’ or ‘NKI’)

skull_strip_fixed_seedbool

Do not use a random seed for skull-stripping - will ensure run-to-run replicability when used with –omp-nthreads 1

freesurferbool

Enable FreeSurfer surface reconstruction (may increase runtime)

hiresbool

Enable sub-millimeter preprocessing in FreeSurfer

templatestr

Name of template targeted by template output space

motion_corr_tostr

Motion correct using the ‘first’ b0 image or use an ‘iterative’ method to motion correct to the midpoint of the b0 images

b0_to_t1w_transform“Rigid” or “Affine”

Use a rigid or full affine transform for b0-T1w registration

intramodal_template_iters: int

Number of iterations for finding the midpoint image from the b0 templates from all groups. Has no effect if there is only one group. If 0, all b0 templates are directly registered to the t1w image.

intramodal_template_transform: str

Transformation used for building the intramodal template.

hmc_model‘none’, ‘3dSHORE’ or ‘MAPMRI’

Model used to generate target images for head motion correction. If ‘none’ the transform from the nearest b0 will be used.

hmc_transform“Rigid” or “Affine”

Type of transform used for head motion correction

impute_slice_thresholdfloat

Impute data in slices that are this many SDs from expected. If 0, no slices will be imputed.

eddy_config: str

Path to a JSON file containing config options for eddy

raw_image_sdc: bool

Use raw (direct from BIDS) images for distortion

prefer_dedicated_fmaps: bool

If a reverse PE fieldmap is available in fmap, use that even if a reverse PE DWI series is available

write_local_bvecsbool

Write out a series of voxelwise bvecs

fmap_bsplinebool

Experimental: Fit B-Spline field using least-squares

fmap_demeanbool

Demean voxel-shift map during unwarp

use_synbool

Experimental: Enable ANTs SyN-based susceptibility distortion correction (SDC). If fieldmaps are present and enabled, this is not run, by default.

force_synbool

Temporary: Always run SyN-based SDC

qsiprep.workflows.base.init_single_subject_wf(subject_id, name, reportlets_dir, output_dir, bids_dir, ignore, debug, write_local_bvecs, low_mem, dwi_only, anat_only, longitudinal, b0_threshold, denoise_before_combining, bids_filters, anatomical_contrast, dwi_denoise_window, denoise_method, unringing_method, b1_biascorrect_stage, no_b0_harmonization, infant_mode, combine_all_dwis, raw_image_sdc, distortion_group_merge, pepolar_method, omp_nthreads, skull_strip_template, force_spatial_normalization, skull_strip_fixed_seed, freesurfer, hires, template, output_resolution, prefer_dedicated_fmaps, motion_corr_to, b0_to_t1w_transform, intramodal_template_iters, intramodal_template_transform, hmc_model, hmc_transform, shoreline_iters, eddy_config, impute_slice_threshold, fmap_bspline, fmap_demean, use_syn, force_syn)[source]

This workflow organizes the preprocessing pipeline for a single subject. It collects and reports information about the subject, and prepares sub-workflows to perform anatomical and diffusion preprocessing.

Anatomical preprocessing is performed in a single workflow, regardless of the number of sessions. Diffusion preprocessing is performed using a separate workflow for each session’s dwi series.

(Source code)

Parameters

subject_idstr

List of subject labels

namestr

Name of workflow

ignorelist

Preprocessing steps to skip (may include “sbref”, “fieldmaps”)

debugbool

Do inaccurate but fast normalization

low_membool

Write uncompressed .nii files in some cases to reduce memory usage

anat_onlybool

Disable dMRI workflows

anatomical_contraststr

Which contrast to use for the anatomical reference

dwi_onlybool

Disable anatomical workflows

longitudinalbool

Treat multiple sessions as longitudinal (may increase runtime) See sub-workflows for specific differences

b0_thresholdint

Images with b-values less than this value will be treated as a b=0 image.

dwi_denoise_windowint

window size in voxels for image-based denoising. Must be odd. If 0, ‘ ‘denoising will not be run’

denoise_methodstr

Either ‘dwidenoise’, ‘patch2self’ or ‘none’

unringing_methodstr

algorithm to use for removing Gibbs ringing. Options: none, mrdegibbs

b1_biascorrect_stagestr

‘final’, ‘none’ or ‘legacy’

no_b0_harmonizationbool

skip rescaling dwi scans to have matching b=0 intensities across scans

denoise_before_combiningbool

‘run dwidenoise before combining dwis. Requires combine_all_dwis

combine_all_dwisBool

Combine all dwi sequences within a session into a single data set

distortion_group_merge: str

How to combine preprocessed scans from different distortion groups. ‘concat’, ‘average’ or ‘none’

pepolar_methodstr

Either ‘DRBUDDI’ or ‘TOPUP’. The method for SDC when EPI fieldmaps are used.

omp_nthreadsint

Maximum number of threads an individual process may use

skull_strip_templatestr

Name of ANTs skull-stripping template (‘OASIS’ or ‘NKI’)

skull_strip_fixed_seedbool

Do not use a random seed for skull-stripping - will ensure run-to-run replicability when used with –omp-nthreads 1

freesurferbool

Enable FreeSurfer surface reconstruction (may increase runtime)

hiresbool

Enable sub-millimeter preprocessing in FreeSurfer

reportlets_dirstr

Directory in which to save reportlets

output_dirstr

Directory in which to save derivatives

bids_dirstr

Root directory of BIDS dataset

templatestr

Name of template targeted by template output space

hmc_model‘none’, ‘3dSHORE’ or ‘eddy’

Model used to generate target images for head motion correction. If ‘none’ the transform from the nearest b0 will be used.

hmc_transform“Rigid” or “Affine”

Type of transform used for head motion correction

impute_slice_thresholdfloat

Impute data in slices that are this many SDs from expected. If 0, no slices will be imputed.

motion_corr_tostr

Motion correct using the ‘first’ b0 image or use an ‘iterative’ method to motion correct to the midpoint of the b0 images

eddy_config: str

Path to a JSON file containing config options for eddy

raw_image_sdc: bool

Use raw (direct from BIDS) images for distortion

fmap_bsplinebool

Experimental: Fit B-Spline field using least-squares

fmap_demeanbool

Demean voxel-shift map during unwarp

use_synbool

Experimental: Enable ANTs SyN-based susceptibility distortion correction (SDC). If fieldmaps are present and enabled, this is not run, by default.

force_synbool

Temporary: Always run SyN-based SDC

eddy_config: str

Path to a JSON file containing config options for eddy

b0_to_t1w_transform“Rigid” or “Affine”

Use a rigid or full affine transform for b0-T1w registration

intramodal_template_iters: int

Number of iterations for finding the midpoint image from the b0 templates from all groups. Has no effect if there is only one group. If 0, all b0 templates are directly registered to the t1w image.

intramodal_template_transform: str

Transformation used for building the intramodal template.

Inputs

subjects_dir

FreeSurfer SUBJECTS_DIR

Pre-processing q-space Image workflows

Orchestrating the dwi-preprocessing workflow

qsiprep.workflows.dwi.base.init_dwi_preproc_wf(dwi_only, scan_groups, output_prefix, ignore, b0_threshold, motion_corr_to, b0_to_t1w_transform, hmc_model, hmc_transform, shoreline_iters, impute_slice_threshold, eddy_config, raw_image_sdc, reportlets_dir, output_dir, dwi_denoise_window, denoise_method, unringing_method, pepolar_method, b1_biascorrect_stage, no_b0_harmonization, denoise_before_combining, template, omp_nthreads, fmap_bspline, fmap_demean, t2w_sdc, low_mem, sloppy, source_file, layout=None)[source]

This workflow controls the dwi preprocessing stages of qsiprep.

(Source code)

Parameters

dwi_groupslist of dicts

List of dicts grouping files by PE-dir

output_prefixstr

beginning of the output file name (eg ‘sub-1_buds-j’)

ignorelist

Preprocessing steps to skip (eg “fieldmaps”)

b0_thresholdint

Images with b-values less than this value will be treated as a b=0 image.

freesurferbool

Enable FreeSurfer functional registration (bbregister) and resampling dwi series to FreeSurfer surface meshes.

motion_corr_tostr

Motion correct using the ‘first’ b0 image or use an ‘iterative’ method to motion correct to the midpoint of the b0 images

b0_to_t1w_transform“Rigid” or “Affine”

Use a rigid or full affine transform for b0-T1w registration

hmc_model‘none’, ‘3dSHORE’, ‘eddy’ or ‘eddy_ingress’

Model used to generate target images for head motion correction. If ‘none’ the transform from the nearest b0 will be used.

hmc_transform“Rigid” or “Affine”

Type of transform used for head motion correction

impute_slice_thresholdfloat

Impute data in slices that are this many SDs from expected. If 0, no slices will be imputed.

eddy_config: str

Path to a JSON file containing config options for eddy

dwi_denoise_windowint

window size in voxels for image-based denoising. Must be odd. If 0, ‘ ‘denoising will not be run’

denoise_methodstr

Either ‘dwidenoise’, ‘patch2self’ or ‘none’

unringing_methodstr

algorithm to use for removing Gibbs ringing. Options: none, mrdegibbs

pepolar_methodstr

Either ‘DRBUDDI’, ‘TOPUP’ or ‘TOPUP+DRBUDDI’. The method for SDC when EPI fieldmaps are used.

b1_biascorrect_stagestr

‘final’, ‘none’ or ‘legacy’

no_b0_harmonizationbool

skip rescaling dwi scans to have matching b=0 intensities across scans

denoise_before_combiningbool

‘run dwidenoise before combining dwis. Requires combine_all_dwis

reportlets_dirstr

Directory in which to save reportlets

templatestr

Name of template targeted by template output space

output_dirstr

Directory in which to save derivatives

omp_nthreadsint

Maximum number of threads an individual process may use

fmap_bsplinebool

Experimental: Fit B-Spline field using least-squares

fmap_demeanbool

Demean voxel-shift map during unwarp

t2w_sdcbool

Are t2w’s going to be included in SDC?

low_membool

Write uncompressed .nii files in some cases to reduce memory usage

layoutBIDSLayout

BIDSLayout structure to enable metadata retrieval

num_dwiint

Total number of dwi files that have been set for preprocessing (default is 1)

sloppybool

Use low-quality settings for motion correction

source_filestr

The file name template used for derivatives

Inputs

t1_preproc

Bias-corrected structural template image

t1_brain

Skull-stripped t1_preproc

t1_mask

Mask of the skull-stripped template image

t1_output_grid

Image to write out DWIs aligned to t1

t1_seg

Segmentation of preprocessed structural image, including gray-matter (GM), white-matter (WM) and cerebrospinal fluid (CSF)

t1_2_mni_forward_transform

ANTs-compatible affine-and-warp transform file

t1_2_mni_reverse_transform

ANTs-compatible affine-and-warp transform file (inverse)

subjects_dir

FreeSurfer SUBJECTS_DIR

subject_id

FreeSurfer subject ID

t1_2_fsnative_forward_transform

LTA-style affine matrix translating from T1w to FreeSurfer-conformed subject space

t1_2_fsnative_reverse_transform

LTA-style affine matrix translating from FreeSurfer-conformed subject space to T1w

dwi_sampling_grid

A NIfTI1 file with the grid spacing and FoV to resample the DWIs

Outputs

dwi_t1

dwi series, resampled to T1w space

dwi_mask_t1

dwi series mask in T1w space

bvals_t1

bvalues of the dwi series

bvecs_t1

bvecs after aligning to the T1w and resampling

local_bvecs_t1

voxelwise bvecs accounting for local displacements

gradient_table_t1

MRTrix-style gradient table

dwi_mni

dwi series, resampled to template space

dwi_mask_mni

dwi series mask in template space

bvals_mni

bvalues of the dwi series

bvecs_mni

bvecs after aligning to the T1w and resampling

local_bvecs_mni

voxelwise bvecs accounting for local displacements

gradient_table_mni

MRTrix-style gradient table

confounds_file

estimated motion parameters and zipper scores

raw_qc_file

DSI Studio QC file for the raw data

raw_concatenated

concatenated raw images for a qc report

carpetplot_data

path to a file containing carpetplot data

Subworkflows

Utility workflows

qsiprep.workflows.dwi.util.init_dwi_reference_wf(omp_nthreads, dwi_file=None, register_t1=False, name='dwi_reference_wf', gen_report=False, source_file=None, desc='initial', sloppy=False)[source]

If register_t2, a skull-stripped T1w image is downsampled to the resolution of the b=0 input image and registered to it.

../_images/index-4.png

(Source code, png, svg, pdf)

Parameters

dwi_filestr

A b=0 image

omp_nthreadsint

Maximum number of threads an individual process may use

namestr

Name of workflow (default: dwi_reference_wf)

gen_reportbool

Whether a mask report node should be appended in the end

Inputs

b0_template

the b0 template used as the motion correction reference

t1_brain

skull-stripped T1w image from the same subject

t1_mask

mask image for t1_brain

wm_seg

white matter segmentation from the T1w image

Outputs

raw_ref_image

Reference image to which dwi series is motion corrected

ref_image

Contrast-enhanced reference image

ref_image_brain

Skull-stripped reference image

dwi_mask

Skull-stripping (rough) mask of reference image

validation_report

HTML reportlet indicating whether dwi_file had a valid affine

Subworkflows

  • init_enhance_and_skullstrip_wf()

Head motion correction

qsiprep.workflows.dwi.hmc.init_dwi_hmc_wf(hmc_transform, hmc_model, hmc_align_to, source_file, num_model_iterations=2, mem_gb=3, omp_nthreads=1, sloppy=False, name='dwi_hmc_wf')[source]

Perform head motion correction and susceptibility distortion correction.

This workflow uses antsRegistration and an iteratively updated signal model to perform motion correction.

Parameters

hmc_transform: ‘Rigid’ or ‘Affine’

How many degrees of freedom to incorporate into motion correction

hmc_model: ‘3dSHORE’, ‘none’ , ‘tensor’ or ‘SH’

Which model to use for generating signal predictions for hmc. ‘3dSHORE’ requires multiple b-values, ‘none’ will only use b0 images for motion correction, ‘tensor’ uses a tensor model for signal predictions for hmc, and ‘SH’ uses spherical harmonics (not implemented yet).

hmc_align_to: ‘first’ or ‘iterative’

Which volume should be used to determine the motion-corrected space?

source_file: str

Path to one of the original dwi files (used for reportlets)

rpe_b0: str

Path to a reverse phase encoding image to be used for 3dQWarp’s TOPUP-style correction

num_model_iterations: int

If hmc_model is '3dSHORE' or 'SH' determines the number of times the model is updated and motion corretion is estimated. Default: 2.

Inputs

dwi_files: list

List of single-volume files across all DWI series

b0_indices: list

Indexes into dwi_files that correspond to b=0 volumes

bvecs: list

List of paths to single-line bvec files

bvals: list

List of paths to single-line bval files

b0_images: list

List of single b=0 volumes

original_files: list

List of the files from which each DWI volume came from.

Outputs

final_template: str

Path to the mean of the coregistered b0 images

forward_transforms: list

List of ITK transforms that motion-correct the images in dwi_files

noise_free_dwis: list

Model-predicted images reverse-transformed into alignment with dwi_files

cnr_image: str

If hmc_model is ‘none’ this is the tsnr of the b=0 images. Otherwise it is the model fit divided by the model error in each voxel.

optimization_data: str

CSV file tracking the motion estimates across shoreline iterations

qsiprep.workflows.dwi.hmc.init_dwi_model_hmc_wf(modelname, transform, mem_gb, omp_nthreads, num_iters=2, name='dwi_model_hmc_wf', metric='Mattes')[source]

Create a model-based hmc workflow.

../_images/index-5.png

(Source code, png, svg, pdf)

Parameters

modelnamestr

one of the models for reconstructing an EAP and producing signal estimates used for motion correction

transformstr

either “Rigid” or “Affine”. Choosing “Affine” may help with Eddy warping

num_itersint

the number of times the model will be updated with transformed data

Inputs

dwi_files

list of 3d dwi files

b0_indices

list of which indices in dwi_files are b0 images

initial_transforms

list of b0-based transforms from dwis to the b0 template

warped_b0_images

list of aligned b0 images

b0_mask

mask of containing brain voxels

bvecs

list of bvec files corresponding to dwi_files

bvals

list of bval files corresponding to dwi_files

Outputs

hmc_transforms

list of transforms, one per file in dwi_files

model_predicted_images: list

Model-predicted images reverse-transformed into alignment with dwi_files

cnr_image: str

If hmc_model is ‘none’ this is the tsnr of the b=0 images. Otherwise it is the model fit divided by the model error in each voxel.

optimization_data: str

CSV file tracking the motion estimates across shoreline iterations

Implementing the FSL preprocessing workflow

qsiprep.workflows.dwi.fsl.init_fsl_hmc_wf(scan_groups, source_file, b0_threshold, impute_slice_threshold, fmap_demean, fmap_bspline, eddy_config, raw_image_sdc, pepolar_method, t2w_sdc, mem_gb=3, omp_nthreads=1, dwi_metadata=None, slice_quality='outlier_n_sqr_stdev_map', sloppy=False, name='fsl_hmc_wf')[source]

This workflow controls the dwi preprocessing stages using FSL tools.

I couldn’t get this to work reliably unless everything was oriented in LAS+ before going to TOPUP and eddy. For this reason, if TOPUP is going to be used (for an epi fieldmap or an RPE series) or there is no fieldmap correction, operations occurring before eddy are done in LAS+. The fieldcoefs are applied during eddy’s run and the corrected series comes out. This is finally converted to LPS+ and sent to the rest of the pipeline.

If a GRE fieldmap is available, the correction is applied to eddy’s outputs after they have been converted back to LPS+.

Finally, if SyN is chosen, it is applied to the LPS+ converted, eddy-resampled data.

Parameters

scan_groups: dict

dictionary with fieldmaps and warp space information for the dwis

impute_slice_threshold: float

threshold for a slice to be replaced with imputed values. Overrides the parameter in eddy_config if set to a number > 0.

pepolar_methodstr

Either ‘DRBUDDI’, ‘TOPUP’ or ‘DRBUDDI+TOPUP’. The method for SDC when EPI fieldmaps are used.

eddy_config: str

Path to a JSON file containing settings for the call to eddy.

Inputs

dwi_file: str

DWI series. Possibly concatenated, denoised, etc

bvec_file: str

bvec file

bval_file: str

bval file

b0_indices: list

Indexes into dwi_files that correspond to b=0 volumes

b0_images: list

List of single b=0 volumes

original_files: list

List of the files from which each DWI volume came. One per original file

t1_brain: str

Skull stripped T1w image

t1_mask: str

mask for t1_brain

Resampling workflows

qsiprep.workflows.dwi.resampling.init_dwi_trans_wf(source_file, template, mem_gb, omp_nthreads, output_resolution, name='dwi_trans_wf', use_compression=True, to_mni=False, write_local_bvecs=False, write_reports=True, concatenate=True)[source]

This workflow samples dwi images to the output_grid in a “single shot” from the original DWI series.

../_images/index-6.png

(Source code, png, svg, pdf)

Parameters

templatestr

Name of template targeted by template output space

mem_gbfloat

Size of DWI file in GB

omp_nthreadsint

Maximum number of threads an individual process may use

namestr

Name of workflow (default: dwi_trans_wf)

use_compressionbool

Save registered DWI series as .nii.gz

use_fieldwarpbool

Include SDC warp in single-shot transform from DWI to MNI

output_resolutionfloat

Voxel size in mm for the output data

to_mnibool

Include warps to MNI

write_local_bvecsbool

if true, local bvec niftis are written

Inputs

itk_b0_to_t1

Affine transform from ref_bold_brain to T1 space (ITK format)

t1_2_mni_forward_transform

ANTs-compatible affine-and-warp transform file

dwi_files

Individual 3D volumes, not motion corrected

cnr_map

Contrast to noise map from model-based hmc

bval_files

individual bval files

bvec_files

one-lined bvec files

b0_ref_image

b0 template for the dwi series

b0_indices

List of indices that contain a b0 image

dwi_mask

Skull-stripping mask of reference image

name_source

DWI series NIfTI file Used to recover original information lost during processing

hmc_xforms

List of affine transforms aligning each volume to ref_image in ITK format

fieldwarps

a DFM in ITK format

output_grid

File defining the output space

t1_mask

Brain mask from the t1w

Outputs

dwi_resampled

DWI series, resampled to template space. One file if concatenate, otherwise a list of files

dwi_ref_resampled

Reference, contrast-enhanced summary of the DWI series, resampled to template space

dwi_mask_resampled

DWI series mask in template space

cnr_map_resampled

Contrast to noise map resampled

bvals

bvals file for the DWI series

rotated_bvecs

bvecs rotated for transforms to output_grid

local_bvecs

NIfTI file containing the bvec rotation matrix (due to transforms) in each voxel. Includes rotations introduced by warpingdenoisin

Calculate dwi confounds

qsiprep.workflows.dwi.confounds.init_dwi_confs_wf(mem_gb, metadata, impute_slice_threshold, name='dwi_confs_wf')[source]

This workflow calculates confounds for a dwi series, and aggregates them into a TSV file, for use as nuisance regressors in a GLM.

The following confounds are calculated, with column headings in parentheses:

  1. Framewise displacement, based on head-motion parameters (framewise_displacement)

  2. Estimated head-motion parameters, in mm and rad (trans_x, trans_y, trans_z, rot_x, rot_y, rot_z)

../_images/index-7.png

(Source code, png, svg, pdf)

Parameters

mem_gbfloat

Size of dwi file in GB - please note that this size should be calculated after resamplings that may extend the FoV

metadatadict

BIDS metadata for dwi file

namestr

Name of workflow (default: dwi_confs_wf)

Inputs

sliceqc_file

dwi image, after the prescribed corrections (STC, HMC and SDC) when available.

motion_params

spm motion params

Outputs

confounds_file

TSV of all aggregated confounds

rois_report

Reportlet visualizing white-matter/CSF mask used for aCompCor, the ROI for tCompCor and the dwi brain mask.

Fieldmap estimation and unwarping workflows

Automatic selection of the appropriate SDC method

If the dataset metadata indicate tha more than one field map acquisition is IntendedFor (see BIDS Specification section 8.9) the following priority will be used:

Table of behavior (fieldmap use-cases):

Fieldmaps found

use_syn

force_syn

Action

True

True

Fieldmaps + SyN

True

False

Fieldmaps

False

True

SyN

False

True

False

SyN

False

False

False

HMC only

qsiprep.workflows.fieldmap.base.init_sdc_wf(fieldmap_info, dwi_meta, omp_nthreads=1, debug=False, fmap_bspline=False, fmap_demean=True)[source]

This workflow implements the heuristics to choose a SDC strategy. When no field map information is present within the BIDS inputs, the EXPERIMENTAL “fieldmap-less SyN” can be performed, using the --use-syn argument. When --force-syn is specified, then the “fieldmap-less SyN” is always executed and reported despite of other fieldmaps available with higher priority. In the latter case (some sort of fieldmap(s) is available and --force-syn is requested), then the SDC method applied is that with the highest priority.

../_images/index-8.png

(Source code, png, svg, pdf)

Parameters

fmapslist of pybids dicts

A list of dictionaries with the available fieldmaps (and their metadata using the key 'metadata' for the case of epi fieldmaps)

dwi_metadict

BIDS metadata dictionary corresponding to the DWI run

omp_nthreadsint

Maximum number of threads an individual process may use

fmap_bsplinebool

Experimental: Fit B-Spline field using least-squares

fmap_demeanbool

Demean voxel-shift map during unwarp

debugbool

Enable debugging outputs

Inputs
b0_ref

A b0 reference calculated at a previous stage

b0_ref_brain

Same as above, but brain-masked

b0_mask

Brain mask for the DWI run

t1_brain

T1w image, brain-masked, for the fieldmap-less SyN method

t1_2_mni_reverse_transform

MNI-to-T1w transform to map prior knowledge to the T1w fo the fieldmap-less SyN method

templatestr

Name of template targeted by template output space

Outputs
b0_ref

An unwarped b0 reference

b0_mask

The corresponding new mask after unwarping

out_warp

The deformation field to unwarp the susceptibility distortions

syn_b0_ref

If --force-syn, an unwarped b0 reference with this method (for reporting purposes)

method

Name of the method used for SDC

fieldmap_hz

The fieldmap in Hz for eddy

Phase Encoding POLARity (PEPOLAR) techniques

qsiprep.workflows.fieldmap.pepolar.init_extended_pepolar_report_wf(segment_t2w, omp_nthreads=1, name='extended_pepolar_report_wf')[source]
qsiprep.workflows.fieldmap.pepolar.init_pepolar_unwarp_wf(dwi_meta, epi_fmaps, omp_nthreads=1, name='pepolar_unwarp_wf')[source]

This workflow takes in a set of EPI files with opposite phase encoding direction than the target file and calculates a displacements field (in other words, an ANTs-compatible warp file).

This procedure works if there is only one ‘_epi’ file is present (as long as it has the opposite phase encoding direction to the target file). The target file will be used to estimate the field distortion. However, if there is another ‘_epi’ file present with a matching phase encoding direction to the target it will be used instead.

Currently, different phase encoding dimension in the target file and the ‘_epi’ file(s) (for example ‘i’ and ‘j’) is not supported.

The warp field correcting for the distortions is estimated using AFNI’s 3dQwarp, with displacement estimation limited to the target file phase encoding direction.

It also calculates a new mask for the input dataset that takes into account the distortions.

../_images/index-9.png

(Source code, png, svg, pdf)

Inputs

in_reference

the reference image

in_reference_brain

the reference image skullstripped

in_mask

a brain mask corresponding to in_reference

Outputs

out_reference

the in_reference after unwarping

out_warp

the corresponding DFM compatible with ANTs

qsiprep.workflows.fieldmap.pepolar.init_prepare_dwi_epi_wf(omp_nthreads, orientation='LPS', name='prepare_epi_wf')[source]

This workflow takes in a set of dwi files with with the same phase encoding direction and returns a single 3D volume ready to be used in field distortion estimation. It removes b>0 volumes.

The procedure involves: estimating a robust template using FreeSurfer’s ‘mri_robust_template’, bias field correction using ANTs N4BiasFieldCorrection and AFNI 3dUnifize, skullstripping using FSL BET and AFNI 3dAutomask, and rigid coregistration to the reference using ANTs.

Phase-difference B0 estimation

The field inhomogeneity inside the scanner (fieldmap) is proportional to the phase drift between two subsequent GRE sequence.

Fieldmap preprocessing workflow for fieldmap data structure 8.9.1 in BIDS 1.0.0: one phase diff and at least one magnitude image 8.9.2 in BIDS 1.0.0: two phases and at least one magnitude image

qsiprep.workflows.fieldmap.phdiff.init_phdiff_wf(omp_nthreads, phasetype='phasediff', name='phdiff_wf')[source]

Estimates the fieldmap using a phase-difference image and one or more magnitude images corresponding to two or more GRE acquisitions. The original code was taken from nipype.

(Source code)

Outputs:

outputnode.fmap_ref - The average magnitude image, skull-stripped
outputnode.fmap_mask - The brain mask applied to the fieldmap
outputnode.fmap - The estimated fieldmap in Hz

Direct B0 mapping sequences

When the fieldmap is directly measured with a prescribed sequence (such as SE), we only need to calculate the corresponding B-Spline coefficients to adapt the fieldmap to the TOPUP tool. This procedure is described with more detail here.

This corresponds to the section 8.9.3 –fieldmap image (and one magnitude image)– of the BIDS specification.

qsiprep.workflows.fieldmap.fmap.init_fmap_wf(omp_nthreads, fmap_bspline, name='fmap_wf')[source]

Fieldmap workflow - when we have a sequence that directly measures the fieldmap we just need to mask it (using the corresponding magnitude image) to remove the noise in the surrounding air region, and ensure that units are Hz.

(Source code)

Fieldmap-less estimation (experimental)

In the absence of direct measurements of fieldmap data, we provide an (experimental) option to estimate the susceptibility distortion based on the ANTs symmetric normalization (SyN) technique. This feature may be enabled, using the --use-syn-sdc flag, and will only be applied if fieldmaps are unavailable.

During the evaluation phase, the --force-syn flag will cause this estimation to be performed in addition to fieldmap-based estimation, to permit the direct comparison of the results of each technique. Note that, even if --force-syn is given, the functional outputs of qsiprep will be corrected using the fieldmap-based estimates.

Feedback will be enthusiastically received.

class qsiprep.workflows.fieldmap.syn.ThreshAndBin(from_file=None, resource_monitor=None, **inputs)[source]

Bases: SimpleInterface

Mandatory Inputs:

in_file (a pathlike object or string representing a file) – File to modify, atlas.

Optional Inputs:

atlas_threshold (an integer) – Threshold value.

Outputs:

out_file (a pathlike object or string representing an existing file) – File with threshold and binarization applied.

qsiprep.workflows.fieldmap.syn.init_syn_sdc_wf(omp_nthreads, bold_pe=None, atlas_threshold=2, name='syn_sdc_wf')[source]

This workflow takes a skull-stripped T1w image and reference b0 image and estimates a susceptibility distortion correction warp, using ANTs symmetric normalization (SyN) and the average fieldmap atlas described in [Treiber2016].

SyN deformation is restricted to the phase-encoding (PE) direction. If no PE direction is specified, anterior-posterior PE is assumed.

SyN deformation is also restricted to regions that are expected to have a >2mm (approximately 1 voxel) warp, based on the fieldmap atlas.

This technique is a variation on those developed in [Huntenburg2014] and [Wang2017].

../_images/index-12.png

(Source code, png, svg, pdf)

Inputs

b0_ref

reference image

templatestr

Name of template targeted by template output space

t1_brain

skull-stripped, bias-corrected structural image

t1_2_mni_reverse_transform

inverse registration transform of T1w image to MNI template

Outputs

out_reference

the bold_ref image after unwarping

out_reference_brain

the bold_ref_brain image after unwarping

out_warp

the corresponding DFM compatible with ANTs

out_mask

mask of the unwarped input file

Unwarping

qsiprep.workflows.fieldmap.unwarp.init_fmap_unwarp_report_wf(name='fmap_unwarp_report_wf', suffix='hmcsdc')[source]

This workflow generates and saves a reportlet showing the effect of fieldmap unwarping a DWI image.

../_images/index-13.png

(Source code, png, svg, pdf)

Parameters

namestr, optional

Workflow name (default: fmap_unwarp_report_wf)

suffixstr, optional

Suffix to be appended to this reportlet

Inputs

in_pre

Reference image, before unwarping

in_post

Reference image, after unwarping

in_seg

Segmentation of preprocessed structural image, including gray-matter (GM), white-matter (WM) and cerebrospinal fluid (CSF)

in_xfm

Affine transform from T1 space to b0 space (ITK format)

qsiprep.workflows.fieldmap.unwarp.init_sdc_unwarp_wf(omp_nthreads, fmap_demean, debug, name='sdc_unwarp_wf')[source]

This workflow takes in a displacements fieldmap and calculates the corresponding displacements field (in other words, an ANTs-compatible warp file).

It also calculates a new mask for the input dataset that takes into account the distortions. The mask is restricted to the field of view of the fieldmap since outside of it corrections could not be performed.

Inputs

in_reference

the reference image

in_reference_brain

the reference image (skull-stripped)

in_mask

a brain mask corresponding to in_reference

metadata

metadata associated to the in_reference EPI input

fmap

the fieldmap in Hz

fmap_ref

the reference (anatomical) image corresponding to fmap

fmap_mask

a brain mask corresponding to fmap

Outputs

out_reference

the in_reference after unwarping

out_reference_brain

the in_reference after unwarping and skullstripping

out_warp

the corresponding DFM compatible with ANTs

out_jacobian

the jacobian of the field (for drop-out alleviation)

out_mask

mask of the unwarped input file

out_hz

fieldmap in Hz that can be sent to eddy

Reconstruction Workflows

..automodule:: qsiprep.workflows.recon