19 #ifndef AMESOS2_CHOLMOD_DECL_HPP
20 #define AMESOS2_CHOLMOD_DECL_HPP
23 #include "Amesos2_SolverCore.hpp"
26 #if defined(KOKKOSKERNELS_ENABLE_SUPERNODAL_SPTRSV) && defined(KOKKOSKERNELS_ENABLE_TPL_CHOLMOD)
27 #include "KokkosKernels_Handle.hpp"
40 template <
class Matrix,
56 using scalar_type = typename super_type::scalar_type;
57 using local_ordinal_type = typename super_type::local_ordinal_type;
58 using global_ordinal_type = typename super_type::global_ordinal_type;
59 using global_size_type = typename super_type::global_size_type;
60 using node_type = typename super_type::node_type;
69 using chol_type = typename type_map::type;
70 using magnitude_type = typename type_map::magnitude_type;
83 Cholmod(Teuchos::RCP<const Matrix> A,
84 Teuchos::RCP<Vector> X,
85 Teuchos::RCP<const Vector> B);
130 int
solve_impl(const Teuchos::Ptr<MultiVecAdapter<Vector> > X,
131 const Teuchos::Ptr<const MultiVecAdapter<Vector> > B) const;
156 const Teuchos::RCP<Teuchos::ParameterList> & parameterList );
180 mutable struct CholData {
183 cholmod_dense *Y, *E;
188 typedef Kokkos::DefaultHostExecutionSpace HostExecSpaceType;
189 typedef typename HostExecSpaceType::memory_space HostMemSpaceType;
195 typedef int size_int_type;
196 typedef int ordinal_int_type;
198 typedef long size_long_type;
199 typedef long ordinal_long_type;
201 typedef Kokkos::View<size_long_type*, HostExecSpaceType> host_size_long_type_array;
202 typedef Kokkos::View<ordinal_long_type*, HostExecSpaceType> host_ordinal_long_type_array;
204 typedef Kokkos::View<size_int_type*, HostExecSpaceType> host_size_int_type_array;
205 typedef Kokkos::View<ordinal_int_type*, HostExecSpaceType> host_ordinal_int_type_array;
207 typedef Kokkos::View<chol_type*, HostExecSpaceType> host_value_type_array;
214 host_size_long_type_array host_rows_long_view_;
217 host_ordinal_long_type_array host_col_ptr_long_view_;
219 typedef typename Kokkos::View<chol_type**, Kokkos::LayoutLeft, HostExecSpaceType>
228 #if defined(KOKKOSKERNELS_ENABLE_SUPERNODAL_SPTRSV) && defined(KOKKOSKERNELS_ENABLE_TPL_CHOLMOD)
230 using DeviceExecSpaceType= Kokkos::DefaultExecutionSpace;
232 #ifdef KOKKOS_ENABLE_CUDA
234 using DeviceMemSpaceType =
typename Kokkos::CudaSpace;
235 #elif KOKKOS_ENABLE_HIP
237 using DeviceMemSpaceType =
typename Kokkos::HIPSpace;
239 using DeviceMemSpaceType =
typename DeviceExecSpaceType::memory_space;
242 typedef Kokkos::View<chol_type**, Kokkos::LayoutLeft, DeviceMemSpaceType>
243 device_solve_array_t;
246 mutable device_solve_array_t device_xValues_;
247 mutable device_solve_array_t device_bValues_;
248 typedef Kokkos::View<int*, HostMemSpaceType> host_int_array;
249 typedef Kokkos::View<int*, DeviceMemSpaceType> device_int_array;
250 host_int_array host_trsv_etree_;
251 host_int_array host_trsv_perm_;
252 device_int_array device_trsv_perm_;
253 mutable device_solve_array_t device_trsv_rhs_;
254 mutable device_solve_array_t device_trsv_sol_;
255 typedef KokkosKernels::Experimental::KokkosKernelsHandle <size_int_type, ordinal_int_type, chol_type,
256 DeviceExecSpaceType, DeviceMemSpaceType, DeviceMemSpaceType> kernel_handle_int_type;
257 typedef KokkosKernels::Experimental::KokkosKernelsHandle <size_long_type, ordinal_long_type, chol_type,
258 DeviceExecSpaceType, DeviceMemSpaceType, DeviceMemSpaceType> kernel_handle_long_type;
259 mutable kernel_handle_int_type device_int_khL_;
260 mutable kernel_handle_int_type device_int_khU_;
261 mutable kernel_handle_long_type device_long_khL_;
262 mutable kernel_handle_long_type device_long_khU_;
270 Teuchos::RCP<const Tpetra::Map<local_ordinal_type,global_ordinal_type,node_type> > map;
273 bool use_triangular_solves_;
274 bool use_cholmod_int_type_;
276 void triangular_solve_symbolic();
277 void triangular_solve_numeric();
280 void triangular_solve()
const;
287 #ifdef HAVE_TEUCHOS_COMPLEX
288 typedef Meta::make_list3<double, std::complex<double>,
289 Kokkos::complex<double>> supported_scalars;
291 typedef Meta::make_list1<double> supported_scalars;
295 template <
typename Scalar,
typename LocalOrdinal,
typename ExecutionSpace>
296 struct solver_supports_matrix<Cholmod,
297 KokkosSparse::CrsMatrix<Scalar, LocalOrdinal, ExecutionSpace>> {
298 static const bool value =
true;
303 #endif // AMESOS2_CHOLMOD_DECL_HPP
int preOrdering_impl()
Performs pre-ordering on the matrix to increase efficiency.
Definition: Amesos2_Cholmod_def.hpp:79
Amesos2::SolverCore: A templated interface for interaction with third-party direct sparse solvers...
Definition: Amesos2_SolverCore_decl.hpp:71
host_solve_array_t host_xValues_
Persisting 1D store for X.
Definition: Amesos2_Cholmod_decl.hpp:223
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters_impl() const
Definition: Amesos2_Cholmod_def.hpp:390
Template for providing a mechanism to map function calls to the correct Solver function based on the ...
Map types to solver-specific data-types and enums.
Definition: Amesos2_TypeMap.hpp:48
bool matrixShapeOK_impl() const
Determines whether the shape of the matrix is OK for this solver.
Definition: Amesos2_Cholmod_def.hpp:319
host_value_type_array host_nzvals_view_
Stores the values of the nonzero entries for CHOLMOD.
Definition: Amesos2_Cholmod_decl.hpp:211
Provides traits about solvers.
Definition: Amesos2_SolverTraits.hpp:37
int solve_impl(const Teuchos::Ptr< MultiVecAdapter< Vector > > X, const Teuchos::Ptr< const MultiVecAdapter< Vector > > B) const
CHOLMOD specific solve.
Definition: Amesos2_Cholmod_def.hpp:208
void setParameters_impl(const Teuchos::RCP< Teuchos::ParameterList > ¶meterList)
Definition: Amesos2_Cholmod_def.hpp:327
host_size_int_type_array host_rows_int_view_
Stores the location in Ai_ and Aval_ that starts row j.
Definition: Amesos2_Cholmod_decl.hpp:213
std::string name() const override
Return the name of this solver.
Definition: Amesos2_SolverCore_def.hpp:725
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_Cholmod_def.hpp:441
int symbolicFactorization_impl()
Perform symbolic factorization of the matrix using CHOLMOD.
Definition: Amesos2_Cholmod_def.hpp:110
host_ordinal_int_type_array host_col_ptr_int_view_
Stores the row indices of the nonzero entries.
Definition: Amesos2_Cholmod_decl.hpp:216
Amesos2 interface to the CHOLMOD package.
Definition: Amesos2_Cholmod_decl.hpp:42
int numericFactorization_impl()
CHOLMOD specific numeric factorization.
Definition: Amesos2_Cholmod_def.hpp:153
host_solve_array_t host_bValues_
Persisting 1D store for B.
Definition: Amesos2_Cholmod_decl.hpp:226