12 #include "../config/config.hpp" 14 #ifdef MFEM_USE_STRUMPACK 25 STRUMPACKRowLocMatrix::STRUMPACKRowLocMatrix(MPI_Comm comm,
26 int num_loc_rows,
int first_loc_row,
27 int glob_nrows,
int glob_ncols,
28 int *I,
int *J,
double *data)
29 : comm_(comm), A_(NULL)
37 MPI_Comm_rank(comm_, &rank);
38 MPI_Comm_size(comm_, &nprocs);
39 int * dist =
new int[nprocs + 1];
40 dist[rank + 1] = first_loc_row + num_loc_rows;
42 MPI_Allgather(MPI_IN_PLACE, 0, MPI_INT, dist + 1, 1, MPI_INT, comm_);
43 A_ =
new CSRMatrixMPI<double,int>(num_loc_rows, I, J, data, dist, comm_,
false);
48 : comm_(hypParMat.GetComm()),
52 hypre_ParCSRMatrix * parcsr_op =
53 (hypre_ParCSRMatrix *)const_cast<HypreParMatrix&>(hypParMat);
55 MFEM_ASSERT(parcsr_op != NULL,
"STRUMPACK: const_cast failed in SetOperator");
60 hypre_CSRMatrix * csr_op = hypre_MergeDiagAndOffd(parcsr_op);
62 hypre_CSRMatrixSetDataOwner(csr_op,0);
63 #if MFEM_HYPRE_VERSION >= 21600 68 hypre_CSRMatrixBigJtoJ(csr_op);
72 width = csr_op->num_rows;
75 MPI_Comm_rank(comm_, &rank);
76 MPI_Comm_size(comm_, &nprocs);
77 int * dist =
new int[nprocs + 1];
78 dist[rank + 1] = parcsr_op->first_row_index + csr_op->num_rows;
80 MPI_Allgather(MPI_IN_PLACE, 0, MPI_INT, dist + 1, 1, MPI_INT, comm_);
81 A_ =
new CSRMatrixMPI<double,int>(csr_op->num_rows, csr_op->i, csr_op->j,
82 csr_op->data, dist, comm_,
false);
86 hypre_CSRMatrixDestroy(csr_op);
92 if ( A_ != NULL ) {
delete A_; }
100 this->
Init(argc, argv);
104 : comm_(A.GetComm()),
119 void STRUMPACKSolver::Init(
int argc,
char* argv[] )
127 solver_ =
new StrumpackSparseSolverMPIDist<double,int>(
comm_, argc, argv,
133 solver_->options().set_from_command_line( );
148 solver_->options().set_Krylov_solver( method );
154 solver_->options().set_reordering_method( method );
159 #if STRUMPACK_VERSION_MAJOR >= 3 160 solver_->options().set_matching( strumpack::MatchingJob::NONE );
162 solver_->options().set_mc64job( strumpack::MC64Job::NONE );
168 #if STRUMPACK_VERSION_MAJOR >= 3 169 solver_->options().set_matching
170 ( strumpack::MatchingJob::MAX_DIAGONAL_PRODUCT_SCALING );
173 ( strumpack::MC64Job::MAX_DIAGONAL_PRODUCT_SCALING );
177 #if STRUMPACK_VERSION_MAJOR >= 3 180 solver_->options().set_matching
181 ( strumpack::MatchingJob::COMBBLAS );
187 solver_->options().set_rel_tol( rtol );
192 solver_->options().set_abs_tol( atol );
198 MFEM_ASSERT(
APtr_ != NULL,
199 "STRUMPACK Error: The operator must be set before" 200 " the system can be solved.");
201 MFEM_ASSERT(x.
Size() ==
Width(),
"invalid x.Size() = " << x.
Size()
202 <<
", expected size = " <<
Width());
203 MFEM_ASSERT(y.
Size() ==
Height(),
"invalid y.Size() = " << y.
Size()
204 <<
", expected size = " <<
Height());
210 ReturnCode ret =
solver_->factor();
213 case ReturnCode::SUCCESS:
break;
214 case ReturnCode::MATRIX_NOT_SET:
216 MFEM_ABORT(
"STRUMPACK: Matrix was not set!");
219 case ReturnCode::REORDERING_ERROR:
221 MFEM_ABORT(
"STRUMPACK: Matrix reordering failed!");
226 MFEM_ABORT(
"STRUMPACK: 'factor()' error code = " << ret);
240 mfem_error(
"STRUMPACKSolver::SetOperator : not STRUMPACKRowLocMatrix!");
253 #endif // MFEM_USE_MPI 254 #endif // MFEM_USE_STRUMPACK const STRUMPACKRowLocMatrix * APtr_
virtual const double * HostRead() const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), false).
STRUMPACKSolver(int argc, char *argv[], MPI_Comm comm)
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
int Size() const
Returns the size of the vector.
virtual double * HostWrite()
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), false).
strumpack::StrumpackSparseSolverMPIDist< double, int > * solver_
void SetFromCommandLine()
void HostRead() const
Update the internal hypre_ParCSRMatrix object, A, to be on host.
void SetPrintFactorStatistics(bool print_stat)
STRUMPACKRowLocMatrix(MPI_Comm comm, int num_loc_rows, int first_loc_row, int glob_nrows, int glob_ncols, int *I, int *J, double *data)
void mfem_error(const char *msg)
Function called when an error is encountered. Used by the macros MFEM_ABORT, MFEM_ASSERT, MFEM_VERIFY.
void SetPrintSolveStatistics(bool print_stat)
strumpack::CSRMatrixMPI< double, int > * getA() const
void SetReorderingStrategy(strumpack::ReorderingStrategy method)
void SetRelTol(double rtol)
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
A class to initialize the size of a Tensor.
int height
Dimension of the output / number of rows in the matrix.
void SetKrylovSolver(strumpack::KrylovSolver method)
void HypreRead() const
Update the internal hypre_ParCSRMatrix object, A, to be in hypre memory space.
void SetAbsTol(double atol)
void SetOperator(const Operator &op)
Set/update the solver for the given operator.
void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Wrapper for hypre's ParCSR matrix class.
int width
Dimension of the input / number of columns in the matrix.
void EnableParallelMatching()