32 #ifndef STRUMPACK_SPARSE_SOLVER_MPI_DIST_HPP 
   33 #define STRUMPACK_SPARSE_SOLVER_MPI_DIST_HPP 
   36 #include "dense/ScaLAPACKWrapper.hpp" 
   37 #include "dense/DistributedVector.hpp" 
   70   template<
typename scalar_t,
typename integer_t>
 
  148     (integer_t N, 
const integer_t* row_ptr, 
const integer_t* col_ind,
 
  149      const scalar_t* values, 
bool symmetric_pattern=
false);
 
  186     (integer_t local_rows, 
const integer_t* row_ptr,
 
  187      const integer_t* col_ind, 
const scalar_t* values,
 
  188      const integer_t* dist, 
bool symmetric_pattern=
false);
 
  196     (integer_t local_rows,
 
  197      const integer_t* d_ptr, 
const integer_t* d_ind, 
const scalar_t* d_val,
 
  198      const integer_t* o_ptr, 
const integer_t* o_ind, 
const scalar_t* o_val,
 
  199      const integer_t* garray);
 
  236     (integer_t local_rows, 
const integer_t* row_ptr,
 
  237      const integer_t* col_ind, 
const scalar_t* values,
 
  238      const integer_t* dist, 
bool symmetric_pattern=
false);
 
  278     (integer_t local_rows, 
const integer_t* d_ptr, 
const integer_t* d_ind,
 
  279      const scalar_t* d_val, 
const integer_t* o_ptr, 
const integer_t* o_ind,
 
  280      const scalar_t* o_val, 
const integer_t* garray);
 
  299     SpMat_t* matrix()
 override { 
return mat_mpi_.get(); }
 
  300     std::unique_ptr<SpMat_t> matrix_nonzero_diag()
 override {
 
  301       using real_t = 
typename RealType<scalar_t>::value_type;
 
  302       return mat_mpi_->add_missing_diagonal
 
  303         (std::sqrt(blas::lamch<real_t>(
'E')) * mat_mpi_->norm1());
 
  305     Reord_t* reordering() 
override;
 
  306     Tree_t* tree()
 override { 
return tree_mpi_dist_.get(); }
 
  307     const SpMat_t* matrix()
 const override { 
return mat_mpi_.get(); }
 
  308     const Reord_t* reordering() 
const override;
 
  309     const Tree_t* tree()
 const override { 
return tree_mpi_dist_.get(); }
 
  311     void setup_tree() 
override;
 
  312     void setup_reordering() 
override;
 
  313     int compute_reordering(
const int* p, 
int base, 
int nx, 
int ny, 
int nz,
 
  314                            int components, 
int width) 
override;
 
  315     void separator_reordering() 
override;
 
  317     void perf_counters_stop(
const std::string& s) 
override;
 
  318     void synchronize()
 override { comm_.
barrier(); }
 
  319     void reduce_flop_counters() 
const override;
 
  321     double max_peak_memory()
 const override {
 
  322       return comm_.
reduce(
double(params::peak_memory), MPI_MAX);
 
  324     double min_peak_memory()
 const override {
 
  325       return comm_.
reduce(
double(params::peak_memory), MPI_MIN);
 
  328     void redistribute_values();
 
  330     void delete_factors_internal() 
override;
 
  333     solve_internal(
const scalar_t* b, scalar_t* x,
 
  334                    bool use_initial_guess=
false) 
override;
 
  336     solve_internal(
const DenseM_t& b, DenseM_t& x,
 
  337                    bool use_initial_guess=
false) 
override;
 
  339     solve_internal(
int nrhs, 
const scalar_t* b, 
int ldb,
 
  340                    scalar_t* x, 
int ldx,
 
  341                    bool use_initial_guess=
false) 
override;
 
  343     std::unique_ptr<CSRMatrixMPI<scalar_t,integer_t>> mat_mpi_;
 
  344     std::unique_ptr<MatrixReorderingMPI<scalar_t,integer_t>> nd_mpi_;
 
  345     std::unique_ptr<EliminationTreeMPIDist<scalar_t,integer_t>> tree_mpi_dist_;
 
  348   template<
typename scalar_t,
typename integer_t>
 
  349   using StrumpackSparseSolverMPIDist =
 
  350     SparseSolverMPIDist<scalar_t,integer_t>;
 
Contains the CSRMatrixMPI class definition, a class representing a block-row distributed compressed s...
Contains the definition of the base (abstract/pure virtual) sparse solver class.
Block-row distributed compressed sparse row storage.
Definition: CSRMatrixMPI.hpp:72
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:1015
This class represents a matrix, stored in column major format, to allow direct use of BLAS/LAPACK rou...
Definition: DenseMatrix.hpp:138
Definition: StrumpackSparseSolverMPIDist.hpp:43
Definition: SparseSolverBase.hpp:53
Wrapper class around an MPI_Comm object.
Definition: MPIWrapper.hpp:194
void barrier() const
Definition: MPIWrapper.hpp:299
T reduce(T t, MPI_Op op) const
Definition: MPIWrapper.hpp:534
Definition: StrumpackSparseSolverMPIDist.hpp:44
Definition: SparseSolverBase.hpp:52
SparseSolverBase is the virtual base for both the sequential/multithreaded and distributed sparse sol...
Definition: SparseSolverBase.hpp:78
This is the fully distributed solver.
Definition: StrumpackSparseSolverMPIDist.hpp:72
SparseSolverMPIDist(MPI_Comm comm, int argc, char *argv[], bool verbose=true)
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 broadcast_matrix(const CSRMatrix< scalar_t, integer_t > &A)
void set_matrix(const CSRMatrixMPI< 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 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)
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)
const MPIComm & Comm() const
Definition: StrumpackSparseSolverMPIDist.hpp:292
SparseSolverMPIDist(MPI_Comm comm, bool verbose=true)
Definition: StrumpackOptions.hpp:43
ReturnCode
Enumeration for the possible return codes.
Definition: StrumpackParameters.hpp:50