Amesos2 - Direct Sparse Solver Interfaces  Version of the Day
Amesos2_MUMPS_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 
18 #ifndef AMESOS2_MUMPS_DECL_HPP
19 #define AMESOS2_MUMPS_DECL_HPP
20 
21 #include "Amesos2_SolverTraits.hpp"
22 #include "Amesos2_SolverCore.hpp"
24 
25 
26 namespace Amesos2 {
27 
28 
49 template <class Matrix,class Vector>
50 class MUMPS : public SolverCore<Amesos2::MUMPS, Matrix, Vector>
51 {
52  friend class SolverCore<Amesos2::MUMPS,Matrix,Vector>; // Give our base access
53  // to our private
54  // implementation funcs
55 public:
56 
58  static const char* name; // declaration. Initialization outside.
59 
60  typedef MUMPS<Matrix,Vector> type;
61  typedef SolverCore<Amesos2::MUMPS,Matrix,Vector> super_type;
62 
63  // Since typedef's are not inheritted, go grab them
64  typedef typename super_type::scalar_type scalar_type;
65  typedef typename super_type::local_ordinal_type local_ordinal_type;
66  typedef typename super_type::global_ordinal_type global_ordinal_type;
67  typedef typename super_type::global_size_type global_size_type;
68 
69  typedef TypeMap<Amesos2::MUMPS,scalar_type> type_map;
70 
71  typedef typename type_map::type mumps_type;
72  typedef typename type_map::magnitude_type magnitude_type;
73  typedef typename type_map::MUMPS_STRUC_C MUMPS_STRUC_C;
74 
75  typedef Kokkos::DefaultHostExecutionSpace HostExecSpaceType;
76  typedef typename HostExecSpaceType::memory_space HostMemSpaceType;
77 
78  typedef Kokkos::View<local_ordinal_type*, HostExecSpaceType> host_ordinal_type_view;
79  typedef Kokkos::View<mumps_type*, HostExecSpaceType> host_value_type_view;
80 
81  MUMPS(Teuchos::RCP<const Matrix> A,
82  Teuchos::RCP<Vector> X,
83  Teuchos::RCP<const Vector> B);
84  ~MUMPS( );
85 
86 
87 private:
88 
94  int preOrdering_impl();
95 
96 
97  int symbolicFactorization_impl();
98 
99 
106 
107 
119  int solve_impl(const Teuchos::Ptr<MultiVecAdapter<Vector> > X,
120  const Teuchos::Ptr<const MultiVecAdapter<Vector> > B) const;
121 
122 
126  bool matrixShapeOK_impl() const;
127 
143  void setParameters_impl(
144  const Teuchos::RCP<Teuchos::ParameterList> & parameterList );
145 
146 
153  Teuchos::RCP<const Teuchos::ParameterList> getValidParameters_impl() const;
154 
155 
164  bool loadA_impl(EPhase current_phase);
165 
166 
167  int ConvertToTriplet();
168 
169  void MUMPS_ERROR() const;
170 
171 #ifdef HAVE_MPI
172  MPI_Comm MUMPSComm;
173 #endif
174 
175 
176  bool MUMPS_MATRIX_LOAD;
177  bool MUMPS_STRUCT;
178  mutable bool MUMPS_MATRIX_LOAD_PREORDERING;
179 
180  // The following Arrays are persisting storage arrays for A, X, and B
182  host_value_type_view host_nzvals_view_;
184  host_ordinal_type_view host_rows_view_;
186  host_ordinal_type_view host_col_ptr_view_;
187 
189  mutable Teuchos::Array<mumps_type> xvals_; local_ordinal_type ldx_;
191  mutable Teuchos::Array<mumps_type> bvals_; local_ordinal_type ldb_;
192 
193  mutable MUMPS_STRUC_C mumps_par;
194 
195  bool is_contiguous_;
196 
197 }; // End class MUMPS
198 
199 
200 template <class Matrix,class Vector>
201 class MUMPSNS : public SolverCore<Amesos2::MUMPS, Matrix, Vector>
202 {
203  friend class SolverCore<Amesos2::MUMPS,Matrix,Vector>; // Give our base access
204  // to our private
205  // implementation funcs
206 public:
207 
208 
209  MUMPSNS(Teuchos::RCP<const Matrix> A,
210  Teuchos::RCP<Vector> X,
211  Teuchos::RCP<const Vector> B);
212  /*
213  {
214  TEUCHOS_TEST_FOR_EXCEPTION(true, std::logic_error,
215  "This solver is not support for these types.");
216 
217  }
218 
219  */
220  ~MUMPSNS( )
221  {}
222 
223 
224 private:
225 
226  int preOrdering_impl()
227  {
228 
229  TEUCHOS_TEST_FOR_EXCEPTION(true, std::logic_error,
230  "This solver is not support for these types.");
231 
232  }
233 
234  int symbolicFactorization_impl()
235  {
236 
237  TEUCHOS_TEST_FOR_EXCEPTION(true, std::logic_error,
238  "This solver is not support for these types.");
239 
240  }
241 
243  {
244 
245  TEUCHOS_TEST_FOR_EXCEPTION(true, std::logic_error,
246  "This solver is not support for these types.");
247 
248  }
249 
250  int solve_impl(const Teuchos::Ptr<MultiVecAdapter<Vector> > X,
251  const Teuchos::Ptr<const MultiVecAdapter<Vector> > B) const
252  {
253 
254  TEUCHOS_TEST_FOR_EXCEPTION(true, std::logic_error,
255  "This solver is not support for these types.");
256 
257  }
258 
259  bool matrixShapeOK_impl() const
260  {
261 
262  TEUCHOS_TEST_FOR_EXCEPTION(true, std::logic_error,
263  "This solver is not support for these types.");
264 
265  }
266 
267  void setParameters_impl(
268  const Teuchos::RCP<Teuchos::ParameterList> & parameterList )
269  {
270 
271  TEUCHOS_TEST_FOR_EXCEPTION(true, std::logic_error,
272  "This solver is not support for these types.");
273 
274  }
275 
276  Teuchos::RCP<const Teuchos::ParameterList> getValidParameters_impl() const
277  {
278 
279  TEUCHOS_TEST_FOR_EXCEPTION(true, std::logic_error,
280  "This solver is not support for these types.");
281 
282  }
283 
284 
285  bool loadA_impl(EPhase current_phase)
286  {
287 
288  TEUCHOS_TEST_FOR_EXCEPTION(true, std::logic_error,
289  "This solver is not support for these types.");
290 
291  }
292 
293 
294  int ConvertToTriplet()
295  {
296 
297  TEUCHOS_TEST_FOR_EXCEPTION(true, std::logic_error,
298  "This solver is not support for these types.");
299 
300  }
301 
302 
303  void MUMPS_ERROR() const
304  {
305 
306  TEUCHOS_TEST_FOR_EXCEPTION(true, std::logic_error,
307  "This solver is not support for these types.");
308 
309  }
310 
311 
312 #ifdef HAVE_MPI
313  MPI_Comm MUMPSComm;
314 #endif
315 
316  bool is_contiguous_;
317 
318 };
319 
320 
321 
322 // Specialize solver_traits struct for MUMPS
323 // TODO
324 template <>
325 struct solver_traits<MUMPS> {
326 #ifdef HAVE_TEUCHOS_COMPLEX
327  typedef Meta::make_list4<float,
328  double,
329  std::complex<float>,
330  std::complex<double> > supported_scalars;
331 #else
332  typedef Meta::make_list2<float, double> supported_scalars;
333 #endif
334 };
335 
336 } // end namespace Amesos2
337 
338 #endif // AMESOS2_MUMPS_DECL_HPP
Provides a mechanism to map function calls to the correct Solver function based on the scalar type of...
host_ordinal_type_view host_col_ptr_view_
Stores the row indices of the nonzero entries.
Definition: Amesos2_MUMPS_decl.hpp:186
Amesos2::SolverCore: A templated interface for interaction with third-party direct sparse solvers...
Definition: Amesos2_SolverCore_decl.hpp:71
int numericFactorization_impl()
MUMPS specific numeric factorization.
Definition: Amesos2_MUMPS_def.hpp:168
EPhase
Used to indicate a phase in the direct solution.
Definition: Amesos2_TypeDecl.hpp:31
Teuchos::Array< mumps_type > xvals_
Persisting 1D store for X.
Definition: Amesos2_MUMPS_decl.hpp:189
Map types to solver-specific data-types and enums.
Definition: Amesos2_TypeMap.hpp:48
int preOrdering_impl()
Performs pre-ordering on the matrix to increase efficiency.
Definition: Amesos2_MUMPS_def.hpp:139
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
MUMPS specific solve.
Definition: Amesos2_MUMPS_def.hpp:191
host_ordinal_type_view host_rows_view_
Stores the location in Ai_ and Aval_ that starts row j.
Definition: Amesos2_MUMPS_decl.hpp:184
host_value_type_view host_nzvals_view_
Stores the values of the nonzero entries for MUMPS.
Definition: Amesos2_MUMPS_decl.hpp:182
bool matrixShapeOK_impl() const
Determines whether the shape of the matrix is OK for this solver.
Definition: Amesos2_MUMPS_def.hpp:254
bool loadA_impl(EPhase current_phase)
Reads matrix data into internal structures.
Definition: Amesos2_MUMPS_def.hpp:334
Amesos2 interface to the MUMPS package.
Definition: Amesos2_MUMPS_decl.hpp:50
Teuchos::Array< mumps_type > bvals_
Persisting 1D store for B.
Definition: Amesos2_MUMPS_decl.hpp:191
std::string name() const override
Return the name of this solver.
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters_impl() const
Definition: Amesos2_MUMPS_def.hpp:306
void setParameters_impl(const Teuchos::RCP< Teuchos::ParameterList > &parameterList)
Definition: Amesos2_MUMPS_def.hpp:263
Provides access to interesting solver traits.
A templated MultiVector class adapter for Amesos2.
Definition: Amesos2_MultiVecAdapter_decl.hpp:142