33#ifndef STRUMPACK_SPARSE_SOLVER_HPP
34#define STRUMPACK_SPARSE_SOLVER_HPP
49 template<
typename scalar_t,
typename integer_t>
class MatrixReordering;
50 template<
typename scalar_t,
typename integer_t>
class EliminationTree;
73 template<
typename scalar_t,
typename integer_t=
int>
97 SparseSolver(
int argc,
char* argv[],
bool verbose=
true,
bool root=
true);
148 const integer_t* row_ptr,
const integer_t* col_ind,
149 const scalar_t* values,
bool symmetric_pattern=
false);
197 const integer_t* row_ptr,
198 const integer_t* col_ind,
199 const scalar_t* values,
200 bool symmetric_pattern=
false);
223 void setup_tree()
override;
224 void setup_reordering()
override;
225 int compute_reordering(
const int* p,
int base,
226 int nx,
int ny,
int nz,
227 int components,
int width)
override;
228 void separator_reordering()
override;
230 SpMat_t* matrix()
override {
return mat_.get(); }
231 Reord_t* reordering()
override {
return nd_.get(); }
232 Tree_t* tree()
override {
return tree_.get(); }
233 const SpMat_t* matrix()
const override {
return mat_.get(); }
234 const Reord_t* reordering()
const override {
return nd_.get(); }
235 const Tree_t* tree()
const override {
return tree_.get(); }
237 void permute_matrix_values();
239 ReturnCode solve_internal(
const scalar_t* b, scalar_t* x,
240 bool use_initial_guess=
false)
override;
241 ReturnCode solve_internal(
const DenseM_t& b, DenseM_t& x,
242 bool use_initial_guess=
false)
override;
244 void delete_factors_internal()
override;
246 void transform_x0(DenseM_t& x, DenseM_t& xtmp);
247 void transform_b(
const DenseM_t& b, DenseM_t& bloc);
248 void transform_x(DenseM_t& x, DenseM_t& xtmp);
250 std::unique_ptr<CSRMatrix<scalar_t,integer_t>> mat_;
251 std::unique_ptr<MatrixReordering<scalar_t,integer_t>> nd_;
252 std::unique_ptr<EliminationTree<scalar_t,integer_t>> tree_;
254 using SPBase_t = SparseSolverBase<scalar_t,integer_t>;
255 using SPBase_t::opts_;
256 using SPBase_t::is_root_;
257 using SPBase_t::matching_;
258 using SPBase_t::equil_;
259 using SPBase_t::factored_;
260 using SPBase_t::reordered_;
261 using SPBase_t::Krylov_its_;
262 using SPBase_t::solve_internal;
265 template<
typename scalar_t,
typename integer_t>
266 using StrumpackSparseSolver = SparseSolver<scalar_t,integer_t>;
Contains the definition of the base (abstract/pure virtual) sparse solver class.
Class for storing a compressed sparse row matrix (single node).
Definition CSRMatrix.hpp:55
Abstract base class for compressed sparse matrix storage.
Definition CompressedSparseMatrix.hpp:149
Like DenseMatrix, this class represents a matrix, stored in column major format, to allow direct use ...
Definition DenseMatrix.hpp:1018
This class represents a matrix, stored in column major format, to allow direct use of BLAS/LAPACK rou...
Definition DenseMatrix.hpp:139
Definition SparseSolverBase.hpp:53
Definition SparseSolverBase.hpp:52
SparseSolverBase is the virtual base for both the sequential/multithreaded and distributed sparse sol...
Definition SparseSolverBase.hpp:78
SparseSolver is the main sequential or multithreaded sparse solver class.
Definition StrumpackSparseSolver.hpp:75
void update_matrix_values(const CSRMatrix< scalar_t, integer_t > &A)
void update_matrix_values(integer_t N, const integer_t *row_ptr, const integer_t *col_ind, const scalar_t *values, bool symmetric_pattern=false)
void set_lower_triangle_matrix(const CSRMatrix< scalar_t, integer_t > &A)
SparseSolver(bool verbose=true, bool root=true)
SparseSolver(int argc, char *argv[], bool verbose=true, bool root=true)
void set_matrix(const CSRMatrix< scalar_t, integer_t > &A)
void set_csr_matrix(integer_t N, const integer_t *row_ptr, const integer_t *col_ind, const scalar_t *values, bool symmetric_pattern=false)
Definition StrumpackOptions.hpp:44
ReturnCode
Enumeration for the possible return codes.
Definition StrumpackParameters.hpp:50