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

This is the fully distributed solver. More...

#include <StrumpackSparseSolverMPIDist.hpp>

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

Public Member Functions

 SparseSolverMPIDist (MPI_Comm comm, int argc, char *argv[], bool verbose=true)
 
 SparseSolverMPIDist (MPI_Comm comm, bool verbose=true)
 
 ~SparseSolverMPIDist ()
 
void broadcast_matrix (const CSRMatrix< scalar_t, integer_t > &A)
 
void broadcast_csr_matrix (integer_t N, const integer_t *row_ptr, const integer_t *col_ind, const scalar_t *values, bool symmetric_pattern=false)
 
void set_matrix (const CSRMatrixMPI< scalar_t, integer_t > &A)
 
void set_distributed_csr_matrix (integer_t local_rows, const integer_t *row_ptr, const integer_t *col_ind, const scalar_t *values, const integer_t *dist, bool symmetric_pattern=false)
 
void set_MPIAIJ_matrix (integer_t local_rows, const integer_t *d_ptr, const integer_t *d_ind, const scalar_t *d_val, const integer_t *o_ptr, const integer_t *o_ind, const scalar_t *o_val, const integer_t *garray)
 
void update_matrix_values (integer_t local_rows, const integer_t *row_ptr, const integer_t *col_ind, const scalar_t *values, const integer_t *dist, bool symmetric_pattern=false)
 
void update_matrix_values (const CSRMatrixMPI< scalar_t, integer_t > &A)
 
void update_MPIAIJ_matrix_values (integer_t local_rows, const integer_t *d_ptr, const integer_t *d_ind, const scalar_t *d_val, const integer_t *o_ptr, const integer_t *o_ind, const scalar_t *o_val, const integer_t *garray)
 
MPI_Comm comm () const
 
const MPICommComm () const
 
- Public Member Functions inherited from strumpack::SparseSolverBase< scalar_t, integer_t >
 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>
class strumpack::SparseSolverMPIDist< scalar_t, integer_t >

This is the fully distributed solver.

All steps are distributed: the symbolic factorization, the factorization and the solve. The sparse factors (fill-in) are distributed over the MPI processes. Only the MC64 phase is not distributed. If MC64 reordering is enabled, the sparse input matrix will be gathered to a single process, where then MC64 is called. This can be a bottleneck for large matrices. Try to disable MC64, since for many matrices, this reordering is not required.

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
SparseSolver

Constructor & Destructor Documentation

◆ SparseSolverMPIDist() [1/2]

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

Constructor taking an MPI communicator and command line arguments. This routine is collective on all ranks in the MPI communicator comm.

Parameters
mpi_commMPI communicator. Can be MPI_COMM_WORLD or a subcommunicator.
argcThe number of arguments, i.e, number of elements in the argv array.
argvCommand line arguments. Add -h or –help to have a description printed.
verbFlag to suppres/enable output. Only the root of comm will print to stdout.

◆ SparseSolverMPIDist() [2/2]

template<typename scalar_t , typename integer_t >
strumpack::SparseSolverMPIDist< scalar_t, integer_t >::SparseSolverMPIDist ( MPI_Comm  comm,
bool  verbose = true 
)

Constructor of the SparseSolver class. This routine is collective on all ranks in the MPI communicator comm.

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

◆ ~SparseSolverMPIDist()

template<typename scalar_t , typename integer_t >
strumpack::SparseSolverMPIDist< scalar_t, integer_t >::~SparseSolverMPIDist ( )

Destructor, virtual.

Member Function Documentation

◆ broadcast_csr_matrix()

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

Set a matrix for this sparse solver. Only the matrix provided by the root process (in comm()) will be referenced. The input matrix will immediately be distributed over all the processes in the communicator associated with the solver. This routine is collective on the MPI communicator from this solver.

Parameters
Nnumber of rows/column in the input matrix
row_ptrrow pointers, points to the start of each row in the col_ind and values arrays. Array of size N+1, with row_ptr[N+1] the total number of nonzeros
col_indfor each nonzeros, the column index
valuesnonzero values
symmetric_patterndenote whether the sparsity pattern (not the numerical values) of the sparse inut matrix is symmetric.
See also
set_matrix, broadcast_matrix

◆ broadcast_matrix()

template<typename scalar_t , typename integer_t >
void strumpack::SparseSolverMPIDist< scalar_t, integer_t >::broadcast_matrix ( const CSRMatrix< scalar_t, integer_t > &  A)

Set a matrix for this sparse solver. Only the matrix provided by the root process (in comm()) will be referenced. The input matrix will immediately be distributed over all the processes in the communicator associated with the solver. This routine is collective on the MPI communicator from this solver.

