1#ifndef SCRAN_PCA_SIMPLE_PCA_HPP
2#define SCRAN_PCA_SIMPLE_PCA_HPP
10#include "tatami_stats/tatami_stats.hpp"
11#include "irlba/irlba.hpp"
13#include "irlba_tatami/irlba_tatami.hpp"
15#include "sanisizer/sanisizer.hpp"
82template<
bool sparse_,
typename Value_,
typename Index_,
class EigenVector_>
83void compute_row_means_and_variances(
const tatami::Matrix<Value_, Index_>& mat,
const int num_threads, EigenVector_& center_v, EigenVector_& scale_v) {
84 const auto ngenes = mat.
nrow();
93 const auto ncells = mat.
ncol();
96 for (Index_ g = start, end = start + length; g < end; ++g) {
97 const auto results = [&]{
98 if constexpr(sparse_) {
99 auto range = ext->fetch(vbuffer.data(), NULL);
100 return tatami_stats::variances::direct(range.value, range.number, ncells,
false);
102 auto ptr = ext->fetch(vbuffer.data());
103 return tatami_stats::variances::direct(ptr, ncells,
false);
106 center_v.coeffRef(g) = results.first;
107 scale_v.coeffRef(g) = results.second;
109 }, ngenes, num_threads);
113 const auto ncells = mat.
ncol();
116 typedef typename EigenVector_::Scalar Scalar;
117 tatami_stats::LocalOutputBuffer<Scalar> cbuffer(t, start, length, center_v.data());
118 tatami_stats::LocalOutputBuffer<Scalar> sbuffer(t, start, length, scale_v.data());
121 if constexpr(sparse_) {
122 return tatami_stats::variances::RunningSparse<Scalar, Value_, Index_>(length, cbuffer.data(), sbuffer.data(),
false, start);
124 return tatami_stats::variances::RunningDense<Scalar, Value_, Index_>(length, cbuffer.data(), sbuffer.data(),
false);
130 if constexpr(sparse_) {
137 for (Index_ c = 0; c < ncells; ++c) {
138 if constexpr(sparse_) {
139 const auto range = ext->fetch(vbuffer.data(), ibuffer.data());
140 running.add(range.value, range.index, range.number);
142 const auto ptr = ext->fetch(vbuffer.data());
150 }, ngenes, num_threads);
154template<
class EigenVector_,
class EigenMatrix_>
155std::unique_ptr<irlba::Matrix<EigenVector_, EigenMatrix_> > prepare_deferred_matrix_for_irlba(
157 const SimplePcaOptions& options,
158 const EigenVector_& center_v,
159 const EigenVector_& scale_v
161 std::unique_ptr<irlba::Matrix<EigenVector_, EigenMatrix_> > alt;
162 alt.reset(
new irlba::CenteredMatrix<EigenVector_, EigenMatrix_, I<
decltype(ptr)>, I<
decltype(¢er_v)> >(std::move(ptr), ¢er_v));
166 alt.reset(
new irlba::ScaledMatrix<EigenVector_, EigenMatrix_, I<
decltype(ptr)>, I<
decltype(&scale_v)> >(std::move(ptr), &scale_v,
true,
true));
173template<
class EigenMatrix_,
typename Value_,
typename Index_,
class EigenVector_>
174std::unique_ptr<irlba::Matrix<EigenVector_, EigenMatrix_> > prepare_sparse_matrix_for_irlba(
176 const SimplePcaOptions& options,
177 EigenVector_& center_v,
178 EigenVector_& scale_v,
179 typename EigenVector_::Scalar& total_var
181 const auto ngenes = mat.
nrow();
182 sanisizer::resize(center_v, ngenes);
183 sanisizer::resize(scale_v, ngenes);
184 std::unique_ptr<irlba::Matrix<EigenVector_, EigenMatrix_> > output;
186 if (options.realize_matrix) {
200 const Index_ ncells = mat.
ncol();
204 I<
decltype(extracted.value)>,
205 I<
decltype(extracted.index)>,
206 I<
decltype(extracted.pointers)>
210 std::move(extracted.value),
211 std::move(extracted.index),
212 std::move(extracted.pointers),
216 output.reset(sparse_ptr);
219 const auto& pointers = sparse_ptr->get_pointers();
220 const auto& values = sparse_ptr->get_values();
221 for (Index_ g = start, end = start + length; g < end; ++g) {
222 const auto offset = pointers[g];
223 const auto next_offset = pointers[g + 1];
224 const Index_ num_nonzero = next_offset - offset;
225 const auto results = tatami_stats::variances::direct(values.data() + offset, num_nonzero, ncells,
false);
226 center_v.coeffRef(g) = results.first;
227 scale_v.coeffRef(g) = results.second;
229 }, ngenes, options.num_threads);
231 total_var = process_scale_vector(options.scale, scale_v);
234 compute_row_means_and_variances<true>(mat, options.num_threads, center_v, scale_v);
235 total_var = process_scale_vector(options.scale, scale_v);
238 new irlba_tatami::Transposed<EigenVector_, EigenMatrix_, Value_, Index_,
decltype(&mat)>(&mat, options.num_threads)
242 return prepare_deferred_matrix_for_irlba(std::move(output), options, center_v, scale_v);
245template<
class EigenMatrix_,
typename Value_,
typename Index_,
class EigenVector_>
246std::unique_ptr<irlba::Matrix<EigenVector_, EigenMatrix_> > prepare_dense_matrix_for_irlba(
248 const SimplePcaOptions& options,
249 EigenVector_& center_v,
250 EigenVector_& scale_v,
251 typename EigenVector_::Scalar& total_var
253 const Index_ ngenes = mat.
nrow();
254 sanisizer::resize(center_v, ngenes);
255 sanisizer::resize(scale_v, ngenes);
257 if (options.realize_matrix) {
259 const Index_ ncells = mat.
ncol();
260 auto emat = std::make_unique<EigenMatrix_>(
261 sanisizer::cast<I<
decltype(std::declval<EigenMatrix_>().rows())> >(ncells),
262 sanisizer::cast<I<
decltype(std::declval<EigenMatrix_>().cols())> >(ngenes)
267 static_assert(!EigenMatrix_::IsRowMajor);
273 tatami::ConvertToDenseOptions opt;
274 opt.num_threads = options.num_threads;
279 center_v.array() = emat->array().colwise().sum();
283 std::fill(center_v.begin(), center_v.end(), std::numeric_limits<typename EigenVector_::Scalar>::quiet_NaN());
285 emat->array().rowwise() -= center_v.adjoint().array();
287 scale_v.array() = emat->array().colwise().squaredNorm();
289 scale_v /= ncells - 1;
291 std::fill(scale_v.begin(), scale_v.end(), std::numeric_limits<typename EigenVector_::Scalar>::quiet_NaN());
294 total_var = process_scale_vector(options.scale, scale_v);
296 emat->array().rowwise() /= scale_v.adjoint().array();
299 return std::unique_ptr<irlba::Matrix<EigenVector_, EigenMatrix_> >(
304 compute_row_means_and_variances<false>(mat, options.num_threads, center_v, scale_v);
305 total_var = process_scale_vector(options.scale, scale_v);
307 std::unique_ptr<irlba::Matrix<EigenVector_, EigenMatrix_> > output(
308 new irlba_tatami::Transposed<EigenVector_, EigenMatrix_, Value_, Index_,
decltype(&mat)>(&mat, options.num_threads)
310 return prepare_deferred_matrix_for_irlba(std::move(output), options, center_v, scale_v);
322template<
typename EigenMatrix_,
typename EigenVector_>
375template<
typename Value_,
typename Index_,
typename EigenMatrix_,
class EigenVector_,
class SubsetFunction_>
376void simple_pca_internal(
380 SubsetFunction_ subset_fun
384 std::unique_ptr<irlba::Matrix<EigenVector_, EigenMatrix_> > ptr;
401 if (!options.
scale) {
402 output.
scale = EigenVector_();
432template<
typename Value_,
typename Index_,
typename EigenMatrix_,
class EigenVector_>
438 [](
const EigenMatrix_&,
const EigenVector_&) ->
void {}
457template<
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
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:433
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:29
bool realize_matrix
Definition simple_pca.hpp:65
irlba::Options< Eigen::VectorXd > irlba_options
Definition simple_pca.hpp:76
bool transpose
Definition simple_pca.hpp:59
int num_threads
Definition simple_pca.hpp:71
bool scale
Definition simple_pca.hpp:53
int number
Definition simple_pca.hpp:46
Results of simple_pca().
Definition simple_pca.hpp:323
EigenMatrix_ components
Definition simple_pca.hpp:332
EigenMatrix_ rotation
Definition simple_pca.hpp:352
bool converged
Definition simple_pca.hpp:369
EigenVector_ center
Definition simple_pca.hpp:358
EigenVector_ variance_explained
Definition simple_pca.hpp:339
EigenVector_ scale
Definition simple_pca.hpp:364
EigenVector_::Scalar total_variance
Definition simple_pca.hpp:345
bool sparse_extract_index