Sisyphe.core.sisypheStatistics
External packages/modules
Cython, static compiler, https://cython.org/
Numpy, scientific computing, https://numpy.org/
Scipy, scientific computing, https://docs.scipy.org/doc/scipy/index.html
SimpleITK, medical image processing, https://simpleitk.org/
- class Sisyphe.core.sisypheStatistics.SisypheDesign
SisypheDesign class
Description
PySisyphe statistical design class.
Model definition class for voxel by voxel statistical parametric mapping analysis.
Reference: Statistical parametric maps in functional imaging: A general linear approach. KJ Friston, AP Holmes, KJ Worsley, JP Poline, CD Frith, RSJ Frackowiak. Human Brain Mapping 1995;2(4):189-210. doi: 10.1002/hbm.460020402.
Analysis of fmri time series revisited. KJ Friston, AP Holmes, JB Poline, PJ Grasby, SCR Williams, RSJ Frackowiak, R Turner. Neuroimage 1995 Mar;2(1):45-53. doi: 10.1006/nimg.1995.1007.
Inheritance
object -> SisypheDesign
Creation: 29/11/2023 Last revision: 06/12/2024
- addCovariableByCondition(name: str, cov: ndarray, estimable: bool = False, zscore: bool = False, logt: bool = False) None
Add a covariable by condition to the model of the current SisypheDesign instance.
Parameters
- namestr
covariable name
- cov: ndarray
covariable vector, a value for each model observation
- estimablebool
covariable estimability, confound variable if False (default False)
- zscorebool
convert covariable values into z-score (i.e. (value - mean) / std, default False)
- logtbool
apply log to covariable values (default False)
- addCovariableByGroup(name: str, cov: ndarray, estimable: bool = False, zscore: bool = False, logt: bool = False) None
Add a covariable by group to the model of the current SisypheDesign instance.
Parameters
- namestr
covariable name
- cov: ndarray
covariable vector, a value for each model observation
- estimablebool
covariable estimability, confound variable if False (default False)
- zscorebool
convert covariable values into z-score (i.e. (value - mean) / std, default False)
- logtbool
apply log to covariable values (default False)
- addCovariableBySubject(name: str, cov: ndarray, estimable: bool = False, zscore: bool = False, logt: bool = False) None
Add a covariable by subject to the model of the current SisypheDesign instance.
Parameters
- namestr
covariable name
- cov: ndarray
covariable vector, a value for each model observation
- estimablebool
covariable estimability, confound variable if False (default False)
- zscorebool
convert covariable values into z-score (i.e. (value - mean) / std, default False)
- logtbool
apply log to covariable values (default False)
- addGlobalCovariable(name: str, cov: ndarray, estimable: bool = False, zscore: bool = False, logt: bool = False) None
Add a global covariable to the model of the current SisypheDesign instance.
Parameters
- namestr
covariable name
- cov: ndarray
covariable vector, a value for each model observation
- estimablebool
covariable estimability, confound variable if False (default False)
- zscorebool
convert covariable values into z-score (i.e. (value - mean) / std, default False)
- logtbool
apply log to covariable values (default False)
- addHRFBoxCarModelToCondition(cond: str) None
Add a BoxCar (with HRF convolution) covariable to the fMRI model of the current SisypheDesign instance.
Parameters
- condstr
fMRI condition name
- addHighPassToCondition(cond: str) None
Add high pass filter covariables to the fMRI model of the current SisypheDesign instance.
Parameters
- condstr
fMRI condition name
- appendFileObsTo(obs: list[str], group: str | None = None, subj: str | None = None, cond: str | None = None) None
Append observations to the model of the current SisypheDesign instance. Observations are SisypheVolume filenames.
Parameters
- obs: list[str]
list of SisypheVolume filenames (observations)
- groupstr | None
group name (default None)
- subjstr | None
subject name (default None)
- condstr | None
condition name (default None)
- classmethod calcMask(obs: list[SisypheVolume]) SisypheVolume
Analysis mask processing
Parameters
- obslist[Sisyphe.core.sisypheVolume.SisypheVolume]
observations
Returns
- Sisyphe.core.sisypheVolume.SisypheVolume
mask
- clearFileObsFrom(group: str | None = None, subj: str | None = None, cond: str | None = None) None
clear observations from the model of the current SisypheDesign instance. Observations are SisypheVolume filenames.
Parameters
- groupstr | None
group name (default None)
- subjstr | None
subject name (default None)
- condstr | None
condition name (default None)
- createXML(doc: Document) None
Write the current SisypheDesign instance attributes to xml instance. This method is called by save() and saveAs() methods, not recommended for use.
Parameters
- docminidom.Document
xml document
- estimate(mask: SisypheROI | SisypheVolume | None = None, roi: SisypheROI | SisypheVolume | None = None, wait: DialogWait | None = None) None
Model estimation of the current SisypheDesign instance. beta, pooled variance and autocorrelations are processed.
Y vector of observations, X design matrix beta = pinv(X) Y = (XtX)-1 Xt Y residuals i.e. model errors = Y - X beta pooled variance = variance(residuals)
Reference: Statistical parametric maps in functional imaging: A general linear approach. KJ Friston, AP Holmes, KJ Worsley, JP Poline, CD Frith, RSJ Frackowiak. Human Brain Mapping 1995;2(4):189-210. doi: 10.1002/hbm.460020402.
Analysis of fmri time series revisited. KJ Friston, AP Holmes, JB Poline, PJ Grasby, SCR Williams, RSJ Frackowiak, R Turner. Neuroimage 1995 Mar;2(1):45-53. doi: 10.1006/nimg.1995.1007.
Parameters
- maskSisyphe.core.SisypheROI.sisypheROI | Sisyphe.core.sisypheVolume.SisypheVolume
statistical analysis mask
- roiSisyphe.core.SisypheROI.sisypheROI | Sisyphe.core.sisypheVolume.SisypheVolume
mask used as reference for signal normalization
- waitSisyphe.gui.dialogWait.DialogWait
progress dialog
- getAgeCovariable() int
Get the age covariable attribute code of the current SisypheDesign instance.
Returns
- int
- Age covariable attribute code:
0: no age covariable,
1: age global covariable,
2: age covariable by group,
3: age covariable by subject,
4: age covariable by condition.
- getAgeCovariableAsString() str
Get the age covariable attribute of the current SisypheDesign instance as string.
Returns
- int
- Age covariable attribute code:
‘no’, no age covariable,
‘global’, age global covariable,
‘by group’, age covariable by group,
‘by subject’, age covariable by subject,
‘by condition’, age covariable by condition.
- getAutocorrelations() list[float] | None
Get autocorrelations of the current SisypheDesign instance.
Reference: Robust smoothness estimation in statistical parametric maps using standardized residuals from the general linear model. SJ Kiebel, JB Poline, KJ Friston, AP Holmes, KJ Worsley. Neuroimage 1999 Dec;10(6):756-66. doi: 10.1006/nimg.1999.0508.
Returns
- list[float]
autocorrelations in x, y and z directions (gaussian point spread function, fwhm in mm)
- getBasename() str
Get the base name attribute of the current SisypheDesign instance. The file name attribute is used to save the current SisypheDesign instance.
Returns
- str
base name (base part of the file name)
- getBeta() SisypheVolume | None
Get the beta volume of the current SisypheDesign instance.
Returns
- Sisyphe.core.sisypheVolume.SisypheVolume
beta
- getConditionCount() int
Get the number of conditions in the model of the current SisypheDesign instance.
Returns
- int
number of conditions
- getConditionName(index: int) str
Get a condition name from the model of the current SisypheDesign instance.
Parameters
- indexint
condition index
Returns
- str
condition name
- getConditionNames() list[str]
Get names of model conditions of the model of the current SisypheDesign instance.
Returns
- list[str]
list of condition names
- getDegreesOfFreedom() int
Get the degrees of freedom of the current SisypheDesign instance.
Returns
- list[float, float, float]
autocorrelations in x, y and directions (mm)
- getDesignMatrix(recalc: bool = False) ndarray
Get the design matrix of the current SisypheDesign instance.
Returns
- ndarray
design matrix X of the statistical model, lines = observations x columns = effects/covariables
- getDirname() str
Get the path name attribute of the current SisypheDesign instance. The file name attribute is used to save the current SisypheDesign instance.
Returns
- str
path name (path part of the file name)
- getEffectInformations() list[tuple[str, int]]
Get the model dictionary of the current SisypheDesign instance.
Returns
- list[tuple[str, int]]
str, column name of the design matrix (covariable/effect name)
- int, covariable/effect estimability
0 not estimable
1 estimable, main effect
2 estimable, global covariable of interest
3 estimable, covariable by group
4 estimable, covariable by subject
5 estimable, covariable by condition
- getFileObsCount(group: str | None = None, subj: str | None = None, cond: str | None = None) int
Get the numer of observations (SisypheVolume) in a given group/subject/condition of the model of the current SisypheDesign instance.
Parameters
- groupstr | None
group name (default None)
- subjstr | None
subject name (default None)
- condstr | None
condition name (default None)
Returns
- int
number of observations in a given group/subject/condition of the model
- getFileObsDict() dict[tuple[str, str, str], list[str]]
Get all model observations of the current SisypheDesign instance. Observations are SisypheVolume filenames.
Returns
- dict[tuple[str, str, str], list[str]]
key, tuple[str, str, str] as (group name, subject name, condition name)
value, list[str], observation filenames
- getFileObsFrom(group: str | None = None, subj: str | None = None, cond: str | None = None) list[str]
Get model observations of the current SisypheDesign instance. Observations are SisypheVolume filenames.
Parameters
- groupstr | None
group name (default None)
- subjstr | None
subject name (default None)
- condstr | None
condition name (default None)
Returns
- list[str]
list of SisypheVolume filenames (observations)
- getFilename() str
Get the file name attribute of the current SisypheDesign instance. The file name attribute is used to save the current SisypheDesign instance.
Returns
- str
file name
- classmethod getFilterExt() str
Get statistical design filter used by QFileDialog.getOpenFileName() and QFileDialog.getSaveFileName().
Returns
- str
‘PySisyphe statistical model (*.xmodel)’
- getGlobalFactorIndex() int | None
Get the column index of the global factor (column of ones) in the design matrix of the current SisypheDesign instance.
Returns
- int | None
column index of the global factor, or None if no global factor in the design matrix
- getGroupCount() int
Get the number of groups in the model of the current SisypheDesign instance.
Returns
- int
number of groups
- getGroupName(index: int) str
Get a group name from the model of the current SisypheDesign instance.
Parameters
- indexint
group index
Returns
- str
group name
- getGroupNames() list[str]
Get names of model groups of the current SisypheDesign instance.
Returns
- list[str]
list of group names
- getMeanObservations() SisypheVolume | None
Get the mean volume of observations of the current SisypheDesign instance.
Returns
- Sisyphe.core.sisypheVolume.SisypheVolume
mean volume of observations
- getNormalizationAsString() str
Get the signal normalization attribute of the current SisypheDesign instance as string.
Returns
- str
- Signal normalization code:
‘No’ signal normalization,
‘Mean scaling’, signal divided by the mean in the analysis mask
‘Median scaling’, signal divided by the median in the analysis mask
‘75th perc. scaling’, signal divided by the 75th percentile in the analysis mask
‘ROI mean scaling’, signal divided by the mean in a reference ROI
‘ROI median scaling’, signal divided by the median in a reference ROI
‘ROI 75th perc. scaling’, scaling, signal divided by the 75th percentile in a reference ROI
‘ANCOVA Mean’, add covariable of the mean in the analysis mask
‘ANCOVA Median’, add covariable of the median in the analysis mask
‘ANCOVA 75th perc.’, add covariable of the 75th percentile in the analysis mask
‘ANCOVA ROI mean’, add covariable of the mean in a reference ROI
‘ANCOVA ROI median’, add covariable of the median in a reference ROI
‘ANCOVA ROI 75th perc.’, add covariable of the 75th percentile in a reference ROI
- getObsCount(group: str | None = None, subj: str | None = None, cond: str | None = None) int
Get the numer of observations (SisypheVolume) required in a given group/subject/condition of the model of the current SisypheDesign instance. Observations are SisypheVolume filenames.
Parameters
- groupstr | None
group name (default None)
- subjstr | None
subject name (default None)
- condstr | None
condition name (default None)
Returns
- int
number of observations in a given group/subject/condition of the model
- getObservations() SisypheVolume | None
Get the observations multicomponent volume of the current SisypheDesign instance.
Returns
- Sisyphe.core.sisypheVolume.SisypheVolume
observations multicomponent volume
- getObservationsDict() dict[tuple[str, str, str], list[int]]
Get the observation dictionary of the current SisypheDesign instance.
Returns
- dict[tuple[str, str, str], list[int]]
- Observation dictionary:
tuple[str, str, str] key as (group name, subject name, condition name)
list[int] value, first element (index 0) = number of observations
- getPooledVariance() SisypheVolume | None
Get the pooled variance volume of the current SisypheDesign instance.
Returns
- Sisyphe.core.sisypheVolume.SisypheVolume
pooled variance
- getSignalCovariable() int
Get the type of ANCOVA covariable for signal normalization of the current SisypheDesign instance.
Returns
- int
type of ANCOVA covariable: 0 ANCOVA global, 1 by group, 2 by subject, 3 by condition
- getSignalCovariableAsString() str
Get the type of ANCOVA covariable for signal normalization of the current SisypheDesign instance.
Returns
- str
type of ANCOVA covariable: ‘Global’, ‘by group’, ‘by subject’,’by condition’
- getSignalNormalization() int
Get the signal normalization attribute of the current SisypheDesign instance.
Returns
- int
- Signal normalization code:
0: No signal normalization,
1: Mean scaling, signal divided by the mean in the analysis mask
2: Median scaling, signal divided by the median in the analysis mask
3: 75th percentile scaling, signal divided by the 75th percentile in the analysis mask
4: ROI mean scaling, signal divided by the mean in a reference ROI
5: ROI median scaling, signal divided by the median in a reference ROI
6: ROI 75th perc. scaling, scaling, signal divided by the 75th percentile in a reference ROI
7: ANCOVA Mean, add covariable of the mean in the analysis mask
8: ANCOVA Median, add covariable of the median in the analysis mask
9: ANCOVA 75th percentile, add covariable of the 75th percentile in the analysis mask
10: ANCOVA ROI mean, add covariable of the mean in a reference ROI
11: ANCOVA ROI median, add covariable of the median in a reference ROI
12: ANCOVA ROI 75th percentile, add covariable of the 75th percentile in a reference ROI
- getSubjectCount() int
Get the number of subjects in the model of the current SisypheDesign instance.
Returns
- int
number of subjects
- getSubjectName(index: int) str
Get a subject name from the model of the current SisypheDesign instance.
Parameters
- indexint
subject index
Returns
- str
subject name
- getSubjectNames() list[str]
Get names of model subjects of the model of the current SisypheDesign instance.
Returns
- list[str]
list of subject names
- getTotalFileObsCount() int
Get the numer of observations (SisypheVolume) in the model of the current SisypheDesign instance.
Returns
- int
total number of observations in the model
- getTotalObsCount() int
Get the number of observations (SisypheVolume) required in the model of the current SisypheDesign instance.
Returns
- int
total number of observations in the model
- hasAllFileObs() bool
Check that there are no missing observation volumes in the model of the current SisypheDesign instance.
Returns
- bool
True if no missing observation volumes
- hasCondition() bool
Check whether the model of the current SisypheDesign instance has condition(s).
Returns
- bool
True if the model has condition(s)
- hasDesignMatrix() bool
Check whether the design matrix of the current SisypheDesign instance is defined.
Returns
- bool
True if the design matrix is defined
- hasFilename() bool
Check whether the filename attribute of the current SisypheDesign instance is defined (i.e. not ‘’)
Returns
- bool
True if the filename attribute is defined, False otherwise
- hasGlobalFactor() bool
Check if the current SisypheDesign instance has a global factor (column of ones in the design matrix).
Returns
- bool
True if design matrix has a global factor
- hasGroup() bool
Check whether the model of the current SisypheDesign instance has group(s).
Returns
- bool
True if the model has group(s)
- hasModel() bool
Check whether the model of the current SisypheDesign instance is defined.
Returns
- bool
True if the model is defined
- hasSubject() bool
Check whether the model of the current SisypheDesign instance has subject(s).
Returns
- bool
True if the model has subject(s)
- isEstimated() bool
Check whether the current SisypheDesign instance is estimated.
Returns
- bool
True if design is estimated
- isfMRIDesign() bool
Check whether the design type of the current SisypheDesign instance is fMRI.
Returns
- bool
True if fMRI design
- load(filename: str, wait: DialogWait | None = None) None
Load the current SisypheDesign instance from a PySisyphe statistical design (.xmodel) file.
Parameters
- filenamestr
PySisyphe statistical design file name
- wait: DialogWait | None
progress dialog (optional)
- parseXML(doc: Document) None
Read the current SisypheDesign instance attributes from xml instance. This method is called by load() method, not recommended for use.
Parameters
- docminidom.Document
xml document
- save(wait: DialogWait | None = None) None
Save the current SisypheDesign instance to a PySisyphe statistical design (.xmodel) file. The file name attribute of the current SisypheDesign instance is used.
Parameters
- wait: DialogWait | None
progress dialog (optional)
- saveAs(filename: str) None
Save the current SisypheDesign instance to a PySisyphe statistical design (.xmodel) file.
Parameters
- filenamestr
PySisyphe statistical design file name
- setAgeCovariableByCondition() None
Add age confound covariable by condition to the model of the current SisypheDesign instance.
- setAgeCovariableByGroup() None
Add age confound covariable by group to the model of the current SisypheDesign instance.
- setAgeCovariableBySubject() None
Add age confound covariable by subject to the model of the current SisypheDesign instance.
- setAgeGlobalCovariable() None
Add age global confound covariable to the model of the current SisypheDesign instance.
- setDictDesign(design: dict) None
Set the observation dictionary of the current SisypheDesign instance.
Parameters
- designdict
Dictionary syntax with 3 levels, key tuple(str = effect name, int = tree level) 3 possible levels: 0=first level (group), 1=second level (subject), 2=third level (condition)
- {(‘Group#1’, 0): {
- (‘Subject#11’, 1): {
(‘Condition#1’, 2): [number of observations] , … , (‘Condition#k’, 2): [number of observations]}
, … , (‘Subject#1j’, 1): {
(‘Condition#1’, 2): [number of observations] , … , (‘Condition#k’, 2): [number of observations]}
}
, … , (‘Group#i’, 0): {
- (‘Subject#i1’, 1): {
(‘Condition#1’, 2): [number of observations] , … , (‘Condition#k’, 2): [number of observations]}
, … , (‘Subject#ij’, 1): {
(‘Condition#1’, 2): [number of observations] , … , (‘Condition#k’, 2): [number of observations]}
}
}
Type of models :
i groups (one-sample t-test, two-sample t-test, ANOVA), no subject, no condition
i conditions, no group, no subject
i subjects, j conditions, no group
i groups, j subjects, no condition
i groups, j subjects, k conditions
- setFileObsTo(obs: list[str], group: str | None = None, subj: str | None = None, cond: str | None = None) None
Set observations to the model of the current SisypheDesign instance. Observations are SisypheVolume filenames.
Parameters
- obs: list[str]
list of SisypheVolume filenames (observations)
- groupstr | None
group name (default None)
- subjstr | None
subject name (default None)
- condstr | None
condition name (default None)
- setFilename(filename: str) None
Set the file name attribute of the current SisypheDesign instance. The file name attribute is used to save the current SisypheDesign instance.
Parameters
- filenamestr
file name
- setNoAgeCovariable() None
Remove age covariable from the model of the current SisypheDesign instance.
- setSignalCovariableByCondition() None
Set the type of ANCOVA covariable for signal normalization of the current SisypheDesign instance. Covariable by condition, one column for each condition in the design matrix.
- setSignalCovariableByGroup() None
Set the type of ANCOVA covariable for signal normalization of the current SisypheDesign instance. Covariable by group, one column for each group in the design matrix.
- setSignalCovariableBySubject() None
Set the type of ANCOVA covariable for signal normalization of the current SisypheDesign instance. Covariable by subject, one column for each subject in the design matrix.
- setSignalGlobalCovariable() None
Set the type of ANCOVA covariable for signal normalization of the current SisypheDesign instance. Global covariable, a single column in the design matrix.
- setSignalNormalization(v: int = 0) None
Set the signal normalization attribute of the current SisypheDesign instance.
Parameters
- vint
- Signal normalization code:
0: No signal normalization,
1: mean scaling, signal divided by the mean in the analysis mask
2: Median scaling, signal divided by the median in the analysis mask
3: 75th percentile scaling, signal divided by the 75th percentile in the analysis mask
4: ROI mean scaling, signal divided by the mean in a reference ROI
5: ROI median scaling, signal divided by the median in a reference ROI
6: ROI 75th perc. scaling, scaling, signal divided by the 75th percentile in a reference ROI
7: ANCOVA Mean, add covariable of the mean in the analysis mask
8: ANCOVA Median, add covariable of the median in the analysis mask
9: ANCOVA 75th percentile, add covariable of the 75th percentile in the analysis mask
10: ANCOVA ROI mean, add covariable of the mean in a reference ROI
11: ANCOVA ROI median, add covariable of the median in a reference ROI
12: ANCOVA ROI 75th percentile, add covariable of the 75th percentile in a reference ROI
- setfMRIDesign() None
Set the design type of the current SisypheDesign instance to fMRI. The current SisypheDesign instance is used to process fMRI.
- Sisyphe.core.sisypheStatistics.autocorrelationsEstimate(error: ndarray, design: ndarray, mask: SisypheROI | SisypheVolume, wait: DialogWait | None = None) tuple[float, float, float]
Estimating autocorrelations of residuals.
Reference: Robust smoothness estimation in statistical parametric maps using standardized residuals from the general linear model. SJ Kiebel, JB Poline, KJ Friston, AP Holmes, KJ Worsley. Neuroimage 1999 Dec;10(6):756-66. doi: 10.1006/nimg.1999.0508.
Parameters
- errornumpy.ndarray
residuals vector / standard deviation of residuals
- designnumpy.ndarray
design matrix X
- masksisyphe.core.sisypheROI.SisypheROI | sisyphe.core.sisypheVolume.SisypheVolume
statistical analysis mask
- waitSisyphe.gui.dialogWait.DialogWait
progress bar
Returns
- tuple[float, float, float]
autocorrelations x, y, z (fwhm in mm)
- Sisyphe.core.sisypheStatistics.clusterCorrectedpvalueToExtent(p: float, ev: float, ec: float) int
Extent processing from cluster corrected p-value.
reference: A unified statistical approach for determining significant signals in images of cerebral activation KJ Worsley, S Marrett, P Neelin, AC Vandal, KJ Friston, AC Evans. Hum Brain Mapp. 1996;4(1):58-73. doi: 10.1002/(SICI)1097-0193(1996)4:1<58::AID-HBM4>3.0.CO;2-O.
Parameters
- pfloat
cluster uncorrected p-value
- evfloat
expected number of voxels
- ecfloat
expected number of clusters
Returns
- int
extent (number of voxels)
- Sisyphe.core.sisypheStatistics.clusterUncorrectedpvalueToExtent(p: float, ev: float, ec: float) int
Extent processing from cluster uncorrected p-value.
reference: A unified statistical approach for determining significant signals in images of cerebral activation KJ Worsley, S Marrett, P Neelin, AC Vandal, KJ Friston, AC Evans. Hum Brain Mapp. 1996;4(1):58-73. doi: 10.1002/(SICI)1097-0193(1996)4:1<58::AID-HBM4>3.0.CO;2-O.
Parameters
- pfloat
cluster uncorrected p-value
- evfloat
expected number of voxels
- ecfloat
expected number of clusters
Returns
- int
extent (number of voxels)
- Sisyphe.core.sisypheStatistics.conjunctionFisher(maps: list[SisypheVolume] | SisypheVolumeCollection) SisypheVolume
Fisher’s method to combine t-maps (t-maps conjunction).
Reference: Combining brains: a survey of methods for statistical pooling of information. Lazar NA, Luna B, Sweeney JA, Eddy WF. Neuroimage 2002 Jun;16(2):538-50. doi: 10.1006/nimg.2002.1107.
Parameters
maps : list[Sisyphe.core.sisypheVolume.SisypheVolume] | Sisyphe.core.sisypheVolume.SisypheVolumeCollection
Returns
- Sisyphe.core.sisypheVolume.SisypheVolume
combined t-map
- Sisyphe.core.sisypheStatistics.conjunctionMudholkar(maps: list[SisypheVolume] | SisypheVolumeCollection) SisypheVolume
Mudholkar’s method to combine t-maps (t-maps conjunction).
Reference: Combining brains: a survey of methods for statistical pooling of information. Lazar NA, Luna B, Sweeney JA, Eddy WF. Neuroimage 2002;Jun;16(2):538-50. doi: 10.1006/nimg.2002.1107.
Parameters
maps : list[Sisyphe.core.sisypheVolume.SisypheVolume] | Sisyphe.core.sisypheVolume.SisypheVolumeCollection
Returns
- Sisyphe.core.sisypheVolume.SisypheVolume
combined t-map
- Sisyphe.core.sisypheStatistics.conjunctionStouffer(maps: list[SisypheVolume] | SisypheVolumeCollection) SisypheVolume
Stouffer’s method to combine t-maps (t-maps conjunction).
Reference: Combining brains: a survey of methods for statistical pooling of information. Lazar NA, Luna B, Sweeney JA, Eddy WF. Neuroimage 2002;Jun;16(2):538-50. doi: 10.1006/nimg.2002.1107.
Parameters
maps : list[Sisyphe.core.sisypheVolume.SisypheVolume] | Sisyphe.core.sisypheVolume.SisypheVolumeCollection
Returns
- Sisyphe.core.sisypheVolume.SisypheVolume
combined t-map
- Sisyphe.core.sisypheStatistics.conjunctionTippett(maps: list[SisypheVolume] | SisypheVolumeCollection) SisypheVolume
Tippett’s method to combine t-maps (t-maps conjunction).
Reference: Combining brains: a survey of methods for statistical pooling of information. Lazar NA, Luna B, Sweeney JA, Eddy WF. Neuroimage 2002;Jun;16(2):538-50. doi: 10.1006/nimg.2002.1107.
Parameters
maps : list[Sisyphe.core.sisypheVolume.SisypheVolume] | Sisyphe.core.sisypheVolume.SisypheVolumeCollection
Returns
- Sisyphe.core.sisypheVolume.SisypheVolume
combined t-map
- Sisyphe.core.sisypheStatistics.conjunctionWorsley(maps: list[SisypheVolume] | SisypheVolumeCollection) SisypheVolume
Worsley’s method to combine t-maps (t-maps conjunction).
Reference: Combining brains: a survey of methods for statistical pooling of information. Lazar NA, Luna B, Sweeney JA, Eddy WF. Neuroimage 2002;Jun;16(2):538-50. doi: 10.1006/nimg.2002.1107.
Parameters
maps : list[Sisyphe.core.sisypheVolume.SisypheVolume] | Sisyphe.core.sisypheVolume.SisypheVolumeCollection
Returns
- Sisyphe.core.sisypheVolume.SisypheVolume
combined t-map
- Sisyphe.core.sisypheStatistics.convolveModelHRF(model: ndarray, iscan: float) ndarray
Convolve the boxcar model with the Hemodynamic Response Function (HRF).
Parameters
- modelnumpy.ndarray
fMRI boxcar model
- iscanfloat
interscan interval (TR in seconds)
Returns
- numpy.ndarray
convolved boxcar model
- Sisyphe.core.sisypheStatistics.expectedVoxelsPerCluster(ev: float, ec: float) float
Expected voxels per cluster processing.
reference: A unified statistical approach for determining significant signals in images of cerebral activation KJ Worsley, S Marrett, P Neelin, AC Vandal, KJ Friston, AC Evans. Hum Brain Mapp. 1996;4(1):58-73. doi: 10.1002/(SICI)1097-0193(1996)4:1<58::AID-HBM4>3.0.CO;2-O.
Parameters
- evfloat
expected number of voxels
- ecfloat
expected number of clusters
Returns
- float
expected voxels per cluster
- Sisyphe.core.sisypheStatistics.extentToClusterCorrectedpvalue(k: int, ev: float, ec: float) float
Cluster corrected p-value processing from extent (number of voxels).
reference: A unified statistical approach for determining significant signals in images of cerebral activation KJ Worsley, S Marrett, P Neelin, AC Vandal, KJ Friston, AC Evans. Hum Brain Mapp. 1996;4(1):58-73. doi: 10.1002/(SICI)1097-0193(1996)4:1<58::AID-HBM4>3.0.CO;2-O.
Parameters
- kint
extent, number of voxels
- evfloat
expected number of voxels
- ecfloat
expected number of clusters
Returns
- float
cluster corrected p-value
- Sisyphe.core.sisypheStatistics.extentToClusterUncorrectedpvalue(k: int, ev: float, ec: float) float
Cluster uncorrected p-value processing from extent (number of voxels).
reference: A unified statistical approach for determining significant signals in images of cerebral activation KJ Worsley, S Marrett, P Neelin, AC Vandal, KJ Friston, AC Evans. Hum Brain Mapp. 1996;4(1):58-73. doi: 10.1002/(SICI)1097-0193(1996)4:1<58::AID-HBM4>3.0.CO;2-O.
Parameters
- kint
extent, number of voxels
- evfloat
expected number of voxels
- ecfloat
expected number of clusters
Returns
- float
cluster uncorrected p-value
- Sisyphe.core.sisypheStatistics.getBoxCarModel(first: int, nactive: int, nrest: int, nblocks: int) ndarray
fMRI boxcar model.
Reference: Analysis of fmri time series revisited. KJ Friston, AP Holmes, JB Poline, PJ Grasby, SCR Williams, RSJ Frackowiak, R Turner. Neuroimage 1995 Mar;2(1):45-53. doi: 10.1006/nimg.1995.1007.
Parameters
- firstint
index of the first scan activity condition
- nactiveint
number of scans in activation blocks
- nrestint
number of scans in rest blocks
- nblocksint
number of blocks (1 block = 1 active block + 1 rest block)
Returns
- numpy.ndarray
fMRI boxcar model
- Sisyphe.core.sisypheStatistics.getDOF(design: ndarray) int
Get the degrees of freedom of a design matrix.
Parameters
- designndarray
design matrix (lines = observations, columns = factors)
Returns
- int
The degrees of freedom
- Sisyphe.core.sisypheStatistics.getHRF(pdelay: float = 6.0, ndelay: float = 16.0, pspread: float = 1.0, nspread: float = 1.0, ratio: float = 6.0, dt: float = 0.1, iscan: float = 2.0) ndarray
Hemodynamic response function. Vector used to convolve fMRI boxcar model.
Reference: Analysis of fmri time series revisited. KJ Friston, AP Holmes, JB Poline, PJ Grasby, SCR Williams, RSJ Frackowiak, R Turner. Neuroimage 1995 Mar;2(1):45-53. doi: 10.1006/nimg.1995.1007.
Parameters
- pdelayfloat
delay in s of positive response (default 6 s)
- ndelayfloat
delay in s of negative response (default 16 s)
- pspreadfloat
spread in s of positive response (default 1 s)
- nspreadfloat
spread in s of negative response (default 1 s)
- ratiofloat
ratio of positive to negative responses (default 6)
- dtfloat
sampling (default 0.1 seconds)
- iscanfloat
interscan interval (TR in seconds)
Returns
- numpy.ndarray
Hemodynamic response function
- Sisyphe.core.sisypheStatistics.getHighPass(nscans: int, nblocks: int) ndarray
High-pass vector processing. These vectors are added to fMRI design matrix as covariable of non-interest.
Reference: Analysis of fmri time series revisited. KJ Friston, AP Holmes, JB Poline, PJ Grasby, SCR Williams, RSJ Frackowiak, R Turner. Neuroimage 1995 Mar;2(1):45-53. doi: 10.1006/nimg.1995.1007.
Parameters
- nscansint
image count
- nblocksint
block count (1 block = 1 active block + 1 rest block) = order
Returns
- numpy.ndarray
High-pass vectors
- Sisyphe.core.sisypheStatistics.isRankDeficient(X: ndarray) bool
Check whether a design matrix is rank-deficient.
Parameters
- Xndarray
design matrix (lines = observations, columns = factors)
Returns
- bool
True if rank-deficient
- Sisyphe.core.sisypheStatistics.modelEstimate(obs: list[SisypheVolume], design: ndarray, mask: SisypheROI | SisypheVolume, scale: ndarray | None = None, wait: DialogWait | None = None) tuple[SisypheVolume, SisypheVolume, SisypheVolume, ndarray]
Model estimation.
Design matrix X (lines = observations, columns = factors): beta = pinv(X) Y = (XtX)-1 Xt Y residuals i.e. model errors = Y - X beta pooled variance = variance(residuals)
Reference: Statistical parametric maps in functional imaging: A general linear approach. KJ Friston, AP Holmes, KJ Worsley, JP Poline, CD Frith, RSJ Frackowiak. Human Brain Mapping 1995;2(4):189-210. doi: 10.1002/hbm.460020402.
Analysis of fmri time series revisited. KJ Friston, AP Holmes, JB Poline, PJ Grasby, SCR Williams, RSJ Frackowiak, R Turner. Neuroimage 1995 Mar;2(1):45-53. doi: 10.1006/nimg.1995.1007.
Parameters
obs : list[sisyphe.core.sisypheVolume.SisypheVolume] mask : sisyphe.core.sisypheROI.SisypheROI | sisyphe.core.sisypheVolume.SisypheVolume
analysis mask
- designnumpy.ndarray
design matrix X
- scalenumpy.ndarray
signal normalization values
- waitSisyphe.gui.dialogWait.DialogWait
progress bar
Returns
- tuple[sisyphe.core.sisypheVolume.SisypheVolume,sisyphe.core.sisypheVolume.SisypheVolume, sisyphe.core.sisypheVolume.SisypheVolume, ndarray]
beta, pooled variance, residuals, autocorrelations
- Sisyphe.core.sisypheStatistics.pCorrectedBonferroni(p: float, n: int) float
Bonferroni multiple-comparison correction.
Parameters
- pfloat
uncorrected p-value
- nint
number of comparisons
Returns
- float
corrected p-value
- Sisyphe.core.sisypheStatistics.pUncorrectedBonferroni(p: float, n: int) float
Get an uncorrected p-value from a Bonferroni corrected p-value.
Parameters
- pfloat
corrected p-value
- nint
number of comparisons
Returns
- float
uncorrected p-value
- Sisyphe.core.sisypheStatistics.pvalueTot(p: float, df: int) float
Get t-value from p-value.
Parameters
- pfloat
p-value
- dfint
degrees of freedom
Returns
- float
t-value
- Sisyphe.core.sisypheStatistics.pvalueToz(p: float) float
Get z-value from p-value.
Parameters
- pfloat
p-value
Returns
- float
z-value
- Sisyphe.core.sisypheStatistics.qFDRTot(q: float, df: int, img: SisypheVolume, independent: bool = False) float
Get t-value from a false discovery rate value.
Reference: Controlling the false discovery rate: a practical and powerful approach to multiple testing. Benjamini Y, Hochberg Y. Journal of the Royal Statistical Society, Series B. 1995;57(1):289–300. doi:10.1111/j.2517-6161.1995.
Parameters
- qfloat
false discovery rate value
- dfint
degrees of freedom
- imgSisyphe.core.sisypheVolume.SisypheVolume
statistical map
independent : bool
Returns
- float
t-value
- Sisyphe.core.sisypheStatistics.qFDRToz(q: float, img: SisypheVolume, independent: bool = False) float
Get z-value from a false discovery rate value.
Reference: Controlling the false discovery rate: a practical and powerful approach to multiple testing. Benjamini Y, Hochberg Y. Journal of the Royal Statistical Society, Series B. 1995;57(1):289–300. doi:10.1111/j.2517-6161.1995.
Parameters
- qfloat
false discovery rate value
- imgSisyphe.core.sisypheVolume.SisypheVolume
statistical map
independent : bool
Returns
- float
z-value
- Sisyphe.core.sisypheStatistics.reselCount(mask: SisypheVolume, autocorrx: float, autocorry: float, autocorrz: float) tuple[float, ...]
Resel counts processing. resel count 0 = Euler characteristic of the statistical map (i.e. gaussian field) resel count 1 = resel diameter resel count 2 = resel surface area resel count 3 = resel volume
reference: A unified statistical approach for determining significant signals in images of cerebral activation KJ Worsley, S Marrett, P Neelin, AC Vandal, KJ Friston, AC Evans. Hum Brain Mapp. 1996;4(1):58-73. doi: 10.1002/(SICI)1097-0193(1996)4:1<58::AID-HBM4>3.0.CO;2-O.
Parameters
- maskSisyphe.core.sisypheVolume.SisypheVolume
mask of analysis
- autocorrxfloat
statistical map gaussian point spread function, fwhm (mm) in x dimension
- autocorryfloat
statistical map gaussian point spread function, fwhm (mm) in y dimension
- autocorrzfloat
statistical map gaussian point spread function, fwhm (mm) in z dimension
Returns
- tuple[float, …]
resel counts
- Sisyphe.core.sisypheStatistics.tFieldEulerCharacteristic(t: float, df: int) tuple[float, ...]
t field Euler characteristic processing.
reference: A unified statistical approach for determining significant signals in images of cerebral activation KJ Worsley, S Marrett, P Neelin, AC Vandal, KJ Friston, AC Evans. Hum Brain Mapp. 1996;4(1):58-73. doi: 10.1002/(SICI)1097-0193(1996)4:1<58::AID-HBM4>3.0.CO;2-O.
Parameters
- tfloat
t-value
- dfint
degrees of freedom
Returns
- tuple[float, …]
Euler characteristics
- Sisyphe.core.sisypheStatistics.tFieldExpectedClusters(t: float, df: int, rc0: float, rc1: float, rc2: float, rc3: float) float
t field expected cluster processing.
reference: A unified statistical approach for determining significant signals in images of cerebral activation KJ Worsley, S Marrett, P Neelin, AC Vandal, KJ Friston, AC Evans. Hum Brain Mapp. 1996;4(1):58-73. doi: 10.1002/(SICI)1097-0193(1996)4:1<58::AID-HBM4>3.0.CO;2-O.
Parameters
- tfloat
t-value
- dfint
degrees of freedom
- rc0float
resel count 0, t-field Euler characteristic
- rc1float
resel count 1, t-field resel diameter
- rc2float
resel count 2, t-field resel surface area
- rc3float
resel count 3, t-field resel volume
Returns
- float
expected clusters
- Sisyphe.core.sisypheStatistics.tFieldExpectedVoxels(t: float, df: int, n: int) float
t field expected voxel processing.
reference: A unified statistical approach for determining significant signals in images of cerebral activation KJ Worsley, S Marrett, P Neelin, AC Vandal, KJ Friston, AC Evans. Hum Brain Mapp. 1996;4(1):58-73. doi: 10.1002/(SICI)1097-0193(1996)4:1<58::AID-HBM4>3.0.CO;2-O.
Parameters
- tfloat
t-value
- dfint
degrees of freedom
- nint
number of voxels in analysis mask
Returns
- float
expected voxels
- Sisyphe.core.sisypheStatistics.tToVoxelCorrectedpvalue(t: float, df: int, rc0: float, rc1: float, rc2: float, rc3: float) float
Voxel corrected p-value processing from t-value.
reference: A unified statistical approach for determining significant signals in images of cerebral activation KJ Worsley, S Marrett, P Neelin, AC Vandal, KJ Friston, AC Evans. Hum Brain Mapp. 1996;4(1):58-73. doi: 10.1002/(SICI)1097-0193(1996)4:1<58::AID-HBM4>3.0.CO;2-O.
Parameters
- tfloat
t-value
- dfint
degrees of freedom
- rc0float
resel count 0, t-field Euler characteristic
- rc1float
resel count 1, t-field resel diameter
- rc2float
resel count 2, t-field resel surface area
- rc3float
resel count 3, t-field resel volume
Returns
- float
t-value
- Sisyphe.core.sisypheStatistics.tTopvalue(t: float, df: int) float
Get p-value from t-value.
Parameters
- tfloat
t-value
- dfint
degrees of freedom
Returns
- float
p-value
- Sisyphe.core.sisypheStatistics.tToz(t: float, df: int) float
Get t-value from z-value.
Parameters
- tfloat
t-value
- dfint
degrees of freedom
Returns
- float
z-value
- Sisyphe.core.sisypheStatistics.tTozmap(tmap: SisypheVolume) SisypheVolume
t to z-map conversion
Parameters
tmap: Sisyphe.core.sisypheVolume.SisypheVolume
Returns
- Sisyphe.core.sisypheVolume.SisypheVolume
z-map
- Sisyphe.core.sisypheStatistics.thresholdMap(stat: float, ext: int, smap: SisypheVolume, lbls: list[SisypheVolume] | None = None) dict
Statistical map thresholding.
Parameters
- statfloat
statistical threshold, t or z value
- extint
extension threshold, number of voxels
- smapSisyphe.core.sisypheVolume.SisypheVolume
statistical map (tmap or zmap)
- lblslist[Sisyphe.core.sisypheVolume.SisypheVolume] | None
label image(s), to process proportion of each label in clusters
Returns
- dict
dict keys: - ‘map’: Sisyphe.core.sisypheVolume.SisypheVolume, thresholded statistical map - ‘c’: list[tuple[float, float, float], …], voxel coordinates of each cluster local maximum - ‘max’: list[float, …], statistical value of each cluster local maximum - ‘extent’: list[int, …], extent i.e. number of voxels of each cluster - label image filenames: list[dict[int, float]], dict: key = label index in the label image, value = proportion of this label in the current cluster (cluster number = list index)
- Sisyphe.core.sisypheStatistics.tmapContrastEstimate(contrast: ndarray, design: ndarray, beta: SisypheVolume, variance: SisypheVolume, df: int, wait: DialogWait | None = None) SisypheVolume
t-test contrast estimation.
Reference: Statistical parametric maps in functional imaging: A general linear approach. KJ Friston, AP Holmes, KJ Worsley, JP Poline, CD Frith, RSJ Frackowiak. Human Brain Mapping 1995;2(4):189-210. doi: 10.1002/hbm.460020402.
Analysis of fmri time series revisited. KJ Friston, AP Holmes, JB Poline, PJ Grasby, SCR Williams, RSJ Frackowiak, R Turner. Neuroimage 1995 Mar;2(1):45-53. doi: 10.1006/nimg.1995.1007.
Parameters
- contrastnumpy.ndarray
contrast vector
- designnumpy.ndarray
design matrix
- betaSisyphe.core.sisypheVolume.SisypheVolume
beta
- varianceSisyphe.core.sisypheVolume.SisypheVolume
pooled variance
- dfint
degrees of freedom (DOF)
- waitSisyphe.gui.dialogWait.DialogWait
progress bar dialog
Returns
- Sisyphe.core.sisypheVolume.SisypheVolume
t-map
- Sisyphe.core.sisypheStatistics.voxelCorrectedpvalueTot(u: float, df: int, rc0: float, rc1: float, rc2: float, rc3: float) float
t-value processing from voxel corrected p-value.
reference: A unified statistical approach for determining significant signals in images of cerebral activation KJ Worsley, S Marrett, P Neelin, AC Vandal, KJ Friston, AC Evans. Hum Brain Mapp. 1996;4(1):58-73. doi: 10.1002/(SICI)1097-0193(1996)4:1<58::AID-HBM4>3.0.CO;2-O.
Parameters
- ufloat
voxel corrected p-value
- dfint
degrees of freedom
- rc0float
resel count 0, t-field Euler characteristic
- rc1float
resel count 1, t-field resel diameter
- rc2float
resel count 2, t-field resel surface area
- rc3float
resel count 3, t-field resel volume
Returns
- float
t-value
- Sisyphe.core.sisypheStatistics.voxelCorrectedpvalueToz(u: float, rc0: float, rc1: float, rc2: float, rc3: float) float
z-value processing from voxel corrected p-value.
reference: A unified statistical approach for determining significant signals in images of cerebral activation KJ Worsley, S Marrett, P Neelin, AC Vandal, KJ Friston, AC Evans. Hum Brain Mapp. 1996;4(1):58-73. doi: 10.1002/(SICI)1097-0193(1996)4:1<58::AID-HBM4>3.0.CO;2-O.
Parameters
- ufloat
voxel corrected p-value
- rc0float
resel count 0, t-field Euler characteristic
- rc1float
resel count 1, t-field resel diameter
- rc2float
resel count 2, t-field resel surface area
- rc3float
resel count 3, t-field resel volume
Returns
- float
z-value
- Sisyphe.core.sisypheStatistics.zFieldEulerCharacteristic(z: float) tuple[float, ...]
Gaussian field Euler characteristic processing.
reference: A unified statistical approach for determining significant signals in images of cerebral activation KJ Worsley, S Marrett, P Neelin, AC Vandal, KJ Friston, AC Evans. Hum Brain Mapp. 1996;4(1):58-73. doi: 10.1002/(SICI)1097-0193(1996)4:1<58::AID-HBM4>3.0.CO;2-O.
Parameters
- zfloat
z-value
Returns
- tuple[float, …]
Euler characteristics
- Sisyphe.core.sisypheStatistics.zFieldExpectedClusters(z: float, rc0: float, rc1: float, rc2: float, rc3: float) float
Gaussian field expected cluster processing.
reference: A unified statistical approach for determining significant signals in images of cerebral activation KJ Worsley, S Marrett, P Neelin, AC Vandal, KJ Friston, AC Evans. Hum Brain Mapp. 1996;4(1):58-73. doi: 10.1002/(SICI)1097-0193(1996)4:1<58::AID-HBM4>3.0.CO;2-O.
Parameters
- zfloat
z-value
- rc0float
resel count 0, t-field Euler characteristic
- rc1float
resel count 1, t-field resel diameter
- rc2float
resel count 2, t-field resel surface area
- rc3float
resel count 3, t-field resel volume
Returns
- float
expected clusters
- Sisyphe.core.sisypheStatistics.zFieldExpectedVoxels(z: float, n: int) float
Gaussian field expected voxel processing.
reference: A unified statistical approach for determining significant signals in images of cerebral activation KJ Worsley, S Marrett, P Neelin, AC Vandal, KJ Friston, AC Evans. Hum Brain Mapp. 1996;4(1):58-73. doi: 10.1002/(SICI)1097-0193(1996)4:1<58::AID-HBM4>3.0.CO;2-O.
Parameters
- zfloat
z-value
- nint
number of voxels in analysis mask
Returns
- float
expected voxels
- Sisyphe.core.sisypheStatistics.zToVoxelCorrectedpvalue(z: float, rc0: float, rc1: float, rc2: float, rc3: float) float
Voxel corrected p-value processing from z-value.
reference: A unified statistical approach for determining significant signals in images of cerebral activation KJ Worsley, S Marrett, P Neelin, AC Vandal, KJ Friston, AC Evans. Hum Brain Mapp. 1996;4(1):58-73. doi: 10.1002/(SICI)1097-0193(1996)4:1<58::AID-HBM4>3.0.CO;2-O.
Parameters
- zfloat
z-value
- rc0float
resel count 0, t-field Euler characteristic
- rc1float
resel count 1, t-field resel diameter
- rc2float
resel count 2, t-field resel surface area
- rc3float
resel count 3, t-field resel volume
Returns
- float
z-value
- Sisyphe.core.sisypheStatistics.zTopvalue(z: float) float
Get p-value from z-value.
Parameters
- zfloat
z-value
Returns
- float
p-value
- Sisyphe.core.sisypheStatistics.zTot(z: float, df: int) float
Get z-value from t-value.
Parameters
- zfloat
z-value
- dfint
degrees of freedom
Returns
- float
t-value
- Sisyphe.core.sisypheStatistics.zmapContrastEstimate(contrast: ndarray, design: ndarray, beta: SisypheVolume, variance: SisypheVolume, df: int, wait: DialogWait | None = None) SisypheVolume
Parameters
- contrastndarray
contrast vector
- designnumpy.ndarray
design matrix
- betaSisyphe.core.sisypheVolume.SisypheVolume
beta
- varianceSisyphe.core.sisypheVolume.SisypheVolume
pooled variance
- dfint
degrees of freedom (DOF)
- waitSisyphe.gui.dialogWait.DialogWait
progress bar dialog
Returns
- Sisyphe.core.sisypheVolume.SisypheVolume
z-map