scran
C++ library for basic single-cell RNA-seq analyses
Loading...
Searching...
No Matches
Classes | Public Types | Public Member Functions | List of all members
scran::ScoreMarkers Class Reference

Score each gene as a candidate marker for each group of cells. More...

#include <ScoreMarkers.hpp>

Classes

struct  Defaults
 Default parameter settings. More...
 
struct  Results
 Results of the marker scoring. More...
 

Public Types

typedef std::array< bool, differential_analysis::n_summaries > ComputeSummaries
 

Public Member Functions

ScoreMarkersset_threshold (double t=Defaults::threshold)
 
ScoreMarkersset_num_threads (int n=Defaults::num_threads)
 
ScoreMarkersset_cache_size (int c=Defaults::cache_size)
 
ScoreMarkersset_block_weight_policy (WeightPolicy w=Defaults::block_weight_policy)
 
ScoreMarkersset_variable_block_weight_parameters (VariableBlockWeightParameters v=Defaults::variable_block_weight_parameters)
 
ScoreMarkersset_compute_cohen (ComputeSummaries s=Defaults::compute_all_summaries())
 
ScoreMarkersset_compute_cohen (bool c)
 
ScoreMarkersset_compute_cohen (differential_analysis::summary s, bool c)
 
ScoreMarkersset_compute_auc (ComputeSummaries s=Defaults::compute_all_summaries())
 
ScoreMarkersset_compute_auc (bool c)
 
ScoreMarkersset_compute_auc (differential_analysis::summary s, bool c)
 
ScoreMarkersset_compute_lfc (ComputeSummaries s=Defaults::compute_all_summaries())
 
ScoreMarkersset_compute_lfc (bool c)
 
ScoreMarkersset_compute_lfc (differential_analysis::summary s, bool c)
 
ScoreMarkersset_compute_delta_detected (ComputeSummaries s=Defaults::compute_all_summaries())
 
ScoreMarkersset_compute_delta_detected (bool c)
 
ScoreMarkersset_compute_delta_detected (differential_analysis::summary s, bool c)
 
ScoreMarkersset_summary_min (bool s)
 
ScoreMarkersset_summary_mean (bool s)
 
ScoreMarkersset_summary_median (bool s)
 
ScoreMarkersset_summary_max (bool s)
 
ScoreMarkersset_summary_min_rank (bool s)
 
template<typename Value_ , typename Index_ , typename Group_ , typename Stat_ >
void run (const tatami::Matrix< Value_, Index_ > *p, const Group_ *group, std::vector< Stat_ * > means, std::vector< Stat_ * > detected, std::vector< std::vector< Stat_ * > > cohen, std::vector< std::vector< Stat_ * > > auc, std::vector< std::vector< Stat_ * > > lfc, std::vector< std::vector< Stat_ * > > delta_detected) const
 
template<typename Value_ , typename Index_ , typename Group_ , typename Block_ , typename Stat_ >
void run_blocked (const tatami::Matrix< Value_, Index_ > *p, const Group_ *group, const Block_ *block, std::vector< Stat_ * > means, std::vector< Stat_ * > detected, std::vector< std::vector< Stat_ * > > cohen, std::vector< std::vector< Stat_ * > > auc, std::vector< std::vector< Stat_ * > > lfc, std::vector< std::vector< Stat_ * > > delta_detected) const
 
template<typename Stat_ = double, typename Value_ , typename Index_ , typename Group_ >
Results< Stat_ > run (const tatami::Matrix< Value_, Index_ > *p, const Group_ *group) const
 
template<typename Stat_ = double, typename Value_ , typename Index_ , typename Group_ , typename Block_ >
Results< Stat_ > run_blocked (const tatami::Matrix< Value_, Index_ > *p, const Group_ *group, const Block_ *block) const
 

Detailed Description

Score each gene as a candidate marker for each group of cells.

Markers are identified by differential expression analyses between pairs of groups of cells (e.g., clusters, cell types). Given n groups, each group is involved in n - 1 pairwise comparisons and thus has n - 1 effect sizes. For each group, we compute summary statistics - e.g., median, mean - of the effect sizes across all of that group's comparisons. Users can then sort by any of these summaries to obtain a ranking of potential marker genes for each group.

