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