|
| 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 MPIComm & | Comm () const |
|
| 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 () |
|
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_t | can be: float, double, std::complex<float> or std::complex<double>. |
integer_t | defaults 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
template<typename scalar_t , typename integer_t >
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
-
A | Sparse matrix, should have the same sparsity pattern as the matrix associated with this solver earlier. |
- See also
- set_csr_matrix, set_matrix
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_rows | Number of rows of the matrix assigned to this process. |
row_ptr | Row 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_ind | Column index array in the typical compressed sparse row representation. This should be the same as used in an earlier call to set_csr_matrix. |
values | Array with numerical nonzero values for the matrix, corresponding to the row_ptr and col_ind compressed sparse row representation. |
dist | Describes the processor distribution, see also set_distributed_csr_matrix |
symmetric_pattern | Denotes 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
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