scran_pca
Principal component analysis for single-cell data
Loading...
Searching...
No Matches
blocked_pca.hpp
Go to the documentation of this file.
1#ifndef SCRAN_PCA_BLOCKED_PCA_HPP
2#define SCRAN_PCA_BLOCKED_PCA_HPP
3
4#include <vector>
5#include <cmath>
6#include <algorithm>
7#include <type_traits>
8#include <cstddef>
9#include <functional>
10
11#include "tatami/tatami.hpp"
12#include "irlba/irlba.hpp"
13#include "irlba/parallel.hpp"
14#include "Eigen/Dense"
16#include "sanisizer/sanisizer.hpp"
17
18#include "utils.hpp"
19
26namespace scran_pca {
27
36 // Avoid throwing an error if too many PCs are requested.
38 }
48 int number = 25;
49
56 bool scale = false;
57
62 bool transpose = true;
63
73 scran_blocks::WeightPolicy block_weight_policy = scran_blocks::WeightPolicy::VARIABLE;
74
80
87
92 bool realize_matrix = true;
93
98 int num_threads = 1;
99
104};
105
109/*****************************************************
110 ************* Blocking data structures **************
111 *****************************************************/
112
113template<typename Index_, class EigenVector_>
114struct BlockingDetails {
115 std::vector<Index_> block_size;
116
117 bool weighted = false;
118 typedef typename EigenVector_::Scalar Weight;
119
120 // The below should only be used if weighted = true.
121 std::vector<Weight> per_element_weight;
122 Weight total_block_weight = 0;
123 EigenVector_ expanded_weights;
124};
125
126template<class EigenVector_, typename Index_, typename Block_>
127BlockingDetails<Index_, EigenVector_> compute_blocking_details(
128 const Index_ ncells,
129 const Block_* block,
130 const scran_blocks::WeightPolicy block_weight_policy,
131 const scran_blocks::VariableWeightParameters& variable_block_weight_parameters)
132{
133 BlockingDetails<Index_, EigenVector_> output;
134 output.block_size = tatami_stats::tabulate_groups(block, ncells);
135 if (block_weight_policy == scran_blocks::WeightPolicy::NONE) {
136 return output;
137 }
138
139 const auto& block_size = output.block_size;
140 const auto nblocks = block_size.size();
141 output.weighted = true;
142 auto& total_weight = output.total_block_weight;
143 auto& element_weight = output.per_element_weight;
144 sanisizer::resize(element_weight, nblocks);
145
146 for (I<decltype(nblocks)> b = 0; b < nblocks; ++b) {
147 const auto bsize = block_size[b];
148
149 // Computing effective block weights that also incorporate division by the
150 // block size. This avoids having to do the division by block size in the
151 // 'compute_blockwise_mean_and_variance*()' functions.
152 if (bsize) {
153 typename EigenVector_::Scalar block_weight = 1;
154 if (block_weight_policy == scran_blocks::WeightPolicy::VARIABLE) {
155 block_weight = scran_blocks::compute_variable_weight(bsize, variable_block_weight_parameters);
156 }
157
158 element_weight[b] = block_weight / bsize;
159 total_weight += block_weight;
160 } else {
161 element_weight[b] = 0;
162 }
163 }
164
165 // Setting a placeholder value to avoid problems with division by zero.
166 if (total_weight == 0) {
167 total_weight = 1;
168 }
169
170 // Expanding them for multiplication in the IRLBA wrappers.
171 auto sqrt_weights = element_weight;
172 for (auto& s : sqrt_weights) {
173 s = std::sqrt(s);
174 }
175
176 auto& expanded = output.expanded_weights;
177 sanisizer::resize(expanded, ncells);
178 for (Index_ c = 0; c < ncells; ++c) {
179 expanded.coeffRef(c) = sqrt_weights[block[c]];
180 }
181
182 return output;
183}
184
185/*****************************************************************
186 ************ Computing the blockwise mean and variance **********
187 *****************************************************************/
188
189template<typename Num_, typename Value_, typename Index_, typename Block_, typename EigenVector_, typename Float_>
190void compute_sparse_mean_and_variance_blocked(
191 const Num_ num_nonzero,
192 const Value_* values,
193 const Index_* indices,
194 const Block_* block,
195 const BlockingDetails<Index_, EigenVector_>& block_details,
196 Float_* centers,
197 Float_& variance,
198 std::vector<Index_>& block_copy,
199 const Num_ num_all)
200{
201 const auto& block_size = block_details.block_size;
202 const auto nblocks = block_size.size();
203
204 std::fill_n(centers, nblocks, 0);
205 for (Num_ i = 0; i < num_nonzero; ++i) {
206 centers[block[indices[i]]] += values[i];
207 }
208 for (I<decltype(nblocks)> b = 0; b < nblocks; ++b) {
209 auto bsize = block_size[b];
210 if (bsize) {
211 centers[b] /= bsize;
212 }
213 }
214
215 // Computing the variance from the sum of squared differences.
216 // This is technically not the correct variance estimate if we
217 // were to consider the loss of residual d.f. from estimating
218 // the block means, but it's what the PCA sees, so whatever.
219 variance = 0;
220 std::copy(block_size.begin(), block_size.end(), block_copy.begin());
221
222 if (block_details.weighted) {
223 for (Num_ i = 0; i < num_nonzero; ++i) {
224 const Block_ curb = block[indices[i]];
225 const auto diff = values[i] - centers[curb];
226 variance += diff * diff * block_details.per_element_weight[curb];
227 --block_copy[curb];
228 }
229 for (I<decltype(nblocks)> b = 0; b < nblocks; ++b) {
230 const auto val = centers[b];
231 variance += val * val * block_copy[b] * block_details.per_element_weight[b];
232 }
233 } else {
234 for (Num_ i = 0; i < num_nonzero; ++i) {
235 const Block_ curb = block[indices[i]];
236 const auto diff = values[i] - centers[curb];
237 variance += diff * diff;
238 --block_copy[curb];
239 }
240 for (I<decltype(nblocks)> b = 0; b < nblocks; ++b) {
241 const auto val = centers[b];
242 variance += val * val * block_copy[b];
243 }
244 }
245
246 // COMMENT ON DENOMINATOR:
247 // If we're not dealing with weights, we compute the actual sample
248 // variance for easy interpretation (and to match up with the
249 // per-PC calculations in clean_up).
250 //
251 // If we're dealing with weights, the concept of the sample variance
252 // becomes somewhat weird, but we just use the same denominator for
253 // consistency in clean_up_projected. Magnitude doesn't matter when
254 // scaling for process_scale_vector anyway.
255 variance /= num_all - 1;
256}
257
258template<class IrlbaSparseMatrix_, typename Block_, class Index_, class EigenVector_, class EigenMatrix_>
259void compute_blockwise_mean_and_variance_realized_sparse(
260 const IrlbaSparseMatrix_& emat, // this should be column-major with genes in the columns.
261 const Block_* block,
262 const BlockingDetails<Index_, EigenVector_>& block_details,
263 EigenMatrix_& centers,
264 EigenVector_& variances,
265 const int nthreads)
266{
267 const auto ngenes = emat.cols();
268 tatami::parallelize([&](const int, const I<decltype(ngenes)> start, const I<decltype(ngenes)> length) -> void {
269 const auto ncells = emat.rows();
270 const auto& values = emat.get_values();
271 const auto& indices = emat.get_indices();
272 const auto& pointers = emat.get_pointers();
273
274 const auto nblocks = block_details.block_size.size();
275 static_assert(!EigenMatrix_::IsRowMajor);
276 auto block_copy = sanisizer::create<std::vector<Index_> >(nblocks);
277
278 for (I<decltype(start)> g = start, end = start + length; g < end; ++g) {
279 const auto offset = pointers[g];
280 const auto next_offset = pointers[g + 1]; // increment won't overflow as 'g < end' and 'end' is of the same type.
281 compute_sparse_mean_and_variance_blocked(
282 static_cast<I<decltype(ncells)> >(next_offset - offset),
283 values.data() + offset,
284 indices.data() + offset,
285 block,
286 block_details,
287 centers.data() + sanisizer::product_unsafe<std::size_t>(g, nblocks),
288 variances[g],
289 block_copy,
290 ncells
291 );
292 }
293 }, ngenes, nthreads);
294}
295
296template<typename Num_, typename Value_, typename Block_, typename Index_, typename EigenVector_, typename Float_>
297void compute_dense_mean_and_variance_blocked(
298 const Num_ number,
299 const Value_* values,
300 const Block_* block,
301 const BlockingDetails<Index_, EigenVector_>& block_details,
302 Float_* centers,
303 Float_& variance)
304{
305 const auto& block_size = block_details.block_size;
306 const auto nblocks = block_size.size();
307 std::fill_n(centers, nblocks, 0);
308 for (Num_ i = 0; i < number; ++i) {
309 centers[block[i]] += values[i];
310 }
311 for (I<decltype(nblocks)> b = 0; b < nblocks; ++b) {
312 const auto& bsize = block_size[b];
313 if (bsize) {
314 centers[b] /= bsize;
315 }
316 }
317
318 variance = 0;
319
320 if (block_details.weighted) {
321 for (Num_ i = 0; i < number; ++i) {
322 const auto curb = block[i];
323 const auto delta = values[i] - centers[curb];
324 variance += delta * delta * block_details.per_element_weight[curb];
325 }
326 } else {
327 for (Num_ i = 0; i < number; ++i) {
328 const auto curb = block[i];
329 const auto delta = values[i] - centers[curb];
330 variance += delta * delta;
331 }
332 }
333
334 variance /= number - 1; // See COMMENT ON DENOMINATOR above.
335}
336
337template<class EigenMatrix_, typename Block_, class Index_, class EigenVector_>
338void compute_blockwise_mean_and_variance_realized_dense(
339 const EigenMatrix_& emat, // this should be column-major with genes in the columns.
340 const Block_* block,
341 const BlockingDetails<Index_, EigenVector_>& block_details,
342 EigenMatrix_& centers,
343 EigenVector_& variances,
344 const int nthreads)
345{
346 const auto ngenes = emat.cols();
347 tatami::parallelize([&](const int, const I<decltype(ngenes)> start, const I<decltype(ngenes)> length) -> void {
348 const auto ncells = emat.rows();
349 static_assert(!EigenMatrix_::IsRowMajor);
350 const auto nblocks = block_details.block_size.size();
351 for (I<decltype(start)> g = start, end = start + length; g < end; ++g) {
352 compute_dense_mean_and_variance_blocked(
353 ncells,
354 emat.data() + sanisizer::product_unsafe<std::size_t>(g, ncells),
355 block,
356 block_details,
357 centers.data() + sanisizer::product_unsafe<std::size_t>(g, nblocks),
358 variances[g]
359 );
360 }
361 }, ngenes, nthreads);
362}
363
364template<typename Value_, typename Index_, typename Block_, class EigenMatrix_, class EigenVector_>
365void compute_blockwise_mean_and_variance_tatami(
366 const tatami::Matrix<Value_, Index_>& mat, // this should have genes in the rows!
367 const Block_* block,
368 const BlockingDetails<Index_, EigenVector_>& block_details,
369 EigenMatrix_& centers,
370 EigenVector_& variances,
371 const int nthreads)
372{
373 const auto& block_size = block_details.block_size;
374 const auto nblocks = block_size.size();
375 const Index_ ngenes = mat.nrow();
376 const Index_ ncells = mat.ncol();
377
378 if (mat.prefer_rows()) {
379 tatami::parallelize([&](const int, const Index_ start, const Index_ length) -> void {
380 static_assert(!EigenMatrix_::IsRowMajor);
381 auto block_copy = sanisizer::create<std::vector<Index_> >(nblocks);
383
384 if (mat.is_sparse()) {
386 auto ext = tatami::consecutive_extractor<true>(mat, true, start, length);
387 for (Index_ g = start, end = start + length; g < end; ++g) {
388 auto range = ext->fetch(vbuffer.data(), ibuffer.data());
389 compute_sparse_mean_and_variance_blocked(
390 range.number,
391 range.value,
392 range.index,
393 block,
394 block_details,
395 centers.data() + sanisizer::product_unsafe<std::size_t>(g, nblocks),
396 variances[g],
397 block_copy,
398 ncells
399 );
400 }
401
402 } else {
403 auto ext = tatami::consecutive_extractor<false>(mat, true, start, length);
404 for (Index_ g = start, end = start + length; g < end; ++g) {
405 auto ptr = ext->fetch(vbuffer.data());
406 compute_dense_mean_and_variance_blocked(
407 ncells,
408 ptr,
409 block,
410 block_details,
411 centers.data() + sanisizer::product_unsafe<std::size_t>(g, nblocks),
412 variances[g]
413 );
414 }
415 }
416 }, ngenes, nthreads);
417
418 } else {
419 typedef typename EigenVector_::Scalar Scalar;
420 std::vector<std::pair<I<decltype(nblocks)>, Scalar> > block_multipliers;
421 block_multipliers.reserve(nblocks);
422
423 for (I<decltype(nblocks)> b = 0; b < nblocks; ++b) {
424 const auto bsize = block_size[b];
425 if (bsize > 1) { // skipping blocks with NaN variances.
426 Scalar mult = bsize - 1; // need to convert variances back into sum of squared differences.
427 if (block_details.weighted) {
428 mult *= block_details.per_element_weight[b];
429 }
430 block_multipliers.emplace_back(b, mult);
431 }
432 }
433
434 tatami::parallelize([&](const int, const Index_ start, const Index_ length) -> void {
435 std::vector<std::vector<Scalar> > re_centers, re_variances;
436 re_centers.reserve(nblocks);
437 re_variances.reserve(nblocks);
438 for (I<decltype(nblocks)> b = 0; b < nblocks; ++b) {
439 re_centers.emplace_back(length);
440 re_variances.emplace_back(length);
441 }
442
444
445 if (mat.is_sparse()) {
446 std::vector<tatami_stats::variances::RunningSparse<Scalar, Value_, Index_> > running;
447 running.reserve(nblocks);
448 for (I<decltype(nblocks)> b = 0; b < nblocks; ++b) {
449 running.emplace_back(length, re_centers[b].data(), re_variances[b].data(), /* skip_nan = */ false, /* subtract = */ start);
450 }
451
453 auto ext = tatami::consecutive_extractor<true>(mat, false, static_cast<Index_>(0), ncells, start, length);
454 for (Index_ c = 0; c < ncells; ++c) {
455 const auto range = ext->fetch(vbuffer.data(), ibuffer.data());
456 running[block[c]].add(range.value, range.index, range.number);
457 }
458
459 for (I<decltype(nblocks)> b = 0; b < nblocks; ++b) {
460 running[b].finish();
461 }
462
463 } else {
464 std::vector<tatami_stats::variances::RunningDense<Scalar, Value_, Index_> > running;
465 running.reserve(nblocks);
466 for (I<decltype(nblocks)> b = 0; b < nblocks; ++b) {
467 running.emplace_back(length, re_centers[b].data(), re_variances[b].data(), /* skip_nan = */ false);
468 }
469
470 auto ext = tatami::consecutive_extractor<false>(mat, false, static_cast<Index_>(0), ncells, start, length);
471 for (Index_ c = 0; c < ncells; ++c) {
472 auto ptr = ext->fetch(vbuffer.data());
473 running[block[c]].add(ptr);
474 }
475
476 for (I<decltype(nblocks)> b = 0; b < nblocks; ++b) {
477 running[b].finish();
478 }
479 }
480
481 static_assert(!EigenMatrix_::IsRowMajor);
482 for (Index_ i = 0; i < length; ++i) {
483 auto mptr = centers.data() + sanisizer::product_unsafe<std::size_t>(start + i, nblocks);
484 for (I<decltype(nblocks)> b = 0; b < nblocks; ++b) {
485 mptr[b] = re_centers[b][i];
486 }
487
488 auto& my_var = variances[start + i];
489 my_var = 0;
490 for (const auto& bm : block_multipliers) {
491 my_var += re_variances[bm.first][i] * bm.second;
492 }
493 my_var /= ncells - 1; // See COMMENT ON DENOMINATOR above.
494 }
495 }, ngenes, nthreads);
496 }
497}
498
499/******************************************************************
500 ************ Project matrices on their rotation vectors **********
501 ******************************************************************/
502
503template<class EigenMatrix_, class EigenVector_>
504const EigenMatrix_& scale_rotation_matrix(const EigenMatrix_& rotation, bool scale, const EigenVector_& scale_v, EigenMatrix_& tmp) {
505 if (scale) {
506 tmp = (rotation.array().colwise() / scale_v.array()).matrix();
507 return tmp;
508 } else {
509 return rotation;
510 }
511}
512
513template<class IrlbaSparseMatrix_, class EigenMatrix_>
514inline void project_matrix_realized_sparse(
515 const IrlbaSparseMatrix_& emat, // cell in rows, genes in the columns, CSC.
516 EigenMatrix_& components, // dims in rows, cells in columns
517 const EigenMatrix_& scaled_rotation, // genes in rows, dims in columns
518 int nthreads)
519{
520 const auto rank = scaled_rotation.cols();
521 const auto ncells = emat.rows();
522 const auto ngenes = emat.cols();
523
524 // Store as transposed for more cache efficiency.
525 components.resize(
526 sanisizer::cast<I<decltype(components.rows())> >(rank),
527 sanisizer::cast<I<decltype(components.cols())> >(ncells)
528 );
529 components.setZero();
530
531 const auto& values = emat.get_values();
532 const auto& indices = emat.get_indices();
533
534 if (nthreads == 1) {
535 const auto& pointers = emat.get_pointers();
536 auto multipliers = sanisizer::create<Eigen::VectorXd>(rank);
537 for (I<decltype(ngenes)> g = 0; g < ngenes; ++g) {
538 multipliers.noalias() = scaled_rotation.row(g);
539 const auto start = pointers[g], end = pointers[g + 1]; // increment is safe as 'g + 1 <= ngenes'.
540 for (auto i = start; i < end; ++i) {
541 components.col(indices[i]).noalias() += values[i] * multipliers;
542 }
543 }
544
545 } else {
546 const auto& row_nonzero_bounds = emat.get_secondary_nonzero_boundaries();
547 irlba::parallelize(nthreads, [&](const int t) -> void {
548 const auto& starts = row_nonzero_bounds[t];
549 const auto& ends = row_nonzero_bounds[t + 1]; // increment is safe as 't + 1 <= nthreads'.
550 auto multipliers = sanisizer::create<Eigen::VectorXd>(rank);
551
552 for (I<decltype(ngenes)> g = 0; g < ngenes; ++g) {
553 multipliers.noalias() = scaled_rotation.row(g);
554 const auto start = starts[g], end = ends[g];
555 for (auto i = start; i < end; ++i) {
556 components.col(indices[i]).noalias() += values[i] * multipliers;
557 }
558 }
559 });
560 }
561}
562
563template<typename Value_, typename Index_, class EigenMatrix_>
564void project_matrix_transposed_tatami(
565 const tatami::Matrix<Value_, Index_>& mat, // genes in rows, cells in columns
566 EigenMatrix_& components,
567 const EigenMatrix_& scaled_rotation, // genes in rows, dims in columns
568 const int nthreads)
569{
570 const auto rank = scaled_rotation.cols();
571 const auto ngenes = mat.nrow();
572 const auto ncells = mat.ncol();
573 typedef typename EigenMatrix_::Scalar Scalar;
574
575 // Store as transposed for more cache efficiency.
576 components.resize(
577 sanisizer::cast<I<decltype(components.rows())> >(rank),
578 sanisizer::cast<I<decltype(components.cols())> >(ncells)
579 );
580
581 if (mat.prefer_rows()) {
582 tatami::parallelize([&](const int, const Index_ start, const Index_ length) -> void {
583 static_assert(!EigenMatrix_::IsRowMajor);
584 const auto vptr = scaled_rotation.data();
586
587 std::vector<std::vector<Scalar> > local_buffers; // create separate buffers to avoid false sharing.
588 local_buffers.reserve(rank);
589 for (I<decltype(rank)> r = 0; r < rank; ++r) {
590 local_buffers.emplace_back(tatami::cast_Index_to_container_size<I<decltype(local_buffers.front())> >(length));
591 }
592
593 if (mat.is_sparse()) {
595 auto ext = tatami::consecutive_extractor<true>(mat, true, static_cast<Index_>(0), ngenes, start, length);
596 for (Index_ g = 0; g < ngenes; ++g) {
597 const auto range = ext->fetch(vbuffer.data(), ibuffer.data());
598 for (I<decltype(rank)> r = 0; r < rank; ++r) {
599 const auto mult = vptr[sanisizer::nd_offset<std::size_t>(g, ngenes, r)];
600 auto& local_buffer = local_buffers[r];
601 for (Index_ i = 0; i < range.number; ++i) {
602 local_buffer[range.index[i] - start] += range.value[i] * mult;
603 }
604 }
605 }
606
607 } else {
608 auto ext = tatami::consecutive_extractor<false>(mat, true, static_cast<Index_>(0), ngenes, start, length);
609 for (Index_ g = 0; g < ngenes; ++g) {
610 const auto ptr = ext->fetch(vbuffer.data());
611 for (I<decltype(rank)> r = 0; r < rank; ++r) {
612 const auto mult = vptr[sanisizer::nd_offset<std::size_t>(g, ngenes, r)];
613 auto& local_buffer = local_buffers[r];
614 for (Index_ i = 0; i < length; ++i) {
615 local_buffer[i] += ptr[i] * mult;
616 }
617 }
618 }
619 }
620
621 for (I<decltype(rank)> r = 0; r < rank; ++r) {
622 for (Index_ c = 0; c < length; ++c) {
623 components.coeffRef(r, c + start) = local_buffers[r][c];
624 }
625 }
626
627 }, ncells, nthreads);
628
629 } else {
630 tatami::parallelize([&](const int, const Index_ start, const Index_ length) -> void {
631 static_assert(!EigenMatrix_::IsRowMajor);
633
634 if (mat.is_sparse()) {
635 std::vector<Index_> ibuffer(ngenes);
636 auto ext = tatami::consecutive_extractor<true>(mat, false, start, length);
637
638 for (Index_ c = start, end = start + length; c < end; ++c) {
639 const auto range = ext->fetch(vbuffer.data(), ibuffer.data());
640 static_assert(!EigenMatrix_::IsRowMajor);
641 for (I<decltype(rank)> r = 0; r < rank; ++r) {
642 auto& output = components.coeffRef(r, c);
643 output = 0;
644 const auto rotptr = scaled_rotation.data() + sanisizer::product_unsafe<std::size_t>(r, ngenes);
645 for (Index_ i = 0; i < range.number; ++i) {
646 output += rotptr[range.index[i]] * range.value[i];
647 }
648 }
649 }
650
651 } else {
652 auto ext = tatami::consecutive_extractor<false>(mat, false, start, length);
653 for (Index_ c = start, end = start + length; c < end; ++c) {
654 const auto ptr = ext->fetch(vbuffer.data());
655 static_assert(!EigenMatrix_::IsRowMajor);
656 for (I<decltype(rank)> r = 0; r < rank; ++r) {
657 const auto rotptr = scaled_rotation.data() + sanisizer::product_unsafe<std::size_t>(r, ngenes);
658 components.coeffRef(r, c) = std::inner_product(rotptr, rotptr + ngenes, ptr, static_cast<Scalar>(0));
659 }
660 }
661 }
662 }, ncells, nthreads);
663 }
664}
665
666template<class EigenMatrix_, class EigenVector_>
667void clean_up_projected(EigenMatrix_& projected, EigenVector_& D) {
668 // Empirically centering to give nice centered PCs, because we can't
669 // guarantee that the projection is centered in this manner.
670 for (I<decltype(projected.rows())> i = 0, prows = projected.rows(); i < prows; ++i) {
671 projected.row(i).array() -= projected.row(i).sum() / projected.cols();
672 }
673
674 // Just dividing by the number of observations - 1 regardless of weighting.
675 const typename EigenMatrix_::Scalar denom = projected.cols() - 1;
676 for (auto& d : D) {
677 d = d * d / denom;
678 }
679}
680
681/*******************************
682 ***** Residual wrapper ********
683 *******************************/
684
685template<class EigenVector_, class IrlbaMatrix_, typename Block_, class CenterMatrix_>
686class ResidualWorkspace final : public irlba::Workspace<EigenVector_> {
687public:
688 ResidualWorkspace(const IrlbaMatrix_& matrix, const Block_* block, const CenterMatrix_& means) :
689 my_work(matrix.new_known_workspace()),
690 my_block(block),
691 my_means(means),
692 my_sub(sanisizer::cast<I<decltype(my_sub.size())> >(my_means.rows()))
693 {}
694
695private:
696 I<decltype(std::declval<IrlbaMatrix_>().new_known_workspace())> my_work;
697 const Block_* my_block;
698 const CenterMatrix_& my_means;
699 EigenVector_ my_sub;
700
701public:
702 void multiply(const EigenVector_& right, EigenVector_& output) {
703 my_work->multiply(right, output);
704
705 my_sub.noalias() = my_means * right;
706 for (I<decltype(output.size())> i = 0, end = output.size(); i < end; ++i) {
707 auto& val = output.coeffRef(i);
708 val -= my_sub.coeff(my_block[i]);
709 }
710 }
711};
712
713template<class EigenVector_, class IrlbaMatrix_, typename Block_, class CenterMatrix_>
714class ResidualAdjointWorkspace final : public irlba::AdjointWorkspace<EigenVector_> {
715public:
716 ResidualAdjointWorkspace(const IrlbaMatrix_& matrix, const Block_* block, const CenterMatrix_& means) :
717 my_work(matrix.new_known_adjoint_workspace()),
718 my_block(block),
719 my_means(means),
720 my_aggr(sanisizer::cast<I<decltype(my_aggr.size())> >(my_means.rows()))
721 {}
722
723private:
724 I<decltype(std::declval<IrlbaMatrix_>().new_known_adjoint_workspace())> my_work;
725 const Block_* my_block;
726 const CenterMatrix_& my_means;
727 EigenVector_ my_aggr;
728
729public:
730 void multiply(const EigenVector_& right, EigenVector_& output) {
731 my_work->multiply(right, output);
732
733 my_aggr.setZero();
734 for (I<decltype(right.size())> i = 0, end = right.size(); i < end; ++i) {
735 my_aggr.coeffRef(my_block[i]) += right.coeff(i);
736 }
737
738 output.noalias() -= my_means.adjoint() * my_aggr;
739 }
740};
741
742template<class EigenMatrix_, class IrlbaMatrix_, typename Block_, class CenterMatrix_>
743class ResidualRealizeWorkspace final : public irlba::RealizeWorkspace<EigenMatrix_> {
744public:
745 ResidualRealizeWorkspace(const IrlbaMatrix_& matrix, const Block_* block, const CenterMatrix_& means) :
746 my_work(matrix.new_known_realize_workspace()),
747 my_block(block),
748 my_means(means)
749 {}
750
751private:
752 I<decltype(std::declval<IrlbaMatrix_>().new_known_realize_workspace())> my_work;
753 const Block_* my_block;
754 const CenterMatrix_& my_means;
755
756public:
757 const EigenMatrix_& realize(EigenMatrix_& buffer) {
758 my_work->realize_copy(buffer);
759 for (I<decltype(buffer.rows())> i = 0, end = buffer.rows(); i < end; ++i) {
760 buffer.row(i) -= my_means.row(my_block[i]);
761 }
762 return buffer;
763 }
764};
765
766// This wrapper class mimics multiplication with the residuals,
767// i.e., after subtracting the per-block mean from each cell.
768template<class EigenVector_, class EigenMatrix_, class IrlbaMatrixPointer_, class Block_, class CenterMatrixPointer_>
769class ResidualMatrix final : public irlba::Matrix<EigenVector_, EigenMatrix_> {
770public:
771 ResidualMatrix(IrlbaMatrixPointer_ mat, const Block_* block, CenterMatrixPointer_ means) :
772 my_matrix(std::move(mat)),
773 my_block(block),
774 my_means(std::move(means))
775 {}
776
777public:
778 Eigen::Index rows() const {
779 return my_matrix->rows();
780 }
781
782 Eigen::Index cols() const {
783 return my_matrix->cols();
784 }
785
786private:
787 IrlbaMatrixPointer_ my_matrix;
788 const Block_* my_block;
789 CenterMatrixPointer_ my_means;
790
791public:
792 std::unique_ptr<irlba::Workspace<EigenVector_> > new_workspace() const {
793 return new_known_workspace();
794 }
795
796 std::unique_ptr<irlba::AdjointWorkspace<EigenVector_> > new_adjoint_workspace() const {
797 return new_known_adjoint_workspace();
798 }
799
800 std::unique_ptr<irlba::RealizeWorkspace<EigenMatrix_> > new_realize_workspace() const {
801 return new_known_realize_workspace();
802 }
803
804public:
805 std::unique_ptr<ResidualWorkspace<EigenVector_, decltype(*my_matrix), Block_, decltype(*my_means)> > new_known_workspace() const {
806 return std::make_unique<ResidualWorkspace<EigenVector_, decltype(*my_matrix), Block_, decltype(*my_means)> >(*my_matrix, my_block, *my_means);
807 }
808
809 std::unique_ptr<ResidualAdjointWorkspace<EigenVector_, decltype(*my_matrix), Block_, decltype(*my_means)> > new_known_adjoint_workspace() const {
810 return std::make_unique<ResidualAdjointWorkspace<EigenVector_, decltype(*my_matrix), Block_, decltype(*my_means)> >(*my_matrix, my_block, *my_means);
811 }
812
813 std::unique_ptr<ResidualRealizeWorkspace<EigenMatrix_, decltype(*my_matrix), Block_, decltype(*my_means)> > new_known_realize_workspace() const {
814 return std::make_unique<ResidualRealizeWorkspace<EigenMatrix_, decltype(*my_matrix), Block_, decltype(*my_means)> >(*my_matrix, my_block, *my_means);
815 }
816};
827template<typename EigenMatrix_, typename EigenVector_>
835 EigenMatrix_ components;
836
842 EigenVector_ variance_explained;
843
848 typename EigenVector_::Scalar total_variance = 0;
849
855 EigenMatrix_ rotation;
856
862 EigenMatrix_ center;
863
868 EigenVector_ scale;
869
873 bool converged = false;
874};
875
923template<typename Value_, typename Index_, typename Block_, typename EigenMatrix_, class EigenVector_>
926 auto block_details = compute_blocking_details<EigenVector_>(mat.ncol(), block, options.block_weight_policy, options.variable_block_weight_parameters);
927
928 const Index_ ngenes = mat.nrow(), ncells = mat.ncol();
929 const auto nblocks = block_details.block_size.size();
930 output.center.resize(
931 sanisizer::cast<I<decltype(output.center.rows())> >(nblocks),
932 sanisizer::cast<I<decltype(output.center.cols())> >(ngenes)
933 );
934 sanisizer::resize(output.scale, ngenes);
935
936 std::unique_ptr<irlba::Matrix<EigenVector_, EigenMatrix_> > ptr;
937 std::function<void(const EigenMatrix_&)> projector;
938
939 if (!options.realize_matrix) {
940 ptr.reset(new TransposedTatamiWrapperMatrix<EigenVector_, EigenMatrix_, Value_, Index_>(mat, options.num_threads));
941 compute_blockwise_mean_and_variance_tatami(mat, block, block_details, output.center, output.scale, options.num_threads);
942
943 projector = [&](const EigenMatrix_& scaled_rotation) -> void {
944 project_matrix_transposed_tatami(mat, output.components, scaled_rotation, options.num_threads);
945 };
946
947 } else if (mat.sparse()) {
948 // 'extracted' contains row-major contents... but we implicitly transpose it to CSC with genes in columns.
950 mat,
951 /* row = */ true,
952 [&]{
954 opt.two_pass = false;
955 opt.num_threads = options.num_threads;
956 return opt;
957 }()
958 );
959
960 // Storing sparse_ptr in the unique pointer should not invalidate the former,
961 // based on a reading of the C++ specification w.r.t. reset();
962 // so we can continue to use it for projection.
963 const auto sparse_ptr = new irlba::ParallelSparseMatrix<
964 EigenVector_,
965 EigenMatrix_,
966 I<decltype(extracted.value)>,
967 I<decltype(extracted.index)>,
968 I<decltype(extracted.pointers)>
969 >(
970 ncells,
971 ngenes,
972 std::move(extracted.value),
973 std::move(extracted.index),
974 std::move(extracted.pointers),
975 true,
976 options.num_threads
977 );
978 ptr.reset(sparse_ptr);
979
980 compute_blockwise_mean_and_variance_realized_sparse(*sparse_ptr, block, block_details, output.center, output.scale, options.num_threads);
981
982 // Make sure to copy sparse_ptr because it doesn't exist outside of this scope.
983 projector = [&,sparse_ptr](const EigenMatrix_& scaled_rotation) -> void {
984 project_matrix_realized_sparse(*sparse_ptr, output.components, scaled_rotation, options.num_threads);
985 };
986
987 } else {
988 // Perform an implicit transposition by performing a row-major extraction into a column-major transposed matrix.
989 auto tmp_ptr = std::make_unique<EigenMatrix_>(
990 sanisizer::cast<I<decltype(std::declval<EigenMatrix_>().rows())> >(ncells),
991 sanisizer::cast<I<decltype(std::declval<EigenMatrix_>().cols())> >(ngenes)
992 );
993 static_assert(!EigenMatrix_::IsRowMajor);
994
996 mat,
997 /* row_major = */ true,
998 tmp_ptr->data(),
999 [&]{
1000 tatami::ConvertToDenseOptions opt;
1001 opt.num_threads = options.num_threads;
1002 return opt;
1003 }()
1004 );
1005
1006 compute_blockwise_mean_and_variance_realized_dense(*tmp_ptr, block, block_details, output.center, output.scale, options.num_threads);
1007 const auto dense_ptr = tmp_ptr.get(); // do this before the move.
1008 ptr.reset(new irlba::SimpleMatrix<EigenVector_, EigenMatrix_, decltype(tmp_ptr)>(std::move(tmp_ptr)));
1009
1010 // Make sure to copy dense_ptr because it doesn't exist outside of this scope.
1011 projector = [&,dense_ptr](const EigenMatrix_& scaled_rotation) -> void {
1012 output.components.noalias() = (*dense_ptr * scaled_rotation).adjoint();
1013 };
1014 }
1015
1016 output.total_variance = process_scale_vector(options.scale, output.scale);
1017
1018 std::unique_ptr<irlba::Matrix<EigenVector_, EigenMatrix_> > alt;
1019 alt.reset(
1020 new ResidualMatrix<
1021 EigenVector_,
1022 EigenMatrix_,
1023 I<decltype(ptr)>,
1024 Block_,
1025 I<decltype(&(output.center))>
1026 >(
1027 std::move(ptr),
1028 block,
1029 &(output.center)
1030 )
1031 );
1032 ptr.swap(alt);
1033
1034 if (options.scale) {
1035 alt.reset(
1037 EigenVector_,
1038 EigenMatrix_,
1039 I<decltype(ptr)>,
1040 I<decltype(&(output.scale))>
1041 >(
1042 std::move(ptr),
1043 &(output.scale),
1044 /* column = */ true,
1045 /* divide = */ true
1046 )
1047 );
1048 ptr.swap(alt);
1049 }
1050
1051 if (block_details.weighted) {
1052 alt.reset(
1054 EigenVector_,
1055 EigenMatrix_,
1056 I<decltype(ptr)>,
1057 I<decltype(&(block_details.expanded_weights))>
1058 >(
1059 std::move(ptr),
1060 &(block_details.expanded_weights),
1061 /* column = */ false,
1062 /* divide = */ false
1063 )
1064 );
1065 ptr.swap(alt);
1066
1067 auto out = irlba::compute(*ptr, options.number, output.components, output.rotation, output.variance_explained, options.irlba_options);
1068 output.converged = out.first;
1069
1070 EigenMatrix_ tmp;
1071 const auto& scaled_rotation = scale_rotation_matrix(output.rotation, options.scale, output.scale, tmp);
1072 projector(scaled_rotation);
1073
1074 // Subtracting each block's mean from the PCs.
1075 if (options.components_from_residuals) {
1076 EigenMatrix_ centering = (output.center * scaled_rotation).adjoint();
1077 for (I<decltype(ncells)> c =0 ; c < ncells; ++c) {
1078 output.components.col(c) -= centering.col(block[c]);
1079 }
1080 }
1081
1082 clean_up_projected(output.components, output.variance_explained);
1083 if (!options.transpose) {
1084 output.components.adjointInPlace();
1085 }
1086
1087 } else {
1088 const auto out = irlba::compute(*ptr, options.number, output.components, output.rotation, output.variance_explained, options.irlba_options);
1089 output.converged = out.first;
1090
1091 if (options.components_from_residuals) {
1092 clean_up(mat.ncol(), output.components, output.variance_explained);
1093 if (options.transpose) {
1094 output.components.adjointInPlace();
1095 }
1096
1097 } else {
1098 EigenMatrix_ tmp;
1099 const auto& scaled_rotation = scale_rotation_matrix(output.rotation, options.scale, output.scale, tmp);
1100 projector(scaled_rotation);
1101
1102 clean_up_projected(output.components, output.variance_explained);
1103 if (!options.transpose) {
1104 output.components.adjointInPlace();
1105 }
1106 }
1107 }
1108
1109 if (!options.scale) {
1110 output.scale = EigenVector_();
1111 }
1112}
1113
1133template<typename EigenMatrix_ = Eigen::MatrixXd, class EigenVector_ = Eigen::VectorXd, typename Value_, typename Index_, typename Block_>
1136 blocked_pca(mat, block, options, output);
1137 return output;
1138}
1139
1140}
1141
1142#endif
virtual Index_ ncol() const=0
virtual Index_ nrow() const=0
virtual bool prefer_rows() const=0
virtual bool is_sparse() const=0
virtual std::unique_ptr< MyopicSparseExtractor< Value_, Index_ > > sparse(bool row, const Options &opt) const=0
void parallelize(Task_ num_tasks, Run_ run_task)
std::pair< bool, int > compute(const Matrix_ &matrix, const Eigen::Index number, EigenMatrix_ &outU, EigenMatrix_ &outV, EigenVector_ &outD, const Options< EigenVector_ > &options)
double compute_variable_weight(const double s, const VariableWeightParameters &params)
Principal component analysis on single-cell data.
void blocked_pca(const tatami::Matrix< Value_, Index_ > &mat, const Block_ *block, const BlockedPcaOptions &options, BlockedPcaResults< EigenMatrix_, EigenVector_ > &output)
Definition blocked_pca.hpp:924
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)
I< decltype(std::declval< Container_ >().size())> cast_Index_to_container_size(const Index_ x)
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 blocked_pca().
Definition blocked_pca.hpp:31
scran_blocks::VariableWeightParameters variable_block_weight_parameters
Definition blocked_pca.hpp:79
int num_threads
Definition blocked_pca.hpp:98
bool components_from_residuals
Definition blocked_pca.hpp:86
bool realize_matrix
Definition blocked_pca.hpp:92
irlba::Options< Eigen::VectorXd > irlba_options
Definition blocked_pca.hpp:103
bool transpose
Definition blocked_pca.hpp:62
scran_blocks::WeightPolicy block_weight_policy
Definition blocked_pca.hpp:73
int number
Definition blocked_pca.hpp:48
bool scale
Definition blocked_pca.hpp:56
Results of blocked_pca().
Definition blocked_pca.hpp:828
EigenVector_::Scalar total_variance
Definition blocked_pca.hpp:848
bool converged
Definition blocked_pca.hpp:873
EigenMatrix_ components
Definition blocked_pca.hpp:835
EigenMatrix_ rotation
Definition blocked_pca.hpp:855
EigenVector_ scale
Definition blocked_pca.hpp:868
EigenMatrix_ center
Definition blocked_pca.hpp:862
EigenVector_ variance_explained
Definition blocked_pca.hpp:842