strumpack::SparseSolverBase< scalar_t, integer_t > Class Template Referenceabstract

SparseSolverBase is the virtual base for both the sequential/multithreaded and distributed sparse solver classes. More...

#include <SparseSolverBase.hpp>

Public Member Functions

 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 (integer_t &neg, integer_t &zero, integer_t &pos)
 
void draw (const std::string &name) const
 
void delete_factors ()
 

Detailed Description

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

SparseSolverBase is the virtual base for both the sequential/multithreaded and distributed sparse solver classes.

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, SparseSolverMixedPrecision

Constructor & Destructor Documentation

◆ SparseSolverBase() [1/2]

template<typename scalar_t , typename integer_t = int>
strumpack::SparseSolverBase< scalar_t, integer_t >::SparseSolverBase ( 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

◆ SparseSolverBase() [2/2]

template<typename scalar_t , typename integer_t = int>
strumpack::SparseSolverBase< scalar_t, integer_t >::SparseSolverBase ( 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

◆ ~SparseSolverBase()

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

(Virtual) destructor of the SparseSolver class.

Member Function Documentation

◆ delete_factors()

template<typename scalar_t , typename integer_t = int>
void strumpack::SparseSolverBase< scalar_t, integer_t >::delete_factors ( )

Free memory held by the factors. This does not delete the symbolic analysis information. So after calling this routine, one can still update the sparse matrix values using for instance update_matrix_values.

◆ draw()

template<typename scalar_t , typename integer_t = int>
void strumpack::SparseSolverBase< scalar_t, integer_t >::draw ( const std::string &  name) const

Create a gnuplot script to draw/plot the sparse factors. Only do this for small matrices! It is very slow!

Parameters
namefilename of the generated gnuplot script. Running
gnuplot plotname.gnuplot 
will generate a pdf file.

◆ factor()

template<typename scalar_t , typename integer_t = int>
ReturnCode strumpack::SparseSolverBase< scalar_t, integer_t >::factor ( )

Perform numerical factorization of the sparse input matrix.

This is the computationally expensive part of the code. However, once the factorization is performed, it can be reused for multiple solves. This requires that a valid matrix has been assigned to this solver. internally this check whether the matrix has been reordered (with a call to reorder), and if not, it will call reorder.

Returns
error code
See also
set_matrix, set_csr_matrix

◆ factor_memory()

template<typename scalar_t , typename integer_t = int>
std::size_t strumpack::SparseSolverBase< scalar_t, integer_t >::factor_memory ( ) const

Return the amount of memory taken by the sparse factorization factors. This is the fill-in. It is simply computed as factor_nonzeros() * sizeof(scalar_t), so it does not include any overhead from the metadata for the datastructures. This should be called after the factorization. For the SparseSolverMPI and SparseSolverMPIDist distributed memory solvers, this routine is collective on the MPI communicator.

◆ factor_nonzeros()

template<typename scalar_t , typename integer_t = int>
std::size_t strumpack::SparseSolverBase< scalar_t, integer_t >::factor_nonzeros ( ) const

Return the number of nonzeros in the (sparse) factors. This is known as the fill-in. This should be called after computing the numerical factorization. For the SparseSolverMPI and SparseSolverMPIDist distributed memory solvers, this routine is collective on the MPI communicator.

◆ inertia()

template<typename scalar_t , typename integer_t = int>
ReturnCode strumpack::SparseSolverBase< scalar_t, integer_t >::inertia ( integer_t &  neg,
integer_t &  zero,
integer_t &  pos 
)

Return the inertia of the matrix. A sparse matrix needs to be set before inertia can be computed. The matrix needs to be factored. If this->factor() was not called already, then it is called inside the inertia routine.

To get accurate inertia the matching needs to be disabled, because the matching applies a non-symmetric permutation. Matching can be disabled using this->options().set_matching(strumpack::MatchingJob::NONE);

The inertia will not be correct if pivoting was performed, in which case the return value will be ReturnCode::INACCURATE_INERTIA. Inertia also cannot be computed when compression is applied (fi, HSS, HODLR, ...).

Parameters
negnumber of negative eigenvalues (if return value is ReturnCode::SUCCESS)
zeronumber of zero eigenvalues (if return value is ReturnCode::SUCCESS)
posnumber of positive eigenvalues (if return value is ReturnCode::SUCCESS)

◆ Krylov_iterations()

template<typename scalar_t , typename integer_t = int>
int strumpack::SparseSolverBase< scalar_t, integer_t >::Krylov_iterations ( ) const

Return the number of iterations performed by the outer (Krylov) iterative solver. Call this after calling the solve routine.

◆ maximum_rank()

template<typename scalar_t , typename integer_t = int>
int strumpack::SparseSolverBase< scalar_t, integer_t >::maximum_rank ( ) const

Return the maximum rank encountered in any of the HSS matrices used to compress the sparse triangular factors. This should be called after the factorization phase. For the SparseSolverMPI and SparseSolverMPIDist distributed memory solvers, this routine is collective on the MPI communicator.

◆ move_to_gpu()

template<typename scalar_t , typename integer_t = int>
void strumpack::SparseSolverBase< scalar_t, integer_t >::move_to_gpu ( )

OBSOLETE if the factors fit in device memory they will already be on the device

◆ options() [1/2]

template<typename scalar_t , typename integer_t = int>
SPOptions< scalar_t > & strumpack::SparseSolverBase< scalar_t, integer_t >::options ( )

Return the object holding the options for this sparse solver.

◆ options() [2/2]

template<typename scalar_t , typename integer_t = int>
const SPOptions< scalar_t > & strumpack::SparseSolverBase< scalar_t, integer_t >::options ( ) const

Return the object holding the options for this sparse solver.

◆ remove_from_gpu()

template<typename scalar_t , typename integer_t = int>
void strumpack::SparseSolverBase< scalar_t, integer_t >::remove_from_gpu ( )

TODO implement this to clear the device memory, but still keep the factors in host memory

◆ reorder() [1/2]

template<typename scalar_t , typename integer_t = int>
ReturnCode strumpack::SparseSolverBase< scalar_t, integer_t >::reorder ( const int *  p,
int  base = 0 
)

Perform sparse matrix reordering, with a user-supplied permutation vector. Using this will ignore the reordering method selected in the options struct.

Parameters
ppermutation vector, should be of size N, the size of the sparse matrix associated with this solver
baseis the permutation 0 or 1 based?

◆ reorder() [2/2]

template<typename scalar_t , typename integer_t = int>
ReturnCode strumpack::SparseSolverBase< scalar_t, integer_t >::reorder ( int  nx = 1,
int  ny = 1,
int  nz = 1,
int  components = 1,
int  width = 1 
)

Compute matrix reorderings for numerical stability and to reduce fill-in.

Start computation of the matrix reorderings. See the relevant options to control the matrix reordering in the manual. A first reordering is the MC64 column permutation for numerical stability. This can be disabled if the matrix has large nonzero diagonal entries. MC64 optionally also performs row and column scaling. Next, a fill-reducing reordering is computed. This is done with the nested dissection algortihms of either (PT-)Scotch, (Par)Metis or a simple geometric nested dissection code which only works on regular meshes.

All parameters here are optional and only used with the geometric ordering, which only makes sense on a regular 1D, 2D or 3D mesh, in the natural ordering.

Parameters
nxthis (optional) parameter is only meaningful when the matrix corresponds to a stencil on a regular mesh. The nx parameter denotes the number of grid points in the first spatial dimension.
nysee parameters nx. Parameter ny denotes the number of gridpoints in the second spatial dimension. This should only be set if the mesh is 2 or 3 dimensional, otherwise it can be 1 (default).
nzSee parameters nx. Parameter nz denotes the number of gridpoints in the third spatial dimension. This should only be set if the mesh is 3 dimensional, otherwise it can be 1 (default).
componentsNumber of degrees of freedom per grid point (default 1)
widthWidth of the stencil, a 1D 3-point stencil needs a separator of width 1, a 1D 5-point wide stencil needs a separator of width 2 (default 1).
Returns
error code
See also
SPOptions

◆ set_from_options() [1/2]

template<typename scalar_t , typename integer_t = int>
void strumpack::SparseSolverBase< scalar_t, integer_t >::set_from_options ( )

Parse the command line options passed in the constructor, and modify the options object accordingly. Run with option -h or –help to see a list of supported options. Or check the SPOptions documentation.

◆ set_from_options() [2/2]

template<typename scalar_t , typename integer_t = int>
void strumpack::SparseSolverBase< scalar_t, integer_t >::set_from_options ( int  argc,
char *  argv[] 
)

Parse the command line options, and modify the options object accordingly. Run with option -h or –help to see a list of supported options. Or check the SPOptions documentation.

Parameters
argcnumber of options in argv
argvlist of options

◆ solve() [1/3]

template<typename scalar_t , typename integer_t = int>
ReturnCode strumpack::SparseSolverBase< scalar_t, integer_t >::solve ( const DenseM_t b,
DenseM_t x,
bool  use_initial_guess = false 
)

Solve a linear system with a single or multiple right-hand sides. Before being able to solve a linear system, the matrix needs to be factored. One can call factor() explicitly, or if this was not yet done, this routine will call factor() internally.

Parameters
binput, will not be modified. DenseMatrix containgin the right-hand side vector/matrix. Should have N rows, with N the dimension of the input matrix for SparseSolver and SparseSolverMPI. For SparseSolverMPIDist, the number or rows of b should be correspond to the partitioning of the block-row distributed input matrix.
xOutput, pointer to the solution vector. Array should be lenght N, the dimension of the input matrix for SparseSolver and SparseSolverMPI. For SparseSolverMPIDist, the length of b should be correspond the partitioning of the block-row distributed input matrix.
use_initial_guessset to true if x contains an intial guess to the solution. This is mainly useful when using an iterative solver. If set to false, x should not be set (but should be allocated).
Returns
error code
See also
DenseMatrix, solve(), factor()

◆ solve() [2/3]

template<typename scalar_t , typename integer_t = int>
ReturnCode strumpack::SparseSolverBase< scalar_t, integer_t >::solve ( const scalar_t *  b,
scalar_t *  x,
bool  use_initial_guess = false 
)

Solve a linear system with a single right-hand side. Before being able to solve a linear system, the matrix needs to be factored. One can call factor() explicitly, or if this was not yet done, this routine will call factor() internally.

Parameters
binput, will not be modified. Pointer to the right-hand side. Array should be lenght N, the dimension of the input matrix for SparseSolver and SparseSolverMPI. For SparseSolverMPIDist, the length of b should be correspond the partitioning of the block-row distributed input matrix.
xOutput, pointer to the solution vector. Array should be lenght N, the dimension of the input matrix for SparseSolver and SparseSolverMPI. For SparseSolverMPIDist, the length of b should be correspond the partitioning of the block-row distributed input matrix.
use_initial_guessset to true if x contains an intial guess to the solution. This is mainly useful when using an iterative solver. If set to false, x should not be set (but should be allocated).
Returns
error code

◆ solve() [3/3]

template<typename scalar_t , typename integer_t = int>
ReturnCode strumpack::SparseSolverBase< scalar_t, integer_t >::solve ( int  nrhs,
const scalar_t *  b,
int  ldb,
scalar_t *  x,
int  ldx,
bool  use_initial_guess = false 
)

Solve a linear system with a single or multiple right-hand sides. Before being able to solve a linear system, the matrix needs to be factored. One can call factor() explicitly, or if this was not yet done, this routine will call factor() internally.

Parameters
nrhsNumber of right hand sides.
binput, will not be modified. DenseMatrix containgin the right-hand side vector/matrix. Should have N rows, with N the dimension of the input matrix for SparseSolver and SparseSolverMPI. For SparseSolverMPIDist, the number or rows of b should be correspond to the partitioning of the block-row distributed input matrix.
ldbleading dimension of b
xOutput, pointer to the solution vector. Array should be lenght N, the dimension of the input matrix for SparseSolver and SparseSolverMPI. For SparseSolverMPIDist, the length of b should be correspond the partitioning of the block-row distributed input matrix.
ldxleading dimension of x
use_initial_guessset to true if x contains an intial guess to the solution. This is mainly useful when using an iterative solver. If set to false, x should not be set (but should be allocated).
Returns
error code
See also
DenseMatrix, solve(), factor()

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