Amesos2 - Direct Sparse Solver Interfaces  Version of the Day
Amesos2_PardisoMKL_def.hpp
Go to the documentation of this file.
1 // @HEADER
2 //
3 // ***********************************************************************
4 //
5 // Amesos2: Templated Direct Sparse Solver Package
6 // Copyright 2011 Sandia Corporation
7 //
8 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9 // the U.S. Government retains certain rights in this software.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact Michael A. Heroux (maherou@sandia.gov)
39 //
40 // ***********************************************************************
41 //
42 // @HEADER
43 
44 
53 #ifndef AMESOS2_PARDISOMKL_DEF_HPP
54 #define AMESOS2_PARDISOMKL_DEF_HPP
55 
56 #include <map>
57 
58 #include <Teuchos_Tuple.hpp>
59 #include <Teuchos_toString.hpp>
60 #include <Teuchos_StandardParameterEntryValidators.hpp>
61 
64 
65 
66 namespace Amesos2 {
67 
68  namespace PMKL {
69 # include <mkl.h>
70 # include <mkl_pardiso.h>
71  }
72 
73  template <class Matrix, class Vector>
74  PardisoMKL<Matrix,Vector>::PardisoMKL(Teuchos::RCP<const Matrix> A,
75  Teuchos::RCP<Vector> X,
76  Teuchos::RCP<const Vector> B)
77  : SolverCore<Amesos2::PardisoMKL,Matrix,Vector>(A, X, B) // instantiate superclass
78  , nzvals_()
79  , colind_()
80  , rowptr_()
81  , n_(Teuchos::as<int_t>(this->globalNumRows_))
82  , perm_(this->globalNumRows_)
83  , nrhs_(0)
84  , is_contiguous_(true)
85  {
86  // set the default matrix type
88 
89  PMKL::_INTEGER_t iparm_temp[64];
90  PMKL::_INTEGER_t mtype_temp = mtype_;
91  PMKL::pardisoinit(pt_, &mtype_temp, iparm_temp);
92 
93  for( int i = 0; i < 64; ++i ){
94  iparm_[i] = iparm_temp[i];
95  }
96 
97  // set single or double precision
98  if( Meta::is_same<solver_magnitude_type, PMKL::_REAL_t>::value ){
99  iparm_[27] = 1; // single-precision
100  } else {
101  iparm_[27] = 0; // double-precision
102  }
103 
104  // Reset some of the default parameters
105  iparm_[34] = 1; // Use zero-based indexing
106 #ifdef HAVE_AMESOS2_DEBUG
107  iparm_[26] = 1; // turn the Pardiso matrix checker on
108 #endif
109  }
110 
111 
112  template <class Matrix, class Vector>
114  {
115  /*
116  * Free any memory allocated by the PardisoMKL library functions
117  */
118  int_t error = 0;
119  void *bdummy, *xdummy;
120 
121  if( this->root_ ){
122  int_t phase = -1; // release all internal solver memory
123  function_map::pardiso( pt_, const_cast<int_t*>(&maxfct_),
124  const_cast<int_t*>(&mnum_), &mtype_, &phase, &n_,
125  nzvals_.getRawPtr(), rowptr_.getRawPtr(),
126  colind_.getRawPtr(), perm_.getRawPtr(), &nrhs_, iparm_,
127  const_cast<int_t*>(&msglvl_), &bdummy, &xdummy, &error );
128  }
129 
130  check_pardiso_mkl_error(Amesos2::CLEAN, error);
131  }
132 
133 
134  template<class Matrix, class Vector>
135  int
137  {
138  // preOrdering done in PardisoMKL during "Analysis" (aka symbolic
139  // factorization) phase
140 
141  return(0);
142  }
143 
144 
145  template <class Matrix, class Vector>
146  int
148  {
149  int_t error = 0;
150 
151  if( this->root_ ){
152 #ifdef HAVE_AMESOS2_TIMERS
153  Teuchos::TimeMonitor symbFactTimer( this->timers_.symFactTime_ );
154 #endif
155 
156  int_t phase = 11;
157  void *bdummy, *xdummy;
158 
159  function_map::pardiso( pt_, const_cast<int_t*>(&maxfct_),
160  const_cast<int_t*>(&mnum_), &mtype_, &phase, &n_,
161  nzvals_.getRawPtr(), rowptr_.getRawPtr(),
162  colind_.getRawPtr(), perm_.getRawPtr(), &nrhs_, iparm_,
163  const_cast<int_t*>(&msglvl_), &bdummy, &xdummy, &error );
164  }
165 
166  check_pardiso_mkl_error(Amesos2::SYMBFACT, error);
167 
168  // Pardiso only lets you retrieve the total number of factor
169  // non-zeros, not for each individually. We should document how
170  // such a situation is reported.
171  this->setNnzLU(iparm_[17]);
172 
173  return(0);
174  }
175 
176 
177  template <class Matrix, class Vector>
178  int
180  {
181  int_t error = 0;
182 
183  if( this->root_ ){
184 #ifdef HAVE_AMESOS2_TIMERS
185  Teuchos::TimeMonitor numFactTimer( this->timers_.numFactTime_ );
186 #endif
187 
188  int_t phase = 22;
189  void *bdummy, *xdummy;
190 
191  function_map::pardiso( pt_, const_cast<int_t*>(&maxfct_),
192  const_cast<int_t*>(&mnum_), &mtype_, &phase, &n_,
193  nzvals_.getRawPtr(), rowptr_.getRawPtr(),
194  colind_.getRawPtr(), perm_.getRawPtr(), &nrhs_, iparm_,
195  const_cast<int_t*>(&msglvl_), &bdummy, &xdummy, &error );
196  }
197 
198  check_pardiso_mkl_error(Amesos2::NUMFACT, error);
199 
200  return( 0 );
201  }
202 
203 
204  template <class Matrix, class Vector>
205  int
207  const Teuchos::Ptr<const MultiVecAdapter<Vector> > B) const
208  {
209  using Teuchos::as;
210 
211  int_t error = 0;
212 
213  // Get B data
214  const global_size_type ld_rhs = this->root_ ? X->getGlobalLength() : 0;
215  nrhs_ = as<int_t>(X->getGlobalNumVectors());
216 
217  const size_t val_store_size = as<size_t>(ld_rhs * nrhs_);
218  xvals_.resize(val_store_size);
219  bvals_.resize(val_store_size);
220 
221  { // Get values from RHS B
222 #ifdef HAVE_AMESOS2_TIMERS
223  Teuchos::TimeMonitor mvConvTimer( this->timers_.vecConvTime_ );
224  Teuchos::TimeMonitor redistTimer( this->timers_.vecRedistTime_ );
225 #endif
226 
227  if ( is_contiguous_ == true ) {
230  solver_scalar_type>::do_get(B, bvals_(),
231  as<size_t>(ld_rhs),
232  ROOTED, this->rowIndexBase_);
233  }
234  else {
237  solver_scalar_type>::do_get(B, bvals_(),
238  as<size_t>(ld_rhs),
239  CONTIGUOUS_AND_ROOTED, this->rowIndexBase_);
240  }
241  }
242 
243  if( this->root_ ){
244 #ifdef HAVE_AMESOS2_TIMERS
245  Teuchos::TimeMonitor solveTimer( this->timers_.solveTime_ );
246 #endif
247 
248  const int_t phase = 33;
249 
250  function_map::pardiso( pt_,
251  const_cast<int_t*>(&maxfct_),
252  const_cast<int_t*>(&mnum_),
253  const_cast<int_t*>(&mtype_),
254  const_cast<int_t*>(&phase),
255  const_cast<int_t*>(&n_),
256  const_cast<solver_scalar_type*>(nzvals_.getRawPtr()),
257  const_cast<int_t*>(rowptr_.getRawPtr()),
258  const_cast<int_t*>(colind_.getRawPtr()),
259  const_cast<int_t*>(perm_.getRawPtr()),
260  &nrhs_,
261  const_cast<int_t*>(iparm_),
262  const_cast<int_t*>(&msglvl_),
263  as<void*>(bvals_.getRawPtr()),
264  as<void*>(xvals_.getRawPtr()), &error );
265  }
266 
267  check_pardiso_mkl_error(Amesos2::SOLVE, error);
268 
269  /* Export X from root to the global space */
270  {
271 #ifdef HAVE_AMESOS2_TIMERS
272  Teuchos::TimeMonitor redistTimer(this->timers_.vecRedistTime_);
273 #endif
274 
275  if ( is_contiguous_ == true ) {
278  solver_scalar_type>::do_put(X, xvals_(),
279  as<size_t>(ld_rhs),
280  ROOTED);
281  }
282  else {
285  solver_scalar_type>::do_put(X, xvals_(),
286  as<size_t>(ld_rhs),
287  CONTIGUOUS_AND_ROOTED);
288  }
289  }
290 
291  return( 0 );
292 }
293 
294 
295  template <class Matrix, class Vector>
296  bool
298  {
299  // PardisoMKL supports square matrices
300  return( this->globalNumRows_ == this->globalNumCols_ );
301  }
302 
303 
304  template <class Matrix, class Vector>
305  void
306  PardisoMKL<Matrix,Vector>::setParameters_impl(const Teuchos::RCP<Teuchos::ParameterList> & parameterList )
307  {
308  using Teuchos::RCP;
309  using Teuchos::getIntegralValue;
310  using Teuchos::ParameterEntryValidator;
311 
312  RCP<const Teuchos::ParameterList> valid_params = getValidParameters_impl();
313 
314  if( parameterList->isParameter("IPARM(2)") )
315  {
316  RCP<const ParameterEntryValidator> fillin_validator = valid_params->getEntry("IPARM(2)").validator();
317  parameterList->getEntry("IPARM(2)").setValidator(fillin_validator);
318  iparm_[1] = getIntegralValue<int>(*parameterList, "IPARM(2)");
319  }
320 
321  if( parameterList->isParameter("IPARM(4)") )
322  {
323  RCP<const ParameterEntryValidator> prec_validator = valid_params->getEntry("IPARM(4)").validator();
324  parameterList->getEntry("IPARM(4)").setValidator(prec_validator);
325  iparm_[3] = getIntegralValue<int>(*parameterList, "IPARM(4)");
326  }
327 
328  if( parameterList->isParameter("IPARM(8)") )
329  {
330  RCP<const ParameterEntryValidator> refine_validator = valid_params->getEntry("IPARM(8)").validator();
331  parameterList->getEntry("IPARM(8)").setValidator(refine_validator);
332  iparm_[7] = getIntegralValue<int>(*parameterList, "IPARM(8)");
333  }
334 
335  if( parameterList->isParameter("IPARM(10)") )
336  {
337  RCP<const ParameterEntryValidator> pivot_perturb_validator = valid_params->getEntry("IPARM(10)").validator();
338  parameterList->getEntry("IPARM(10)").setValidator(pivot_perturb_validator);
339  iparm_[9] = getIntegralValue<int>(*parameterList, "IPARM(10)");
340  }
341 
342  // First check if the control object requests a transpose solve.
343  // Then solver specific options can override this.
344  iparm_[11] = this->control_.useTranspose_ ? 2 : 0;
345 
346  if( parameterList->isParameter("IPARM(12)") )
347  {
348  RCP<const ParameterEntryValidator> trans_validator = valid_params->getEntry("IPARM(12)").validator();
349  parameterList->getEntry("IPARM(12)").setValidator(trans_validator);
350  iparm_[11] = getIntegralValue<int>(*parameterList, "IPARM(12)");
351  }
352 
353  if( parameterList->isParameter("IPARM(18)") )
354  {
355  RCP<const ParameterEntryValidator> report_validator = valid_params->getEntry("IPARM(18)").validator();
356  parameterList->getEntry("IPARM(18)").setValidator(report_validator);
357  iparm_[17] = getIntegralValue<int>(*parameterList, "IPARM(18)");
358  }
359 
360  if( parameterList->isParameter("IPARM(24)") )
361  {
362  RCP<const ParameterEntryValidator> par_fact_validator = valid_params->getEntry("IPARM(24)").validator();
363  parameterList->getEntry("IPARM(24)").setValidator(par_fact_validator);
364  iparm_[23] = getIntegralValue<int>(*parameterList, "IPARM(24)");
365  }
366 
367  if( parameterList->isParameter("IPARM(25)") )
368  {
369  RCP<const ParameterEntryValidator> par_fbsolve_validator = valid_params->getEntry("IPARM(25)").validator();
370  parameterList->getEntry("IPARM(25)").setValidator(par_fbsolve_validator);
371  iparm_[24] = getIntegralValue<int>(*parameterList, "IPARM(25)");
372  }
373 
374  if( parameterList->isParameter("IPARM(60)") )
375  {
376  RCP<const ParameterEntryValidator> ooc_validator = valid_params->getEntry("IPARM(60)").validator();
377  parameterList->getEntry("IPARM(60)").setValidator(ooc_validator);
378  iparm_[59] = getIntegralValue<int>(*parameterList, "IPARM(60)");
379  }
380 
381  if( parameterList->isParameter("IsContiguous") ){
382  is_contiguous_ = parameterList->get<bool>("IsContiguous");
383  }
384  }
385 
386 
387 /*
388  * TODO: It would be nice if the parameters could be expressed as
389  * either all string or as all integers. I see no way of doing this
390  * at present with the standard validators. However, we could create
391  * our own validators or kindly ask the Teuchos team to add some
392  * features for use.
393  *
394  * The issue is that with the current validators we cannot specify
395  * arbitrary sets of numbers that are the only allowed parameters.
396  * For example the IPARM(2) parameter can take only the values 0, 2,
397  * and 3. The EnhancedNumberValidator can take a min value, and max
398  * value, and a step size, but with those options there is no way to
399  * specify the needed set.
400  *
401  * Another missing feature is the ability to give docstrings for such
402  * numbers. For example IPARM(25) can take on the values 0 and 1.
403  * This would be easy enough to accomplish with just a number
404  * validator, but then have no way to document the effect of each
405  * value.
406  */
407 template <class Matrix, class Vector>
408 Teuchos::RCP<const Teuchos::ParameterList>
410 {
411  using std::string;
412  using Teuchos::as;
413  using Teuchos::RCP;
414  using Teuchos::tuple;
415  using Teuchos::toString;
416  using Teuchos::EnhancedNumberValidator;
417  using Teuchos::setStringToIntegralParameter;
418  using Teuchos::anyNumberParameterEntryValidator;
419 
420  static Teuchos::RCP<const Teuchos::ParameterList> valid_params;
421 
422  if( is_null(valid_params) ){
423  Teuchos::RCP<Teuchos::ParameterList> pl = Teuchos::parameterList();
424 
425  // Use pardisoinit to get some default values;
426  void *pt_dummy[64];
427  PMKL::_INTEGER_t mtype_temp = mtype_;
428  PMKL::_INTEGER_t iparm_temp[64];
429  PMKL::pardisoinit(pt_dummy,
430  const_cast<PMKL::_INTEGER_t*>(&mtype_temp),
431  const_cast<PMKL::_INTEGER_t*>(iparm_temp));
432 
433  setStringToIntegralParameter<int>("IPARM(2)", toString(iparm_temp[1]),
434  "Fill-in reducing ordering for the input matrix",
435  tuple<string>("0", "2", "3"),
436  tuple<string>("The minimum degree algorithm",
437  "Nested dissection algorithm from METIS",
438  "OpenMP parallel nested dissection algorithm"),
439  tuple<int>(0, 2, 3),
440  pl.getRawPtr());
441 
442  Teuchos::RCP<EnhancedNumberValidator<int> > iparm_4_validator
443  = Teuchos::rcp( new EnhancedNumberValidator<int>() );
444  iparm_4_validator->setMin(0);
445  pl->set("IPARM(4)" , as<int>(iparm_temp[3]) , "Preconditioned CGS/CG",
446  iparm_4_validator);
447 
448  setStringToIntegralParameter<int>("IPARM(12)", toString(iparm_temp[11]),
449  "Solve with transposed or conjugate transposed matrix A",
450  tuple<string>("0", "1", "2"),
451  tuple<string>("Non-transposed",
452  "Conjugate-transposed",
453  "Transposed"),
454  tuple<int>(0, 1, 2),
455  pl.getRawPtr());
456 
457  setStringToIntegralParameter<int>("IPARM(24)", toString(iparm_temp[23]),
458  "Parallel factorization control",
459  tuple<string>("0", "1"),
460  tuple<string>("PARDISO uses the previous algorithm for factorization",
461  "PARDISO uses the new two-level factorization algorithm"),
462  tuple<int>(0, 1),
463  pl.getRawPtr());
464 
465  setStringToIntegralParameter<int>("IPARM(25)", toString(iparm_temp[24]),
466  "Parallel forward/backward solve control",
467  tuple<string>("0", "1"),
468  tuple<string>("PARDISO uses the parallel algorithm for the solve step",
469  "PARDISO uses the sequential forward and backward solve"),
470  tuple<int>(0, 1),
471  pl.getRawPtr());
472 
473  setStringToIntegralParameter<int>("IPARM(60)", toString(iparm_temp[59]),
474  "PARDISO mode (OOC mode)",
475  tuple<string>("0", "2"),
476  tuple<string>("In-core PARDISO",
477  "Out-of-core PARDISO. The OOC PARDISO can solve very "
478  "large problems by holding the matrix factors in files "
479  "on the disk. Hence the amount of RAM required by OOC "
480  "PARDISO is significantly reduced."),
481  tuple<int>(0, 2),
482  pl.getRawPtr());
483 
484  Teuchos::AnyNumberParameterEntryValidator::EPreferredType preferred_int =
485  Teuchos::AnyNumberParameterEntryValidator::PREFER_INT;
486 
487  Teuchos::AnyNumberParameterEntryValidator::AcceptedTypes accept_int( false );
488  accept_int.allowInt( true );
489 
490  pl->set("IPARM(8)" , as<int>(iparm_temp[8]) , "Iterative refinement step",
491  anyNumberParameterEntryValidator(preferred_int, accept_int));
492 
493  pl->set("IPARM(10)", as<int>(iparm_temp[9]) , "Pivoting perturbation",
494  anyNumberParameterEntryValidator(preferred_int, accept_int));
495 
496  pl->set("IPARM(18)", as<int>(iparm_temp[17]), "Report the number of non-zero elements in the factors",
497  anyNumberParameterEntryValidator(preferred_int, accept_int));
498 
499  pl->set("IsContiguous", true, "Whether GIDs contiguous");
500 
501  valid_params = pl;
502  }
503 
504  return valid_params;
505 }
506 
507 
508 
509 template <class Matrix, class Vector>
510 bool
512 {
513 #ifdef HAVE_AMESOS2_TIMERS
514  Teuchos::TimeMonitor convTimer(this->timers_.mtxConvTime_);
515 #endif
516 
517  // PardisoMKL does not need matrix data in the pre-ordering phase
518  if( current_phase == PREORDERING ) return( false );
519 
520  if( this->root_ ){
521  nzvals_.resize(this->globalNumNonZeros_);
522  colind_.resize(this->globalNumNonZeros_);
523  rowptr_.resize(this->globalNumRows_ + 1);
524  }
525 
526  int_t nnz_ret = 0;
527  {
528 #ifdef HAVE_AMESOS2_TIMERS
529  Teuchos::TimeMonitor mtxRedistTimer( this->timers_.mtxRedistTime_ );
530 #endif
531 
532  if ( is_contiguous_ == true ) {
535  solver_scalar_type,
536  int_t,int_t>::do_get(this->matrixA_.ptr(),
537  nzvals_(), colind_(), rowptr_(),
538  nnz_ret, ROOTED, SORTED_INDICES, this->rowIndexBase_);
539  }
540  else {
543  solver_scalar_type,
544  int_t,int_t>::do_get(this->matrixA_.ptr(),
545  nzvals_(), colind_(), rowptr_(),
546  nnz_ret, CONTIGUOUS_AND_ROOTED, SORTED_INDICES, this->rowIndexBase_);
547  }
548 }
549 
550  return( true );
551 }
552 
553 
554 template <class Matrix, class Vector>
555 void
557  int_t error) const
558 {
559  int error_i = error;
560  Teuchos::broadcast(*(this->getComm()), 0, &error_i); // We only care about root's value
561 
562  if( error == 0 ) return; // No error
563 
564  std::string errmsg = "Other error";
565  switch( error ){
566  case -1:
567  errmsg = "PardisoMKL reported error: 'Input inconsistent'";
568  break;
569  case -2:
570  errmsg = "PardisoMKL reported error: 'Not enough memory'";
571  break;
572  case -3:
573  errmsg = "PardisoMKL reported error: 'Reordering problem'";
574  break;
575  case -4:
576  errmsg =
577  "PardisoMKL reported error: 'Zero pivot, numerical "
578  "factorization or iterative refinement problem'";
579  break;
580  case -5:
581  errmsg = "PardisoMKL reported error: 'Unclassified (internal) error'";
582  break;
583  case -6:
584  errmsg = "PardisoMKL reported error: 'Reordering failed'";
585  break;
586  case -7:
587  errmsg = "PardisoMKL reported error: 'Diagonal matrix is singular'";
588  break;
589  case -8:
590  errmsg = "PardisoMKL reported error: '32-bit integer overflow problem'";
591  break;
592  case -9:
593  errmsg = "PardisoMKL reported error: 'Not enough memory for OOC'";
594  break;
595  case -10:
596  errmsg = "PardisoMKL reported error: 'Problems with opening OOC temporary files'";
597  break;
598  case -11:
599  errmsg = "PardisoMKL reported error: 'Read/write problem with OOC data file'";
600  break;
601  }
602 
603  TEUCHOS_TEST_FOR_EXCEPTION( true, std::runtime_error, errmsg );
604 }
605 
606 
607 template <class Matrix, class Vector>
608 void
610 {
611  if( mtype == 0 ){
612  if( complex_ ){
613  mtype_ = 13; // complex, unsymmetric
614  } else {
615  mtype_ = 11; // real, unsymmetric
616  }
617  } else {
618  switch( mtype ){
619  case 11:
620  TEUCHOS_TEST_FOR_EXCEPTION( complex_,
621  std::invalid_argument,
622  "Cannot set a real Pardiso matrix type with scalar type complex" );
623  mtype_ = 11; break;
624  case 13:
625  TEUCHOS_TEST_FOR_EXCEPTION( !complex_,
626  std::invalid_argument,
627  "Cannot set a complex Pardiso matrix type with non-complex scalars" );
628  mtype_ = 13; break;
629  default:
630  TEUCHOS_TEST_FOR_EXCEPTION( true,
631  std::invalid_argument,
632  "Symmetric matrices are not yet supported by the Amesos2 interface" );
633  }
634  }
635 }
636 
637 
638 template <class Matrix, class Vector>
639 const char* PardisoMKL<Matrix,Vector>::name = "PARDISOMKL";
640 
641 template <class Matrix, class Vector>
642 const typename PardisoMKL<Matrix,Vector>::int_t
644 
645 template <class Matrix, class Vector>
646 const typename PardisoMKL<Matrix,Vector>::int_t
648 
649 template <class Matrix, class Vector>
650 const typename PardisoMKL<Matrix,Vector>::int_t
652 
653 
654 } // end namespace Amesos
655 
656 #endif // AMESOS2_PARDISOMKL_DEF_HPP
Amesos2::SolverCore: A templated interface for interaction with third-party direct sparse solvers...
Definition: Amesos2_SolverCore_decl.hpp:105
~PardisoMKL()
Destructor.
Definition: Amesos2_PardisoMKL_def.hpp:113
PardisoMKL(Teuchos::RCP< const Matrix > A, Teuchos::RCP< Vector > X, Teuchos::RCP< const Vector > B)
Initialize from Teuchos::RCP.
Definition: Amesos2_PardisoMKL_def.hpp:74
EPhase
Used to indicate a phase in the direct solution.
Definition: Amesos2_TypeDecl.hpp:65
Similar to get_ccs_helper , but used to get a CRS representation of the given matrix.
Definition: Amesos2_Util.hpp:780
void set_pardiso_mkl_matrix_type(int_t mtype=0)
Definition: Amesos2_PardisoMKL_def.hpp:609
Helper class for getting 1-D copies of multivectors.
Definition: Amesos2_MultiVecAdapter_decl.hpp:266
void * pt_[64]
PardisoMKL internal data address pointer.
Definition: Amesos2_PardisoMKL_decl.hpp:285
A template class that does nothing useful besides show developers what, in general, needs to be done to add a new solver interface to the Amesos2 collection.
Amesos2 interface to the PardisoMKL package.
Definition: Amesos2_PardisoMKL_decl.hpp:83
int numericFactorization_impl()
PardisoMKL specific numeric factorization.
Definition: Amesos2_PardisoMKL_def.hpp:179
int symbolicFactorization_impl()
Perform symbolic factorization of the matrix using PardisoMKL.
Definition: Amesos2_PardisoMKL_def.hpp:147
A Matrix adapter interface for Amesos2.
Definition: Amesos2_MatrixAdapter_decl.hpp:76
bool loadA_impl(EPhase current_phase)
Reads matrix data into internal structures.
Definition: Amesos2_PardisoMKL_def.hpp:511
void setParameters_impl(const Teuchos::RCP< Teuchos::ParameterList > &parameterList)
Definition: Amesos2_PardisoMKL_def.hpp:306
int solve_impl(const Teuchos::Ptr< MultiVecAdapter< Vector > > X, const Teuchos::Ptr< const MultiVecAdapter< Vector > > B) const
PardisoMKL specific solve.
Definition: Amesos2_PardisoMKL_def.hpp:206
void check_pardiso_mkl_error(EPhase phase, int_t error) const
Throws an appropriate runtime error in the event that error &lt; 0 .
Definition: Amesos2_PardisoMKL_def.hpp:556
Helper class for putting 1-D data arrays into multivectors.
Definition: Amesos2_MultiVecAdapter_decl.hpp:344
A templated MultiVector class adapter for Amesos2.
Definition: Amesos2_MultiVecAdapter_decl.hpp:176
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters_impl() const
Definition: Amesos2_PardisoMKL_def.hpp:409
bool matrixShapeOK_impl() const
Determines whether the shape of the matrix is OK for this solver.
Definition: Amesos2_PardisoMKL_def.hpp:297
int_t mtype_
The matrix type. We deal only with unsymmetrix matrices.
Definition: Amesos2_PardisoMKL_decl.hpp:287
int_t iparm_[64]
Definition: Amesos2_PardisoMKL_decl.hpp:297
int preOrdering_impl()
Performs pre-ordering on the matrix to increase efficiency.
Definition: Amesos2_PardisoMKL_def.hpp:136