Amesos2 - Direct Sparse Solver Interfaces  Version of the Day
Amesos2_Superlumt_decl.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Amesos2: Templated Direct Sparse Solver Package
4 //
5 // Copyright 2011 NTESS and the Amesos2 contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
19 #ifndef AMESOS2_SUPERLUMT_DECL_HPP
20 #define AMESOS2_SUPERLUMT_DECL_HPP
21 
22 #include "Amesos2_SolverTraits.hpp"
25 
26 
27 namespace Amesos2 {
28 
29 
41 template <class Matrix,
42  class Vector>
43 class Superlumt : public SolverCore<Amesos2::Superlumt, Matrix, Vector>
44 {
45  friend class SolverCore<Amesos2::Superlumt,Matrix,Vector>; // Give our base access
46  // to our private
47  // implementation funcs
48 public:
49 
51  static const char* name; // declaration. Initialization outside.
52 
53  typedef Superlumt<Matrix,Vector> type;
54  typedef SolverCore<Amesos2::Superlumt,Matrix,Vector> super_type;
55 
56  // Since typedef's are not inheritted, go grab them
57  typedef typename super_type::scalar_type scalar_type;
58  typedef typename super_type::local_ordinal_type local_ordinal_type;
59  typedef typename super_type::global_ordinal_type global_ordinal_type;
60  typedef typename super_type::global_size_type global_size_type;
61 
62  typedef TypeMap<Amesos2::Superlumt,scalar_type> type_map;
63 
64  typedef typename type_map::type slu_type;
65  typedef typename type_map::magnitude_type magnitude_type;
66 
67  typedef FunctionMap<Amesos2::Superlumt,slu_type> function_map;
68 
69 
71 
72 
79  Superlumt(Teuchos::RCP<const Matrix> A,
80  Teuchos::RCP<Vector> X,
81  Teuchos::RCP<const Vector> B);
82 
83 
85  ~Superlumt( );
86 
88 
89 private:
90 
97  int preOrdering_impl();
98 
99 
108 
109 
120 
121 
133  int solve_impl(const Teuchos::Ptr<MultiVecAdapter<Vector> > X,
134  const Teuchos::Ptr<const MultiVecAdapter<Vector> > B) const;
135 
136 
142  bool matrixShapeOK_impl() const;
143 
144 
182  void setParameters_impl(
183  const Teuchos::RCP<Teuchos::ParameterList> & parameterList );
184  /* Parameters to support in the future:
185  *
186  * <li> \c "IterRefine" : { \c "NO" | \c "SINGLE" | \c "DOUBLE" | \c "EXTRA"
187  * }. Specifies whether to perform iterative refinement, and in
188  * what precision to compute the residual. (Not currently supported)</li>
189  */
190 
197  Teuchos::RCP<const Teuchos::ParameterList> getValidParameters_impl() const;
198 
199 
208  bool loadA_impl(EPhase current_phase);
209 
210 
211  // struct holds all data necessary to make a superlu factorization or solve call
212  mutable struct SLUData {
213  SLUMT::SuperMatrix A, BX, L, U;
214  SLUMT::SuperMatrix AC;
215 
216  SLUMT::superlumt_options_t options;
217  SLUMT::superlu_memusage_t mem_usage;
218  SLUMT::Gstat_t stat;
219 
220  Teuchos::Array<magnitude_type> berr;
221  Teuchos::Array<magnitude_type> ferr;
222  Teuchos::Array<int> perm_r;
223  Teuchos::Array<int> perm_c;
224  Teuchos::Array<magnitude_type> R;
225  Teuchos::Array<magnitude_type> C;
226 
227  // in contrast to SuperLU, memory for etree will be allocated by
228  // pxgssvx and the pointer will be stored in `options'
229 
230  SLUMT::equed_t equed;
231  bool rowequ, colequ;
232  } data_;
233 
234  // The following Arrays are persisting storage arrays for A, X, and B
236  Teuchos::Array<typename TypeMap<Amesos2::Superlumt,scalar_type>::type> nzvals_;
238  Teuchos::Array<int> rowind_;
240  Teuchos::Array<int> colptr_;
241 
242  /* Note: In the above, must use "Amesos2::Superlumt" rather than
243  * "Superlumt" because otherwise the compiler references the
244  * specialized type of the class, and not the templated type that is
245  * required for Amesos2::TypeMap
246  */
247 
248  /* SuperLU can accept input in either compressed-row or
249  * compressed-column storage. We will store and pass matrices in
250  * *compressed-row* format because that is the format Amesos used.
251  */
252 
253  bool is_contiguous_;
254 
255 }; // End class Superlumt
256 
257 
258 // Specialize the solver_traits struct for SuperLU_MT
259 template <>
260 struct solver_traits<Superlumt> {
261 #ifdef HAVE_TEUCHOS_COMPLEX
262  typedef Meta::make_list6<float,
263  double,
264  std::complex<float>,
265  std::complex<double>,
266  SLUMT::C::complex,
267  SLUMT::Z::doublecomplex> supported_scalars;
268 #else
269  typedef Meta::make_list2<float, double> supported_scalars;
270 #endif
271 };
272 
273 } // end namespace Amesos2
274 
275 #endif // AMESOS2_SUPERLUMT_DECL_HPP
Amesos2::SolverCore: A templated interface for interaction with third-party direct sparse solvers...
Definition: Amesos2_SolverCore_decl.hpp:71
bool matrixShapeOK_impl() const
Determines whether the shape of the matrix is OK for this solver.
Definition: Amesos2_Superlumt_def.hpp:411
Amesos2 interface to the Multi-threaded version of SuperLU.
Definition: Amesos2_Superlumt_decl.hpp:43
int solve_impl(const Teuchos::Ptr< MultiVecAdapter< Vector > > X, const Teuchos::Ptr< const MultiVecAdapter< Vector > > B) const
SuperLU_MT specific solve.
Definition: Amesos2_Superlumt_def.hpp:297
Map types to solver-specific data-types and enums.
Definition: Amesos2_TypeMap.hpp:48
Teuchos::Array< int > colptr_
Stores the location in Ai_ and Aval_ that starts row j.
Definition: Amesos2_Superlumt_decl.hpp:240
int symbolicFactorization_impl()
Perform symbolic factorization of the matrix using SuperLU_MT.
Definition: Amesos2_Superlumt_def.hpp:168
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters_impl() const
Definition: Amesos2_Superlumt_def.hpp:476
Templated core-functionality class for Amesos2 solvers.
Provides traits about solvers.
Definition: Amesos2_SolverTraits.hpp:37
int preOrdering_impl()
Performs pre-ordering on the matrix to increase efficiency.
Definition: Amesos2_Superlumt_def.hpp:148
Provides a mechanism to map function calls to the correct Solver function based on the scalar type of...
std::string name() const override
Return the name of this solver.
Definition: Amesos2_SolverCore_def.hpp:725
int numericFactorization_impl()
SuperLU_MT specific numeric factorization.
Definition: Amesos2_Superlumt_def.hpp:193
Interface to Amesos2 solver objects.
Definition: Amesos2_Solver_decl.hpp:44
void setParameters_impl(const Teuchos::RCP< Teuchos::ParameterList > &parameterList)
Definition: Amesos2_Superlumt_def.hpp:422
Teuchos::Array< typename TypeMap< Amesos2::Superlumt, scalar_type >::type > nzvals_
Stores the values of the nonzero entries for SuperLU.
Definition: Amesos2_Superlumt_decl.hpp:236
Passes functions to TPL functions based on type.
Definition: Amesos2_FunctionMap.hpp:42
Provides access to interesting solver traits.
bool loadA_impl(EPhase current_phase)
Reads matrix data into internal structures.
Definition: Amesos2_Superlumt_def.hpp:565
Teuchos::Array< int > rowind_
Stores the row indices of the nonzero entries.
Definition: Amesos2_Superlumt_decl.hpp:238