MFEM v4.8.0
Finite element discretization library
Loading...
Searching...
No Matches
strumpack.hpp
Go to the documentation of this file.
1// Copyright (c) 2010-2025, 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 override
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 override;
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,
104 Array<Vector *> &Y) const override;
105
106 /** @brief Set the operator/matrix.
107 \note @a A must be a STRUMPACKRowLocMatrix. */
108 void SetOperator(const Operator &op) override;
109
110 /** @brief Set options that were captured from the command line.
111
112 These were captured in the constructer STRUMPACKSolverBase. Refer
113 to the STRUMPACK documentation for details. */
114 void SetFromCommandLine();
115
116 /// Set up verbose printing during the factor step
117 void SetPrintFactorStatistics(bool print_stat);
118
119 /// Set up verbose printing during the solve step
120 void SetPrintSolveStatistics(bool print_stat);
121
122 /// Set the relative tolerance for interative solvers
123 void SetRelTol(double rtol);
124
125 /// Set the absolute tolerance for iterative solvers
126 void SetAbsTol(double atol);
127
128 /// Set the maximum number of iterations for iterative solvers
129 void SetMaxIter(int max_it);
130
131 /** @brief Set the flag controlling reuse of the symbolic factorization for
132 multiple operators.
133
134 This method must be called before repeated calls to SetOperator. */
135 void SetReorderingReuse(bool reuse);
136
137 /** @brief Enable GPU off-loading available if STRUMPACK was compiled with
138 CUDA.
139 @note Input/Output from MFEM to STRUMPACK is all still through host
140 memory. */
141 void EnableGPU();
142
143 /** @brief Disable GPU off-loading available if STRUMPACK was compiled with
144 CUDA.
145 @note Input/Output from MFEM to STRUMPACK is all still through host
146 memory. */
147 void DisableGPU();
148
149 /** @brief Set the Krylov solver method to use
150 *
151 * STRUMPACK is an (approximate) direct solver. It can be used as a direct
152 * solver or as a preconditioner. To use STRUMPACK as only a preconditioner,
153 * set the Krylov solver to DIRECT. STRUMPACK also provides iterative solvers
154 * which can use the preconditioner, and these iterative solvers can also be
155 * used without preconditioner.
156 *
157 * Supported values are:
158 * - AUTO: Use iterative refinement if no HSS compression is
159 * used, otherwise use GMRES
160 * - DIRECT: No outer iterative solver, just a single application
161 * of the multifrontal solver
162 * - REFINE: Iterative refinement
163 * - PREC_GMRES: Preconditioned GMRES
164 * The preconditioner is the (approx) multifrontal solver
165 * - GMRES: UN-preconditioned GMRES (for testing mainly)
166 * - PREC_BICGSTAB: Preconditioned BiCGStab
167 * The preconditioner is the (approx) multifrontal solver
168 * - BICGSTAB: UN-preconditioned BiCGStab. (for testing mainly)
169 */
170 void SetKrylovSolver(strumpack::KrylovSolver method);
171
172 /** @brief Set matrix reordering strategy
173 *
174 * Supported reorderings are:
175 * - NATURAL: Do not reorder the system
176 * - METIS: Use Metis nested-dissection reordering (default)
177 * - PARMETIS: Use ParMetis nested-dissection reordering
178 * - SCOTCH: Use Scotch nested-dissection reordering
179 * - PTSCOTCH: Use PT-Scotch nested-dissection reordering
180 * - RCM: Use RCM reordering
181 * - GEOMETRIC: A simple geometric nested dissection code that
182 * only works for regular meshes
183 * - AMD: Approximate minimum degree
184 * - MMD: Multiple minimum degree
185 * - AND: Nested dissection
186 * - MLF: Minimum local fill
187 * - SPECTRAL: Spectral nested dissection
188 */
189 void SetReorderingStrategy(strumpack::ReorderingStrategy method);
190
191 /** @brief Configure static pivoting for stability.
192 *
193 * The static pivoting in STRUMPACK
194 * permutes the sparse input matrix in order to get large (nonzero) elements
195 * on the diagonal. If the input matrix is already diagonally dominant, this
196 * reordering can be disabled.
197 *
198 * Supported matching algorithms are:
199 * - NONE: Don't do anything
200 * - MAX_CARDINALITY: Maximum cardinality
201 * - MAX_SMALLEST_DIAGONAL: Maximum smallest diagonal value
202 * - MAX_SMALLEST_DIAGONAL_2: Same as MAX_SMALLEST_DIAGONAL
203 * but different algorithm
204 * - MAX_DIAGONAL_SUM: Maximum sum of diagonal values
205 * - MAX_DIAGONAL_PRODUCT_SCALING: Maximum product of diagonal values
206 * and row and column scaling (default)
207 * - COMBBLAS: Use AWPM from CombBLAS (only with
208 * version >= 3)
209 */
210 void SetMatching(strumpack::MatchingJob job);
211
212 /** @brief Select compression for sparse data types
213 *
214 * Enable support for rank-structured data formats, which can be used
215 * for compression within the sparse solver.
216 *
217 * Supported compression types are:
218 * - NONE: No compression, purely direct solver (default)
219 * - HSS: HSS compression of frontal matrices
220 * - BLR: Block low-rank compression of fronts
221 * - HODLR: Hierarchically Off-diagonal Low-Rank
222 * compression of frontal matrices
223 * - BLR_HODLR: Block low-rank compression of medium
224 * fronts and Hierarchically Off-diagonal
225 * Low-Rank compression of large fronts
226 * - ZFP_BLR_HODLR: ZFP compression for small fronts,
227 * Block low-rank compression of medium
228 * fronts and Hierarchically Off-diagonal
229 * Low-Rank compression of large fronts
230 * - LOSSLESS: Lossless compression
231 * - LOSSY: Lossy compression
232 *
233 * For versions of STRUMPACK < 5, we support only NONE, HSS, and BLR.
234 * BLR_HODLR and ZPR_BLR_HODLR are supported in STRUMPACK >= 6.
235 */
236 void SetCompression(strumpack::CompressionType type);
237
238 /** @brief Set the relative tolerance for low rank compression methods
239 *
240 * This currently affects BLR, HSS, and HODLR. Use
241 * STRUMPACKSolverBase::SetCompression to set the proper compression type.
242 */
243 void SetCompressionRelTol(double rtol);
244
245 /** @brief Set the absolute tolerance for low rank compression methods
246 *
247 * This currently affects BLR, HSS, and HODLR. Use
248 * STRUMPACKSolverBase::SetCompression to set the proper compression type.
249 */
250 void SetCompressionAbsTol(double atol);
251
252#if STRUMPACK_VERSION_MAJOR >= 5
253 /** @brief Set the precision for the lossy compression option
254 *
255 * Use STRUMPACKSolverBase::SetCompression to set the proper compression
256 * type. */
257 void SetCompressionLossyPrecision(int precision);
258
259 /** @brief Set the number of butterfly levels for the HODLR compression
260 * option.
261 *
262 * Use STRUMPACKSolverBase::SetCompression to set the proper compression
263 * type. */
264 void SetCompressionButterflyLevels(int levels);
265#endif
266
267private:
268 /// Helper method for calling the STRUMPACK factoriation routine.
269 void FactorInternal() const;
270
271protected:
273 STRUMPACKSolverType *solver_;
274
278
279 mutable Vector rhs_, sol_;
280 mutable int nrhs_;
281};
282
284 public STRUMPACKSolverBase<strumpack::
285 SparseSolverMPIDist<double, HYPRE_BigInt>>
286{
287public:
288 /// Constructor with the MPI Comm parameter
289 STRUMPACKSolver(MPI_Comm comm);
290
291 /// Constructor with STRUMPACK matrix object.
293
294 /** @brief Constructor with MPI_Comm parameter and command line arguments.
295
296 SetFromCommandLine must be called for the command line arguments
297 to be used.
298 */
299 STRUMPACKSolver(MPI_Comm comm, int argc, char *argv[]);
300 MFEM_DEPRECATED STRUMPACKSolver(int argc, char *argv[], MPI_Comm comm)
301 : STRUMPACKSolver(comm, argc, argv) {}
302
303 /** @brief Constructor with STRUMPACK matrix object and command line
304 arguments.
305
306 SetFromCommandLine must be called for the command line arguments
307 to be used.
308 */
309 STRUMPACKSolver(STRUMPACKRowLocMatrix &A, int argc, char *argv[]);
310
311 /// Destructor.
313};
314
315#if STRUMPACK_VERSION_MAJOR >= 7
317 public STRUMPACKSolverBase<strumpack::
318 SparseSolverMixedPrecisionMPIDist<float, double, HYPRE_BigInt>>
319{
320public:
321 /// Constructor with MPI_Comm parameter.
322 STRUMPACKMixedPrecisionSolver(MPI_Comm comm);
323
324 /// Constructor with STRUMPACK matrix object.
326
327 /** @brief Constructor with MPI_Comm parameter and command line arguments.
328
329 SetFromCommandLine must be called for the command line arguments
330 to be used.
331 */
332 STRUMPACKMixedPrecisionSolver(MPI_Comm comm, int argc, char *argv[]);
333
334 /** @brief Constructor with STRUMPACK matrix object and command line
335 arguments.
336
337 SetFromCommandLine must be called for the command line arguments
338 to be used.
339 */
341 int argc, char *argv[]);
342
343 /// Destructor.
345};
346#endif
347
348} // namespace mfem
349
350#endif // MFEM_USE_MPI
351#endif // MFEM_USE_STRUMPACK
352#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 override
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 SetAbsTol(double atol)
Set the absolute tolerance for iterative solvers.
void SetMaxIter(int max_it)
Set the maximum number of iterations for iterative solvers.
STRUMPACKSolverType * solver_
void ArrayMult(const Array< const Vector * > &X, Array< Vector * > &Y) const override
Factor and solve the linear systems across the array of vectors.
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 override
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 SetOperator(const Operator &op) override
Set the operator/matrix.
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:780
Vector data type.
Definition vector.hpp:82
HYPRE_Int HYPRE_BigInt