Title: | Bridging Simultaneous Confidence Corridors and PET Neuroimaging |
---|---|
Description: | Tools for the structured processing of PET neuroimaging data in preparation for the estimation of Simultaneous Confidence Corridors (SCCs) for one-group, two-group, or single-patient vs group comparisons. The package facilitates PET image loading, data restructuring, integration into a Functional Data Analysis framework, contour extraction, identification of significant results, and performance evaluation. It bridges established packages (e.g., 'oro.nifti') with novel statistical methodologies (e.g., 'ImageSCC') and enables reproducible analysis pipelines, including comparison with Statistical Parametric Mapping ('SPM'). |
Authors: | Juan A. Arias Lopez [aut, cre, cph]
|
Maintainer: | Juan A. Arias Lopez <[email protected]> |
License: | MIT + file LICENSE |
Version: | 1.0.0 |
Built: | 2025-03-28 08:34:13 UTC |
Source: | https://github.com/iguanamarina/neuroscc |
The neuroSCC
package provides tools to preprocess and structure neuroimaging data
for functional data analysis using Simultaneous Confidence Corridors (SCCs). It wraps external packages
to prepare data from PET images, extract contours, generate meshes, and evaluate regions of statistical significance.
The methods implemented support both group comparisons and single-subject vs. group inference, following the methodology described in Wang et al. (2020) and the author's PhD thesis.
This package serves as a bridge between neuroimaging file formats (e.g., NIfTI) and advanced
statistical tools like ImageSCC::scc.image
. It includes the following key components.
Loading and cleaning PET image data.
Extracting ROIs and constructing functional data matrices.
Generating synthetic Poisson clones for 1-vs-group settings.
Extracting SCC-detected points and evaluating detection metrics.
Maintainer: Juan A. Arias Lopez [email protected] (ORCID) [copyright holder]
Other contributors:
Virgilio Gomez Rubio [email protected] (ORCID) [reviewer]
Pablo Aguiar Fernandez [email protected] (ORCID) [thesis advisor]
Andrew Haddon Kemp [email protected] (ORCID) [thesis advisor]
neuroCleaner
, databaseCreator
, getPoints
Computes Sensitivity, Specificity, Positive Predictive Value (PPV), and Negative Predictive Value (NPV) by comparing detected points with ground truth ROI points. This function is used to assess the accuracy of SCC- or SPM-based detection in neuroimaging analysis.
calculateMetrics(detectedPoints, truePoints, totalCoords, regionName)
calculateMetrics(detectedPoints, truePoints, totalCoords, regionName)
detectedPoints |
A data frame containing detected coordinates ( |
truePoints |
A data frame with ground truth ROI coordinates ( |
totalCoords |
A list with the full voxel grid dimensions, created by |
regionName |
A character string used to label the output region. |
This function requires three precomputed objects
detectedPoints
: SCC-detected points from getPoints
or SPM-detected points from getSPMbinary
.
truePoints
: Ground truth ROIs extracted using processROIs
.
totalCoords
: Full voxel coordinate grid from getDimensions
.
A data frame with the following evaluation metrics
region
: Name of the analyzed region.
sensitivity
: True positive rate (TP / (TP + FN) * 100).
specificity
: True negative rate (TN / (TN + FP) * 100).
PPV
: Positive predictive value (TP / (TP + FP) * 100).
NPV
: Negative predictive value (TN / (TN + FN) * 100).
getPoints
for SCC-detected regions. getSPMbinary
for binary SPM-detected points. processROIs
for defining ground truth ROIs. getDimensions
for generating the coordinate grid.
# Load precomputed inputs for the example data("calculateMetricsExample", package = "neuroSCC") # Evaluate SCC and SPM detection performance with(calculateMetricsExample, { metricsSCC <- calculateMetrics(detectedSCC, trueROI, totalCoords, "Region2_SCC") metricsSPM <- calculateMetrics(detectedSPM, trueROI, totalCoords, "Region2_SPM") print(metricsSCC) print(metricsSPM) })
# Load precomputed inputs for the example data("calculateMetricsExample", package = "neuroSCC") # Evaluate SCC and SPM detection performance with(calculateMetricsExample, { metricsSCC <- calculateMetrics(detectedSCC, trueROI, totalCoords, "Region2_SCC") metricsSPM <- calculateMetrics(detectedSPM, trueROI, totalCoords, "Region2_SPM") print(metricsSCC) print(metricsSPM) })
A dataset containing all necessary inputs for demonstrating calculateMetrics
.
It enables reproducible and fast example code that compares SCC-detected and SPM-detected
points against a known ground truth ROI.
These inputs were generated using sample PET and ROI files included in the neuroSCC
package.
data("calculateMetricsExample")
data("calculateMetricsExample")
A set of four objects:
detectedSCC
Data frame of SCC-detected coordinates (from getPoints
).
detectedSPM
Data frame of SPM-detected coordinates (from getSPMbinary
).
trueROI
Ground truth ROI voxel data (from processROIs
).
totalCoords
List with full image grid dimensions (from getDimensions
).
Simulated PET neuroimaging study for testing SCC and SPM detection accuracy.
calculateMetrics
, getPoints
, getSPMbinary
,
processROIs
, getDimensions
Processes multiple PET image files matching a specified filename pattern.
Each file is processed using neuroCleaner
, and the results are aggregated
into a unified data frame for functional data analysis. This function serves as a key step
in the neuroSCC
workflow, bridging raw image data and Simultaneous Confidence Corridors (SCC) computation.
databaseCreator( pattern, control = TRUE, useSequentialNumbering = FALSE, demo = NULL, quiet = FALSE )
databaseCreator( pattern, control = TRUE, useSequentialNumbering = FALSE, demo = NULL, quiet = FALSE )
pattern |
|
control |
|
useSequentialNumbering |
|
demo |
|
quiet |
|
The function performs the following steps
Identifies image files matching the given pattern.
Processes each file using neuroCleaner
, optionally merging demographic data.
Adds a subject identifier column (CN_number
or AD_number
).
Aggregates all results into a single data frame.
If no files are successfully processed, an empty data frame is returned with a warning.
This function is typically followed by matrixCreator
, which converts the output
into a matrix format for functional analysis.
A data.frame
combining processed voxel-level data from all matched files.
Each row represents a voxel (3D pixel). The column structure depends on input
For the control group: CN_number
, z
, x
, y
, pet
For the pathological group: AD_number
, z
, x
, y
, pet
If demographics are included: additional columns PPT
, Group
, Sex
, Age
neuroCleaner
for the underlying image processing function. matrixCreator
for the next step in the workflow that converts
the database to a matrix format for SCC analysis.
# NOTE: To keep runtime below CRAN limits, this example processes only 1 subject. # You can expand the pattern to include all subjects for real use. # Example: Create a database from a single synthetic PET image (control group) controlPattern <- "^syntheticControl1\\.nii\\.gz$" databaseControls <- databaseCreator(pattern = controlPattern, control = TRUE, quiet = TRUE) head(databaseControls)
# NOTE: To keep runtime below CRAN limits, this example processes only 1 subject. # You can expand the pattern to include all subjects for real use. # Example: Create a database from a single synthetic PET image (control group) controlPattern <- "^syntheticControl1\\.nii\\.gz$" databaseControls <- databaseCreator(pattern = controlPattern, control = TRUE, quiet = TRUE) head(databaseControls)
Generates synthetic clones of a PET data matrix by adding Poisson-distributed noise to each non-zero voxel. This approach helps address the limitations of functional data analysis (FDA) in single-subject versus group (1 vs. Group) setups, where a single subject lacks sufficient variability to reliably estimate Simultaneous Confidence Corridors (SCCs).
generatePoissonClones(originalMatrix, numClones, lambdaFactor)
generatePoissonClones(originalMatrix, numClones, lambdaFactor)
originalMatrix |
A numeric matrix where each row represents a flattened PET image. |
numClones |
An integer specifying the number of synthetic clones to generate. |
lambdaFactor |
A positive numeric value that scales the magnitude of Poisson noise. |
Values equal to 0
remain unchanged to preserve background regions.
NA
values are replaced with 0
before adding noise.
Poisson noise is applied only to positive values, scaled by lambdaFactor
.
Enables valid SCC estimation in single-subject settings by artificially increasing sample size.
A numeric matrix with numClones
rows, each representing a noisy version
of originalMatrix
with Poisson noise added.
# Load example input matrix for Poisson cloning data("generatePoissonClonesExample", package = "neuroSCC") # Select 10 random voxel positions for display set.seed(123) sampledCols <- sample(ncol(generatePoissonClonesExample), 10) # Generate 1 synthetic clone clones <- generatePoissonClones(generatePoissonClonesExample, numClones = 1, lambdaFactor = 0.25) # Show voxel intensity values after cloning clones[, sampledCols]
# Load example input matrix for Poisson cloning data("generatePoissonClonesExample", package = "neuroSCC") # Select 10 random voxel positions for display set.seed(123) sampledCols <- sample(ncol(generatePoissonClonesExample), 10) # Generate 1 synthetic clone clones <- generatePoissonClones(generatePoissonClonesExample, numClones = 1, lambdaFactor = 0.25) # Show voxel intensity values after cloning clones[, sampledCols]
A full single-subject PET matrix used to demonstrate generatePoissonClones
.
This matrix was extracted from simulated neuroimaging data included in the neuroSCC
package.
The example avoids long runtime by generating only one synthetic clone.
data("generatePoissonClonesExample")
data("generatePoissonClonesExample")
A numeric matrix named generatePoissonClonesExample
, with 1 row and all voxel columns.
Simulated PET neuroimaging dataset included with neuroSCC
.
Extracts voxel dimension information from a NIfTI or similar neuroimaging file.
This function is designed to work with neuroCleaner
, but it can also be used
independently to inspect image dimensions.
getDimensions(file)
getDimensions(file)
file |
A NIfTI file object or a file path pointing to a NIfTI image. |
The function accepts either a file path or a preloaded nifti
object.
If a file path is provided, it uses oro.nifti::readNIfTI()
to load the image.
This function ensures consistent dimension extraction across the neuroSCC
pipeline.
A named list with the following elements
xDim
– Number of voxels along the X axis.
yDim
– Number of voxels along the Y axis.
zDim
– Number of slices along the Z axis.
dim
– Total number of voxels in a 2D slice (calculated as xDim * yDim
).
# Get the file path for a sample NIfTI file niftiFile <- system.file("extdata", "syntheticControl1.nii.gz", package = "neuroSCC") # Extract dimensions from the NIfTI file dimensions <- getDimensions(niftiFile) # Display the extracted dimensions print(dimensions)
# Get the file path for a sample NIfTI file niftiFile <- system.file("extdata", "syntheticControl1.nii.gz", package = "neuroSCC") # Extract dimensions from the NIfTI file dimensions <- getDimensions(niftiFile) # Display the extracted dimensions print(dimensions)
Identifies and extracts coordinates where differences fall outside the simultaneous confidence corridors (SCCs),
indicating statistically significant regions. This function processes the results from
ImageSCC::scc.image()
and returns the voxel locations that represent either hypo- or hyperactivity.
Interpretation depends on the order of inputs in the SCC computation.
If SCC was computed as scc.image(Ya = Y_AD, Yb = Y_CN, ...)
(i.e., the Control group is the second argument).
positivePoints
— Regions where Control minus Pathological is significantly above the SCC.
These correspond to areas where the Pathological group (AD) is hypoactive relative to Controls.
negativePoints
— Regions where Control minus Pathological is significantly below the SCC.
These correspond to areas where the Pathological group is hyperactive relative to Controls.
Always confirm the order of Ya
and Yb
in the SCC computation
to interpret the directionality correctly.
getPoints(sccResult)
getPoints(sccResult)
sccResult |
A list of SCC computation results, as returned by
|
A named list with two elements
positivePoints
— A data frame with coordinates where the first group (Ya) shows
significantly lower activity than the second group (Yb).
negativePoints
— A data frame with coordinates where the first group (Ya) shows
significantly higher activity than the second group (Yb).
ImageSCC::scc.image
for SCC computation.
# Load precomputed SCC example data("SCCcomp", package = "neuroSCC") # Extract significant SCC points significantPoints <- getPoints(SCCcomp) # Show extracted points (interpretation depends on SCC setup; see description) head(significantPoints$positivePoints) # Pathological hypoactive vs. Control head(significantPoints$negativePoints) # Pathological hyperactive vs. Control
# Load precomputed SCC example data("SCCcomp", package = "neuroSCC") # Extract significant SCC points significantPoints <- getPoints(SCCcomp) # Show extracted points (interpretation depends on SCC setup; see description) head(significantPoints$positivePoints) # Pathological hypoactive vs. Control head(significantPoints$negativePoints) # Pathological hyperactive vs. Control
Extracts voxel coordinates where pet = 1
(i.e., statistically significant points)
from a binary NIfTI file produced by an external SPM analysis.
Only voxels from a specific brain slice (z = paramZ
) are retained.
The output data frame is structured identically to that of getPoints
,
allowing direct comparison between SCC- and SPM-detected regions via calculateMetrics
.
getSPMbinary(niftiFile, paramZ = 35)
getSPMbinary(niftiFile, paramZ = 35)
niftiFile |
|
paramZ |
|
This function converts externally generated SPM results into a format compatible
with SCC analysis tools in neuroSCC
.
Use getDimensions
to inspect the full coordinate space if needed.
A data frame with the following columns:
x
, y
– Coordinates of significant voxels at the specified slice.
getPoints
for SCC-based detection. getDimensions
for obtaining full coordinate grids. calculateMetrics
for evaluating SCC vs. SPM detection performance.
# Load a sample binary NIfTI file (SPM result) niftiFile <- system.file("extdata", "binary.nii.gz", package = "neuroSCC") detectedSPM <- getSPMbinary(niftiFile, paramZ = 35) # Show detected points head(detectedSPM)
# Load a sample binary NIfTI file (SPM result) niftiFile <- system.file("extdata", "binary.nii.gz", package = "neuroSCC") detectedSPM <- getSPMbinary(niftiFile, paramZ = 35) # Show detected points head(detectedSPM)
Converts a PET image database (created via databaseCreator
) into
a matrix format suitable for functional data analysis.
Each row of the resulting matrix corresponds to a subject, and each column to a voxel's PET intensity
values at a specified brain slice.
matrixCreator( database, paramZ = 35, useSequentialNumbering = FALSE, quiet = FALSE )
matrixCreator( database, paramZ = 35, useSequentialNumbering = FALSE, quiet = FALSE )
database |
A data frame created by |
paramZ |
An integer specifying the z-coordinate (slice) to extract. Default is |
useSequentialNumbering |
|
quiet |
|
This function performs the following steps
Verifies that the specified z-slice exists in the database.
Identifies the correct subject grouping column (CN_number
or AD_number
).
Determines the matrix dimensions using x
and y
coordinates.
Extracts PET intensities per subject at the given slice.
Replaces any NaN
values with 0
to ensure numerical stability.
This function typically follows databaseCreator
and precedes
meanNormalization
in the neuroSCC workflow.
A numeric matrix where
Each row represents one subject's PET values at the selected z-slice.
Each column corresponds to a voxel (flattened as a 1D row).
databaseCreator
for generating the input database. meanNormalization
for scaling matrix data prior to SCC computation.
# NOTE: To keep example runtime short, only one synthetic PET file is used. # For full analysis, expand the filename pattern accordingly. # Step 1: Generate a database for a single subject controlPattern <- "^syntheticControl1\\.nii\\.gz$" databaseControls <- databaseCreator(pattern = controlPattern, control = TRUE, quiet = TRUE) # Step 2: Convert the database into a matrix format matrixControls <- matrixCreator(databaseControls, paramZ = 35, quiet = TRUE) # Display dimensions of the matrix dim(matrixControls)
# NOTE: To keep example runtime short, only one synthetic PET file is used. # For full analysis, expand the filename pattern accordingly. # Step 1: Generate a database for a single subject controlPattern <- "^syntheticControl1\\.nii\\.gz$" databaseControls <- databaseCreator(pattern = controlPattern, control = TRUE, quiet = TRUE) # Step 2: Convert the database into a matrix format matrixControls <- matrixCreator(databaseControls, paramZ = 35, quiet = TRUE) # Display dimensions of the matrix dim(matrixControls)
Normalizes each row of a matrix by dividing its elements by the row mean, ignoring NA
values.
This step is commonly used to adjust for global intensity differences across subjects before
applying statistical comparisons or functional data analysis.
meanNormalization( matrixData, handleInvalidRows = c("warn", "error", "omit"), returnDetails = FALSE, quiet = FALSE )
meanNormalization( matrixData, handleInvalidRows = c("warn", "error", "omit"), returnDetails = FALSE, quiet = FALSE )
matrixData |
A |
handleInvalidRows |
|
returnDetails |
|
quiet |
|
The function performs the following steps
Computes the row means of the input matrix, ignoring NA
s.
Divides each row by its corresponding mean.
Replaces NaN
values (from division by 0) with 0
if applicable.
Handles problematic rows according to the selected handleInvalidRows
option:
"warn"
(default) issues a warning, "error"
stops execution,
and "omit"
removes the affected rows from the result.
This step is often used prior to applying SCC methods to ensure comparability across subjects.
A normalized matrix, or a list if returnDetails = TRUE
.
normalizedMatrix
– The normalized matrix.
problemRows
– Indices of rows that had zero or NA
means.
matrixCreator
for building the matrix input to normalize.
# Generate a minimal database and create a matrix (1 control subject) dataDir <- system.file("extdata", package = "neuroSCC") controlPattern <- "^syntheticControl1\\.nii\\.gz$" databaseControls <- databaseCreator(pattern = controlPattern, control = TRUE, quiet = TRUE) matrixControls <- matrixCreator(databaseControls, paramZ = 35, quiet = TRUE) # Normalize the matrix (with diagnostics) normalizationResult <- meanNormalization(matrixControls, returnDetails = TRUE, quiet = FALSE)
# Generate a minimal database and create a matrix (1 control subject) dataDir <- system.file("extdata", package = "neuroSCC") controlPattern <- "^syntheticControl1\\.nii\\.gz$" databaseControls <- databaseCreator(pattern = controlPattern, control = TRUE, quiet = TRUE) matrixControls <- matrixCreator(databaseControls, paramZ = 35, quiet = TRUE) # Normalize the matrix (with diagnostics) normalizationResult <- meanNormalization(matrixControls, returnDetails = TRUE, quiet = FALSE)
Loads a NIfTI-format neuroimaging file and transforms it into a structured data frame,
organizing voxel-level information for downstream analysis. This function is the first step
in the neuroimaging processing pipeline in neuroSCC
, converting raw PET data into
a format suitable for functional data analysis. SCCs are later computed using functions
from the ImageSCC
package, such as ImageSCC::scc.image()
.
neuroCleaner(name, demo = NULL, demoRow = 1)
neuroCleaner(name, demo = NULL, demoRow = 1)
name |
|
demo |
Optional |
demoRow |
|
The function performs the following steps
Loads the NIfTI file using oro.nifti::readNIfTI()
.
Converts the 3D image into a tidy data frame.
Adds z
, x
, and y
voxel coordinates.
If demographic data is provided, attempts to match based on PPT
(case-insensitive). If no match is found, demoRow
is used.
The resulting data frame serves as input for databaseCreator
, matrixCreator
,
and other core functions in the neuroSCC
pipeline.
A data frame where each row represents a voxel (3D pixel).
If demographics are provided: the columns include PPT
, Group
, Sex
, Age
, z
, x
, y
, and pet
.
If demographics are not provided: the columns include z
, x
, y
, and pet
.
The pet
column contains the PET intensity value at each voxel location.
databaseCreator
for batch image processing. readNIfTI
for reading NIfTI-format files.
# Load a sample Control NIfTI file niftiFile <- system.file("extdata", "syntheticControl1.nii.gz", package = "neuroSCC") # Example Without demographic data petData <- neuroCleaner(niftiFile) petData[sample(nrow(petData), 10), ] # Show 10 random voxels
# Load a sample Control NIfTI file niftiFile <- system.file("extdata", "syntheticControl1.nii.gz", package = "neuroSCC") # Example Without demographic data petData <- neuroCleaner(niftiFile) petData[sample(nrow(petData), 10), ] # Show 10 random voxels
This function extracts contours from a neuroimaging NIFTI file where values change according to specified levels.
It processes the NIFTI file with neuroCleaner
to extract structured neuroimaging data, then extracts contours
using contoureR::getContourLines
. These contours serve as input for Triangulation::TriMesh
,
which is used in Simultaneous Confidence Corridors (SCCs) calculations.
While not mandatory, it is highly recommended that the input NIFTI file be pre-processed
such that zero values represent the background and non-zero values represent regions of interest.
The function's default behavior extracts contours at level 0
, which is ideal for well-masked data.
neuroContour(niftiFile, paramZ = 35, levels = c(0), plotResult = FALSE)
neuroContour(niftiFile, paramZ = 35, levels = c(0), plotResult = FALSE)
niftiFile |
|
paramZ |
|
levels |
|
plotResult |
|
This function extracts contours from a NIFTI file, typically a masked image where background values
are set to zero, and regions of interest contain non-zero values. While users can specify a different
boundary level, the recommended approach is to use levels = 0
for masked data.
The extracted contours are typically used as input to Triangulation::TriMesh
to create
a triangular mesh of the region, which is then used for Simultaneous Confidence Corridors calculations.
A list
of data frames, where each data frame contains the x and y coordinates of a contour.
The first element typically represents the external boundary, while subsequent elements (if present)
represent internal contours or holes. Each data frame has two columns:
x
– x-coordinates of the contour points.
y
– y-coordinates of the contour points.
getContourLines for the underlying contour extraction. Triangulation::TriMesh
for the next step in the SCC calculation process.
# Get the file path for a sample NIfTI file niftiFile <- system.file("extdata", "syntheticControl1.nii.gz", package = "neuroSCC") # Extract contours at level 0 contours <- neuroContour(niftiFile, paramZ = 35, levels = 0, plotResult = TRUE) # Display the first few points of the main contour head(contours[[1]])
# Get the file path for a sample NIfTI file niftiFile <- system.file("extdata", "syntheticControl1.nii.gz", package = "neuroSCC") # Extract contours at level 0 contours <- neuroContour(niftiFile, paramZ = 35, levels = 0, plotResult = TRUE) # Display the first few points of the main contour head(contours[[1]])
Processes Regions of Interest (ROIs) from a binary NIfTI file by extracting voxel-level
coordinates and labeling each voxel as part of the ROI or not. The function preserves the
spatial structure and is typically used to prepare ground truth ROIs for comparison with
SCC-detected regions via calculateMetrics
.
processROIs( roiFile, region, number, save = TRUE, outputDir = tempdir(), verbose = TRUE )
processROIs( roiFile, region, number, save = TRUE, outputDir = tempdir(), verbose = TRUE )
roiFile |
|
region |
|
number |
|
save |
|
outputDir |
|
verbose |
|
The function uses neuroCleaner
to load and flatten the NIfTI file into a structured
data frame. All voxels are retained, with the pet
column indicating which ones are part
of the ROI (1
) versus background (0
). An ROI label is added in the group
column.
This output is used as ground truth for evaluating detection performance in SCC analyses.
A data frame with voxel-level ROI information.
group
– Combined identifier built from region
and number
.
z
, x
, y
– Voxel coordinates.
pet
– Binary value indicating ROI membership (1
= ROI, 0
= non-ROI).
If save = TRUE
, the data frame is saved as an .RDS
file and not returned to the console.
calculateMetrics
for evaluating SCC detection performance. neuroCleaner
for reading and structuring voxel data.
# Load and process a sample ROI NIfTI file (console output)
# Load and process a sample ROI NIfTI file (console output)
A precomputed example of a Simultaneous Confidence Corridor (SCC) analysis
comparing a group of pathological subjects against controls. This object was generated
using the ImageSCC::scc.image
function and represents a realistic output
from SCC-based neuroimaging group comparisons.
This dataset is used in the examples of getPoints
and calculateMetrics
,
allowing users to explore SCC outputs without needing to recompute them.
data("SCCcomp")
data("SCCcomp")
A named list of class "image"
with the following elements
scc
3D array of SCC confidence bands, dimensions [n, 2, alpha]
.
Z.band
Matrix of grid coordinates corresponding to evaluated locations.
ind.inside.cover
Integer vector of indices for grid points inside the SCC band.
V.est.a
, V.est.b
Vertex matrices for triangulated domains (pathological and control groups).
Tr.est.a
, Tr.est.b
Triangle index matrices corresponding to the domain meshes.
alpha
Vector of confidence levels used (e.g., 0.1, 0.05, 0.01).
d.est
Spline degree used in mean function estimation.
r
Smoothing parameter used during fitting.
Simulated PET neuroimaging study for evaluating SCC methodology.
getPoints
, calculateMetrics
, ImageSCC::scc.image