MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
strumpack.hpp
Go to the documentation of this file.
1// Copyright (c) 2010-2024, Lawrence Livermore National Security, LLC. Produced
2// at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3// LICENSE and NOTICE for details. LLNL-CODE-806117.
4//
5// This file is part of the MFEM library. For more information and source code
6// availability visit https://mfem.org.
7//
8// MFEM is free software; you can redistribute it and/or modify it under the
9// terms of the BSD-3 license. We welcome feedback and contributions, see file
10// CONTRIBUTING.md for details.
11
12#ifndef MFEM_STRUMPACK
13#define MFEM_STRUMPACK
14
15#include "../config/config.hpp"
16
17#ifdef MFEM_USE_STRUMPACK
18#ifdef MFEM_USE_MPI
19
20#include "operator.hpp"
21#include "hypre.hpp"
22#include <mpi.h>
23
24// STRUMPACK headers
25#include "StrumpackSparseSolverMPIDist.hpp"
26#include "StrumpackSparseSolverMixedPrecisionMPIDist.hpp"
27
28namespace mfem
29{
30
32{
33public:
34 /** @brief Creates a general parallel matrix from a local CSR matrix on each
35 processor.
36
37 The CSR matrix is described by the I, J and data arrays. The local matrix
38 should be of size (local) nrows by (global) glob_ncols. The new parallel
39 matrix contains copies of all input arrays (so they can be deleted). */
40 STRUMPACKRowLocMatrix(MPI_Comm comm,
41 int num_loc_rows, HYPRE_BigInt first_loc_row,
42 HYPRE_BigInt glob_nrows, HYPRE_BigInt glob_ncols,
43 int *I, HYPRE_BigInt *J, double *data,
44 bool sym_sparse = false);
45
46 /** @brief Creates a copy of the parallel matrix hypParMat in STRUMPACK's
47 RowLoc format.
48
49 All data is copied so the original matrix may be deleted. */
50 STRUMPACKRowLocMatrix(const Operator &op, bool sym_sparse = false);
51
53
54 /// Matrix vector products are not supported for this type of matrix.
55 void Mult(const Vector &x, Vector &y) const
56 {
57 MFEM_ABORT("STRUMPACKRowLocMatrix::Mult: Matrix vector products are not "
58 "supported!");
59 }
60
61 /// Get the MPI Comm being used by the parallel matrix.
62 MPI_Comm GetComm() const { return A_->comm(); }
63
64 /// Get access to the internal CSR matrix.
65 strumpack::CSRMatrixMPI<double, HYPRE_BigInt> *GetA() const { return A_; }
66
67private:
68 strumpack::CSRMatrixMPI<double, HYPRE_BigInt> *A_;
69};
70
71/** The MFEM STRUMPACK Direct Solver class.
72
73 The mfem::STRUMPACKSolver class uses the STRUMPACK library to perform LU
74 factorization of a parallel sparse matrix. The solver is capable of handling
75 double precision types. See
76 http://portal.nersc.gov/project/sparse/strumpack/. */
77template <typename STRUMPACKSolverType>
79{
80protected:
81 /** @brief Constructor with MPI_Comm parameter and command line arguments.
82
83 STRUMPACKSolverBase::SetFromCommandLine must be called for the command
84 line arguments to be used. */
85 STRUMPACKSolverBase(MPI_Comm comm, int argc, char *argv[]);
86
87 /** @brief Constructor with STRUMPACK matrix object and command line
88 arguments.
89
90 STRUMPACKSolverBase::SetFromCommandLine must be called for the command
91 line arguments to be used. */
92 STRUMPACKSolverBase(STRUMPACKRowLocMatrix &A, int argc, char *argv[]);
93
94public:
95 /// Default destructor.
96 virtual ~STRUMPACKSolverBase();
97
98 /// Factor and solve the linear system $y = Op^{-1} x $.
99 void Mult(const Vector &x, Vector &y) const;
100
101 /** @brief Factor and solve the linear systems $ Y_i = Op^{-1} X_i $
102 across the array of vectors. */
103 void ArrayMult(const Array<const Vector *> &X, Array<Vector *> &Y) const;
104
105 /** @brief Set the operator/matrix.
106 \note @a A must be a STRUMPACKRowLocMatrix. */
107 void SetOperator(const Operator &op);
108
109 /** @brief Set options that were captured from the command line.
110
111 These were captured in the constructer STRUMPACKSolverBase. Refer
112 to the STRUMPACK documentation for details. */
113 void SetFromCommandLine();
114
115 /// Set up verbose printing during the factor step
116 void SetPrintFactorStatistics(bool print_stat);
117
118 /// Set up verbose printing during the solve step
119 void SetPrintSolveStatistics(bool print_stat);
120
121 /// Set the relative tolerance for interative solvers
122 void SetRelTol(double rtol);
123
124 /// Set the absolute tolerance for iterative solvers
125 void SetAbsTol(double atol);
126
127 /// Set the maximum number of iterations for iterative solvers
128 void SetMaxIter(int max_it);
129
130 /** @brief Set the flag controlling reuse of the symbolic factorization for
131 multiple operators.
132
133 This method must be called before repeated calls to SetOperator. */
134 void SetReorderingReuse(bool reuse);
135
136 /** @brief Enable GPU off-loading available if STRUMPACK was compiled with
137 CUDA.
138 @note Input/Output from MFEM to STRUMPACK is all still through host
139 memory. */
140 void EnableGPU();
141
142 /** @brief Disable GPU off-loading available if STRUMPACK was compiled with
143 CUDA.
144 @note Input/Output from MFEM to STRUMPACK is all still through host
145 memory. */
146 void DisableGPU();
147
148 /** @brief Set the Krylov solver method to use
149 *
150 * STRUMPACK is an (approximate) direct solver. It can be used as a direct
151 * solver or as a preconditioner. To use STRUMPACK as only a preconditioner,
152 * set the Krylov solver to DIRECT. STRUMPACK also provides iterative solvers
153 * which can use the preconditioner, and these iterative solvers can also be
154 * used without preconditioner.
155 *
156 * Supported values are:
157 * - AUTO: Use iterative refinement if no HSS compression is
158 * used, otherwise use GMRES
159 * - DIRECT: No outer iterative solver, just a single application
160 * of the multifrontal solver
161 * - REFINE: Iterative refinement
162 * - PREC_GMRES: Preconditioned GMRES
163 * The preconditioner is the (approx) multifrontal solver
164 * - GMRES: UN-preconditioned GMRES (for testing mainly)
165 * - PREC_BICGSTAB: Preconditioned BiCGStab
166 * The preconditioner is the (approx) multifrontal solver
167 * - BICGSTAB: UN-preconditioned BiCGStab. (for testing mainly)
168 */
169 void SetKrylovSolver(strumpack::KrylovSolver method);
170
171 /** @brief Set matrix reordering strategy
172 *
173 * Supported reorderings are:
174 * - NATURAL: Do not reorder the system
175 * - METIS: Use Metis nested-dissection reordering (default)
176 * - PARMETIS: Use ParMetis nested-dissection reordering
177 * - SCOTCH: Use Scotch nested-dissection reordering
178 * - PTSCOTCH: Use PT-Scotch nested-dissection reordering
179 * - RCM: Use RCM reordering
180 * - GEOMETRIC: A simple geometric nested dissection code that
181 * only works for regular meshes
182 * - AMD: Approximate minimum degree
183 * - MMD: Multiple minimum degree
184 * - AND: Nested dissection
185 * - MLF: Minimum local fill
186 * - SPECTRAL: Spectral nested dissection
187 */
188 void SetReorderingStrategy(strumpack::ReorderingStrategy method);
189
190 /** @brief Configure static pivoting for stability.
191 *
192 * The static pivoting in STRUMPACK
193 * permutes the sparse input matrix in order to get large (nonzero) elements
194 * on the diagonal. If the input matrix is already diagonally dominant, this
195 * reordering can be disabled.
196 *
197 * Supported matching algorithms are:
198 * - NONE: Don't do anything
199 * - MAX_CARDINALITY: Maximum cardinality
200 * - MAX_SMALLEST_DIAGONAL: Maximum smallest diagonal value
201 * - MAX_SMALLEST_DIAGONAL_2: Same as MAX_SMALLEST_DIAGONAL
202 * but different algorithm
203 * - MAX_DIAGONAL_SUM: Maximum sum of diagonal values
204 * - MAX_DIAGONAL_PRODUCT_SCALING: Maximum product of diagonal values
205 * and row and column scaling (default)
206 * - COMBBLAS: Use AWPM from CombBLAS (only with
207 * version >= 3)
208 */
209 void SetMatching(strumpack::MatchingJob job);
210
211 /** @brief Select compression for sparse data types
212 *
213 * Enable support for rank-structured data formats, which can be used
214 * for compression within the sparse solver.
215 *
216 * Supported compression types are:
217 * - NONE: No compression, purely direct solver (default)
218 * - HSS: HSS compression of frontal matrices
219 * - BLR: Block low-rank compression of fronts
220 * - HODLR: Hierarchically Off-diagonal Low-Rank
221 * compression of frontal matrices
222 * - BLR_HODLR: Block low-rank compression of medium
223 * fronts and Hierarchically Off-diagonal
224 * Low-Rank compression of large fronts
225 * - ZFP_BLR_HODLR: ZFP compression for small fronts,
226 * Block low-rank compression of medium
227 * fronts and Hierarchically Off-diagonal
228 * Low-Rank compression of large fronts
229 * - LOSSLESS: Lossless compression
230 * - LOSSY: Lossy compression
231 *
232 * For versions of STRUMPACK < 5, we support only NONE, HSS, and BLR.
233 * BLR_HODLR and ZPR_BLR_HODLR are supported in STRUMPACK >= 6.
234 */
235 void SetCompression(strumpack::CompressionType type);
236
237 /** @brief Set the relative tolerance for low rank compression methods
238 *
239 * This currently affects BLR, HSS, and HODLR. Use
240 * STRUMPACKSolverBase::SetCompression to set the proper compression type.
241 */
242 void SetCompressionRelTol(double rtol);
243
244 /** @brief Set the absolute tolerance for low rank compression methods
245 *
246 * This currently affects BLR, HSS, and HODLR. Use
247 * STRUMPACKSolverBase::SetCompression to set the proper compression type.
248 */
249 void SetCompressionAbsTol(double atol);
250
251#if STRUMPACK_VERSION_MAJOR >= 5
252 /** @brief Set the precision for the lossy compression option
253 *
254 * Use STRUMPACKSolverBase::SetCompression to set the proper compression
255 * type. */
256 void SetCompressionLossyPrecision(int precision);
257
258 /** @brief Set the number of butterfly levels for the HODLR compression
259 * option.
260 *
261 * Use STRUMPACKSolverBase::SetCompression to set the proper compression
262 * type. */
263 void SetCompressionButterflyLevels(int levels);
264#endif
265
266private:
267 /// Helper method for calling the STRUMPACK factoriation routine.
268 void FactorInternal() const;
269
270protected:
272 STRUMPACKSolverType *solver_;
273
277
278 mutable Vector rhs_, sol_;
279 mutable int nrhs_;
280};
281
283 public STRUMPACKSolverBase<strumpack::
284 SparseSolverMPIDist<double, HYPRE_BigInt>>
285{
286public:
287 /// Constructor with the MPI Comm parameter
288 STRUMPACKSolver(MPI_Comm comm);
289
290 /// Constructor with STRUMPACK matrix object.
292
293 /** @brief Constructor with MPI_Comm parameter and command line arguments.
294
295 SetFromCommandLine must be called for the command line arguments
296 to be used.
297 */
298 STRUMPACKSolver(MPI_Comm comm, int argc, char *argv[]);
299 MFEM_DEPRECATED STRUMPACKSolver(int argc, char *argv[], MPI_Comm comm)
300 : STRUMPACKSolver(comm, argc, argv) {}
301
302 /** @brief Constructor with STRUMPACK matrix object and command line
303 arguments.
304
305 SetFromCommandLine must be called for the command line arguments
306 to be used.
307 */
308 STRUMPACKSolver(STRUMPACKRowLocMatrix &A, int argc, char *argv[]);
309
310 /// Destructor.
312};
313
314#if STRUMPACK_VERSION_MAJOR >= 7
316 public STRUMPACKSolverBase<strumpack::
317 SparseSolverMixedPrecisionMPIDist<float, double, HYPRE_BigInt>>
318{
319public:
320 /// Constructor with MPI_Comm parameter.
321 STRUMPACKMixedPrecisionSolver(MPI_Comm comm);
322
323 /// Constructor with STRUMPACK matrix object.
325
326 /** @brief Constructor with MPI_Comm parameter and command line arguments.
327
328 SetFromCommandLine must be called for the command line arguments
329 to be used.
330 */
331 STRUMPACKMixedPrecisionSolver(MPI_Comm comm, int argc, char *argv[]);
332
333 /** @brief Constructor with STRUMPACK matrix object and command line
334 arguments.
335
336 SetFromCommandLine must be called for the command line arguments
337 to be used.
338 */
340 int argc, char *argv[]);
341
342 /// Destructor.
344};
345#endif
346
347} // namespace mfem
348
349#endif // MFEM_USE_MPI
350#endif // MFEM_USE_STRUMPACK
351#endif // MFEM_STRUMPACK
Abstract operator.
Definition operator.hpp:25
STRUMPACKMixedPrecisionSolver(MPI_Comm comm)
Constructor with MPI_Comm parameter.
STRUMPACKRowLocMatrix(MPI_Comm comm, int num_loc_rows, HYPRE_BigInt first_loc_row, HYPRE_BigInt glob_nrows, HYPRE_BigInt glob_ncols, int *I, HYPRE_BigInt *J, double *data, bool sym_sparse=false)
Creates a general parallel matrix from a local CSR matrix on each processor.
Definition strumpack.cpp:22
void Mult(const Vector &x, Vector &y) const
Matrix vector products are not supported for this type of matrix.
Definition strumpack.hpp:55
strumpack::CSRMatrixMPI< double, HYPRE_BigInt > * GetA() const
Get access to the internal CSR matrix.
Definition strumpack.hpp:65
MPI_Comm GetComm() const
Get the MPI Comm being used by the parallel matrix.
Definition strumpack.hpp:62
void SetOperator(const Operator &op)
Set the operator/matrix.
void SetAbsTol(double atol)
Set the absolute tolerance for iterative solvers.
void ArrayMult(const Array< const Vector * > &X, Array< Vector * > &Y) const
Factor and solve the linear systems across the array of vectors.
void SetMaxIter(int max_it)
Set the maximum number of iterations for iterative solvers.
STRUMPACKSolverType * solver_
void SetMatching(strumpack::MatchingJob job)
Configure static pivoting for stability.
void SetCompression(strumpack::CompressionType type)
Select compression for sparse data types.
STRUMPACKSolverBase(MPI_Comm comm, int argc, char *argv[])
Constructor with MPI_Comm parameter and command line arguments.
void SetPrintFactorStatistics(bool print_stat)
Set up verbose printing during the factor step.
void SetKrylovSolver(strumpack::KrylovSolver method)
Set the Krylov solver method to use.
void Mult(const Vector &x, Vector &y) const
Factor and solve the linear system .
void SetCompressionButterflyLevels(int levels)
Set the number of butterfly levels for the HODLR compression option.
void SetCompressionAbsTol(double atol)
Set the absolute tolerance for low rank compression methods.
void EnableGPU()
Enable GPU off-loading available if STRUMPACK was compiled with CUDA.
virtual ~STRUMPACKSolverBase()
Default destructor.
void DisableGPU()
Disable GPU off-loading available if STRUMPACK was compiled with CUDA.
const STRUMPACKRowLocMatrix * APtr_
void SetFromCommandLine()
Set options that were captured from the command line.
void SetRelTol(double rtol)
Set the relative tolerance for interative solvers.
void SetCompressionLossyPrecision(int precision)
Set the precision for the lossy compression option.
void SetReorderingReuse(bool reuse)
Set the flag controlling reuse of the symbolic factorization for multiple operators.
void SetPrintSolveStatistics(bool print_stat)
Set up verbose printing during the solve step.
void SetCompressionRelTol(double rtol)
Set the relative tolerance for low rank compression methods.
void SetReorderingStrategy(strumpack::ReorderingStrategy method)
Set matrix reordering strategy.
STRUMPACKSolver(MPI_Comm comm)
Constructor with the MPI Comm parameter.
~STRUMPACKSolver()
Destructor.
MFEM_DEPRECATED STRUMPACKSolver(int argc, char *argv[], MPI_Comm comm)
Base class for solvers.
Definition operator.hpp:683
Vector data type.
Definition vector.hpp:80
HYPRE_Int HYPRE_BigInt