The choice of effect size and summary statistic determines the characteristics of the resulting marker set. For the effect sizes: we compute Cohen's d, the area under the curve (AUC), the log-fold change and the delta-detected, which are described in more detail in the documentation for PairwiseEffects, For the summary statistics: we compute the minimum, mean, median, maximum and min-rank of the effect sizes across each group's pairwise comparisons, which are described in SummarizeEffects.

If the dataset contains blocking factors such as batch or sample, we compute the effect size within each level of the blocking factor. This avoids interference from batch effects or sample-to-sample variation. Users can also adjust the effect size to account for a minimum log-fold change threshold, in order to focus on markers with larger changes in expression. See PairwiseEffects for more details.

As a courtesy, we also compute the mean expression a

Member Typedef Documentation

◆ ComputeSummaries

typedef std::array<bool, differential_analysis::n_summaries> scran::ScoreMarkers::ComputeSummaries

Array type indicating whether each summary statistic should be computed. Each entry corresponds to a summary statistic enumerated in differential_analysis::summary.

Member Function Documentation

◆ set_threshold()

ScoreMarkers & scran::ScoreMarkers::set_threshold ( double  t = Defaults::threshold)
inline
Parameters
tThreshold on the log-fold change. This should be non-negative.
Returns
A reference to this ScoreMarkers object.

◆ set_num_threads()

ScoreMarkers & scran::ScoreMarkers::set_num_threads ( int  n = Defaults::num_threads)
inline
Parameters
nNumber of threads to use.
Returns
A reference to this ScoreMarkers object.

◆ set_cache_size()

ScoreMarkers & scran::ScoreMarkers::set_cache_size ( int  c = Defaults::cache_size)
inline
Parameters
cSize of the cache, in terms of the number of pairwise comparisons. Larger values speed up the comparisons at the cost of higher memory consumption.
Returns
A reference to this ScoreMarkers object.

◆ set_block_weight_policy()

ScoreMarkers & scran::ScoreMarkers::set_block_weight_policy ( WeightPolicy  w = Defaults::block_weight_policy)
inline
Parameters
wPolicy to use for weighting blocks when computing average statistics/effect sizes across blocks.
Returns
A reference to this ScoreMarkers instance.

◆ set_variable_block_weight_parameters()

ScoreMarkers & scran::ScoreMarkers::set_variable_block_weight_parameters ( VariableBlockWeightParameters  v = Defaults::variable_block_weight_parameters)
inline
Parameters
vParameters for the variable block weights, see variable_block_weight() for more details. Only used when the block weight policy is set to WeightPolicy::VARIABLE.
Returns
A reference to this ScoreMarkers instance.

◆ set_compute_cohen() [1/3]

ScoreMarkers & scran::ScoreMarkers::set_compute_cohen ( ComputeSummaries  s = Defaults::compute_all_summaries())
inline
Parameters
sWhich summary statistics to compute for Cohen's d.

This only has an effect for run() methods that return Results. Otherwise, we make this decision based on the validity of the input pointers.

Returns
A reference to this ScoreMarkers object.

◆ set_compute_cohen() [2/3]

ScoreMarkers & scran::ScoreMarkers::set_compute_cohen ( bool  c)
inline
Parameters
cWhether to compute Cohen's d at all.

This is an alias for set_compute_cohen() where c = true is equivalent to s = Defaults::compute_all_summaries() and c = false is equivalent to s = Defaults::compute_no_summaries().

This only has an effect for run() methods that return Results. Otherwise, we make this decision based on the validity of the input pointers.

Returns
A reference to this ScoreMarkers object.

◆ set_compute_cohen() [3/3]

ScoreMarkers & scran::ScoreMarkers::set_compute_cohen ( differential_analysis::summary  s,
bool  c 
)
inline
Parameters
sA summary statistic of interest.
cWhether to compute the summary statistic s for Cohen's d.

This only has an effect for run() methods that return Results. Otherwise, we make this decision based on the validity of the input pointers.

Returns
A reference to this ScoreMarkers object.

◆ set_compute_auc() [1/3]

ScoreMarkers & scran::ScoreMarkers::set_compute_auc ( ComputeSummaries  s = Defaults::compute_all_summaries())
inline
Parameters
sWhich summary statistics to compute for the AUC.

This only has an effect for run() methods that return Results. Otherwise, we make this decision based on the validity of the input pointers.

