MFEM  v4.6.0
Finite element discretization library
superlu.hpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2023, 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_SUPERLU
13 #define MFEM_SUPERLU
14 
15 #include "../config/config.hpp"
16 
17 #ifdef MFEM_USE_SUPERLU
18 #ifdef MFEM_USE_MPI
19 
20 #include "operator.hpp"
21 #include "hypre.hpp"
22 #include <mpi.h>
23 
24 namespace mfem
25 {
26 
27 namespace superlu
28 {
29 
30 // Copy selected enumerations from SuperLU (from superlu_enum_consts.h)
31 #ifdef MFEM_USE_SUPERLU5
33 #else
35 #endif
38  } ColPerm;
41 
42 } // namespace superlu
43 
45 {
46 public:
47  /** Creates a general parallel matrix from a local CSR matrix on each
48  processor described by the I, J and data arrays. The local matrix should
49  be of size (local) nrows by (global) glob_ncols. The new parallel matrix
50  contains copies of all input arrays (so they can be deleted). */
51  SuperLURowLocMatrix(MPI_Comm comm,
52  int num_loc_rows, HYPRE_BigInt first_loc_row,
53  HYPRE_BigInt glob_nrows, HYPRE_BigInt glob_ncols,
54  int *I, HYPRE_BigInt *J, double *data);
55 
56  /** Creates a copy of the parallel matrix hypParMat in SuperLU's RowLoc
57  format. All data is copied so the original matrix may be deleted. */
58  SuperLURowLocMatrix(const Operator &op);
59 
61 
62  void Mult(const Vector &x, Vector &y) const
63  {
64  MFEM_ABORT("SuperLURowLocMatrix::Mult: Matrix vector products are not "
65  "supported!");
66  }
67 
68  void *InternalData() const { return rowLocPtr_; }
69 
70  MPI_Comm GetComm() const { return comm_; }
71 
72  HYPRE_BigInt GetGlobalNumRows() const { return num_global_rows_; }
73 
74  HYPRE_BigInt GetGlobalNumColumns() const { return num_global_cols_; }
75 
76 private:
77  MPI_Comm comm_;
78  void *rowLocPtr_;
79  HYPRE_BigInt num_global_rows_, num_global_cols_;
80 };
81 
82 /** The MFEM SuperLU Direct Solver class.
83 
84  The mfem::SuperLUSolver class uses the SuperLU_DIST library to perform LU
85  factorization of a parallel sparse matrix. The solver is capable of handling
86  double precision types. It is currently maintained by Xiaoye Sherry Li at
87  NERSC, see http://crd-legacy.lbl.gov/~xiaoye/SuperLU/.
88 */
89 class SuperLUSolver : public Solver
90 {
91 public:
92  // Constructor with MPI_Comm parameter.
93  SuperLUSolver(MPI_Comm comm, int npdep = 1);
94 
95  // Constructor with SuperLU matrix object.
96  SuperLUSolver(SuperLURowLocMatrix &A, int npdep = 1);
97 
98  // Default destructor.
100 
101  // Set the operator.
102  void SetOperator(const Operator &op);
103 
104  // Factor and solve the linear system y = Op^{-1} x.
105  // Note: Factorization modifies the operator matrix.
106  void Mult(const Vector &x, Vector &y) const;
107  void ArrayMult(const Array<const Vector *> &X, Array<Vector *> &Y) const;
108 
109  // Factor and solve the linear system y = Op^{-T} x.
110  // Note: Factorization modifies the operator matrix.
111  void MultTranspose(const Vector &x, Vector &y) const;
113  Array<Vector *> &Y) const;
114 
115  // Set various solver options. Refer to SuperLU_DIST documentation for
116  // details.
117  void SetPrintStatistics(bool print_stat);
118  void SetEquilibriate(bool equil);
119  void SetColumnPermutation(superlu::ColPerm col_perm);
120  void SetRowPermutation(superlu::RowPerm row_perm);
122  void SetReplaceTinyPivot(bool rtp);
123  void SetNumLookAheads(int num_lookaheads);
124  void SetLookAheadElimTree(bool etree);
125  void SetSymmetricPattern(bool sym);
126  void SetParSymbFact(bool par);
127  void SetFact(superlu::Fact fact);
128 
129  // Processor grid for SuperLU_DIST.
130  const int nprow_, npcol_, npdep_;
131 
132 private:
133  // Initialize the solver.
134  void Init(MPI_Comm comm);
135 
136  // Handle error message from call to SuperLU solver.
137  void HandleError(int info) const;
138 
139 protected:
141  mutable Vector sol_;
142  mutable int nrhs_;
143 
144  /** The actual types of the following pointers are hidden to avoid exposing
145  the SuperLU header files to the entire library. Their types are given in
146  the trailing comments. The reason that this is necessary is that SuperLU
147  defines these structs differently for use with its real and complex
148  solvers. If we want to add support for SuperLU's complex solvers one day
149  we will need to hide these types to avoid name conflicts. */
150  void *optionsPtr_; // superlu_options_t *
151  void *ScalePermstructPtr_; // ScalePermsruct_t *
152  void *LUstructPtr_; // LUstruct_t *
153  void *SOLVEstructPtr_; // SOLVEstruct_t *
154  void *gridPtr_; // gridinfo_t * or gridinfo3d_t *
155 };
156 
157 } // namespace mfem
158 
159 #endif // MFEM_USE_MPI
160 #endif // MFEM_USE_SUPERLU
161 #endif // MFEM_SUPERLU
void SetRowPermutation(superlu::RowPerm row_perm)
Definition: superlu.cpp:415
void SetFact(superlu::Fact fact)
Definition: superlu.cpp:468
MPI_Comm GetComm() const
Definition: superlu.hpp:70
void ArrayMultTranspose(const Array< const Vector *> &X, Array< Vector *> &Y) const
Action of the transpose operator on a matrix: Y=A^t(X).
Definition: superlu.cpp:711
void * ScalePermstructPtr_
Definition: superlu.hpp:151
const int nprow_
Definition: superlu.hpp:130
HYPRE_BigInt GetGlobalNumRows() const
Definition: superlu.hpp:72
void SetSymmetricPattern(bool sym)
Definition: superlu.cpp:454
HYPRE_BigInt GetGlobalNumColumns() const
Definition: superlu.hpp:74
void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Definition: superlu.cpp:587
void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition: superlu.cpp:475
void * InternalData() const
Definition: superlu.hpp:68
void MultTranspose(const Vector &x, Vector &y) const
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
Definition: superlu.cpp:700
void SetIterativeRefine(superlu::IterRefine iter_ref)
Definition: superlu.cpp:427
void SetColumnPermutation(superlu::ColPerm col_perm)
Definition: superlu.cpp:399
const int npdep_
Definition: superlu.hpp:130
void SetNumLookAheads(int num_lookaheads)
Definition: superlu.cpp:441
void SetPrintStatistics(bool print_stat)
Definition: superlu.cpp:385
HYPRE_Int HYPRE_BigInt
A class to initialize the size of a Tensor.
Definition: dtensor.hpp:54
void * SOLVEstructPtr_
Definition: superlu.hpp:153
void SetLookAheadElimTree(bool etree)
Definition: superlu.cpp:447
void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Definition: superlu.hpp:62
void ArrayMult(const Array< const Vector *> &X, Array< Vector *> &Y) const
Operator application on a matrix: Y=A(X).
Definition: superlu.cpp:596
void SetReplaceTinyPivot(bool rtp)
Definition: superlu.cpp:434
Vector data type.
Definition: vector.hpp:58
SuperLURowLocMatrix(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)
Definition: superlu.cpp:107
const SuperLURowLocMatrix * APtr_
Definition: superlu.hpp:140
void SetParSymbFact(bool par)
Definition: superlu.cpp:461
Base class for solvers.
Definition: operator.hpp:682
Abstract operator.
Definition: operator.hpp:24
const int npcol_
Definition: superlu.hpp:130
SuperLUSolver(MPI_Comm comm, int npdep=1)
Definition: superlu.cpp:262
void SetEquilibriate(bool equil)
Definition: superlu.cpp:392