1#ifndef SCRAN_PCA_SIMPLE_PCA_HPP
2#define SCRAN_PCA_SIMPLE_PCA_HPP
9#include "tatami_stats/tatami_stats.hpp"
13#include "sanisizer/sanisizer.hpp"
80template<
bool sparse_,
typename Value_,
typename Index_,
class EigenVector_>
81void compute_row_means_and_variances(
const tatami::Matrix<Value_, Index_>& mat,
const int num_threads, EigenVector_& center_v, EigenVector_& scale_v) {
82 const auto ngenes = mat.
nrow();
91 const auto ncells = mat.
ncol();
94 for (Index_ g = start, end = start + length; g < end; ++g) {
95 const auto results = [&]{
96 if constexpr(sparse_) {
97 auto range = ext->fetch(vbuffer.data(), NULL);
98 return tatami_stats::variances::direct(range.value, range.number, ncells,
false);
100 auto ptr = ext->fetch(vbuffer.data());
101 return tatami_stats::variances::direct(ptr, ncells,
false);
104 center_v.coeffRef(g) = results.first;
105 scale_v.coeffRef(g) = results.second;
107 }, ngenes, num_threads);
111 const auto ncells = mat.
ncol();
114 typedef typename EigenVector_::Scalar Scalar;
115 tatami_stats::LocalOutputBuffer<Scalar> cbuffer(t, start, length, center_v.data());
116 tatami_stats::LocalOutputBuffer<Scalar> sbuffer(t, start, length, scale_v.data());
119 if constexpr(sparse_) {
120 return tatami_stats::variances::RunningSparse<Scalar, Value_, Index_>(length, cbuffer.data(), sbuffer.data(),
false, start);
122 return tatami_stats::variances::RunningDense<Scalar, Value_, Index_>(length, cbuffer.data(), sbuffer.data(),
false);
128 if constexpr(sparse_) {
135 for (Index_ c = 0; c < ncells; ++c) {
136 if constexpr(sparse_) {
137 const auto range = ext->fetch(vbuffer.data(), ibuffer.data());
138 running.add(range.value, range.index, range.number);
140 const auto ptr = ext->fetch(vbuffer.data());
148 }, ngenes, num_threads);
152template<
class EigenVector_,
class EigenMatrix_>
153std::unique_ptr<irlba::Matrix<EigenVector_, EigenMatrix_> > prepare_deferred_matrix_for_irlba(
155 const SimplePcaOptions& options,
156 const EigenVector_& center_v,
157 const EigenVector_& scale_v
159 std::unique_ptr<irlba::Matrix<EigenVector_, EigenMatrix_> > alt;
160 alt.reset(
new irlba::CenteredMatrix<EigenVector_, EigenMatrix_, I<
decltype(ptr)>, I<
decltype(¢er_v)>>(std::move(ptr), ¢er_v));
164 alt.reset(
new irlba::ScaledMatrix<EigenVector_, EigenMatrix_, I<
decltype(ptr)>, I<
decltype(&scale_v)>>(std::move(ptr), &scale_v,
true,
true));
171template<
class EigenMatrix_,
typename Value_,
typename Index_,
class EigenVector_>
172std::unique_ptr<irlba::Matrix<EigenVector_, EigenMatrix_> > prepare_sparse_matrix_for_irlba(
174 const SimplePcaOptions& options,
175 EigenVector_& center_v,
176 EigenVector_& scale_v,
177 typename EigenVector_::Scalar& total_var
179 const auto ngenes = mat.
nrow();
180 sanisizer::resize(center_v, ngenes);
181 sanisizer::resize(scale_v, ngenes);
182 std::unique_ptr<irlba::Matrix<EigenVector_, EigenMatrix_> > output;
184 if (options.realize_matrix) {
198 const Index_ ncells = mat.
ncol();
202 I<
decltype(extracted.value)>,
203 I<
decltype(extracted.index)>,
204 I<
decltype(extracted.pointers)>
208 std::move(extracted.value),
209 std::move(extracted.index),
210 std::move(extracted.pointers),
214 output.reset(sparse_ptr);
217 const auto& pointers = sparse_ptr->get_pointers();
218 const auto& values = sparse_ptr->get_values();
219 for (Index_ g = start, end = start + length; g < end; ++g) {
220 const auto offset = pointers[g];
221 const auto next_offset = pointers[g + 1];
222 const Index_ num_nonzero = next_offset - offset;
223 const auto results = tatami_stats::variances::direct(values.data() + offset, num_nonzero, ncells,
false);
224 center_v.coeffRef(g) = results.first;
225 scale_v.coeffRef(g) = results.second;
227 }, ngenes, options.num_threads);
229 total_var = process_scale_vector(options.scale, scale_v);
232 compute_row_means_and_variances<true>(mat, options.num_threads, center_v, scale_v);
233 total_var = process_scale_vector(options.scale, scale_v);
235 output.reset(
new TransposedTatamiWrapperMatrix<EigenVector_, EigenMatrix_, Value_, Index_>(mat, options.num_threads));
238 return prepare_deferred_matrix_for_irlba(std::move(output), options, center_v, scale_v);
241template<
class EigenMatrix_,
typename Value_,
typename Index_,
class EigenVector_>
242std::unique_ptr<irlba::Matrix<EigenVector_, EigenMatrix_> > prepare_dense_matrix_for_irlba(
244 const SimplePcaOptions& options,
245 EigenVector_& center_v,
246 EigenVector_& scale_v,
247 typename EigenVector_::Scalar& total_var
249 const Index_ ngenes = mat.
nrow();
250 sanisizer::resize(center_v, ngenes);
251 sanisizer::resize(scale_v, ngenes);
253 if (options.realize_matrix) {
255 const Index_ ncells = mat.
ncol();
256 auto emat = std::make_unique<Eigen::MatrixXd>(
257 sanisizer::cast<I<
decltype(std::declval<EigenMatrix_>().rows())> >(ncells),
258 sanisizer::cast<I<
decltype(std::declval<EigenMatrix_>().cols())> >(ngenes)
268 tatami::ConvertToDenseOptions opt;
269 opt.num_threads = options.num_threads;
274 center_v.array() = emat->array().colwise().sum();
278 std::fill(center_v.begin(), center_v.end(), std::numeric_limits<typename EigenVector_::Scalar>::quiet_NaN());
280 emat->array().rowwise() -= center_v.adjoint().array();
282 scale_v.array() = emat->array().colwise().squaredNorm();
284 scale_v /= ncells - 1;
286 std::fill(scale_v.begin(), scale_v.end(), std::numeric_limits<typename EigenVector_::Scalar>::quiet_NaN());
289 total_var = process_scale_vector(options.scale, scale_v);
291 emat->array().rowwise() /= scale_v.adjoint().array();
294 return std::unique_ptr<irlba::Matrix<EigenVector_, EigenMatrix_> >(
299 compute_row_means_and_variances<false>(mat, options.num_threads, center_v, scale_v);
300 total_var = process_scale_vector(options.scale, scale_v);
302 std::unique_ptr<irlba::Matrix<EigenVector_, EigenMatrix_> > output(
303 new TransposedTatamiWrapperMatrix<EigenVector_, EigenMatrix_, Value_, Index_>(mat, options.num_threads)
305 return prepare_deferred_matrix_for_irlba(std::move(output), options, center_v, scale_v);
317template<
typename EigenMatrix_,
typename EigenVector_>
387template<
typename Value_,
typename Index_,
typename EigenMatrix_,
class EigenVector_>
391 std::unique_ptr<irlba::Matrix<EigenVector_, EigenMatrix_> > ptr;
406 if (!options.
scale) {
407 output.
scale = EigenVector_();
426template<
typename EigenMatrix_ = Eigen::MatrixXd,
class EigenVector_ = Eigen::VectorXd,
typename Value_,
typename Index_>
virtual Index_ ncol() const=0
virtual Index_ nrow() const=0
virtual bool prefer_rows() const=0
virtual std::unique_ptr< MyopicSparseExtractor< Value_, Index_ > > sparse(bool row, const Options &opt) const=0
std::pair< bool, int > compute(const Matrix_ &matrix, const Eigen::Index number, EigenMatrix_ &outU, EigenMatrix_ &outV, EigenVector_ &outD, const Options< EigenVector_ > &options)
Principal component analysis on single-cell data.
void simple_pca(const tatami::Matrix< Value_, Index_ > &mat, const SimplePcaOptions &options, SimplePcaResults< EigenMatrix_, EigenVector_ > &output)
Definition simple_pca.hpp:388
void parallelize(Function_ fun, const Index_ tasks, const int threads)
CompressedSparseContents< StoredValue_, StoredIndex_, StoredPointer_ > retrieve_compressed_sparse_contents(const Matrix< InputValue_, InputIndex_ > &matrix, const bool row, const RetrieveCompressedSparseContentsOptions &options)
void convert_to_dense(const Matrix< InputValue_, InputIndex_ > &matrix, const bool row_major, StoredValue_ *const store, const ConvertToDenseOptions &options)
Container_ create_container_of_Index_size(const Index_ x, Args_ &&... args)
auto consecutive_extractor(const Matrix< Value_, Index_ > &matrix, const bool row, const Index_ iter_start, const Index_ iter_length, Args_ &&... args)
Options for simple_pca().
Definition simple_pca.hpp:27
bool realize_matrix
Definition simple_pca.hpp:63
irlba::Options< Eigen::VectorXd > irlba_options
Definition simple_pca.hpp:74
bool transpose
Definition simple_pca.hpp:57
int num_threads
Definition simple_pca.hpp:69
bool scale
Definition simple_pca.hpp:51
int number
Definition simple_pca.hpp:44
Results of simple_pca().
Definition simple_pca.hpp:318
EigenMatrix_ components
Definition simple_pca.hpp:325
EigenMatrix_ rotation
Definition simple_pca.hpp:344
bool converged
Definition simple_pca.hpp:361
EigenVector_ center
Definition simple_pca.hpp:350
EigenVector_ variance_explained
Definition simple_pca.hpp:331
EigenVector_ scale
Definition simple_pca.hpp:356
EigenVector_::Scalar total_variance
Definition simple_pca.hpp:337
bool sparse_extract_index