Returns
A reference to this ScoreMarkers object.

◆ set_compute_auc() [2/3]

ScoreMarkers & scran::ScoreMarkers::set_compute_auc ( bool  c)
inline
Parameters
cWhether to compute the AUC at all.

This is an alias for set_compute_auc() where c = true is equivalent to s = Defaults::compute_all_summaries() and c = false is equivalent to s = Defaults::compute_no_summaries().

This only has an effect for run() methods that return Results. Otherwise, we make this decision based on the validity of the input pointers.

Returns
A reference to this ScoreMarkers object.

◆ set_compute_auc() [3/3]

ScoreMarkers & scran::ScoreMarkers::set_compute_auc ( differential_analysis::summary  s,
bool  c 
)
inline
Parameters
sA summary statistic of interest.
cWhether to compute the summary statistic s for the AUC.

This only has an effect for run() methods that return Results. Otherwise, we make this decision based on the validity of the input pointers.

Returns
A reference to this ScoreMarkers object.

◆ set_compute_lfc() [1/3]

ScoreMarkers & scran::ScoreMarkers::set_compute_lfc ( ComputeSummaries  s = Defaults::compute_all_summaries())
inline
Parameters
sWhich summary statistics to compute for the LFC.

This only has an effect for run() methods that return Results. Otherwise, we make this decision based on the validity of the input pointers.

Returns
A reference to this ScoreMarkers object.

◆ set_compute_lfc() [2/3]

ScoreMarkers & scran::ScoreMarkers::set_compute_lfc ( bool  c)
inline
Parameters
cWhether to compute the LFC at all.

This is an alias for set_compute_lfc() where c = true is equivalent to s = Defaults::compute_all_summaries() and c = false is equivalent to s = Defaults::compute_no_summaries().

This only has an effect for run() methods that return Results. Otherwise, we make this decision based on the validity of the input pointers.

Returns
A reference to this ScoreMarkers object.

◆ set_compute_lfc() [3/3]

ScoreMarkers & scran::ScoreMarkers::set_compute_lfc ( differential_analysis::summary  s,
bool  c 
)
inline
Parameters
sA summary statistic of interest.
cWhether to compute the summary statistic s for the LFC.

This only has an effect for run() methods that return Results. Otherwise, we make this decision based on the validity of the input pointers.

Returns
A reference to this ScoreMarkers object.

◆ set_compute_delta_detected() [1/3]

ScoreMarkers & scran::ScoreMarkers::set_compute_delta_detected ( ComputeSummaries  s = Defaults::compute_all_summaries())
inline
Parameters
sWhich summary statistics to compute for the delta detected.

This only has an effect for run() methods that return Results. Otherwise, we make this decision based on the validity of the input pointers.

Returns
A reference to this ScoreMarkers object.

◆ set_compute_delta_detected() [2/3]

ScoreMarkers & scran::ScoreMarkers::set_compute_delta_detected ( bool  c)
inline
Parameters
cWhether to compute the delta detected at all.

This is an alias for set_compute_delta_detected() where c = true is equivalent to s = Defaults::compute_all_summaries() and c = false is equivalent to s = Defaults::compute_no_summaries().

This only has an effect for run() methods that return Results. Otherwise, we make this decision based on the validity of the input pointers.

Returns
A reference to this ScoreMarkers object.

◆ set_compute_delta_detected() [3/3]

ScoreMarkers & scran::ScoreMarkers::set_compute_delta_detected ( differential_analysis::summary  s,
bool  c 
)
inline
Parameters
sA summary statistic of interest.
cWhether to compute the summary statistic s for the delta detected.

This only has an effect for run() methods that return Results. Otherwise, we make this decision based on the validity of the input pointers.

Returns
A reference to this ScoreMarkers object.

◆ set_summary_min()

ScoreMarkers & scran::ScoreMarkers::set_summary_min ( bool  s)
inline
Parameters
sWhether to compute the minimum summary statistic for any effect size.

This overrides any previous settings for the minimum from the effect-size-specific setters, e.g., set_compute_cohen(). However, it can also be overridden by later calls to those setters.

This only has an effect for run() methods that return Results. Otherwise, we make this decision based on the validity of the input pointers.

Returns
A reference to this ScoreMarkers object.

◆ set_summary_mean()

