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_>
82 auto ngenes = mat.
nrow();
91 auto ncells = mat.
ncol();
94 for (Index_ g = start, end = start + length; g < end; ++g) {
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 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 auto range = ext->fetch(vbuffer.data(), ibuffer.data());
138 running.add(range.value, range.index, range.number);
140 auto ptr = ext->fetch(vbuffer.data());
148 }, ngenes, num_threads);
152template<
class IrlbaMatrix_,
class EigenMatrix_,
class EigenVector_>
153auto run_irlba_deferred(
154 const IrlbaMatrix_& mat,
155 const SimplePcaOptions& options,
156 EigenMatrix_& components,
157 EigenMatrix_& rotation,
158 EigenVector_& variance_explained,
159 EigenVector_& center_v,
160 EigenVector_& scale_v)
164 irlba::Scaled<
true,
decltype(centered), EigenVector_> scaled(centered, scale_v,
true);
165 return irlba::compute(scaled, options.number, components, rotation, variance_explained, options.irlba_options);
167 return irlba::compute(centered, options.number, components, rotation, variance_explained, options.irlba_options);
171template<
typename Value_,
typename Index_,
class EigenMatrix_,
class EigenVector_>
174 const SimplePcaOptions& options,
175 EigenMatrix_& components,
176 EigenMatrix_& rotation,
177 EigenVector_& variance_explained,
178 EigenVector_& center_v,
179 EigenVector_& scale_v,
180 typename EigenVector_::Scalar& total_var,
183 auto ngenes = mat.
nrow();
187 if (options.realize_matrix) {
201 Index_ ncells = mat.
ncol();
205 std::move(extracted.value),
206 std::move(extracted.index),
207 std::move(extracted.pointers),
213 const auto& pointers = emat.get_pointers();
214 const auto& values = emat.get_values();
215 for (Index_ g = start, end = start + length; g < end; ++g) {
216 auto offset = pointers[g];
217 auto next_offset = pointers[g + 1];
218 Index_ num_nonzero = next_offset - offset;
219 auto results = tatami_stats::variances::direct(values.data() + offset, num_nonzero, ncells,
false);
220 center_v.coeffRef(g) = results.first;
221 scale_v.coeffRef(g) = results.second;
223 }, ngenes, options.num_threads);
225 total_var = internal::process_scale_vector(options.scale, scale_v);
226 auto out = run_irlba_deferred(emat, options, components, rotation, variance_explained, center_v, scale_v);
227 converged = out.first;
230 compute_row_means_and_variances<true>(mat, options.num_threads, center_v, scale_v);
231 total_var = internal::process_scale_vector(options.scale, scale_v);
232 auto out = run_irlba_deferred(
233 internal::TransposedTatamiWrapper<EigenVector_, Value_, Index_>(mat, options.num_threads),
241 converged = out.first;
245template<
typename Value_,
typename Index_,
class EigenMatrix_,
class EigenVector_>
248 const SimplePcaOptions& options,
249 EigenMatrix_& components,
250 EigenMatrix_& rotation,
251 EigenVector_& variance_explained,
252 EigenVector_& center_v,
253 EigenVector_& scale_v,
254 typename EigenVector_::Scalar& total_var,
257 Index_ ngenes = mat.
nrow();
261 if (options.realize_matrix) {
263 Index_ ncells = mat.
ncol();
265 sanisizer::cast<
decltype(std::declval<EigenMatrix_>().rows())>(ncells),
266 sanisizer::cast<
decltype(std::declval<EigenMatrix_>().cols())>(ngenes)
276 tatami::ConvertToDenseOptions opt;
277 opt.num_threads = options.num_threads;
282 center_v.array() = emat.array().colwise().sum();
286 std::fill(center_v.begin(), center_v.end(), std::numeric_limits<typename EigenVector_::Scalar>::quiet_NaN());
288 emat.array().rowwise() -= center_v.adjoint().array();
290 scale_v.array() = emat.array().colwise().squaredNorm();
292 scale_v /= ncells - 1;
294 std::fill(scale_v.begin(), scale_v.end(), std::numeric_limits<typename EigenVector_::Scalar>::quiet_NaN());
297 total_var = internal::process_scale_vector(options.scale, scale_v);
299 emat.array().rowwise() /= scale_v.adjoint().array();
302 auto out =
irlba::compute(emat, options.number, components, rotation, variance_explained, options.irlba_options);
303 converged = out.first;
306 compute_row_means_and_variances<false>(mat, options.num_threads, center_v, scale_v);
307 total_var = internal::process_scale_vector(options.scale, scale_v);
308 auto out = run_irlba_deferred(
309 internal::TransposedTatamiWrapper<EigenVector_, Value_, Index_>(mat, options.num_threads),
317 converged = out.first;
331template<
typename EigenMatrix_,
typename EigenVector_>
401template<
typename Value_,
typename Index_,
typename EigenMatrix_,
class EigenVector_>
416 if (!options.
scale) {
417 output.
scale = EigenVector_();
436template<
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, Eigen::Index number, EigenMatrix_ &outU, EigenMatrix_ &outV, EigenVector_ &outD, const Options &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:402
CompressedSparseContents< StoredValue_, StoredIndex_, StoredPointer_ > retrieve_compressed_sparse_contents(const Matrix< InputValue_, InputIndex_ > &matrix, bool row, const RetrieveCompressedSparseContentsOptions &options)
decltype(std::declval< Container_ >().size()) cast_Index_to_container_size(Index_ x)
void parallelize(Function_ fun, Index_ tasks, int threads)
Container_ create_container_of_Index_size(Index_ x, Args_ &&... args)
void convert_to_dense(const Matrix< InputValue_, InputIndex_ > &matrix, bool row_major, StoredValue_ *store, const ConvertToDenseOptions &options)
auto consecutive_extractor(const Matrix< Value_, Index_ > &matrix, bool row, Index_ iter_start, Index_ iter_length, Args_ &&... args)
Options for simple_pca().
Definition simple_pca.hpp:27
bool realize_matrix
Definition simple_pca.hpp:67
bool transpose
Definition simple_pca.hpp:55
int num_threads
Definition simple_pca.hpp:61
irlba::Options irlba_options
Definition simple_pca.hpp:72
bool scale
Definition simple_pca.hpp:49
int number
Definition simple_pca.hpp:43
Results of simple_pca().
Definition simple_pca.hpp:332
EigenMatrix_ components
Definition simple_pca.hpp:339
EigenMatrix_ rotation
Definition simple_pca.hpp:358
bool converged
Definition simple_pca.hpp:375
EigenVector_ center
Definition simple_pca.hpp:364
EigenVector_ variance_explained
Definition simple_pca.hpp:345
EigenVector_ scale
Definition simple_pca.hpp:370
EigenVector_::Scalar total_variance
Definition simple_pca.hpp:351
bool sparse_extract_index