Amesos2 - Direct Sparse Solver Interfaces Version of the Day
Amesos2_CssMKL_def.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
10
19#ifndef AMESOS2_CSSMKL_DEF_HPP
20#define AMESOS2_CSSMKL_DEF_HPP
21
22#include <map>
23
24#include <Teuchos_Tuple.hpp>
25#include <Teuchos_toString.hpp>
26#include <Teuchos_StandardParameterEntryValidators.hpp>
27
30
31
32namespace Amesos2 {
33
34 namespace PMKL {
35# include <mkl.h>
36# include <mkl_pardiso.h>
37 }
38
39 template <class Matrix, class Vector>
40 CssMKL<Matrix,Vector>::CssMKL(Teuchos::RCP<const Matrix> A,
41 Teuchos::RCP<Vector> X,
42 Teuchos::RCP<const Vector> B)
43 : SolverCore<Amesos2::CssMKL,Matrix,Vector>(A, X, B) // instantiate superclass
44 , n_(Teuchos::as<int_t>(this->globalNumRows_))
45 , perm_(this->globalNumRows_)
46 , nrhs_(0)
47 , css_initialized_(false)
48 , is_contiguous_(true)
49 , msglvl_(0)
50 {
51 // Matrix info
52 Teuchos::RCP<const Teuchos::Comm<int> > matComm = this->matrixA_->getComm ();
53 const global_ordinal_type indexBase = this->matrixA_->getRowMap ()->getIndexBase ();
54 const local_ordinal_type nrows = this->matrixA_->getLocalNumRows();
55
56 // rowmap for loadA (to have locally contiguous)
57 css_rowmap_ =
58 Teuchos::rcp (new map_type (this->globalNumRows_, nrows, indexBase, matComm));
59 css_contig_rowmap_ = Teuchos::rcp (new map_type (0, 0, indexBase, matComm));
60 css_contig_colmap_ = Teuchos::rcp (new map_type (0, 0, indexBase, matComm));
61
62 // set the default matrix type
64 set_css_mkl_default_parameters(pt_, iparm_);
65
66 // index base
67 iparm_[34] = (indexBase == 0 ? 1 : 0); /* Use one or zero-based indexing */
68 // 1D block-row distribution (using Contiguous map)
69 auto frow = css_rowmap_->getMinGlobalIndex();
70 iparm_[39] = 2; /* Matrix input format. */
71 iparm_[40] = frow; /* > Beginning of input domain. */
72 iparm_[41] = frow+nrows-1; /* > End of input domain. */
73
74 // get MPI Comm
75 TEUCHOS_TEST_FOR_EXCEPTION(
76 matComm.is_null (), std::logic_error, "Amesos2::CssMKL "
77 "constructor: The matrix's communicator is null!");
78 Teuchos::RCP<const Teuchos::MpiComm<int> > matMpiComm =
79 Teuchos::rcp_dynamic_cast<const Teuchos::MpiComm<int> > (matComm);
80 TEUCHOS_TEST_FOR_EXCEPTION(
81 matMpiComm.is_null (), std::logic_error, "Amesos2::CssMKL "
82 "constructor: The matrix's communicator is not an MpiComm!");
83 TEUCHOS_TEST_FOR_EXCEPTION(
84 matMpiComm->getRawMpiComm ().is_null (), std::logic_error, "Amesos2::"
85 "CssMKL constructor: The matrix's communicator claims to be a "
86 "Teuchos::MpiComm<int>, but its getRawPtrComm() method returns "
87 "Teuchos::null! This means that the underlying MPI_Comm doesn't even "
88 "exist, which likely implies that the Teuchos::MpiComm was constructed "
89 "incorrectly. It means something different than if the MPI_Comm were "
90 "MPI_COMM_NULL.");
91 MPI_Comm CssComm = *(matMpiComm->getRawMpiComm ());
92 CssComm_ = MPI_Comm_c2f(CssComm);
93 }
94
95
96 template <class Matrix, class Vector>
98 {
99 /*
100 * Free any memory allocated by the CssMKL library functions
101 */
102 int_t error = 0;
103 if (css_initialized_)
104 {
105 int_t phase = -1; // release all internal solver memory
106 void *bdummy, *xdummy;
107 const MPI_Fint CssComm = CssComm_;
108 function_map::cluster_sparse_solver( pt_, const_cast<int_t*>(&maxfct_),
109 const_cast<int_t*>(&mnum_), &mtype_, &phase, &n_,
110 nzvals_view_.data(), rowptr_view_.data(),
111 colind_view_.data(), perm_.getRawPtr(), &nrhs_, iparm_,
112 const_cast<int_t*>(&msglvl_), &bdummy, &xdummy, &CssComm, &error );
113 css_initialized_ = false;
114 }
115 check_css_mkl_error(Amesos2::CLEAN, error);
116 }
117
118
119 template<class Matrix, class Vector>
120 int
122 {
123 // preOrdering done during "Analysis" (aka symbolic
124 // factorization) phase
125 return(0);
126 }
127
128
129 template <class Matrix, class Vector>
130 int
132 {
133 if (msglvl_ > 0 && this->matrixA_->getComm()->getRank() == 0) {
134 std::cout << " CssMKL::symbolicFactorization:\n" << std::endl;
135 for (int i=0; i < 64; i++) std::cout << " * IPARM[" << i << "] = " << iparm_[i] << std::endl;
136 }
137 int_t error = 0;
138 {
139#ifdef HAVE_AMESOS2_TIMERS
140 Teuchos::TimeMonitor symbFactTimer( this->timers_.symFactTime_ );
141#endif
142
143 int_t phase = 11; // Analysis
144 void *bdummy, *xdummy;
145 const MPI_Fint CssComm = CssComm_;
146 function_map::cluster_sparse_solver( pt_, const_cast<int_t*>(&maxfct_),
147 const_cast<int_t*>(&mnum_), &mtype_, &phase, &n_,
148 nzvals_view_.data(), rowptr_view_.data(),
149 colind_view_.data(), perm_.getRawPtr(), &nrhs_, iparm_,
150 const_cast<int_t*>(&msglvl_), &bdummy, &xdummy, &CssComm, &error );
151 }
152 check_css_mkl_error(Amesos2::SYMBFACT, error);
153 if (msglvl_ > 0 && this->matrixA_->getComm()->getRank() == 0) {
154 std::cout << " CssMKL::symbolicFactorization done:" << std::endl;
155 std::cout << " * Time : " << this->timers_.symFactTime_.totalElapsedTime() << std::endl;
156 }
157
158 // Pardiso only lets you retrieve the total number of factor
159 // non-zeros, not for each individually. We should document how
160 // such a situation is reported.
161 this->setNnzLU(iparm_[17]);
162 css_initialized_ = true;
163 return(0);
164 }
165
166
167 template <class Matrix, class Vector>
168 int
170 {
171 if (msglvl_ > 0 && this->matrixA_->getComm()->getRank() == 0) {
172 std::cout << " CssMKL::numericFactorization:\n" << std::endl;
173 }
174 int_t error = 0;
175 {
176#ifdef HAVE_AMESOS2_TIMERS
177 Teuchos::TimeMonitor numFactTimer( this->timers_.numFactTime_ );
178#endif
179
180 //int_t phase = 12; // Analysis, numerical factorization
181 int_t phase = 22; // Numerical factorization
182 void *bdummy, *xdummy;
183 const MPI_Fint CssComm = CssComm_;
184 function_map::cluster_sparse_solver( pt_, const_cast<int_t*>(&maxfct_),
185 const_cast<int_t*>(&mnum_), &mtype_, &phase, &n_,
186 nzvals_view_.data(), rowptr_view_.data(),
187 colind_view_.data(), perm_.getRawPtr(), &nrhs_, iparm_,
188 const_cast<int_t*>(&msglvl_), &bdummy, &xdummy, &CssComm, &error );
189 }
190 check_css_mkl_error(Amesos2::NUMFACT, error);
191 if (msglvl_ > 0 && this->matrixA_->getComm()->getRank() == 0) {
192 std::cout << " CssMKL::numericFactorization done:" << std::endl;
193 std::cout << " Time : " << this->timers_.numFactTime_.totalElapsedTime() << std::endl;
194 }
195
196 return( 0 );
197 }
198
199
200 template <class Matrix, class Vector>
201 int
203 const Teuchos::Ptr<const MultiVecAdapter<Vector> > B) const
204 {
205 using Teuchos::as;
206
207 // Get B data
208 const local_ordinal_type ld_rhs = this->matrixA_->getLocalNumRows();
209 nrhs_ = as<int_t>(X->getGlobalNumVectors());
210
211 const size_t val_store_size = as<size_t>(ld_rhs * nrhs_);
212 xvals_.resize(val_store_size);
213 bvals_.resize(val_store_size);
214 {
215#ifdef HAVE_AMESOS2_TIMERS
216 Teuchos::TimeMonitor mvConvTimer( this->timers_.vecConvTime_ );
217 Teuchos::TimeMonitor redistTimer( this->timers_.vecRedistTime_ );
218#endif
219
222 solver_scalar_type>::do_get(B, bvals_(),
223 as<size_t>(ld_rhs),
224 Teuchos::ptrInArg(*css_rowmap_));
225 }
226
227 int_t error = 0;
228 {
229#ifdef HAVE_AMESOS2_TIMERS
230 Teuchos::TimeMonitor solveTimer( this->timers_.solveTime_ );
231#endif
232
233 const int_t phase = 33; // Solve, iterative refinement
234 const MPI_Fint CssComm = CssComm_;
235 function_map::cluster_sparse_solver( pt_,
236 const_cast<int_t*>(&maxfct_),
237 const_cast<int_t*>(&mnum_),
238 const_cast<int_t*>(&mtype_),
239 const_cast<int_t*>(&phase),
240 const_cast<int_t*>(&n_),
241 const_cast<solver_scalar_type*>(nzvals_view_.data()),
242 const_cast<int_t*>(rowptr_view_.data()),
243 const_cast<int_t*>(colind_view_.data()),
244 const_cast<int_t*>(perm_.getRawPtr()),
245 &nrhs_,
246 const_cast<int_t*>(iparm_),
247 const_cast<int_t*>(&msglvl_),
248 as<void*>(bvals_.getRawPtr()),
249 as<void*>(xvals_.getRawPtr()), &CssComm, &error );
250 }
251 check_css_mkl_error(Amesos2::SOLVE, error);
252
253 /* Get values to X */
254 {
255#ifdef HAVE_AMESOS2_TIMERS
256 Teuchos::TimeMonitor redistTimer(this->timers_.vecRedistTime_);
257#endif
258
261 solver_scalar_type>::do_put(X, xvals_(),
262 as<size_t>(ld_rhs),
263 Teuchos::ptrInArg(*css_rowmap_));
264 }
265
266 return( 0 );
267 }
268
269
270 template <class Matrix, class Vector>
271 bool
273 {
274 // CssMKL supports square matrices
275 return( this->globalNumRows_ == this->globalNumCols_ );
276 }
277
278
279 template <class Matrix, class Vector>
280 void
281 CssMKL<Matrix,Vector>::setParameters_impl(const Teuchos::RCP<Teuchos::ParameterList> & parameterList )
282 {
283 using Teuchos::RCP;
284 using Teuchos::getIntegralValue;
285 using Teuchos::ParameterEntryValidator;
286
287 RCP<const Teuchos::ParameterList> valid_params = getValidParameters_impl();
288
289 // 2: Fill-in reordering from METIS, 3: thread dissection, 10: MPI version of the nested dissection
290 if( parameterList->isParameter("IPARM(2)") )
291 {
292 RCP<const ParameterEntryValidator> fillin_validator = valid_params->getEntry("IPARM(2)").validator();
293 parameterList->getEntry("IPARM(2)").setValidator(fillin_validator);
294 iparm_[1] = getIntegralValue<int>(*parameterList, "IPARM(2)");
295 }
296
297 // Max numbers of iterative refinement steps
298 if( parameterList->isParameter("IPARM(8)") )
299 {
300 RCP<const ParameterEntryValidator> refine_validator = valid_params->getEntry("IPARM(8)").validator();
301 parameterList->getEntry("IPARM(8)").setValidator(refine_validator);
302 iparm_[7] = getIntegralValue<int>(*parameterList, "IPARM(8)");
303 }
304
305 // Perturb the pivot elements
306 if( parameterList->isParameter("IPARM(10)") )
307 {
308 RCP<const ParameterEntryValidator> pivot_perturb_validator = valid_params->getEntry("IPARM(10)").validator();
309 parameterList->getEntry("IPARM(10)").setValidator(pivot_perturb_validator);
310 iparm_[9] = getIntegralValue<int>(*parameterList, "IPARM(10)");
311 }
312
313 // First check if the control object requests a transpose solve.
314 // Then solver specific options can override this.
315 iparm_[11] = this->control_.useTranspose_ ? 2 : 0;
316 // Normal solve (0), or a transpose solve (1)
317 if( parameterList->isParameter("IPARM(12)") )
318 {
319 RCP<const ParameterEntryValidator> trans_validator = valid_params->getEntry("IPARM(12)").validator();
320 parameterList->getEntry("IPARM(12)").setValidator(trans_validator);
321 iparm_[11] = getIntegralValue<int>(*parameterList, "IPARM(12)");
322 }
323
324 // (Non-)symmetric matchings : detault 1 for nonsymmetric and 0 for symmetric matrix (default is nonsymmetric)
325 if( parameterList->isParameter("IPARM(13)") )
326 {
327 RCP<const ParameterEntryValidator> trans_validator = valid_params->getEntry("IPARM(13)").validator();
328 parameterList->getEntry("IPARM(13)").setValidator(trans_validator);
329 iparm_[12] = getIntegralValue<int>(*parameterList, "IPARM(13)");
330 }
331
332 // Output: Number of nonzeros in the factor LU
333 if( parameterList->isParameter("IPARM(18)") )
334 {
335 RCP<const ParameterEntryValidator> report_validator = valid_params->getEntry("IPARM(18)").validator();
336 parameterList->getEntry("IPARM(18)").setValidator(report_validator);
337 iparm_[17] = getIntegralValue<int>(*parameterList, "IPARM(18)");
338 }
339
340 // Check input matrix is sorted
341 if( parameterList->isParameter("IPARM(28)") )
342 {
343 RCP<const ParameterEntryValidator> report_validator = valid_params->getEntry("IPARM(28)").validator();
344 parameterList->getEntry("IPARM(28)").setValidator(report_validator);
345 iparm_[27] = getIntegralValue<int>(*parameterList, "IPARM(28)");
346 }
347
348 if( parameterList->isParameter("IsContiguous") ){
349 is_contiguous_ = parameterList->get<bool>("IsContiguous");
350 }
351
352 if( parameterList->isParameter("verbose") ){
353 msglvl_ = parameterList->get<int>("verbose");
354 }
355 }
356
357
358/*
359 * TODO: It would be nice if the parameters could be expressed as
360 * either all string or as all integers. I see no way of doing this
361 * at present with the standard validators. However, we could create
362 * our own validators or kindly ask the Teuchos team to add some
363 * features for use.
364 *
365 * The issue is that with the current validators we cannot specify
366 * arbitrary sets of numbers that are the only allowed parameters.
367 * For example the IPARM(2) parameter can take only the values 0, 2,
368 * and 3. The EnhancedNumberValidator can take a min value, and max
369 * value, and a step size, but with those options there is no way to
370 * specify the needed set.
371 *
372 * Another missing feature is the ability to give docstrings for such
373 * numbers. For example IPARM(25) can take on the values 0 and 1.
374 * This would be easy enough to accomplish with just a number
375 * validator, but then have no way to document the effect of each
376 * value.
377 */
378template <class Matrix, class Vector>
379Teuchos::RCP<const Teuchos::ParameterList>
381{
382 using std::string;
383 using Teuchos::as;
384 using Teuchos::RCP;
385 using Teuchos::tuple;
386 using Teuchos::toString;
387 using Teuchos::EnhancedNumberValidator;
388 using Teuchos::setStringToIntegralParameter;
389 using Teuchos::anyNumberParameterEntryValidator;
390
391 static Teuchos::RCP<const Teuchos::ParameterList> valid_params;
392
393 if( is_null(valid_params) ){
394 Teuchos::RCP<Teuchos::ParameterList> pl = Teuchos::parameterList();
395
396 void* pt_temp[64];
397 int_t iparm_temp[64];
398 set_css_mkl_default_parameters(pt_temp, iparm_temp);
399 setStringToIntegralParameter<int>("IPARM(2)", toString(iparm_temp[1]),
400 "Fill-in reducing ordering for the input matrix",
401 tuple<string>("2", "3", "10"),
402 tuple<string>("Nested dissection algorithm from METIS",
403 "Parallel version of the nested dissection algorithm",
404 "MPI version of the nested dissection and symbolic factorization algorithms"),
405 tuple<int>(2, 3, 10),
406 pl.getRawPtr());
407
408 setStringToIntegralParameter<int>("IPARM(12)", toString(iparm_temp[11]),
409 "Solve with transposed or conjugate transposed matrix A",
410 tuple<string>("0", "1", "2"),
411 tuple<string>("Non-transposed",
412 "Conjugate-transposed",
413 "Transposed"),
414 tuple<int>(0, 1, 2),
415 pl.getRawPtr());
416
417 setStringToIntegralParameter<int>("IPARM(13)", toString(iparm_temp[12]),
418 "Use weighted matching",
419 tuple<string>("0", "1"),
420 tuple<string>("No matching", "Use matching"),
421 tuple<int>(0, 1),
422 pl.getRawPtr());
423
424 Teuchos::AnyNumberParameterEntryValidator::EPreferredType preferred_int =
425 Teuchos::AnyNumberParameterEntryValidator::PREFER_INT;
426
427 Teuchos::AnyNumberParameterEntryValidator::AcceptedTypes accept_int( false );
428 accept_int.allowInt( true );
429
430 pl->set("IPARM(8)" , as<int>(iparm_temp[7]) , "Iterative refinement step",
431 anyNumberParameterEntryValidator(preferred_int, accept_int));
432
433 pl->set("IPARM(10)", as<int>(iparm_temp[9]) , "Pivoting perturbation",
434 anyNumberParameterEntryValidator(preferred_int, accept_int));
435
436 pl->set("IPARM(18)", as<int>(iparm_temp[17]), "Report the number of non-zero elements in the factors",
437 anyNumberParameterEntryValidator(preferred_int, accept_int));
438
439 pl->set("IPARM(28)", as<int>(iparm_temp[27]), "Check input matrix is sorted",
440 anyNumberParameterEntryValidator(preferred_int, accept_int));
441
442 pl->set("IsContiguous", true, "Whether GIDs contiguous");
443
444 pl->set("verbose", 0, "Verbosity Message Level");
445
446 valid_params = pl;
447 }
448
449 return valid_params;
450}
451
452
453
454template <class Matrix, class Vector>
455bool
457{
458#ifdef HAVE_AMESOS2_TIMERS
459 Teuchos::TimeMonitor convTimer(this->timers_.mtxConvTime_);
460#endif
461
462 // CssMKL does not need matrix data in the pre-ordering phase
463 if( current_phase == PREORDERING ) return( false );
464
465 // is_contiguous : input is contiguous
466 // CONTIGUOUS_AND_ROOTED : input is not contiguous, so make output contiguous
467 EDistribution dist_option = (iparm_[39] != 0 ? DISTRIBUTED_NO_OVERLAP : ((is_contiguous_ == true) ? ROOTED : CONTIGUOUS_AND_ROOTED));
468 if (dist_option == DISTRIBUTED_NO_OVERLAP && !is_contiguous_) {
469 // Neeed to form contiguous matrix
470 #if 1
471 // Only reinex GIDs
472 css_rowmap_ = this->matrixA_->getRowMap(); // use original map to redistribute vectors in solve
473 Teuchos::RCP<const MatrixAdapter<Matrix> > contig_mat = this->matrixA_->reindex(css_contig_rowmap_, css_contig_colmap_, current_phase);
474 #else
475 // Redistribued matrixA into contiguous GIDs
476 Teuchos::RCP<const MatrixAdapter<Matrix> > contig_mat = this->matrixA_->get(ptrInArg(*css_rowmap_));
477 //css_rowmap_ = contig_mat->getRowMap(); // use new map to redistribute vectors in solve
478 #endif
479 // Copy into local views
480 if (current_phase == SYMBFACT) {
481 Kokkos::resize(nzvals_temp_, contig_mat->getLocalNNZ());
482 Kokkos::resize(nzvals_view_, contig_mat->getLocalNNZ());
483 Kokkos::resize(colind_view_, contig_mat->getLocalNNZ());
484 Kokkos::resize(rowptr_view_, contig_mat->getLocalNumRows() + 1);
485 }
486 int_t nnz_ret = 0;
487 {
488#ifdef HAVE_AMESOS2_TIMERS
489 Teuchos::TimeMonitor mtxRedistTimer( this->timers_.mtxRedistTime_ );
490#endif
492 host_value_type_array,host_ordinal_type_array, host_size_type_array >::do_get(
493 contig_mat.ptr(),
494 nzvals_temp_, colind_view_, rowptr_view_,
495 nnz_ret,
496 ptrInArg(*(contig_mat->getRowMap())),
497 #if 1
499 #else
500 ROOTED,
501 #endif
503 Kokkos::deep_copy(nzvals_view_, nzvals_temp_);
504 }
505 } else {
506 if (current_phase == SYMBFACT) {
507 if (dist_option == DISTRIBUTED_NO_OVERLAP) {
508 Kokkos::resize(nzvals_temp_, this->matrixA_->getLocalNNZ());
509 Kokkos::resize(nzvals_view_, this->matrixA_->getLocalNNZ());
510 Kokkos::resize(colind_view_, this->matrixA_->getLocalNNZ());
511 Kokkos::resize(rowptr_view_, this->matrixA_->getLocalNumRows() + 1);
512 } else {
513 if( this->root_ ) {
514 Kokkos::resize(nzvals_temp_, this->matrixA_->getGlobalNNZ());
515 Kokkos::resize(nzvals_view_, this->matrixA_->getGlobalNNZ());
516 Kokkos::resize(colind_view_, this->matrixA_->getGlobalNNZ());
517 Kokkos::resize(rowptr_view_, this->matrixA_->getGlobalNumRows() + 1);
518 } else {
519 Kokkos::resize(nzvals_temp_, 0);
520 Kokkos::resize(nzvals_view_, 0);
521 Kokkos::resize(colind_view_, 0);
522 Kokkos::resize(rowptr_view_, 0);
523 }
524 }
525 }
526 int_t nnz_ret = 0;
527 {
528#ifdef HAVE_AMESOS2_TIMERS
529 Teuchos::TimeMonitor mtxRedistTimer( this->timers_.mtxRedistTime_ );
530#endif
532 host_value_type_array,host_ordinal_type_array, host_size_type_array >::do_get(
533 this->matrixA_.ptr(),
534 nzvals_temp_, colind_view_, rowptr_view_,
535 nnz_ret,
536 ptrInArg(*(this->matrixA_->getRowMap())),
537 dist_option,
539 Kokkos::deep_copy(nzvals_view_, nzvals_temp_);
540 }
541 }
542 return( true );
543}
544
545
546template <class Matrix, class Vector>
547void
549 int_t error) const
550{
551 int error_i = error;
552 Teuchos::broadcast(*(this->getComm()), 0, &error_i); // We only care about root's value
553
554 if( error == 0 ) return; // No error
555
556 std::string errmsg = "Other error";
557 switch( error ){
558 case -1:
559 errmsg = "CssMKL reported error: 'Input inconsistent'";
560 break;
561 case -2:
562 errmsg = "CssMKL reported error: 'Not enough memory'";
563 break;
564 case -3:
565 errmsg = "CssMKL reported error: 'Reordering problem'";
566 break;
567 case -4:
568 errmsg =
569 "CssMKL reported error: 'Zero pivot, numerical "
570 "factorization or iterative refinement problem'";
571 break;
572 case -5:
573 errmsg = "CssMKL reported error: 'Unclassified (internal) error'";
574 break;
575 case -6:
576 errmsg = "CssMKL reported error: 'Reordering failed'";
577 break;
578 case -7:
579 errmsg = "CssMKL reported error: 'Diagonal matrix is singular'";
580 break;
581 case -8:
582 errmsg = "CssMKL reported error: '32-bit integer overflow problem'";
583 break;
584 case -9:
585 errmsg = "CssMKL reported error: 'Not enough memory for OOC'";
586 break;
587 case -10:
588 errmsg = "CssMKL reported error: 'Problems with opening OOC temporary files'";
589 break;
590 case -11:
591 errmsg = "CssMKL reported error: 'Read/write problem with OOC data file'";
592 break;
593 }
594 errmsg += (" at phase = "+std::to_string(phase));
595
596 TEUCHOS_TEST_FOR_EXCEPTION( true, std::runtime_error, errmsg );
597}
598
599
600template <class Matrix, class Vector>
601void
603{
604 if( mtype == 0 ){
605 if( complex_ ){
606 mtype_ = 13; // complex, unsymmetric
607 } else {
608 mtype_ = 11; // real, unsymmetric
609 }
610 } else {
611 switch( mtype ){
612 case 11:
613 TEUCHOS_TEST_FOR_EXCEPTION( complex_,
614 std::invalid_argument,
615 "Cannot set a real Pardiso matrix type with scalar type complex" );
616 mtype_ = 11; break;
617 case 13:
618 TEUCHOS_TEST_FOR_EXCEPTION( !complex_,
619 std::invalid_argument,
620 "Cannot set a complex Pardiso matrix type with non-complex scalars" );
621 mtype_ = 13; break;
622 default:
623 TEUCHOS_TEST_FOR_EXCEPTION( true,
624 std::invalid_argument,
625 "Symmetric matrices are not yet supported by the Amesos2 interface" );
626 }
627 }
628}
629
630template <class Matrix, class Vector>
631void
632CssMKL<Matrix,Vector>::set_css_mkl_default_parameters(void* pt[], int_t iparm[]) const
633{
634 for( int i = 0; i < 64; ++i ){
635 pt[i] = nullptr;
636 iparm[i] = 0;
637 }
638 iparm[0] = 1; /* No solver default */
639 // Reset some of the default parameters
640 iparm[1] = 10; /* 2: Fill-in reordering from METIS, 3: thread dissection, 10: MPI version of the nested dissection and symbolic factorization*/
641 iparm[7] = 0; /* Max numbers of iterative refinement steps */
642 iparm[10] = 0; /* Disable nonsymmetric permutation and scaling MPS */
643 iparm[11] = 0; /* Normal solve (0), or a transpose solve (1) */
644 iparm[12] = 0; /* Do not use (non-)symmetric matchings */
645 iparm[17] = -1; /* Output: Number of nonzeros in the factor LU */
646 iparm[20] = 1; /* Pivoting for symmetric indefinite matrices */
647 iparm[26] = 1; /* Check input matrix is sorted */
648
649 // diagonal pertubation
650 if (mtype_ == -2 || mtype_ == -4) {
651 // symmetric indefinite
652 iparm[9] = 8; /* Perturb the pivot elements with 1E-8 */
653 } else {
654 // non-symmetric
655 iparm[9] = 13; /* Perturb the pivot elements with 1E-13 */
656 }
657
658 // set single or double precision
659 if constexpr ( std::is_same_v<solver_magnitude_type, PMKL::_REAL_t> ) {
660 iparm[27] = 1; // single-precision
661 } else {
662 iparm[27] = 0; // double-precision
663 }
664 iparm[34] = 1; /* Use zero-based indexing */
665}
666
667
668template <class Matrix, class Vector>
669const char* CssMKL<Matrix,Vector>::name = "CSSMKL";
670
671template <class Matrix, class Vector>
672const typename CssMKL<Matrix,Vector>::int_t
673CssMKL<Matrix,Vector>::maxfct_ = 1;
674
675template <class Matrix, class Vector>
676const typename CssMKL<Matrix,Vector>::int_t
677CssMKL<Matrix,Vector>::mnum_ = 1;
678
679
680} // end namespace Amesos
681
682#endif // AMESOS2_CSSMKL_DEF_HPP
A template class that does nothing useful besides show developers what, in general,...
EDistribution
Definition Amesos2_TypeDecl.hpp:89
@ DISTRIBUTED_NO_OVERLAP
Definition Amesos2_TypeDecl.hpp:91
@ ROOTED
Definition Amesos2_TypeDecl.hpp:93
@ CONTIGUOUS_AND_ROOTED
Definition Amesos2_TypeDecl.hpp:94
@ SORTED_INDICES
Definition Amesos2_TypeDecl.hpp:108
Amesos2 interface to the CssMKL package.
Definition Amesos2_CssMKL_decl.hpp:50
void * pt_[64]
CssMKL internal data address pointer.
Definition Amesos2_CssMKL_decl.hpp:261
bool loadA_impl(EPhase current_phase)
Reads matrix data into internal structures.
Definition Amesos2_CssMKL_def.hpp:456
~CssMKL()
Destructor.
Definition Amesos2_CssMKL_def.hpp:97
void check_css_mkl_error(EPhase phase, int_t error) const
Throws an appropriate runtime error in the event that error < 0 .
Definition Amesos2_CssMKL_def.hpp:548
static const char * name
The name of this solver interface.
Definition Amesos2_CssMKL_decl.hpp:57
int preOrdering_impl()
Performs pre-ordering on the matrix to increase efficiency.
Definition Amesos2_CssMKL_def.hpp:121
int numericFactorization_impl()
CssMKL specific numeric factorization.
Definition Amesos2_CssMKL_def.hpp:169
int symbolicFactorization_impl()
Perform symbolic factorization of the matrix using CssMKL.
Definition Amesos2_CssMKL_def.hpp:131
bool matrixShapeOK_impl() const
Determines whether the shape of the matrix is OK for this solver.
Definition Amesos2_CssMKL_def.hpp:272
int solve_impl(const Teuchos::Ptr< MultiVecAdapter< Vector > > X, const Teuchos::Ptr< const MultiVecAdapter< Vector > > B) const
CssMKL specific solve.
Definition Amesos2_CssMKL_def.hpp:202
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters_impl() const
Definition Amesos2_CssMKL_def.hpp:380
int_t iparm_[64]
Definition Amesos2_CssMKL_decl.hpp:279
void set_css_mkl_matrix_type(int_t mtype=0)
Definition Amesos2_CssMKL_def.hpp:602
void setParameters_impl(const Teuchos::RCP< Teuchos::ParameterList > &parameterList)
Definition Amesos2_CssMKL_def.hpp:281
CssMKL(Teuchos::RCP< const Matrix > A, Teuchos::RCP< Vector > X, Teuchos::RCP< const Vector > B)
Initialize from Teuchos::RCP.
Definition Amesos2_CssMKL_def.hpp:40
Amesos2::SolverCore: A templated interface for interaction with third-party direct sparse solvers.
Definition Amesos2_SolverCore_decl.hpp:72
Teuchos::RCP< const MatrixAdapter< Matrix > > matrixA_
The LHS operator.
Definition Amesos2_SolverCore_decl.hpp:421
Teuchos::RCP< const Teuchos::Comm< int > > getComm() const override
Returns a pointer to the Teuchos::Comm communicator with this operator.
Definition Amesos2_SolverCore_decl.hpp:329
global_size_type globalNumRows_
Number of global rows in matrixA_.
Definition Amesos2_SolverCore_decl.hpp:442
EPhase
Used to indicate a phase in the direct solution.
Definition Amesos2_TypeDecl.hpp:31
A templated MultiVector class adapter for Amesos2.
Definition Amesos2_MultiVecAdapter_decl.hpp:142
Helper class for getting 1-D copies of multivectors.
Definition Amesos2_MultiVecAdapter_decl.hpp:233
Similar to get_ccs_helper , but used to get a CRS representation of the given matrix.
Definition Amesos2_Util.hpp:644
Helper class for putting 1-D data arrays into multivectors.
Definition Amesos2_MultiVecAdapter_decl.hpp:339