ScoreMarkers & scran::ScoreMarkers::set_summary_mean ( bool  s)
inline
Parameters
sWhether to compute the mean summary statistic for any effect size.

This overrides any previous settings for the mean from the effect-size-specific setters, e.g., set_compute_cohen(). However, it can also be overridden by later calls to those setters.

This only has an effect for run() methods that return Results. Otherwise, we make this decision based on the validity of the input pointers.

Returns
A reference to this ScoreMarkers object.

◆ set_summary_median()

ScoreMarkers & scran::ScoreMarkers::set_summary_median ( bool  s)
inline
Parameters
sWhether to compute the median summary statistic for any effect size.

This overrides any previous settings for the median from the effect-size-specific setters, e.g., set_compute_cohen(). However, it can also be overridden by later calls to those setters.

This only has an effect for run() methods that return Results. Otherwise, we make this decision based on the validity of the input pointers.

Returns
A reference to this ScoreMarkers object.

◆ set_summary_max()

ScoreMarkers & scran::ScoreMarkers::set_summary_max ( bool  s)
inline
Parameters
sWhether to compute the maximum summary statistic for any effect size.

This overrides any previous settings for the maximum from the effect-size-specific setters, e.g., set_compute_cohen(). However, it can also be overridden by later calls to those setters.

This only has an effect for run() methods that return Results. Otherwise, we make this decision based on the validity of the input pointers.

Returns
A reference to this ScoreMarkers object.

◆ set_summary_min_rank()

ScoreMarkers & scran::ScoreMarkers::set_summary_min_rank ( bool  s)
inline
Parameters
sWhether to compute the minimum rank summary statistic for any effect size.

This overrides any previous settings for the minimum rank from the effect-size-specific setters, e.g., set_compute_cohen(). However, it can also be overridden by later calls to those setters.

This only has an effect for run() methods that return Results. Otherwise, we make this decision based on the validity of the input pointers.

Returns
A reference to this ScoreMarkers object.

◆ run() [1/2]

template<typename Value_ , typename Index_ , typename Group_ , typename Stat_ >
void scran::ScoreMarkers::run ( const tatami::Matrix< Value_, Index_ > *  p,
const Group_ *  group,
std::vector< Stat_ * >  means,
std::vector< Stat_ * >  detected,
std::vector< std::vector< Stat_ * > >  cohen,
std::vector< std::vector< Stat_ * > >  auc,
std::vector< std::vector< Stat_ * > >  lfc,
std::vector< std::vector< Stat_ * > >  delta_detected 
) const
inline

Score potential marker genes by computing summary statistics across pairwise comparisons between groups. On completion, means, detected, cohen, auc, lfc and delta_detected are filled with their corresponding statistics.

If cohen is of length 0, Cohen's d is not computed. If any of the inner vectors are of length 0, the corresponding summary statistic is not computed. The same applies to auc, lfc and delta_detected. (set_compute_cohen() and related functions have no effect here.)

Template Parameters
Data_Matrix data type.
Index_Matrix index type.
Group_Integer type for the group assignments.
Stat_Floating-point type to store the statistics.
Parameters
pPointer to a tatami matrix instance.
[in]groupPointer to an array of length equal to the number of columns in p, containing the group assignments. These should be 0-based and consecutive.
[out]meansPointers to arrays of length equal to the number of rows in p, used to store the mean expression of each group.
[out]detectedPointers to arrays of length equal to the number of rows in p, used to store the proportion of detected expression in each group.
[out]cohenVector of vector of pointers to arrays of length equal to the number of rows in p. Each inner vector corresponds to a summary statistic for Cohen's d, ordered as in differential_analysis::summary. Each pointer corresponds to a group and is filled with the relevant summary statistic for that group.
[out]aucVector of vector of pointers as described for cohen, but instead storing summary statistics for the AUC.
[out]lfcVector of vector of pointers as described for cohen, but instead storing summary statistics for the log-fold change instead of Cohen's d.
[out]delta_detectedVector of vector of pointers as described for cohen, but instead storing summary statistics for the delta in the detected proportions.

◆ run_blocked() [1/2]

