Amesos2 - Direct Sparse Solver Interfaces Version of the Day
Amesos2_Superlu_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_SUPERLU_DECL_HPP
20#define AMESOS2_SUPERLU_DECL_HPP
21
23#include "Amesos2_SolverCore.hpp"
25
26#if defined(KOKKOSKERNELS_ENABLE_SUPERNODAL_SPTRSV) && defined(KOKKOSKERNELS_ENABLE_TPL_SUPERLU)
27#include "KokkosKernels_Handle.hpp"
28#endif
29
30namespace Amesos2 {
31
32
40template <class Matrix,
41 class Vector>
42class Superlu : public SolverCore<Amesos2::Superlu, Matrix, Vector>
43{
44 friend class SolverCore<Amesos2::Superlu,Matrix,Vector>; // Give our base access
45 // to our private
46 // implementation funcs
47public:
48
50 static const char* name; // declaration. Initialization outside.
51
52 typedef Superlu<Matrix,Vector> type;
53 typedef SolverCore<Amesos2::Superlu,Matrix,Vector> super_type;
54
55 // Since typedef's are not inheritted, go grab them
56 typedef typename super_type::scalar_type scalar_type;
57 typedef typename super_type::local_ordinal_type local_ordinal_type;
58 typedef typename super_type::global_ordinal_type global_ordinal_type;
59 typedef typename super_type::global_size_type global_size_type;
60
61 typedef TypeMap<Amesos2::Superlu,scalar_type> type_map;
62
63 /*
64 * The SuperLU interface will need two other typedef's, which are:
65 * - the superlu type that corresponds to scalar_type and
66 * - the corresponding type to use for magnitude
67 */
68 typedef typename type_map::type slu_type;
69 typedef typename type_map::convert_type slu_convert_type;
70 typedef typename type_map::magnitude_type magnitude_type;
71
72 typedef FunctionMap<Amesos2::Superlu,slu_type> function_map;
73
75
76
83 Superlu(Teuchos::RCP<const Matrix> A,
84 Teuchos::RCP<Vector> X,
85 Teuchos::RCP<const Vector> B);
86
87
89 ~Superlu( );
90
92
94 std::string description() const override;
95
96private:
97
103 int preOrdering_impl();
104
105
114
115
122
123
135 int solve_impl(const Teuchos::Ptr<MultiVecAdapter<Vector> > X,
136 const Teuchos::Ptr<const MultiVecAdapter<Vector> > B) const;
137
138
142 bool matrixShapeOK_impl() const;
143
144
179 const Teuchos::RCP<Teuchos::ParameterList> & parameterList );
180
181
188 Teuchos::RCP<const Teuchos::ParameterList> getValidParameters_impl() const;
189
190
199 bool loadA_impl(EPhase current_phase);
200
201 typedef Kokkos::DefaultHostExecutionSpace HostExecSpaceType;
202
203 // struct holds all data necessary to make a superlu factorization or solve call
204 mutable struct SLUData {
205 SLU::SuperMatrix A, B, X, L, U; // matrix A in NCformat
206 SLU::SuperMatrix AC; // permuted matrix A in NCPformat
207
208 SLU::superlu_options_t options;
209 SLU::mem_usage_t mem_usage;
210#ifdef HAVE_AMESOS2_SUPERLU5_API
211 SLU::GlobalLU_t lu; // Use for gssvx and gsisx in SuperLU 5.0
212#endif
213 SLU::SuperLUStat_t stat;
214
215
216
217 typedef Kokkos::View<magnitude_type*, HostExecSpaceType> host_mag_array;
218 typedef Kokkos::View<int*, HostExecSpaceType> host_int_array;
219 host_mag_array berr;
220 host_mag_array ferr;
221 host_int_array perm_r;
222 host_int_array perm_c;
223 host_int_array etree;
224 host_mag_array R;
225 host_mag_array C;
226
227#if defined(KOKKOSKERNELS_ENABLE_SUPERNODAL_SPTRSV) && defined(KOKKOSKERNELS_ENABLE_TPL_SUPERLU)
228 host_int_array parents;
229#endif
230
231 char equed;
232 bool rowequ, colequ; // flags what type of equilibration
233 // has been performed
234 magnitude_type anorm, rcond; // condition number estimate
235
236 int relax;
237 int panel_size;
238 } data_;
239
240 typedef int size_type;
241 typedef int ordinal_type;
242 typedef Kokkos::View<size_type*, HostExecSpaceType> host_size_type_array;
243 typedef Kokkos::View<ordinal_type*, HostExecSpaceType> host_ordinal_type_array;
244 typedef Kokkos::View<slu_type*, HostExecSpaceType> host_value_type_array;
245
246 // The following Arrays are persisting storage arrays for A, X, and B
248 host_value_type_array host_nzvals_view_;
249 Teuchos::Array<slu_convert_type> convert_nzvals_; // copy to SuperLU native array before calling SuperLU
250
252 host_size_type_array host_rows_view_;
254 host_ordinal_type_array host_col_ptr_view_;
255
256 typedef typename Kokkos::View<slu_type**, Kokkos::LayoutLeft, HostExecSpaceType>
257 host_solve_array_t;
258
260 mutable host_solve_array_t host_xValues_;
261 mutable Teuchos::Array<slu_convert_type> convert_xValues_; // copy to SuperLU native array before calling SuperLU
262
264 mutable host_solve_array_t host_bValues_;
265 mutable Teuchos::Array<slu_convert_type> convert_bValues_; // copy to SuperLU native array before calling SuperLU
266
267#if defined(KOKKOSKERNELS_ENABLE_SUPERNODAL_SPTRSV) && defined(KOKKOSKERNELS_ENABLE_TPL_SUPERLU)
268 typedef Kokkos::DefaultExecutionSpace DeviceExecSpaceType;
269
270 #ifdef KOKKOS_ENABLE_CUDA
271 // solver will be UVM off even though Tpetra is CudaUVMSpace
272 typedef typename Kokkos::CudaSpace DeviceMemSpaceType;
273 #else
274 typedef typename DeviceExecSpaceType::memory_space DeviceMemSpaceType;
275 #endif
276
277 typedef Kokkos::View<slu_type**, Kokkos::LayoutLeft, DeviceMemSpaceType>
278 device_solve_array_t;
279 // For triangular solves we have both host and device versions of xValues and
280 // bValues because a parameter can turn it on or off.
281 mutable device_solve_array_t device_xValues_;
282 mutable device_solve_array_t device_bValues_;
283 typedef Kokkos::View<int*, DeviceMemSpaceType> device_int_array;
284 typedef Kokkos::View<magnitude_type*, DeviceMemSpaceType> device_mag_array;
285 device_int_array device_trsv_perm_r_;
286 device_int_array device_trsv_perm_c_;
287 device_mag_array device_trsv_R_;
288 device_mag_array device_trsv_C_;
289 mutable device_solve_array_t device_trsv_rhs_;
290 mutable device_solve_array_t device_trsv_sol_;
291 typedef KokkosKernels::Experimental::KokkosKernelsHandle <size_type, ordinal_type, slu_type,
292 DeviceExecSpaceType, DeviceMemSpaceType, DeviceMemSpaceType> kernel_handle_type;
293 mutable kernel_handle_type device_khL_;
294 mutable kernel_handle_type device_khU_;
295 /* parameters for SpTRSV */
296 bool sptrsv_invert_diag_;
297 bool sptrsv_invert_offdiag_;
298 bool sptrsv_u_in_csr_;
299 bool sptrsv_merge_supernodes_;
300 bool sptrsv_use_spmv_;
301#endif
302
303 /* Note: In the above, must use "Amesos2::Superlu" rather than
304 * "Superlu" because otherwise the compiler references the
305 * specialized type of the class, and not the templated type that is
306 * required for Amesos2::TypeMap
307 */
308
309 /* SuperLU can accept input in either compressed-row or
310 * compressed-column storage. We will store and pass matrices in
311 * *compressed-column* format.
312 */
313
314 /*
315 * Internal flag that is used for the numericFactorization_impl
316 * routine. If true, then the superlu gstrf routine should have
317 * SamePattern_SameRowPerm in its options. Otherwise, it should
318 * factor from scratch.
319 *
320 * This is somewhat of a kludge to get around the fact that the
321 * superlu routines all expect something different from the options
322 * struct. The big issue is that we don't want gstrf doing the
323 * symbolic factorization if it doesn't need to. On the other hand,
324 * we can't leave options.Fact set to SamePattern_SameRowPerm
325 * because the solver driver needs it to be set at FACTORED. But
326 * having it set at FACTORED upon re-entrance into
327 * numericFactorization prompts gstrf to redo the symbolic
328 * factorization.
329 */
330 bool same_symbolic_;
331 bool ILU_Flag_;
332
333 bool is_contiguous_;
334 bool use_triangular_solves_;
335
336 void triangular_solve_factor();
337
338 /* call metis before SuperLU */
339 bool use_metis_;
340 bool symmetrize_metis_;
341
342 public: // for GPU
343 void triangular_solve() const; // Only for internal use - public to support kernels
344}; // End class Superlu
345
346
347// Specialize solver_traits struct for SuperLU
348template <>
349struct solver_traits<Superlu> {
350#ifdef HAVE_TEUCHOS_COMPLEX
351 typedef Meta::make_list6<float, double,
352 std::complex<float>, std::complex<double>,
353 Kokkos::complex<float>, Kokkos::complex<double>>
354 supported_scalars;
355#else
356 typedef Meta::make_list2<float, double> supported_scalars;
357#endif
358};
359
360template <typename Scalar, typename LocalOrdinal, typename ExecutionSpace>
361struct solver_supports_matrix<Superlu,
362 KokkosSparse::CrsMatrix<Scalar, LocalOrdinal, ExecutionSpace>> {
363 static const bool value = true;
364};
365
366} // end namespace Amesos2
367
368#endif // AMESOS2_SUPERLU_DECL_HPP
Provides access to interesting solver traits.
Provides a mechanism to map function calls to the correct Solver function based on the scalar type of...
Amesos2::SolverCore: A templated interface for interaction with third-party direct sparse solvers.
Definition Amesos2_SolverCore_decl.hpp:72
Amesos2 interface to the SuperLU package.
Definition Amesos2_Superlu_decl.hpp:43
bool matrixShapeOK_impl() const
Determines whether the shape of the matrix is OK for this solver.
Definition Amesos2_Superlu_def.hpp:674
int numericFactorization_impl()
Superlu specific numeric factorization.
Definition Amesos2_Superlu_def.hpp:283
host_solve_array_t host_xValues_
Persisting 1D store for X.
Definition Amesos2_Superlu_decl.hpp:260
std::string description() const override
Returns a short description of this Solver.
Definition Amesos2_Superlu_def.hpp:117
host_ordinal_type_array host_col_ptr_view_
Stores the row indices of the nonzero entries.
Definition Amesos2_Superlu_decl.hpp:254
int preOrdering_impl()
Performs pre-ordering on the matrix to increase efficiency.
Definition Amesos2_Superlu_def.hpp:171
static const char * name
Name of this solver interface.
Definition Amesos2_Superlu_decl.hpp:50
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters_impl() const
Definition Amesos2_Superlu_def.hpp:793
int solve_impl(const Teuchos::Ptr< MultiVecAdapter< Vector > > X, const Teuchos::Ptr< const MultiVecAdapter< Vector > > B) const
Superlu specific solve.
Definition Amesos2_Superlu_def.hpp:473
host_size_type_array host_rows_view_
Stores the location in Ai_ and Aval_ that starts row j.
Definition Amesos2_Superlu_decl.hpp:252
int symbolicFactorization_impl()
Perform symbolic factorization of the matrix using Superlu.
Definition Amesos2_Superlu_def.hpp:196
void setParameters_impl(const Teuchos::RCP< Teuchos::ParameterList > &parameterList)
Definition Amesos2_Superlu_def.hpp:685
host_value_type_array host_nzvals_view_
Stores the values of the nonzero entries for SuperLU.
Definition Amesos2_Superlu_decl.hpp:248
host_solve_array_t host_bValues_
Persisting 1D store for B.
Definition Amesos2_Superlu_decl.hpp:264
bool loadA_impl(EPhase current_phase)
Reads matrix data into internal structures.
Definition Amesos2_Superlu_def.hpp:937
Passes functions to TPL functions based on type.
Definition Amesos2_FunctionMap.hpp:43
Map types to solver-specific data-types and enums.
Definition Amesos2_TypeMap.hpp:48
Provides traits about solvers.
Definition Amesos2_SolverTraits.hpp:37