1#ifndef SCRAN_PCA_SIMPLE_PCA_HPP
2#define SCRAN_PCA_SIMPLE_PCA_HPP
10#include "tatami_stats/tatami_stats.hpp"
13#include "irlba_tatami/irlba_tatami.hpp"
15#include "sanisizer/sanisizer.hpp"
31template<
typename EigenVector_ = Eigen::VectorXd>
87template<
bool sparse_,
typename Value_,
typename Index_,
class EigenVector_>
88void compute_row_means_and_variances(
const tatami::Matrix<Value_, Index_>& mat,
const int num_threads, EigenVector_& center_v, EigenVector_& scale_v) {
89 const auto ngenes = mat.
nrow();
98 const auto ncells = mat.
ncol();
101 for (Index_ g = start, end = start + length; g < end; ++g) {
102 const auto results = [&]{
103 if constexpr(sparse_) {
104 auto range = ext->fetch(vbuffer.data(), NULL);
105 return tatami_stats::variances::direct(range.value, range.number, ncells,
false);
107 auto ptr = ext->fetch(vbuffer.data());
108 return tatami_stats::variances::direct(ptr, ncells,
false);
111 center_v.coeffRef(g) = results.first;
112 scale_v.coeffRef(g) = results.second;
114 }, ngenes, num_threads);
118 const auto ncells = mat.
ncol();
121 typedef typename EigenVector_::Scalar Scalar;
122 tatami_stats::LocalOutputBuffer<Scalar> cbuffer(t, start, length, center_v.data());
123 tatami_stats::LocalOutputBuffer<Scalar> sbuffer(t, start, length, scale_v.data());
126 if constexpr(sparse_) {
127 return tatami_stats::variances::RunningSparse<Scalar, Value_, Index_>(length, cbuffer.data(), sbuffer.data(),
false, start);
129 return tatami_stats::variances::RunningDense<Scalar, Value_, Index_>(length, cbuffer.data(), sbuffer.data(),
false);
135 if constexpr(sparse_) {
142 for (Index_ c = 0; c < ncells; ++c) {
143 if constexpr(sparse_) {
144 const auto range = ext->fetch(vbuffer.data(), ibuffer.data());
145 running.add(range.value, range.index, range.number);
147 const auto ptr = ext->fetch(vbuffer.data());
155 }, ngenes, num_threads);
159template<
class EigenVector_,
class EigenMatrix_>
160std::unique_ptr<irlba::Matrix<EigenVector_, EigenMatrix_> > prepare_deferred_matrix_for_irlba(
162 const SimplePcaOptions<EigenVector_>& options,
163 const EigenVector_& center_v,
164 const EigenVector_& scale_v
166 std::unique_ptr<irlba::Matrix<EigenVector_, EigenMatrix_> > alt;
167 alt.reset(
new irlba::CenteredMatrix<EigenVector_, EigenMatrix_, I<
decltype(ptr)>, I<
decltype(¢er_v)> >(std::move(ptr), ¢er_v));
171 alt.reset(
new irlba::ScaledMatrix<EigenVector_, EigenMatrix_, I<
decltype(ptr)>, I<
decltype(&scale_v)> >(std::move(ptr), &scale_v,
true,
true));
178template<
class EigenMatrix_,
typename Value_,
typename Index_,
class EigenVector_>
179std::unique_ptr<irlba::Matrix<EigenVector_, EigenMatrix_> > prepare_sparse_matrix_for_irlba(
181 const SimplePcaOptions<EigenVector_>& options,
182 EigenVector_& center_v,
183 EigenVector_& scale_v,
184 typename EigenVector_::Scalar& total_var
186 const auto ngenes = mat.
nrow();
187 sanisizer::resize(center_v, ngenes);
188 sanisizer::resize(scale_v, ngenes);
189 std::unique_ptr<irlba::Matrix<EigenVector_, EigenMatrix_> > output;
191 if (options.realize_matrix) {
205 const Index_ ncells = mat.
ncol();
209 I<
decltype(extracted.value)>,
210 I<
decltype(extracted.index)>,
211 I<
decltype(extracted.pointers)>
215 std::move(extracted.value),
216 std::move(extracted.index),
217 std::move(extracted.pointers),
221 output.reset(sparse_ptr);
224 const auto& pointers = sparse_ptr->get_pointers();
225 const auto& values = sparse_ptr->get_values();
226 for (Index_ g = start, end = start + length; g < end; ++g) {
227 const auto offset = pointers[g];
228 const auto next_offset = pointers[g + 1];
229 const Index_ num_nonzero = next_offset - offset;
230 const auto results = tatami_stats::variances::direct(values.data() + offset, num_nonzero, ncells,
false);
231 center_v.coeffRef(g) = results.first;
232 scale_v.coeffRef(g) = results.second;
234 }, ngenes, options.num_threads);
236 total_var = process_scale_vector(options.scale, scale_v);
239 compute_row_means_and_variances<true>(mat, options.num_threads, center_v, scale_v);
240 total_var = process_scale_vector(options.scale, scale_v);
243 new irlba_tatami::Transposed<EigenVector_, EigenMatrix_, Value_, Index_,
decltype(&mat)>(&mat, options.num_threads)
247 return prepare_deferred_matrix_for_irlba(std::move(output), options, center_v, scale_v);
250template<
class EigenMatrix_,
typename Value_,
typename Index_,
class EigenVector_>
251std::unique_ptr<irlba::Matrix<EigenVector_, EigenMatrix_> > prepare_dense_matrix_for_irlba(
253 const SimplePcaOptions<EigenVector_>& options,
254 EigenVector_& center_v,
255 EigenVector_& scale_v,
256 typename EigenVector_::Scalar& total_var
258 const Index_ ngenes = mat.
nrow();
259 sanisizer::resize(center_v, ngenes);
260 sanisizer::resize(scale_v, ngenes);
262 if (options.realize_matrix) {
264 const Index_ ncells = mat.
ncol();
265 auto emat = std::make_unique<EigenMatrix_>(
266 sanisizer::cast<I<
decltype(std::declval<EigenMatrix_>().rows())> >(ncells),
267 sanisizer::cast<I<
decltype(std::declval<EigenMatrix_>().cols())> >(ngenes)
272 static_assert(!EigenMatrix_::IsRowMajor);
278 tatami::ConvertToDenseOptions opt;
279 opt.num_threads = options.num_threads;
284 center_v.array() = emat->array().colwise().sum();
288 std::fill(center_v.begin(), center_v.end(), std::numeric_limits<typename EigenVector_::Scalar>::quiet_NaN());
290 emat->array().rowwise() -= center_v.adjoint().array();
292 scale_v.array() = emat->array().colwise().squaredNorm();
294 scale_v /= ncells - 1;
296 std::fill(scale_v.begin(), scale_v.end(), std::numeric_limits<typename EigenVector_::Scalar>::quiet_NaN());
299 total_var = process_scale_vector(options.scale, scale_v);
301 emat->array().rowwise() /= scale_v.adjoint().array();
304 return std::unique_ptr<irlba::Matrix<EigenVector_, EigenMatrix_> >(
309 compute_row_means_and_variances<false>(mat, options.num_threads, center_v, scale_v);
310 total_var = process_scale_vector(options.scale, scale_v);
312 std::unique_ptr<irlba::Matrix<EigenVector_, EigenMatrix_> > output(
313 new irlba_tatami::Transposed<EigenVector_, EigenMatrix_, Value_, Index_,
decltype(&mat)>(&mat, options.num_threads)
315 return prepare_deferred_matrix_for_irlba(std::move(output), options, center_v, scale_v);
327template<
typename EigenMatrix_,
typename EigenVector_>
380template<
typename Value_,
typename Index_,
typename EigenMatrix_,
class EigenVector_,
class SubsetFunction_>
381void simple_pca_internal(
385 SubsetFunction_ subset_fun
389 std::unique_ptr<irlba::Matrix<EigenVector_, EigenMatrix_> > ptr;
405 if (!options.
scale) {
406 output.
scale = EigenVector_();
436template<
typename Value_,
typename Index_,
typename EigenMatrix_,
class EigenVector_>
442 [](
const EigenMatrix_&,
const EigenVector_&) ->
void {}
461template<
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
Metrics 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< EigenVector_ > &options, SimplePcaResults< EigenMatrix_, EigenVector_ > &output)
Definition simple_pca.hpp:437
CompressedSparseContents< StoredValue_, StoredIndex_, StoredPointer_ > retrieve_compressed_sparse_contents(const Matrix< InputValue_, InputIndex_ > &matrix, const bool row, const RetrieveCompressedSparseContentsOptions &options)
int parallelize(Function_ fun, const Index_ tasks, const int workers)
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:32
bool scale
Definition simple_pca.hpp:56
int number
Definition simple_pca.hpp:49
bool realize_matrix
Definition simple_pca.hpp:68
int num_threads
Definition simple_pca.hpp:76
bool transpose
Definition simple_pca.hpp:62
irlba::Options< EigenVector_ > irlba_options
Definition simple_pca.hpp:81
Results of simple_pca().
Definition simple_pca.hpp:328
EigenMatrix_ components
Definition simple_pca.hpp:337
EigenMatrix_ rotation
Definition simple_pca.hpp:357
EigenVector_ center
Definition simple_pca.hpp:363
EigenVector_ variance_explained
Definition simple_pca.hpp:344
EigenVector_ scale
Definition simple_pca.hpp:369
EigenVector_::Scalar total_variance
Definition simple_pca.hpp:350
irlba::Metrics metrics
Definition simple_pca.hpp:374
bool sparse_extract_index