View source: R/MsExperiment-db.R
dbWriteSampleData | R Documentation |
For MsExperiment
objects with their MS data represented by a Spectra
object that use a MsBackendSql
backend, its sample annotations can be
written to the backend's SQL database with the dbWriteSampleData()
function.
The content of the object's [sampleData()]
(as well as eventually present
linking between samples and spectra) will be stored in two separate
database tables sample_data and sample_to_msms_spectrum in the same
database.
This requires that the MS data of the experiment is represented by a
MsBackendSql
backend (see help on the createMsBackendSqlDatabase
or the
MsBackendSql package vignette for more information on how to create or
use such SQL databases).
dbWriteSampleData(x)
x |
|
Johannes Rainer, Laurent Gatto
library(MsExperiment)
## Create a MsBackendSql database from two mzML files.
## Connect first to an empty SQLite database (for the example we create
## a database in a temporary file).
library(RSQLite)
sqlite_db <- tempfile()
con <- dbConnect(SQLite(), sqlite_db)
## Define the files from which we import the data
fls <- dir(system.file("sciex", package = "msdata"), pattern = "mzML",
full.names = TRUE)
## Create a MsBackendSql database containing the full MS data
library(MsBackendSql)
createMsBackendSqlDatabase(con, fls)
## Note: alternatively it would be possible to first import the MS data
## to a `Spectra` object and then change the backend to a `MsBackendSql`
## using the `setBackend` function.
## Load this data as a `Spectra` object (using a `MsBackendOfflineSql`
## backend)
library(Spectra)
sps <- Spectra(sqlite_db, source = MsBackendOfflineSql(),
drv = SQLite())
sps
## Define sample annotations for the two data files. Adding one column
## `"file"` that contains the file name of the data files.
df <- data.frame(sample = c("QC1", "QC2"), file = basename(fls))
## Add a spectra variable `"file"` to the `Spectra` object with
## the raw data files' file names to simplify the linking between
## samples and spectra performed later.
sps$file <- basename(dataOrigin(sps))
## Create a MsExperiment with the spectra and sample data.
mse <- MsExperiment(spectra = sps, sampleData = df)
## Establish the link (mapping) between samples and spectra
## using the column `"file"` in the `sampleData` and the spectra
## variable `"file"`.
mse <- linkSampleData(mse, with = "sampleData.file = spectra.file")
mse
## Write sample data (and the sample to spectra mapping) to the
## *MsBackendSql* database.
dbWriteSampleData(mse)
## List the tables in the database
dbListTables(con)
## Sample data was thus stored to the database.
dbGetQuery(con, "select * from sample_data;")
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.