Amesos2 - Direct Sparse Solver Interfaces Version of the Day
Amesos2_STRUMPACK_decl.hpp
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
10
11#ifndef AMESOS2_STRUMPACK_DECL_HPP
12#define AMESOS2_STRUMPACK_DECL_HPP
13
15#include "Amesos2_SolverCore.hpp"
16
17#ifdef HAVE_MPI
18#include "StrumpackSparseSolverMPIDist.hpp"
19#else
20#include "StrumpackSparseSolver.hpp"
21#endif
22
23namespace Amesos2 {
24
25
35template <class Matrix,
36 class Vector>
37class STRUMPACK : public SolverCore<Amesos2::STRUMPACK, Matrix, Vector>
38{
39 friend class SolverCore<Amesos2::STRUMPACK,Matrix,Vector>; // Give our base access
40 // to our private
41 // implementation funcs
42public:
43
45 static const char* name; // declaration. Initialization outside.
46
47 typedef STRUMPACK<Matrix,Vector> type;
48 typedef SolverCore<Amesos2::STRUMPACK,Matrix,Vector> strum_type;
49
50 typedef Matrix matrix_type;
51 typedef Vector vector_type;
52
53 // Since typedef's are not inheritted, go grab them
54 typedef typename strum_type::scalar_type scalar_type;
55 typedef typename strum_type::local_ordinal_type local_ordinal_type;
56 typedef typename strum_type::global_ordinal_type global_ordinal_type;
57 typedef typename strum_type::global_size_type global_size_type;
58 typedef typename strum_type::node_type node_type;
59
60 typedef Kokkos::DefaultHostExecutionSpace HostExecSpaceType;
61 typedef Kokkos::View<global_ordinal_type*, HostExecSpaceType> host_ordinal_type_array;
62 typedef Kokkos::View<scalar_type*, HostExecSpaceType> host_value_type_array;
63
65
66
73 STRUMPACK(Teuchos::RCP<const Matrix> A,
74 Teuchos::RCP<Vector> X,
75 Teuchos::RCP<const Vector> B);
76
77
79 ~STRUMPACK( );
80
82
83private:
84
91 int preOrdering_impl();
92
93
102
103
110
111
123 int solve_impl(const Teuchos::Ptr<MultiVecAdapter<Vector> > X,
124 const Teuchos::Ptr<const MultiVecAdapter<Vector> > B) const;
125
126
132 bool matrixShapeOK_impl() const;
133
134
169 const Teuchos::RCP<Teuchos::ParameterList> & parameterList );
170
171
178 Teuchos::RCP<const Teuchos::ParameterList> getValidParameters_impl() const;
179
180
195 bool loadA_impl(EPhase current_phase);
196
197#ifdef HAVE_MPI
198 Teuchos::RCP<strumpack::StrumpackSparseSolverMPIDist<scalar_type,global_ordinal_type>> sp_;
199#else
200 Teuchos::RCP<strumpack::StrumpackSparseSolver<scalar_type,global_ordinal_type>> sp_;
201#endif
202
203
204 // The following Arrays are persisting storage arrays for A, X, and B
206 host_value_type_array nzvals_view_;
208 host_ordinal_type_array colind_view_;
210 host_ordinal_type_array rowptr_view_;
211 // /// 1D store for B values
212 mutable Teuchos::Array<scalar_type> bvals_;
213 // /// 1D store for X values
214 mutable Teuchos::Array<scalar_type> xvals_;
215
216#ifdef HAVE_MPI
218 Teuchos::RCP<const Tpetra::Map<local_ordinal_type,
219 global_ordinal_type,
220 node_type> > strumpack_rowmap_;
221#endif
222}; // End class STRUMPACK
223
224
225// Specialize the solver_traits template for STRUMPACK
226template <>
227struct solver_traits<STRUMPACK> {
228#if defined(HAVE_TEUCHOS_COMPLEX) && !defined(__clang__)
229 typedef Meta::make_list4<float, double, std::complex<float>, std::complex<double>> supported_scalars;
230#else
231 typedef Meta::make_list2<float, double> supported_scalars;
232#endif
233};
234
235} // end namespace Amesos2
236
237#endif // AMESOS2_STRUMPACK_DECL_HPP
Provides access to interesting solver traits.
Amesos2 interface to STRUMPACK direct solver and preconditioner.
Definition Amesos2_STRUMPACK_decl.hpp:38
void setParameters_impl(const Teuchos::RCP< Teuchos::ParameterList > &parameterList)
Definition Amesos2_STRUMPACK_def.hpp:286
static const char * name
Name of this solver interface.
Definition Amesos2_STRUMPACK_decl.hpp:45
bool loadA_impl(EPhase current_phase)
Reads matrix data into internal solver structures.
Definition Amesos2_STRUMPACK_def.hpp:456
host_value_type_array nzvals_view_
Stores the values of the nonzero entries for STRUMPACK.
Definition Amesos2_STRUMPACK_decl.hpp:206
bool matrixShapeOK_impl() const
Determines whether the shape of the matrix is OK for this solver.
Definition Amesos2_STRUMPACK_def.hpp:268
host_ordinal_type_array rowptr_view_
Stores the location in Ai_ and Aval_ that starts row j.
Definition Amesos2_STRUMPACK_decl.hpp:210
int solve_impl(const Teuchos::Ptr< MultiVecAdapter< Vector > > X, const Teuchos::Ptr< const MultiVecAdapter< Vector > > B) const
STRUMPACK specific solve.
Definition Amesos2_STRUMPACK_def.hpp:169
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters_impl() const
Definition Amesos2_STRUMPACK_def.hpp:343
int preOrdering_impl()
Performs pre-ordering on the matrix to increase efficiency.
Definition Amesos2_STRUMPACK_def.hpp:102
int numericFactorization_impl()
STRUMPACK specific numeric factorization.
Definition Amesos2_STRUMPACK_def.hpp:147
int symbolicFactorization_impl()
Perform symbolic factorization of the matrix using STRUMPACK.
Definition Amesos2_STRUMPACK_def.hpp:122
host_ordinal_type_array colind_view_
Stores the row indices of the nonzero entries.
Definition Amesos2_STRUMPACK_decl.hpp:208
Amesos2::SolverCore: A templated interface for interaction with third-party direct sparse solvers.
Definition Amesos2_SolverCore_decl.hpp:72
Provides traits about solvers.
Definition Amesos2_SolverTraits.hpp:37