Parameters
Ainput sparse matrix, should only be provided on the root process. The matrix will be copied internally, so it can be safely modified/deleted after calling this function.
See also
set_csr_matrix

◆ comm()

template<typename scalar_t , typename integer_t >
MPI_Comm strumpack::SparseSolverMPIDist< scalar_t, integer_t >::comm ( ) const

Return the MPI_Comm object associated with this solver.

Returns
MPI_Comm object for this solver.

◆ Comm()

template<typename scalar_t , typename integer_t >
const MPIComm & strumpack::SparseSolverMPIDist< scalar_t, integer_t >::Comm ( ) const
inline

Return the const MPIComm object associated with this solver.

Returns
MPIComm object for this solver.

◆ set_distributed_csr_matrix()

template<typename scalar_t , typename integer_t >
void strumpack::SparseSolverMPIDist< scalar_t, integer_t >::set_distributed_csr_matrix ( integer_t  local_rows,
const integer_t *  row_ptr,
const integer_t *  col_ind,
const scalar_t *  values,
const integer_t *  dist,
bool  symmetric_pattern = false 
)

Associate a block-row distributed CSR matrix with the solver object.

Parameters
local_rowsthe number of rows of the input matrix assigned to this MPI process. This should equal to dist[rank+1]-dist[rank].
row_ptrindices in col_ind and values for the start of each row. Nonzeros for row r+dist[rank] are in [row_ptr[r],row_ptr[r+1]).
col_indcolumn indices of each nonzero.
valuesnonzero values. Should have at least (row_ptr[dist[p+1]-dist[p]]-row_ptr[0]) elements.
distspecifies the block-row distribution. A process with rank p owns rows [dist[p],dist[p+1]).
symmetric_patterndenotes whether the sparsity pattern (not the numerical values) of the input matrix is symmetric.
See also
set_matrix, set_csr_matrix

◆ set_matrix()

template<typename scalar_t , typename integer_t >
void strumpack::SparseSolverMPIDist< scalar_t, integer_t >::set_matrix ( const CSRMatrixMPI< scalar_t, integer_t > &  A)

Set a (distributed) matrix for this sparse 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. This routine is collective on the MPI communicator associated with the solver.

Parameters
Ainput sparse matrix, should be provided on all ranks. The matrix will be copied internally, so it can be safely modified/deleted after calling this function.
See also
set_csr_matrix

◆ set_MPIAIJ_matrix()

template<typename scalar_t , typename integer_t >
void strumpack::SparseSolverMPIDist< scalar_t, integer_t >::set_MPIAIJ_matrix ( integer_t  local_rows,
const integer_t *  d_ptr,
const integer_t *  d_ind,
const scalar_t *  d_val,
const integer_t *  o_ptr,
const integer_t *  o_ind,
const scalar_t *  o_val,
const integer_t *  garray 
)

Associate a (PETSc) MPIAIJ block-row distributed CSR matrix with the solver object. See the PETSc manual for a description of MPIAIJ matrix format.

◆ update_matrix_values() [1/2]

template<typename scalar_t , typename integer_t >
void strumpack::SparseSolverMPIDist< scalar_t, integer_t >::update_matrix_values ( const CSRMatrixMPI< 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 (row_ptr and col_ind) and distribution as used to set the initial matrix (using set_matrix or set_distributed_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 >
void strumpack::SparseSolverMPIDist< scalar_t, integer_t >::update_matrix_values ( integer_t  local_rows,
const integer_t *  row_ptr,
const integer_t *  col_ind,
const scalar_t *  values,
const integer_t *  dist,
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) and distribution as used to set the initial matrix (using set_matrix or set_distributed_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
local_rowsNumber of rows of the matrix assigned to this process.
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.
distDescribes the processor distribution, see also set_distributed_csr_matrix
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

◆ update_MPIAIJ_matrix_values()

template<typename scalar_t , typename integer_t >
void strumpack::SparseSolverMPIDist< scalar_t, integer_t >::update_MPIAIJ_matrix_values ( integer_t  local_rows,
const integer_t *  d_ptr,
const integer_t *  d_ind,
const scalar_t *  d_val,
const integer_t *  o_ptr,
const integer_t *  o_ind,
const scalar_t *  o_val,
const integer_t *  garray 
)

This can only be used to UPDATE the nonzero values of the matrix. So it should be called with exactly the same sparsity pattern (d_ptr, d_ind, o_ptr and o_ind) and distribution (garray) as used to set the initial matrix (using set_matrix or set_distributed_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.

See also
set_MPIAIJ_matrix

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