template<typename Value_ , typename Index_ , typename Group_ , typename Block_ , typename Stat_ >
void scran::ScoreMarkers::run_blocked ( const tatami::Matrix< Value_, Index_ > *  p,
const Group_ *  group,
const Block_ *  block,
std::vector< Stat_ * >  means,
std::vector< Stat_ * >  detected,
std::vector< std::vector< Stat_ * > >  cohen,
std::vector< std::vector< Stat_ * > >  auc,
std::vector< std::vector< Stat_ * > >  lfc,
std::vector< std::vector< Stat_ * > >  delta_detected 
) const
inline

Score potential marker genes by computing summary statistics across pairwise comparisons between groups in multiple blocks. On completion, means, detected, cohen, auc, lfc and delta_detected are filled with their corresponding statistics.

If cohen is of length 0, Cohen's d is not computed. If any of the inner vectors are of length 0, the corresponding summary statistic is not computed. The same applies to auc, lfc and delta_detected. (set_compute_cohen() and related functions have no effect here.)

Template Parameters
Data_Matrix data type.
Index_Matrix index type.
Group_Integer type for the group assignments.
Block_Integer type for the block assignments.
Stat_Floating-point type to store the statistics.
Parameters
pPointer to a tatami matrix instance.
[in]groupPointer to an array of length equal to the number of columns in p, containing the group assignments. These should be 0-based and consecutive.
[in]blockPointer to an array of length equal to the number of columns in p, containing the blocking factor. Levels should be 0-based and consecutive. If NULL, this is ignored and the result is the same as calling run().
[out]meansVector of vectors of pointers to arrays of length equal to the number of rows in p. Each inner vector corresponds to a group and each pointer therein contains the mean expression in a blocking level.
[out]detectedPointers to arrays of length equal to the number of rows in p. Each inner vector corresponds to a group and each pointer therein contains the proportion of detected expression in a blocking level.
[out]cohenVector of vector of pointers to arrays of length equal to the number of rows in p. Each inner vector corresponds to a summary statistic for Cohen's d, ordered as in differential_analysis::summary. Each pointer corresponds to a group and is filled with the relevant summary statistic for that group.
[out]aucVector of vector of pointers as described for cohen, but instead storing summary statistics for the AUC.
[out]lfcVector of vector of pointers as described for cohen, but instead storing summary statistics for the log-fold change instead of Cohen's d.
[out]delta_detectedVector of vector of pointers as described for cohen, but instead storing summary statistics for the delta in the detected proportions.

◆ run() [2/2]

template<typename Stat_ = double, typename Value_ , typename Index_ , typename Group_ >
Results< Stat_ > scran::ScoreMarkers::run ( const tatami::Matrix< Value_, Index_ > *  p,
const Group_ *  group 
) const
inline

Score potential marker genes by computing summary statistics across pairwise comparisons between groups.

Template Parameters
Stat_Floating-point type to store the statistics.
Data_Matrix data type.
Index_Matrix index type.
Group_Integer type for the group assignments.
Parameters
pPointer to a tatami matrix instance.
[in]groupPointer to an array of length equal to the number of columns in p, containing the group assignments. These should be 0-based and consecutive.
Returns
A Results object containing the summary statistics and the other per-group statistics. Whether particular statistics are computed depends on the configuration from set_compute_cohen() and related setters.

◆ run_blocked() [2/2]

template<typename Stat_ = double, typename Value_ , typename Index_ , typename Group_ , typename Block_ >
Results< Stat_ > scran::ScoreMarkers::run_blocked ( const tatami::Matrix< Value_, Index_ > *  p,
const Group_ *  group,
const Block_ *  block 
) const
inline

Score potential marker genes by computing summary statistics across pairwise comparisons between groups in multiple blocks.

Template Parameters
Stat_Floating-point type to store the statistics.
Data_Matrix data type.
Index_Matrix index type.
Group_Integer type for the group assignments.
Block_Integer type for the block assignments.
Parameters
pPointer to a tatami matrix instance.
[in]groupPointer to an array of length equal to the number of columns in p, containing the group assignments. These should be 0-based and consecutive.
[in]blockPointer to an array of length equal to the number of columns in p, containing the blocking factor. Levels should be 0-based and consecutive. If NULL, this is ignored and the result is the same as calling run().
Returns
A Results object containing the summary statistics and the other per-group statistics. Whether particular statistics are computed depends on the configuration from set_compute_cohen() and related setters.

The documentation for this class was generated from the following file: