CompressedSparseMatrix.hpp
Go to the documentation of this file.
1/*
2 * STRUMPACK -- STRUctured Matrices PACKage, Copyright (c) 2014, The
3 * Regents of the University of California, through Lawrence Berkeley
4 * National Laboratory (subject to receipt of any required approvals
5 * from the U.S. Dept. of Energy). All rights reserved.
6 *
7 * If you have questions about your rights to use or distribute this
8 * software, please contact Berkeley Lab's Technology Transfer
9 * Department at TTD@lbl.gov.
10 *
11 * NOTICE. This software is owned by the U.S. Department of Energy. As
12 * such, the U.S. Government has been granted for itself and others
13 * acting on its behalf a paid-up, nonexclusive, irrevocable,
14 * worldwide license in the Software to reproduce, prepare derivative
15 * works, and perform publicly and display publicly. Beginning five
16 * (5) years after the date permission to assert copyright is obtained
17 * from the U.S. Department of Energy, and subject to any subsequent
18 * five (5) year renewals, the U.S. Government is granted for itself
19 * and others acting on its behalf a paid-up, nonexclusive,
20 * irrevocable, worldwide license in the Software to reproduce,
21 * prepare derivative works, distribute copies to the public, perform
22 * publicly and display publicly, and to permit others to do so.
23 *
24 * Developers: Pieter Ghysels, Francois-Henry Rouet, Xiaoye S. Li.
25 * (Lawrence Berkeley National Lab, Computational Research
26 * Division).
27 *
28 */
34#ifndef COMPRESSED_SPARSE_MATRIX_HPP
35#define COMPRESSED_SPARSE_MATRIX_HPP
36
37#include <vector>
38#include <string>
39#include <tuple>
40
41#include "misc/Tools.hpp"
42#include "misc/Triplet.hpp"
43#include "dense/DenseMatrix.hpp"
44#include "StrumpackOptions.hpp"
45
46namespace strumpack {
47
48 // forward declarations
49 template<typename integer_t> class CSRGraph;
50 template<typename scalar_t> class DenseMatrix;
51 template<typename scalar_t> class DistributedMatrix;
52
53
54 template<typename scalar_t, typename integer_t,
55 typename real_t = typename RealType<scalar_t>::value_type>
57 public:
58 MatchingData() {}
59 MatchingData(MatchingJob j, std::size_t n) : job(j) {
60 if (job != MatchingJob::NONE)
61 Q.resize(n);
63 R.resize(n);
64 C.resize(n);
65 }
66 }
67
69 std::vector<integer_t> Q;
70 std::vector<real_t> R, C;
71
72 integer_t mc64_work_int(std::size_t n, std::size_t nnz) const {
73 switch (job) {
75 case MatchingJob::MAX_SMALLEST_DIAGONAL_2: return 10*n + nnz;
79 default: return 5*n;
80 }
81 }
82
83 integer_t mc64_work_double(std::size_t n, std::size_t nnz) const {
84 switch (job) {
85 case MatchingJob::MAX_CARDINALITY: return 0;
88 case MatchingJob::MAX_DIAGONAL_SUM: return 2*n + nnz;
90 default: return 3*n + nnz;
91 }
92 }
93 };
94
95 template<typename scalar_t,
96 typename real_t = typename RealType<scalar_t>::value_type>
98 public:
99 Equilibration() {}
100 Equilibration(std::size_t N) : R(N), C(N) {}
101 // use this one for block-row distributed sparse matrix, with rows
102 // local rows, and cols global columns
103 Equilibration(std::size_t rows, std::size_t cols) : R(rows), C(cols) {}
104
105 int info = 0;
106 EquilibrationType type = EquilibrationType::NONE;
107 real_t rcond = 1, ccond = 1, Amax = 0;
108 std::vector<real_t> R, C;
109
110 void set_type() {
111 const real_t thres = 0.1;
112 const real_t small = blas::lamch<real_t>('S') / blas::lamch<real_t>('P');
113 const real_t large = 1. / small;
114 if (rcond >= thres && Amax >= small && Amax <= large) {
115 R.clear();
116 if (ccond >= thres) {
117 type = EquilibrationType::NONE;
118 C.clear();
119 } else type = EquilibrationType::COLUMN;
120 } else {
121 if (ccond >= thres) {
122 type = EquilibrationType::ROW;
123 C.clear();
124 } else type = EquilibrationType::BOTH;
125 }
126 }
127 };
128
148 template<typename scalar_t,typename integer_t>
151 using real_t = typename RealType<scalar_t>::value_type;
152#if defined(STRUMPACK_USE_MPI)
154#endif
157
158 public:
163
169 integer_t size() const { return n_; }
170
176 integer_t nnz() const { return nnz_; }
177
187 const integer_t* ptr() const { return ptr_.data(); }
188
196 const integer_t* ind() const { return ind_.data(); }
197
205 const scalar_t* val() const { return val_.data(); }
206
215 integer_t* ptr() { return ptr_.data(); }
216
224 integer_t* ind() { return ind_.data(); }
225
232 scalar_t* val() { return val_.data(); }
233
239 const integer_t& ptr(integer_t i) const { assert(i >= 0 && i <= size()); return ptr_[i]; }
240
247 const integer_t& ind(integer_t i) const { assert(i >= 0 && i < nnz()); return ind_[i]; }
248
254 const scalar_t& val(integer_t i) const { assert(i >= 0 && i < nnz()); return val_[i]; }
255
261 integer_t& ptr(integer_t i) { assert(i <= size()); return ptr_[i]; }
262
269 integer_t& ind(integer_t i) { assert(i < nnz()); return ind_[i]; }
270
276 scalar_t& val(integer_t i) { assert(i < nnz()); return val_[i]; }
277
278 virtual real_t norm1() const = 0; //{ assert(false); return -1.; };
279
284 bool symm_sparse() const { return symm_sparse_; }
285
293 void set_symm_sparse(bool symm_sparse=true) { symm_sparse_ = symm_sparse; }
294
295
309 virtual void spmv(const DenseM_t& x, DenseM_t& y) const = 0;
310
324 virtual void spmv(const scalar_t* x, scalar_t* y) const = 0;
325
330 virtual void permute(const integer_t* iorder, const integer_t* order);
331
332 virtual void permute(const std::vector<integer_t>& iorder,
333 const std::vector<integer_t>& order) {
334 permute(iorder.data(), order.data());
335 }
336
337 virtual void permute_columns(const std::vector<integer_t>& perm) = 0;
338
339 virtual Equil_t equilibration() const { return Equil_t(this->size()); }
340
341 virtual void equilibrate(const Equil_t&) {}
342
343 virtual Match_t matching(MatchingJob, bool apply=true);
344
345 virtual void apply_matching(const Match_t&);
346
347 virtual void symmetrize_sparsity();
348
349 virtual void print() const;
350 virtual void print_dense(const std::string& name) const {
351 std::cerr << "print_dense not implemented for this matrix type"
352 << std::endl;
353 }
354 virtual void print_matrix_market(const std::string& filename) const {
355 std::cerr << "print_matrix_market not implemented for this matrix type"
356 << std::endl;
357 }
358
359 virtual int read_matrix_market(const std::string& filename) = 0;
360
361 virtual real_t max_scaled_residual(const scalar_t* x,
362 const scalar_t* b) const = 0;
363 virtual real_t max_scaled_residual(const DenseM_t& x,
364 const DenseM_t& b) const = 0;
365
366#ifndef DOXYGEN_SHOULD_SKIP_THIS
367 virtual CSRGraph<integer_t>
368 extract_graph(int ordering_level, integer_t lo, integer_t hi) const = 0;
369 virtual CSRGraph<integer_t>
370 extract_graph_sep_CB(int ordering_level, integer_t lo, integer_t hi,
371 const std::vector<integer_t>& upd) const = 0;
372 virtual CSRGraph<integer_t>
373 extract_graph_CB_sep(int ordering_level, integer_t lo, integer_t hi,
374 const std::vector<integer_t>& upd) const = 0;
375 virtual CSRGraph<integer_t>
376 extract_graph_CB(int ordering_level,
377 const std::vector<integer_t>& upd) const = 0;
378
379 virtual void
380 extract_separator(integer_t sep_end, const std::vector<std::size_t>& I,
381 const std::vector<std::size_t>& J, DenseM_t& B,
382 int depth) const = 0;
383 virtual void
384 extract_front(DenseM_t& F11, DenseM_t& F12, DenseM_t& F21,
385 integer_t slo, integer_t shi,
386 const std::vector<integer_t>& upd,
387 int depth) const = 0;
388 virtual void
389 push_front_elements(integer_t, integer_t, const std::vector<integer_t>&,
390 std::vector<Triplet<scalar_t>>&,
391 std::vector<Triplet<scalar_t>>&,
392 std::vector<Triplet<scalar_t>>&) const = 0;
393 virtual void
394 set_front_elements(integer_t, integer_t, const std::vector<integer_t>&,
395 Triplet<scalar_t>*, Triplet<scalar_t>*,
396 Triplet<scalar_t>*) const = 0;
397 virtual void
398 count_front_elements(integer_t, integer_t, const std::vector<integer_t>&,
399 std::size_t&, std::size_t&, std::size_t&) const = 0;
400
401 virtual void
402 front_multiply(integer_t slo, integer_t shi,
403 const std::vector<integer_t>& upd,
404 const DenseM_t& R, DenseM_t& Sr, DenseM_t& Sc,
405 int depth) const = 0;
406
407 virtual void
408 front_multiply_F11(Trans op, integer_t slo, integer_t shi,
409 const DenseM_t& R, DenseM_t& S, int depth) const = 0;
410 virtual void
411 front_multiply_F12(Trans op, integer_t slo, integer_t shi,
412 const std::vector<integer_t>& upd,
413 const DenseM_t& R, DenseM_t& S, int depth) const = 0;
414 virtual void
415 front_multiply_F21(Trans op, integer_t slo, integer_t shi,
416 const std::vector<integer_t>& upd,
417 const DenseM_t& R, DenseM_t& S, int depth) const = 0;
418
419#if defined(STRUMPACK_USE_MPI)
420 virtual void
421 extract_F11_block(scalar_t* F, integer_t ldF,
422 integer_t row, integer_t nr_rows,
423 integer_t col, integer_t nr_cols) const = 0;
424 virtual void
425 extract_F12_block(scalar_t* F, integer_t ldF,
426 integer_t row, integer_t nr_rows,
427 integer_t col, integer_t nr_cols,
428 const integer_t* upd) const = 0;
429 virtual void
430 extract_F21_block(scalar_t* F, integer_t ldF,
431 integer_t row, integer_t nr_rows,
432 integer_t col, integer_t nr_cols,
433 const integer_t* upd) const = 0;
434 virtual void
435 extract_separator_2d(integer_t sep_end,
436 const std::vector<std::size_t>& I,
437 const std::vector<std::size_t>& J,
438 DistM_t& B) const = 0;
439 virtual void
440 front_multiply_2d(integer_t sep_begin, integer_t sep_end,
441 const std::vector<integer_t>& upd,
442 const DistM_t& R, DistM_t& Srow, DistM_t& Scol,
443 int depth) const = 0;
444 virtual void
445 front_multiply_2d(Trans op, integer_t sep_begin, integer_t sep_end,
446 const std::vector<integer_t>& upd, const DistM_t& R,
447 DistM_t& S, int depth) const = 0;
448#endif //STRUMPACK_USE_MPI
449#endif //DOXYGEN_SHOULD_SKIP_THIS
450
451 protected:
452 integer_t n_, nnz_;
453 std::vector<integer_t> ptr_, ind_;
454 std::vector<scalar_t> val_;
455 bool symm_sparse_;
456
457 enum MMsym {GENERAL, SYMMETRIC, SKEWSYMMETRIC, HERMITIAN};
458
459 CompressedSparseMatrix();
460 CompressedSparseMatrix(integer_t n, integer_t nnz,
461 bool symm_sparse=false);
462 CompressedSparseMatrix(integer_t n,
463 const integer_t* row_ptr,
464 const integer_t* col_ind,
465 const scalar_t* values, bool symm_sparsity);
466
467 std::vector<std::tuple<integer_t,integer_t,scalar_t>>
468 read_matrix_market_entries(const std::string& filename);
469
470 virtual int strumpack_mc64(MatchingJob, Match_t&) { return 0; }
471
472 virtual void scale(const std::vector<scalar_t>&,
473 const std::vector<scalar_t>&) = 0;
474 virtual void scale_real(const std::vector<real_t>&,
475 const std::vector<real_t>&) = 0;
476
477 long long spmv_flops() const;
478 long long spmv_bytes() const;
479 };
480
481} //end namespace strumpack
482
483#endif // COMPRESSED_SPARSE_MATRIX_HPP
Contains the DenseMatrix and DenseMatrixWrapper classes, simple wrappers around BLAS/LAPACK style den...
Holds options for the sparse solver.
Definition: CompressedSparseMatrix.hpp:49
Abstract base class for compressed sparse matrix storage.
Definition: CompressedSparseMatrix.hpp:149
const integer_t & ptr(integer_t i) const
Definition: CompressedSparseMatrix.hpp:239
bool symm_sparse() const
Definition: CompressedSparseMatrix.hpp:284
integer_t * ind()
Definition: CompressedSparseMatrix.hpp:224
scalar_t & val(integer_t i)
Definition: CompressedSparseMatrix.hpp:276
scalar_t * val()
Definition: CompressedSparseMatrix.hpp:232
integer_t nnz() const
Definition: CompressedSparseMatrix.hpp:176
virtual void spmv(const scalar_t *x, scalar_t *y) const =0
const integer_t * ptr() const
Definition: CompressedSparseMatrix.hpp:187
const scalar_t & val(integer_t i) const
Definition: CompressedSparseMatrix.hpp:254
const scalar_t * val() const
Definition: CompressedSparseMatrix.hpp:205
integer_t size() const
Definition: CompressedSparseMatrix.hpp:169
integer_t & ind(integer_t i)
Definition: CompressedSparseMatrix.hpp:269
integer_t & ptr(integer_t i)
Definition: CompressedSparseMatrix.hpp:261
virtual ~CompressedSparseMatrix()
Definition: CompressedSparseMatrix.hpp:162
virtual void permute(const integer_t *iorder, const integer_t *order)
virtual void spmv(const DenseM_t &x, DenseM_t &y) const =0
void set_symm_sparse(bool symm_sparse=true)
Definition: CompressedSparseMatrix.hpp:293
const integer_t & ind(integer_t i) const
Definition: CompressedSparseMatrix.hpp:247
const integer_t * ind() const
Definition: CompressedSparseMatrix.hpp:196
integer_t * ptr()
Definition: CompressedSparseMatrix.hpp:215
This class represents a matrix, stored in column major format, to allow direct use of BLAS/LAPACK rou...
Definition: DenseMatrix.hpp:138
2D block cyclicly distributed matrix, as used by ScaLAPACK.
Definition: DistributedMatrix.hpp:84
Definition: CompressedSparseMatrix.hpp:97
Definition: CompressedSparseMatrix.hpp:56
Definition: StrumpackOptions.hpp:43
MatchingJob
Definition: StrumpackOptions.hpp:119
Trans
Definition: DenseMatrix.hpp:51