strumpack::SparseSolver< scalar_t, integer_t > Class Template Reference

SparseSolver is the main sequential or multithreaded sparse solver class. More...

#include <StrumpackSparseSolver.hpp>

Inheritance diagram for strumpack::SparseSolver< scalar_t, integer_t >:
Collaboration diagram for strumpack::SparseSolver< scalar_t, integer_t >:

Public Member Functions

 SparseSolver (int argc, char *argv[], bool verbose=true, bool root=true)
 
 SparseSolver (bool verbose=true, bool root=true)
 
 ~SparseSolver ()
 
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)
 
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 update_matrix_values (const CSRMatrix< scalar_t, integer_t > &A)
 
- Public Member Functions inherited from strumpack::SparseSolverBase< scalar_t, int >
 SparseSolverBase (int argc, char *argv[], bool verbose=true, bool root=true)
 
 SparseSolverBase (bool verbose=true, bool root=true)
 
virtual ~SparseSolverBase ()
 
ReturnCode reorder (int nx=1, int ny=1, int nz=1, int components=1, int width=1)
 
ReturnCode reorder (const int *p, int base=0)
 
ReturnCode factor ()
 
void move_to_gpu ()
 
void remove_from_gpu ()
 
ReturnCode solve (const scalar_t *b, scalar_t *x, bool use_initial_guess=false)
 
ReturnCode solve (const DenseM_t &b, DenseM_t &x, bool use_initial_guess=false)
 
ReturnCode solve (int nrhs, const scalar_t *b, int ldb, scalar_t *x, int ldx, bool use_initial_guess=false)
 
SPOptions< scalar_t > & options ()
 
const SPOptions< scalar_t > & options () const
 
void set_from_options ()
 
void set_from_options (int argc, char *argv[])
 
int maximum_rank () const
 
std::size_t factor_nonzeros () const
 
std::size_t factor_memory () const
 
int Krylov_iterations () const
 
ReturnCode inertia (int &neg, int &zero, int &pos)
 
void draw (const std::string &name) const
 
void delete_factors ()
 

Detailed Description

template<typename scalar_t, typename integer_t = int>
class strumpack::SparseSolver< scalar_t, integer_t >

SparseSolver is the main sequential or multithreaded sparse solver class.

This is the main interface to STRUMPACK's sparse solver. Use this for a sequential or multithreaded sparse solver. For the fully distributed solver, see SparseSolverMPIDist.

Template Parameters
scalar_tcan be: float, double, std::complex<float> or std::complex<double>.
integer_tdefaults to a regular int. If regular int causes 32 bit integer overflows, you should switch to integer_t=int64_t instead. This should be a signed integer type.
See also
SparseSolverMPIDist

Constructor & Destructor Documentation

◆ SparseSolver() [1/2]

template<typename scalar_t , typename integer_t = int>
strumpack::SparseSolver< scalar_t, integer_t >::SparseSolver ( int  argc,
char *  argv[],
bool  verbose = true,
bool  root = true 
)

Constructor of the SparseSolver class, taking command line arguments.

Parameters
argcnumber of arguments, i.e, number of elements in the argv array
argvcommand line arguments. Add -h or –help to have a description printed
verboseflag to enable/disable output to cout
rootflag to denote whether this process is the root MPI process, only the root will print certain messages to cout

◆ SparseSolver() [2/2]

template<typename scalar_t , typename integer_t = int>
strumpack::SparseSolver< scalar_t, integer_t >::SparseSolver ( bool  verbose = true,
bool  root = true 
)

Constructor of the SparseSolver class.

Parameters
verboseflag to enable/disable output to cout
rootflag to denote whether this process is the root MPI process. Only the root will print certain messages
See also
set_from_options

◆ ~SparseSolver()

template<typename scalar_t , typename integer_t = int>
strumpack::SparseSolver< scalar_t, integer_t >::~SparseSolver ( )

(Virtual) destructor of the SparseSolver class.

Member Function Documentation

◆ set_csr_matrix()

template<typename scalar_t , typename integer_t = int>
void strumpack::SparseSolver< scalar_t, integer_t >::set_csr_matrix ( integer_t  N,
const integer_t *  row_ptr,
const integer_t *  col_ind,
const scalar_t *  values,
bool  symmetric_pattern = false 
)

Associate a (sequential) NxN CSR matrix with this solver.

This matrix will not be modified. An internal copy will be made, so it is safe to delete the data immediately after calling this function. See the manual for a description of the CSR format. You can also use the CSRMatrix class.

Parameters
Nnumber of rows and columns of the CSR input matrix.
row_ptrindices in col_ind and values for the start of each row. Nonzeros for row r are in [row_ptr[r],row_ptr[r+1])
col_indcolumn indices of each nonzero
valuesnonzero values
symmetric_patterndenotes whether the sparsity pattern of the input matrix is symmetric, does not require the matrix values to be symmetric
See also
set_matrix

◆ set_matrix()

template<typename scalar_t , typename integer_t = int>
void strumpack::SparseSolver< scalar_t, integer_t >::set_matrix ( const CSRMatrix< scalar_t, integer_t > &  A)

Associate a (sequential) CSRMatrix with this solver.

This matrix will not be modified. An internal copy will be made, so it is safe to delete the data immediately after calling this function.

Parameters
AA CSRMatrix<scalar_t,integer_t> object, will internally be duplicated
See also
set_csr_matrix

◆ update_matrix_values() [1/2]

template<typename scalar_t , typename integer_t = int>
void strumpack::SparseSolver< scalar_t, integer_t >::update_matrix_values ( const CSRMatrix< scalar_t, integer_t > &  A)

This can only be used to UPDATE the nonzero values of the matrix. So it should be called with exactly the same sparsity pattern as used to set the initial matrix (using set_matrix or set_csr_matrix). This routine can be called after having performed a factorization of a different matrix with the same sparsity pattern. In that case, when this solver is used for another solve, with the updated matrix values, the permutation vector previously computed will be reused to permute the updated matrix values, instead of recomputing the permutation. The numerical factorization will automatically be redone.

Parameters
ASparse matrix, should have the same sparsity pattern as the matrix associated with this solver earlier.
See also
set_csr_matrix, set_matrix

◆ update_matrix_values() [2/2]

template<typename scalar_t , typename integer_t = int>
void strumpack::SparseSolver< scalar_t, integer_t >::update_matrix_values ( integer_t  N,
const integer_t *  row_ptr,
const integer_t *  col_ind,
const scalar_t *  values,
bool  symmetric_pattern = false 
)

This can only be used to UPDATE the nonzero values of the matrix. So it should be called with exactly the same sparsity pattern (row_ptr and col_ind) as used to set the initial matrix (using set_matrix or set_csr_matrix). This routine can be called after having performed a factorization of a different matrix with the same sparsity pattern. In that case, when this solver is used for another solve, with the updated matrix values, the permutation vector previously computed will be reused to permute the updated matrix values, instead of recomputing the permutation. The numerical factorization will automatically be redone.

Parameters
NNumber of rows in the matrix.
row_ptrRow pointer array in the typical compressed sparse row representation. This should be the same as used in an earlier call to set_csr_matrix.
col_indColumn index array in the typical compressed sparse row representation. This should be the same as used in an earlier call to set_csr_matrix.
valuesArray with numerical nonzero values for the matrix, corresponding to the row_ptr and col_ind compressed sparse row representation.
symmetric_patternDenotes whether the sparsity pattern of the input matrix is symmetric, does not require the matrix values to be symmetric
See also
set_csr_matrix, set_matrix

The documentation for this class was generated from the following file: