| Title: | Tools and Operations for Reading, Writing, Viewing, and Manipulating SKIFTI Files |
|---|---|
| Description: | SKIFTI files contain brain imaging data in coordinates across Tract Based Spatial Statistics (TBSS) skeleton, which represent the brain white matter intensity values. 'skiftiTools' provides a unified environment for reading, writing, visualizing and manipulating SKIFTI-format data. It supports the "subsetting", "concatenating", and using data as data.frame for R statistical functions. The SKIFTI data is structured for convenient access to the data and metadata, and includes support for visualizations. For more information see Merisaari et al. (2024) <doi:10.57736/87d2-0608>. |
| Authors: | Harri Merisaari [aut] (ORCID: <https://orcid.org/0000-0002-8515-5399>), Ilkka Suuronen [cre] (ORCID: <https://orcid.org/0009-0001-6516-4116>) |
| Maintainer: | Ilkka Suuronen <[email protected]> |
| License: | GPL-3 |
| Version: | 0.1.0 |
| Built: | 2026-05-14 07:54:13 UTC |
| Source: | https://github.com/haanme/skiftitools |
Concatenate Skifti data
concat(Skifti_data1, Skifti_data2)concat(Skifti_data1, Skifti_data2)
Skifti_data1 |
Skifti data object1 |
Skifti_data2 |
Skifti data object2 |
concatenated Skifti data object
library(RNifti) data<-array(0,dim=list(10,10,10,10)) for(t in 1:10) { for(x in 1:10) { for(y in 1:10) { for(z in 1:10) { data[x,y,z,t]<-t+x } } } } data_Nifti<-RNifti::retrieveNifti(data) RNifti::writeNifti(data_Nifti, "data_Nifti.nii.gz", template = NULL, datatype = "auto") data_skeleton<-array(0,dim=list(10,10,10)) data_skeleton[5,5,5]<-1 data_skeleton[6,6,6]<-1 data_skeleton[7,7,7]<-1 data_skeleton_Nifti<-RNifti::retrieveNifti(data_skeleton) RNifti::writeNifti(data_skeleton_Nifti, "data_skeleton_Nifti.nii.gz", datatype = "auto") data_Skifti<-Nifti2Skifti(Nifti_data="data_Nifti.nii.gz", Nifti_skeleton="data_skeleton_Nifti.nii.gz", selected_volumes=1:10, Nifti_labels=NULL, write_coordinates=TRUE, verbose=FALSE) data_Skifti_subset<-subset(data_Skifti, c(1,5,10)) m<-matrix(c(6,10,15,7,11,16,8,12,17), nrow=3, ncol=3) rownames(m)<-c("vol1", "vol5", "vol10") data_Skifti_subset1<-subset(data_Skifti, c(1,5)) data_Skifti_subset2<-subset(data_Skifti, c(10)) data_Skifti_concat<-concat(data_Skifti_subset1, data_Skifti_subset2) m<-matrix(c(6,10,15,7,11,16,8,12,17), nrow=3, ncol=3) rownames(m)<-c("vol1", "vol5", "")library(RNifti) data<-array(0,dim=list(10,10,10,10)) for(t in 1:10) { for(x in 1:10) { for(y in 1:10) { for(z in 1:10) { data[x,y,z,t]<-t+x } } } } data_Nifti<-RNifti::retrieveNifti(data) RNifti::writeNifti(data_Nifti, "data_Nifti.nii.gz", template = NULL, datatype = "auto") data_skeleton<-array(0,dim=list(10,10,10)) data_skeleton[5,5,5]<-1 data_skeleton[6,6,6]<-1 data_skeleton[7,7,7]<-1 data_skeleton_Nifti<-RNifti::retrieveNifti(data_skeleton) RNifti::writeNifti(data_skeleton_Nifti, "data_skeleton_Nifti.nii.gz", datatype = "auto") data_Skifti<-Nifti2Skifti(Nifti_data="data_Nifti.nii.gz", Nifti_skeleton="data_skeleton_Nifti.nii.gz", selected_volumes=1:10, Nifti_labels=NULL, write_coordinates=TRUE, verbose=FALSE) data_Skifti_subset<-subset(data_Skifti, c(1,5,10)) m<-matrix(c(6,10,15,7,11,16,8,12,17), nrow=3, ncol=3) rownames(m)<-c("vol1", "vol5", "vol10") data_Skifti_subset1<-subset(data_Skifti, c(1,5)) data_Skifti_subset2<-subset(data_Skifti, c(10)) data_Skifti_concat<-concat(data_Skifti_subset1, data_Skifti_subset2) m<-matrix(c(6,10,15,7,11,16,8,12,17), nrow=3, ncol=3) rownames(m)<-c("vol1", "vol5", "")
Skeleton mask and corresponding image intensity data must be in Nifti format. The skeleton mask is used to determine the coordinates of intensity data. If optional label file is given, that is used to label the voxels.
Nifti2Skifti( Nifti_data = NULL, Nifti_skeleton = NULL, selected_volumes = NULL, Nifti_labels = NULL, write_coordinates = FALSE, verbose = FALSE )Nifti2Skifti( Nifti_data = NULL, Nifti_skeleton = NULL, selected_volumes = NULL, Nifti_labels = NULL, write_coordinates = FALSE, verbose = FALSE )
Nifti_data |
Intensity data in Nifti format (required) |
Nifti_skeleton |
Skeleton at same imaging space as the data, in Nifti format (required) |
selected_volumes |
Selected volume indexes starting from 1 (default==NULL, selecting all) |
Nifti_labels |
Labeling data to be used inside mask, writing extra line to the output about labels (default==NULL, no extra labeling) |
write_coordinates |
TRUE/FALSE(default) write coordinates of voxels in x,y,z ASCII format, in the same order as they appear in the skifti |
verbose |
TRUE/FALSE(default) for verbose messages |
skifti object with default rownames as vol1, vol2 .... volN as indexes from the nifti data
#source('../../R/Skifti2Nifti.R') #source('../../R/Nifti2Skifti.R') library(RNifti) data<-array(0,dim=list(10,10,10,10)) for(t in 1:10) { for(x in 1:10) { for(y in 1:10) { for(z in 1:10) { data[x,y,z,t]<-t+x } } } } data_Nifti<-RNifti::retrieveNifti(data) RNifti::writeNifti(data_Nifti, "data_Nifti.nii.gz", template = NULL, datatype = "auto") data_skeleton<-array(0,dim=list(10,10,10)) data_skeleton[5,5,5]<-1 data_skeleton[6,6,6]<-1 data_skeleton[7,7,7]<-1 data_skeleton_Nifti<-RNifti::retrieveNifti(data_skeleton) RNifti::writeNifti(data_skeleton_Nifti, "data_skeleton_Nifti.nii.gz", datatype = "auto") data_Skifti<-Nifti2Skifti(Nifti_data="data_Nifti.nii.gz", Nifti_skeleton="data_skeleton_Nifti.nii.gz", selected_volumes=c(1), Nifti_labels=NULL, write_coordinates=TRUE, verbose=FALSE) # Create Skifti data_Nifti2<-Skifti2Nifti(data_Skifti) RNifti::writeNifti(data_Nifti2[[1]], "data_Nifti.nii.gz", datatype = "auto") data_Nifti2<-RNifti::readNifti("data_Nifti.nii.gz", internal = TRUE, volumes = NULL)#source('../../R/Skifti2Nifti.R') #source('../../R/Nifti2Skifti.R') library(RNifti) data<-array(0,dim=list(10,10,10,10)) for(t in 1:10) { for(x in 1:10) { for(y in 1:10) { for(z in 1:10) { data[x,y,z,t]<-t+x } } } } data_Nifti<-RNifti::retrieveNifti(data) RNifti::writeNifti(data_Nifti, "data_Nifti.nii.gz", template = NULL, datatype = "auto") data_skeleton<-array(0,dim=list(10,10,10)) data_skeleton[5,5,5]<-1 data_skeleton[6,6,6]<-1 data_skeleton[7,7,7]<-1 data_skeleton_Nifti<-RNifti::retrieveNifti(data_skeleton) RNifti::writeNifti(data_skeleton_Nifti, "data_skeleton_Nifti.nii.gz", datatype = "auto") data_Skifti<-Nifti2Skifti(Nifti_data="data_Nifti.nii.gz", Nifti_skeleton="data_skeleton_Nifti.nii.gz", selected_volumes=c(1), Nifti_labels=NULL, write_coordinates=TRUE, verbose=FALSE) # Create Skifti data_Nifti2<-Skifti2Nifti(data_Skifti) RNifti::writeNifti(data_Nifti2[[1]], "data_Nifti.nii.gz", datatype = "auto") data_Nifti2<-RNifti::readNifti("data_Nifti.nii.gz", internal = TRUE, volumes = NULL)
Read Skifti data
readSkifti(filename, verbose = FALSE)readSkifti(filename, verbose = FALSE)
filename |
file to read |
verbose |
TRUE/FALSE(default), for verbosity |
Skifti data object
#source('../../R/Skifti2Nifti.R') #source('../../R/Nifti2Skifti.R') library(RNifti) data<-array(0,dim=list(10,10,10,10)) for(t in 1:10) { for(x in 1:10) { for(y in 1:10) { for(z in 1:10) { data[x,y,z,t]<-t+x } } } } data_Nifti<-RNifti::retrieveNifti(data) RNifti::writeNifti(data_Nifti, "data_Nifti.nii.gz", template = NULL, datatype = "auto") data_skeleton<-array(0,dim=list(10,10,10)) data_skeleton[5,5,5]<-1 data_skeleton[6,6,6]<-1 data_skeleton[7,7,7]<-1 data_skeleton_Nifti<-RNifti::retrieveNifti(data_skeleton) RNifti::writeNifti(data_skeleton_Nifti, "data_skeleton_Nifti.nii.gz", datatype = "auto") data_Skifti<-Nifti2Skifti(Nifti_data="data_Nifti.nii.gz", Nifti_skeleton="data_skeleton_Nifti.nii.gz", selected_volumes=c(1), Nifti_labels=NULL, write_coordinates=TRUE, verbose=FALSE) # Create Skifti data_Nifti2<-Skifti2Nifti(data_Skifti) RNifti::writeNifti(data_Nifti2[[1]], "data_Nifti.nii.gz", datatype = "auto") data_Nifti2<-RNifti::readNifti("data_Nifti.nii.gz", internal = TRUE, volumes = NULL)#source('../../R/Skifti2Nifti.R') #source('../../R/Nifti2Skifti.R') library(RNifti) data<-array(0,dim=list(10,10,10,10)) for(t in 1:10) { for(x in 1:10) { for(y in 1:10) { for(z in 1:10) { data[x,y,z,t]<-t+x } } } } data_Nifti<-RNifti::retrieveNifti(data) RNifti::writeNifti(data_Nifti, "data_Nifti.nii.gz", template = NULL, datatype = "auto") data_skeleton<-array(0,dim=list(10,10,10)) data_skeleton[5,5,5]<-1 data_skeleton[6,6,6]<-1 data_skeleton[7,7,7]<-1 data_skeleton_Nifti<-RNifti::retrieveNifti(data_skeleton) RNifti::writeNifti(data_skeleton_Nifti, "data_skeleton_Nifti.nii.gz", datatype = "auto") data_Skifti<-Nifti2Skifti(Nifti_data="data_Nifti.nii.gz", Nifti_skeleton="data_skeleton_Nifti.nii.gz", selected_volumes=c(1), Nifti_labels=NULL, write_coordinates=TRUE, verbose=FALSE) # Create Skifti data_Nifti2<-Skifti2Nifti(data_Skifti) RNifti::writeNifti(data_Nifti2[[1]], "data_Nifti.nii.gz", datatype = "auto") data_Nifti2<-RNifti::readNifti("data_Nifti.nii.gz", internal = TRUE, volumes = NULL)
Skeleton mask and corresponding image intensity data must be in Nifti format. The skeleton mask is used to determine the coordinates of intensity data.
save_skeleton( mask, data, img_hdr, output, legend_title, scale, keep_temp = FALSE, palette = "lajolla", verbose = FALSE )save_skeleton( mask, data, img_hdr, output, legend_title, scale, keep_temp = FALSE, palette = "lajolla", verbose = FALSE )
mask |
Intensity data in Nifti object format |
data |
Skeleton at same imaging space as the data, in Nifti format |
img_hdr |
Nifti header object |
output |
Output PNG filename |
legend_title |
Title to be shown |
scale |
scaling for intensity values, tune for better color depth |
keep_temp |
TRUE/FALSE(default) to keep temporary png images |
palette |
color palette |
verbose |
TRUE/FALSE(default), for verbosity |
No output, as resuts are saved to a png file
Skeleton mask and corresponding image intensity data in Nifti format. The skeleton mask is used to determine the coordinates of intensity data. If optional label file is given, that is used to label the voxels.
Skifti2CSV(Skifti_data, filename, overwrite = FALSE, sep = ";")Skifti2CSV(Skifti_data, filename, overwrite = FALSE, sep = ";")
Skifti_data |
Intensity data in Nifti format |
filename |
file to read' |
overwrite |
TRUE/FALSE(default) to overwrite existing data |
sep |
file separator to be written default ';' |
CSV filename
library(RNifti) data<-array(0,dim=list(10,10,10,10)) for(t in 1:10) { for(x in 1:10) { for(y in 1:10) { for(z in 1:10) { data[x,y,z,t]<-t+x } } } } data_Nifti<-RNifti::retrieveNifti(data) RNifti::writeNifti(data_Nifti, "data_Nifti.nii.gz", template = NULL, datatype = "auto") data_skeleton<-array(0,dim=list(10,10,10)) data_skeleton[5,5,5]<-1 data_skeleton[6,6,6]<-1 data_skeleton[7,7,7]<-1 data_skeleton_Nifti<-RNifti::retrieveNifti(data_skeleton) RNifti::writeNifti(data_skeleton_Nifti, "data_skeleton_Nifti.nii.gz", datatype = "auto") data_Skifti<-Nifti2Skifti(Nifti_data="data_Nifti.nii.gz", Nifti_skeleton="data_skeleton_Nifti.nii.gz", selected_volumes=1:10, Nifti_labels=NULL, write_coordinates=TRUE, verbose=FALSE) Skifti2CSV(data_Skifti, "data_Skifti.csv", overwrite=TRUE, sep=';') data_csv<-read.csv2("data_Skifti.csv", ';', header = FALSE, row.names = NULL)library(RNifti) data<-array(0,dim=list(10,10,10,10)) for(t in 1:10) { for(x in 1:10) { for(y in 1:10) { for(z in 1:10) { data[x,y,z,t]<-t+x } } } } data_Nifti<-RNifti::retrieveNifti(data) RNifti::writeNifti(data_Nifti, "data_Nifti.nii.gz", template = NULL, datatype = "auto") data_skeleton<-array(0,dim=list(10,10,10)) data_skeleton[5,5,5]<-1 data_skeleton[6,6,6]<-1 data_skeleton[7,7,7]<-1 data_skeleton_Nifti<-RNifti::retrieveNifti(data_skeleton) RNifti::writeNifti(data_skeleton_Nifti, "data_skeleton_Nifti.nii.gz", datatype = "auto") data_Skifti<-Nifti2Skifti(Nifti_data="data_Nifti.nii.gz", Nifti_skeleton="data_skeleton_Nifti.nii.gz", selected_volumes=1:10, Nifti_labels=NULL, write_coordinates=TRUE, verbose=FALSE) Skifti2CSV(data_Skifti, "data_Skifti.csv", overwrite=TRUE, sep=';') data_csv<-read.csv2("data_Skifti.csv", ';', header = FALSE, row.names = NULL)
Skeleton mask and corresponding image intensity data in Nifti format. The skeleton mask is used to determine the coordinates of intensity data. If optional label file is given, that is used to label the voxels.
Skifti2Nifti(Skifti_data)Skifti2Nifti(Skifti_data)
Skifti_data |
Intensity data in Nifti format |
Nifti skeleton file for Skifti data
#source('../../R/Skifti2Nifti.R') #source('../../R/Nifti2Skifti.R') library(RNifti) data<-array(0,dim=list(10,10,10,10)) for(t in 1:10) { for(x in 1:10) { for(y in 1:10) { for(z in 1:10) { data[x,y,z,t]<-t+x } } } } data_Nifti<-RNifti::retrieveNifti(data) RNifti::writeNifti(data_Nifti, "data_Nifti.nii.gz", template = NULL, datatype = "auto") data_skeleton<-array(0,dim=list(10,10,10)) data_skeleton[5,5,5]<-1 data_skeleton[6,6,6]<-1 data_skeleton[7,7,7]<-1 data_skeleton_Nifti<-RNifti::retrieveNifti(data_skeleton) RNifti::writeNifti(data_skeleton_Nifti, "data_skeleton_Nifti.nii.gz", datatype = "auto") data_Skifti<-Nifti2Skifti(Nifti_data="data_Nifti.nii.gz", Nifti_skeleton="data_skeleton_Nifti.nii.gz", selected_volumes=c(1), Nifti_labels=NULL, write_coordinates=TRUE, verbose=FALSE) # Create Skifti data_Nifti2<-Skifti2Nifti(data_Skifti) RNifti::writeNifti(data_Nifti2[[1]], "data_Nifti.nii.gz", datatype = "auto") data_Nifti2<-RNifti::readNifti("data_Nifti.nii.gz", internal = TRUE, volumes = NULL)#source('../../R/Skifti2Nifti.R') #source('../../R/Nifti2Skifti.R') library(RNifti) data<-array(0,dim=list(10,10,10,10)) for(t in 1:10) { for(x in 1:10) { for(y in 1:10) { for(z in 1:10) { data[x,y,z,t]<-t+x } } } } data_Nifti<-RNifti::retrieveNifti(data) RNifti::writeNifti(data_Nifti, "data_Nifti.nii.gz", template = NULL, datatype = "auto") data_skeleton<-array(0,dim=list(10,10,10)) data_skeleton[5,5,5]<-1 data_skeleton[6,6,6]<-1 data_skeleton[7,7,7]<-1 data_skeleton_Nifti<-RNifti::retrieveNifti(data_skeleton) RNifti::writeNifti(data_skeleton_Nifti, "data_skeleton_Nifti.nii.gz", datatype = "auto") data_Skifti<-Nifti2Skifti(Nifti_data="data_Nifti.nii.gz", Nifti_skeleton="data_skeleton_Nifti.nii.gz", selected_volumes=c(1), Nifti_labels=NULL, write_coordinates=TRUE, verbose=FALSE) # Create Skifti data_Nifti2<-Skifti2Nifti(data_Skifti) RNifti::writeNifti(data_Nifti2[[1]], "data_Nifti.nii.gz", datatype = "auto") data_Nifti2<-RNifti::readNifti("data_Nifti.nii.gz", internal = TRUE, volumes = NULL)
Get subset of Skifti data
subset(Skifti_data, volumes)subset(Skifti_data, volumes)
Skifti_data |
Skifti data object |
volumes |
selection |
Skifti data object of subset
Write Skifti data
writeSkifti( Skifti_data, basename, overwrite = FALSE, compress = "none", verbose = FALSE )writeSkifti( Skifti_data, basename, overwrite = FALSE, compress = "none", verbose = FALSE )
Skifti_data |
Skifti data object |
basename |
basename to write without suffix |
overwrite |
TRUE/FALSE(default) to overwrite existing data |
compress |
bz2/zip/none(default) to select compression method |
verbose |
TRUE/FALSE(default), for verbosity |
filename where Skifti data was written
#source('../../R/Skifti2Nifti.R') #source('../../R/Nifti2Skifti.R') library(RNifti) data<-array(0,dim=list(10,10,10,10)) for(t in 1:10) { for(x in 1:10) { for(y in 1:10) { for(z in 1:10) { data[x,y,z,t]<-t+x } } } } data_Nifti<-RNifti::retrieveNifti(data) RNifti::writeNifti(data_Nifti, "data_Nifti.nii.gz", template = NULL, datatype = "auto") data_skeleton<-array(0,dim=list(10,10,10)) data_skeleton[5,5,5]<-1 data_skeleton[6,6,6]<-1 data_skeleton[7,7,7]<-1 data_skeleton_Nifti<-RNifti::retrieveNifti(data_skeleton) RNifti::writeNifti(data_skeleton_Nifti, "data_skeleton_Nifti.nii.gz", datatype = "auto") data_Skifti<-Nifti2Skifti(Nifti_data="data_Nifti.nii.gz", Nifti_skeleton="data_skeleton_Nifti.nii.gz", selected_volumes=c(1), Nifti_labels=NULL, write_coordinates=TRUE, verbose=FALSE) # Create Skifti data_Nifti2<-Skifti2Nifti(data_Skifti) RNifti::writeNifti(data_Nifti2[[1]], "data_Nifti.nii.gz", datatype = "auto") data_Nifti2<-RNifti::readNifti("data_Nifti.nii.gz", internal = TRUE, volumes = NULL)#source('../../R/Skifti2Nifti.R') #source('../../R/Nifti2Skifti.R') library(RNifti) data<-array(0,dim=list(10,10,10,10)) for(t in 1:10) { for(x in 1:10) { for(y in 1:10) { for(z in 1:10) { data[x,y,z,t]<-t+x } } } } data_Nifti<-RNifti::retrieveNifti(data) RNifti::writeNifti(data_Nifti, "data_Nifti.nii.gz", template = NULL, datatype = "auto") data_skeleton<-array(0,dim=list(10,10,10)) data_skeleton[5,5,5]<-1 data_skeleton[6,6,6]<-1 data_skeleton[7,7,7]<-1 data_skeleton_Nifti<-RNifti::retrieveNifti(data_skeleton) RNifti::writeNifti(data_skeleton_Nifti, "data_skeleton_Nifti.nii.gz", datatype = "auto") data_Skifti<-Nifti2Skifti(Nifti_data="data_Nifti.nii.gz", Nifti_skeleton="data_skeleton_Nifti.nii.gz", selected_volumes=c(1), Nifti_labels=NULL, write_coordinates=TRUE, verbose=FALSE) # Create Skifti data_Nifti2<-Skifti2Nifti(data_Skifti) RNifti::writeNifti(data_Nifti2[[1]], "data_Nifti.nii.gz", datatype = "auto") data_Nifti2<-RNifti::readNifti("data_Nifti.nii.gz", internal = TRUE, volumes = NULL)