updateTables | R Documentation |
metabCombiner
ObjectsThis method updates the feature list (featData) and aligned
table (combinedTable
) within a metabCombiner
object. Manual
changes to the (combinedTable
) as well as unmatched X & Y dataset
features can be incorporated into the object and the corresponding results.
This function is typically paired with link{reduceTable}
or other
forms of table reduction performed by the user.
updateTables(object, xdata = NULL, ydata = NULL, combinedTable = NULL)
object |
|
xdata |
|
ydata |
|
combinedTable |
merged table which may be altered by the user. This must
have the |
There are two points where features can be removed from the
combinedTable
report: during m/z grouping and during the table
reduction step. It is also possible for user-specified changes to the report
to remove certain features entirely. This function allows for the missed
features to be brought back into the table as non-matched entities. For xdata
features, the Y columns will be entirely missing values, and ydata features
will have missing X information. The feature data (featData) will also be
updated for use in subsequent alignments, but only features present in the
representative dataset will be retained by default.
metabCombiner
object with updates to combinedTable
to
include features that have been missed or changes by the user.
Duplicated sample & extra column names cannot be copied from the original data they feature in, therefore they are left as missing values.
data(plasma30)
data(plasma20)
p30 <- metabData(plasma30, samples = "CHEAR")
p20 <- metabData(plasma20, samples = "Red")
p.comb <- metabCombiner(xdata = p30, ydata = p20, xid = "p30", yid = "p20",
binGap = 0.0075)
##extracting, modifying, and updating combinedTable
cTable <- combinedTable(p.comb)
cTable <- dplyr::filter(cTable, rty < 17.25)
p.comb <- updateTables(p.comb, combinedTable = cTable)
p.comb <- selectAnchors(p.comb, tolmz = 0.003, tolQ = 0.3, windy = 0.02)
p.comb <- fit_gam(p.comb, k = 20, iterFilter = 1)
p.comb <- calcScores(p.comb, A = 90, B = 14, C = 0.5)
p.comb <- reduceTable(p.comb, delta = 0.2, maxRTerr = 0.5)
##updating to include features removed from xdata & ydata
p.comb <- updateTables(p.comb, xdata = p30, ydata = p20)
#view results
cTable <- combinedTable(p.comb)
fdata <- featData(p.comb)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.