52 #ifndef AMESOS2_SUPERLUDIST_DEF_HPP
53 #define AMESOS2_SUPERLUDIST_DEF_HPP
55 #include <Teuchos_Tuple.hpp>
56 #include <Teuchos_StandardParameterEntryValidators.hpp>
57 #include <Teuchos_DefaultMpiComm.hpp>
67 template <
class Matrix,
class Vector>
69 Teuchos::RCP<Vector> X,
70 Teuchos::RCP<const Vector> B)
78 , is_contiguous_(true)
83 using Teuchos::MpiComm;
84 using Teuchos::outArg;
85 using Teuchos::ParameterList;
86 using Teuchos::parameterList;
89 using Teuchos::rcp_dynamic_cast;
90 using Teuchos::REDUCE_SUM;
91 using Teuchos::reduceAll;
92 typedef global_ordinal_type GO;
93 typedef Tpetra::Map<local_ordinal_type, GO, node_type> map_type;
99 RCP<const Comm<int> > comm = this->
getComm ();
100 const int myRank = comm->getRank ();
101 const int numProcs = comm->getSize ();
103 SLUD::int_t nprow, npcol;
110 RCP<const Comm<int> > matComm = this->
matrixA_->getComm ();
111 TEUCHOS_TEST_FOR_EXCEPTION(
112 matComm.is_null (), std::logic_error,
"Amesos2::Superlustdist "
113 "constructor: The matrix's communicator is null!");
114 RCP<const MpiComm<int> > matMpiComm =
115 rcp_dynamic_cast<
const MpiComm<int> > (matComm);
120 TEUCHOS_TEST_FOR_EXCEPTION(
121 matMpiComm.is_null (), std::logic_error,
"Amesos2::Superlustdist "
122 "constructor: The matrix's communicator is not an MpiComm!");
123 TEUCHOS_TEST_FOR_EXCEPTION(
124 matMpiComm->getRawMpiComm ().is_null (), std::logic_error,
"Amesos2::"
125 "Superlustdist constructor: The matrix's communicator claims to be a "
126 "Teuchos::MpiComm<int>, but its getRawPtrComm() method returns "
127 "Teuchos::null! This means that the underlying MPI_Comm doesn't even "
128 "exist, which likely implies that the Teuchos::MpiComm was constructed "
129 "incorrectly. It means something different than if the MPI_Comm were "
131 MPI_Comm rawMpiComm = (* (matMpiComm->getRawMpiComm ())) ();
132 data_.mat_comm = rawMpiComm;
138 SLUD::superlu_gridinit(data_.mat_comm, nprow, npcol, &(data_.grid));
147 SLUD::set_default_options_dist(&data_.options);
149 RCP<ParameterList> default_params =
154 data_.options.Fact = SLUD::DOFACT;
155 data_.equed = SLUD::NOEQUIL;
156 data_.options.SolveInitialized = SLUD::NO;
157 data_.options.RefineInitialized = SLUD::NO;
158 data_.rowequ =
false;
159 data_.colequ =
false;
167 data_.symb_comm = MPI_COMM_NULL;
172 data_.domains = (int) ( pow(2.0, floor(log10((
double)nprow*npcol)/log10(2.0))) );
174 const int color = (myRank < data_.domains) ? 0 : MPI_UNDEFINED;
175 MPI_Comm_split (data_.mat_comm, color, myRank, &(data_.symb_comm));
187 int myProcParticipates = 0;
188 if (myRank < nprow * npcol) {
190 myProcParticipates = 1;
195 int numParticipatingProcs = 0;
196 reduceAll<int, int> (*comm, REDUCE_SUM, myProcParticipates,
197 outArg (numParticipatingProcs));
198 TEUCHOS_TEST_FOR_EXCEPTION(
200 std::logic_error,
"Amesos2::Superludist constructor: The matrix has "
201 << this->
globalNumRows_ <<
" > 0 global row(s), but no processes in the "
202 "communicator want to participate in its factorization! nprow = "
203 << nprow <<
" and npcol = " << npcol <<
".");
206 size_t myNumRows = 0;
209 const GO quotient = (numParticipatingProcs == 0) ? static_cast<GO> (0) :
210 GNR /
static_cast<GO
> (numParticipatingProcs);
212 GNR - quotient *
static_cast<GO
> (numParticipatingProcs);
213 const GO lclNumRows = (
static_cast<GO
> (myRank) < remainder) ?
214 (quotient +
static_cast<GO
> (1)) : quotient;
215 myNumRows =
static_cast<size_t> (lclNumRows);
219 const GO indexBase = this->
matrixA_->getRowMap ()->getIndexBase ();
221 rcp (
new map_type (this->
globalNumRows_, myNumRows, indexBase, comm));
227 data_.A.Store = NULL;
229 SLUD::PStatInit(&(data_.stat));
232 data_.scale_perm.perm_r = data_.perm_r.getRawPtr();
233 data_.scale_perm.perm_c = data_.perm_c.getRawPtr();
237 template <
class Matrix,
class Vector>
249 if ( this->status_.getNumPreOrder() > 0 ){
251 free( data_.fstVtxSep );
256 if( data_.A.Store != NULL ){
257 SLUD::Destroy_SuperMatrix_Store_dist( &(data_.A) );
261 if ( this->status_.getNumNumericFact() > 0 ){
262 function_map::Destroy_LU(this->globalNumRows_, &(data_.grid), &(data_.lu));
264 function_map::LUstructFree(&(data_.lu));
270 if ( this->status_.symbolicFactorizationDone() &&
271 !this->status_.numericFactorizationDone() ){
272 if ( data_.pslu_freeable.xlsub != NULL ){
273 free( data_.pslu_freeable.xlsub );
274 free( data_.pslu_freeable.lsub );
276 if ( data_.pslu_freeable.xusub != NULL ){
277 free( data_.pslu_freeable.xusub );
278 free( data_.pslu_freeable.usub );
280 if ( data_.pslu_freeable.supno_loc != NULL ){
281 free( data_.pslu_freeable.supno_loc );
282 free( data_.pslu_freeable.xsup_beg_loc );
283 free( data_.pslu_freeable.xsup_end_loc );
285 free( data_.pslu_freeable.globToLoc );
288 SLUD::PStatFree( &(data_.stat) ) ;
293 if ( data_.options.SolveInitialized == SLUD::YES )
294 function_map::SolveFinalize(&(data_.options), &(data_.solve_struct));
296 SLUD::superlu_gridexit(&(data_.grid));
300 if ( data_.symb_comm != MPI_COMM_NULL ) MPI_Comm_free(&(data_.symb_comm));
303 template<
class Matrix,
class Vector>
312 SLUD::int_t slu_rows_ub = Teuchos::as<SLUD::int_t>(this->globalNumRows_);
313 for( SLUD::int_t i = 0; i < slu_rows_ub; ++i ) data_.perm_r[i] = i;
323 if( this->status_.getNumPreOrder() > 0 ){
325 free( data_.fstVtxSep );
327 #ifdef HAVE_AMESOS2_TIMERS
328 Teuchos::TimeMonitor preOrderTime( this->timers_.preOrderTime_ );
332 info = SLUD::get_perm_c_parmetis( &(data_.A),
333 data_.perm_r.getRawPtr(), data_.perm_c.getRawPtr(),
334 data_.grid.nprow * data_.grid.npcol, data_.domains,
335 &(data_.sizes), &(data_.fstVtxSep),
336 &(data_.grid), &(data_.symb_comm) );
338 TEUCHOS_TEST_FOR_EXCEPTION( info > 0.0,
340 "SuperLU_DIST pre-ordering ran out of memory after allocating "
341 << info <<
" bytes of memory" );
353 template <
class Matrix,
class Vector>
361 #ifdef HAVE_AMESOS2_TIMERS
362 Teuchos::TimeMonitor symFactTime( this->timers_.symFactTime_ );
366 info = SLUD::symbfact_dist((data_.grid.nprow) * (data_.grid.npcol),
367 data_.domains, &(data_.A), data_.perm_c.getRawPtr(),
368 data_.perm_r.getRawPtr(), data_.sizes,
369 data_.fstVtxSep, &(data_.pslu_freeable),
370 &(data_.grid.comm), &(data_.symb_comm),
373 TEUCHOS_TEST_FOR_EXCEPTION( info > 0.0,
375 "SuperLU_DIST symbolic factorization ran out of memory after"
376 " allocating " << info <<
" bytes of memory" );
378 same_symbolic_ =
false;
379 same_solve_struct_ =
false;
385 template <
class Matrix,
class Vector>
404 size_t nnz_loc = ((SLUD::NRformat_loc*)data_.A.Store)->nnz_loc;
405 for(
size_t i = 0; i < nnz_loc; ++i ) colind_[i] = data_.perm_c[colind_[i]];
408 if( same_symbolic_ ){
413 function_map::pdistribute(SLUD::SamePattern_SameRowPerm,
414 as<SLUD::int_t>(this->globalNumRows_),
415 &(data_.A), &(data_.scale_perm),
416 &(data_.glu_freeable), &(data_.lu),
419 function_map::dist_psymbtonum(SLUD::DOFACT,
420 as<SLUD::int_t>(this->globalNumRows_),
421 &(data_.A), &(data_.scale_perm),
422 &(data_.pslu_freeable), &(data_.lu),
427 double anorm = function_map::plangs((
char *)
"I", &(data_.A), &(data_.grid));
431 #ifdef HAVE_AMESOS2_TIMERS
432 Teuchos::TimeMonitor numFactTimer(this->timers_.numFactTime_);
435 function_map::gstrf(&(data_.options), this->globalNumRows_,
436 this->globalNumCols_, anorm, &(data_.lu),
437 &(data_.grid), &(data_.stat), &info);
441 TEUCHOS_TEST_FOR_EXCEPTION( info > 0,
443 "L and U factors have been computed but U("
444 << info <<
"," << info <<
") is exactly zero "
445 "(i.e. U is singular)");
452 data_.options.Fact = SLUD::FACTORED;
453 same_symbolic_ =
true;
459 template <
class Matrix,
class Vector>
468 const size_t local_len_rhs = superlu_rowmap_->getNodeNumElements();
469 const global_size_type nrhs = X->getGlobalNumVectors();
470 const global_ordinal_type first_global_row_b = superlu_rowmap_->getMinGlobalIndex();
473 bvals_.resize(nrhs * local_len_rhs);
474 xvals_.resize(nrhs * local_len_rhs);
480 #ifdef HAVE_AMESOS2_TIMERS
481 Teuchos::TimeMonitor convTimer(this->timers_.vecConvTime_);
491 #ifdef HAVE_AMESOS2_TIMERS
492 Teuchos::TimeMonitor redistTimer(this->timers_.vecRedistTime_);
498 copy_helper::do_get(B,
501 Teuchos::ptrInArg(*superlu_rowmap_));
522 if( !same_solve_struct_ ){
523 if( data_.options.SolveInitialized == SLUD::YES ){
524 function_map::SolveFinalize(&(data_.options), &(data_.solve_struct));
526 function_map::SolveInit(&(data_.options), &(data_.A), data_.perm_r.getRawPtr(),
527 data_.perm_c.getRawPtr(), as<SLUD::int_t>(nrhs), &(data_.lu),
528 &(data_.grid), &(data_.solve_struct));
532 same_solve_struct_ =
true;
537 #ifdef HAVE_AMESOS2_TIMERS
538 Teuchos::TimeMonitor solveTimer(this->timers_.solveTime_);
541 function_map::gstrs(as<SLUD::int_t>(this->globalNumRows_), &(data_.lu),
542 &(data_.scale_perm), &(data_.grid), bvals_.getRawPtr(),
543 as<SLUD::int_t>(local_len_rhs), as<SLUD::int_t>(first_global_row_b),
544 as<SLUD::int_t>(local_len_rhs), as<int>(nrhs),
545 &(data_.solve_struct), &(data_.stat), &ierr);
548 TEUCHOS_TEST_FOR_EXCEPTION( ierr < 0,
550 "Argument " << -ierr <<
" to gstrs had an illegal value" );
565 #ifdef HAVE_AMESOS2_TIMERS
566 Teuchos::TimeMonitor redistTimer(this->timers_.vecRedistTime_);
568 SLUD::int_t ld = as<SLUD::int_t>(local_len_rhs);
569 function_map::permute_Dense_Matrix(as<SLUD::int_t>(first_global_row_b),
570 as<SLUD::int_t>(local_len_rhs),
571 data_.solve_struct.row_to_proc,
572 data_.solve_struct.inv_perm_c,
573 bvals_.getRawPtr(), ld,
574 xvals_.getRawPtr(), ld,
582 #ifdef HAVE_AMESOS2_TIMERS
583 Teuchos::TimeMonitor redistTimer(this->timers_.vecRedistTime_);
587 put_helper::do_put(X,
590 Teuchos::ptrInArg(*superlu_rowmap_));
597 template <
class Matrix,
class Vector>
602 return( this->globalNumRows_ == this->globalNumCols_ );
606 template <
class Matrix,
class Vector>
612 using Teuchos::getIntegralValue;
613 using Teuchos::ParameterEntryValidator;
615 RCP<const Teuchos::ParameterList> valid_params = getValidParameters_impl();
617 if( parameterList->isParameter(
"npcol") || parameterList->isParameter(
"nprow") ){
618 TEUCHOS_TEST_FOR_EXCEPTION( !(parameterList->isParameter(
"nprow") &&
619 parameterList->isParameter(
"npcol")),
620 std::invalid_argument,
621 "nprow and npcol must be set together" );
623 SLUD::int_t nprow = parameterList->template get<SLUD::int_t>(
"nprow");
624 SLUD::int_t npcol = parameterList->template get<SLUD::int_t>(
"npcol");
626 TEUCHOS_TEST_FOR_EXCEPTION( nprow * npcol > this->getComm()->getSize(),
627 std::invalid_argument,
628 "nprow and npcol combination invalid" );
630 if( (npcol != data_.grid.npcol) || (nprow != data_.grid.nprow) ){
632 SLUD::superlu_gridexit(&(data_.grid));
634 SLUD::superlu_gridinit(data_.mat_comm, nprow, npcol, &(data_.grid));
638 TEUCHOS_TEST_FOR_EXCEPTION( this->control_.useTranspose_,
639 std::invalid_argument,
640 "SuperLU_DIST does not support solving the tranpose system" );
642 data_.options.Trans = SLUD::NOTRANS;
647 data_.options.Equil = SLUD::NO;
649 if( parameterList->isParameter(
"ColPerm") ){
650 RCP<const ParameterEntryValidator> colperm_validator = valid_params->getEntry(
"ColPerm").validator();
651 parameterList->getEntry(
"ColPerm").setValidator(colperm_validator);
653 data_.options.ColPerm = getIntegralValue<SLUD::colperm_t>(*parameterList,
"ColPerm");
660 data_.options.RowPerm = SLUD::NOROWPERM;
669 data_.options.IterRefine = SLUD::NOREFINE;
671 bool replace_tiny = parameterList->get<
bool>(
"ReplaceTinyPivot",
true);
672 data_.options.ReplaceTinyPivot = replace_tiny ? SLUD::YES : SLUD::NO;
674 if( parameterList->isParameter(
"IsContiguous") ){
675 is_contiguous_ = parameterList->get<
bool>(
"IsContiguous");
680 template <
class Matrix,
class Vector>
681 Teuchos::RCP<const Teuchos::ParameterList>
685 using Teuchos::tuple;
686 using Teuchos::ParameterList;
687 using Teuchos::EnhancedNumberValidator;
688 using Teuchos::setStringToIntegralParameter;
689 using Teuchos::stringToIntegralParameterEntryValidator;
691 static Teuchos::RCP<const Teuchos::ParameterList> valid_params;
693 if( is_null(valid_params) ){
694 Teuchos::RCP<Teuchos::ParameterList> pl = Teuchos::parameterList();
696 Teuchos::RCP<EnhancedNumberValidator<SLUD::int_t> > col_row_validator
697 = Teuchos::rcp(
new EnhancedNumberValidator<SLUD::int_t>() );
698 col_row_validator->setMin(1);
700 pl->set(
"npcol", data_.grid.npcol,
701 "Number of columns in the processor grid. "
702 "Must be set with nprow", col_row_validator);
703 pl->set(
"nprow", data_.grid.nprow,
704 "Number of rows in the SuperLU_DIST processor grid. "
705 "Must be set together with npcol", col_row_validator);
708 setStringToIntegralParameter<SLUD::trans_t>(
"Trans",
"NOTRANS",
709 "Solve for the transpose system or not",
710 tuple<string>(
"NOTRANS"),
711 tuple<string>(
"Do not solve with transpose"),
712 tuple<SLUD::trans_t>(SLUD::NOTRANS),
728 pl->set(
"ReplaceTinyPivot",
true,
729 "Specifies whether to replace tiny diagonals during LU factorization");
731 setStringToIntegralParameter<SLUD::colperm_t>(
"ColPerm",
"PARMETIS",
732 "Specifies how to permute the columns of the "
733 "matrix for sparsity preservation",
734 tuple<string>(
"NATURAL",
"PARMETIS"),
735 tuple<string>(
"Natural ordering",
736 "ParMETIS ordering on A^T + A"),
737 tuple<SLUD::colperm_t>(SLUD::NATURAL,
741 pl->set(
"IsContiguous",
true,
"Whether GIDs contiguous");
750 template <
class Matrix,
class Vector>
754 SLUD::int_t& npcol)
const {
755 TEUCHOS_TEST_FOR_EXCEPTION( nprocs < 1,
756 std::invalid_argument,
757 "Number of MPI processes must be at least 1" );
758 SLUD::int_t c, r = 1;
759 while( r*r <= nprocs ) r++;
762 while( (r--)*c != nprocs ){
768 if( r > 1 || nprocs < 9){
775 template <
class Matrix,
class Vector>
779 using Teuchos::Array;
780 using Teuchos::ArrayView;
781 using Teuchos::ptrInArg;
786 #ifdef HAVE_AMESOS2_TIMERS
787 Teuchos::TimeMonitor convTimer(this->timers_.mtxConvTime_);
791 if( data_.A.Store != NULL ){
792 SLUD::Destroy_SuperMatrix_Store_dist( &(data_.A) );
793 data_.A.Store = NULL;
796 Teuchos::RCP<const MatrixAdapter<Matrix> > redist_mat
797 = this->matrixA_->get(ptrInArg(*superlu_rowmap_));
799 int_t l_nnz, l_rows, g_rows, g_cols, fst_global_row;
800 l_nnz = as<int_t>(redist_mat->getLocalNNZ());
801 l_rows = as<int_t>(redist_mat->getLocalNumRows());
802 g_rows = as<int_t>(redist_mat->getGlobalNumRows());
804 fst_global_row = as<int_t>(superlu_rowmap_->getMinGlobalIndex());
806 nzvals_.resize(l_nnz);
807 colind_.resize(l_nnz);
808 rowptr_.resize(l_rows + 1);
812 #ifdef HAVE_AMESOS2_TIMERS
813 Teuchos::TimeMonitor mtxRedistTimer( this->timers_.mtxRedistTime_ );
818 slu_type, int_t, int_t >::do_get(redist_mat.ptr(),
819 nzvals_(), colind_(), rowptr_(),
821 ptrInArg(*superlu_rowmap_),
826 TEUCHOS_TEST_FOR_EXCEPTION( nnz_ret != l_nnz,
828 "Did not get the expected number of non-zero vals");
831 SLUD::Dtype_t dtype = type_map::dtype;
834 function_map::create_CompRowLoc_Matrix(&(data_.A),
836 l_nnz, l_rows, fst_global_row,
841 dtype, SLUD::SLU_GE);
848 template<
class Matrix,
class Vector>
854 #endif // AMESOS2_SUPERLUDIST_DEF_HPP
Amesos2::SolverCore: A templated interface for interaction with third-party direct sparse solvers...
Definition: Amesos2_SolverCore_decl.hpp:105
Amesos2 interface to the distributed memory version of SuperLU.
Definition: Amesos2_Superludist_decl.hpp:90
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:665
global_size_type globalNumCols_
Number of global columns in matrixA_.
Definition: Amesos2_SolverCore_decl.hpp:478
Superludist(Teuchos::RCP< const Matrix > A, Teuchos::RCP< Vector > X, Teuchos::RCP< const Vector > B)
Initialize from Teuchos::RCP.
Definition: Amesos2_Superludist_def.hpp:68
void setParameters_impl(const Teuchos::RCP< Teuchos::ParameterList > ¶meterList)
Definition: Amesos2_Superludist_def.hpp:608
global_size_type globalNumRows_
Number of global rows in matrixA_.
Definition: Amesos2_SolverCore_decl.hpp:475
Helper class for getting 1-D copies of multivectors.
Definition: Amesos2_MultiVecAdapter_decl.hpp:266
Utility functions for Amesos2.
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters_impl() const
Definition: Amesos2_Superludist_def.hpp:682
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
Return a const parameter list of all of the valid parameters that this->setParameterList(...) will accept.
Definition: Amesos2_SolverCore_def.hpp:314
Provides definition of SuperLU_DIST types as well as conversions and type traits. ...
int preOrdering_impl()
Performs pre-ordering on the matrix to increase efficiency.
Definition: Amesos2_Superludist_def.hpp:305
Teuchos::RCP< const Teuchos::Comm< int > > getComm() const
Returns a pointer to the Teuchos::Comm communicator with this operator.
Definition: Amesos2_SolverCore_decl.hpp:362
A Matrix adapter interface for Amesos2.
Definition: Amesos2_MatrixAdapter_decl.hpp:76
bool matrixShapeOK_impl() const
Determines whether the shape of the matrix is OK for this solver.
Definition: Amesos2_Superludist_def.hpp:599
Teuchos::RCP< const Tpetra::Map< local_ordinal_type, global_ordinal_type, node_type > > superlu_rowmap_
Maps rows of the matrix to processors in the SuperLU_DIST processor grid.
Definition: Amesos2_Superludist_decl.hpp:327
~Superludist()
Destructor.
Definition: Amesos2_Superludist_def.hpp:238
super_type & setParameters(const Teuchos::RCP< Teuchos::ParameterList > ¶meterList)
Set/update internal variables and solver options.
Definition: Amesos2_SolverCore_def.hpp:282
int symbolicFactorization_impl()
Perform symbolic factorization of the matrix using SuperLU_DIST.
Definition: Amesos2_Superludist_def.hpp:355
void get_default_grid_size(int nprocs, SLUD::int_t &nprow, SLUD::int_t &npcol) const
Definition: Amesos2_Superludist_def.hpp:752
bool in_grid_
true if this processor is in SuperLU_DISTS's 2D process grid
Definition: Amesos2_Superludist_decl.hpp:320
Helper class for putting 1-D data arrays into multivectors.
Definition: Amesos2_MultiVecAdapter_decl.hpp:322
int solve_impl(const Teuchos::Ptr< MultiVecAdapter< Vector > > X, const Teuchos::Ptr< const MultiVecAdapter< Vector > > B) const
SuperLU_DIST specific solve.
Definition: Amesos2_Superludist_def.hpp:461
A templated MultiVector class adapter for Amesos2.
Definition: Amesos2_MultiVecAdapter_decl.hpp:176
bool loadA_impl(EPhase current_phase)
Reads matrix data into internal solver structures.
Definition: Amesos2_Superludist_def.hpp:777
Teuchos::RCP< const MatrixAdapter< Matrix > > matrixA_
The LHS operator.
Definition: Amesos2_SolverCore_decl.hpp:454
int numericFactorization_impl()
SuperLU_DIST specific numeric factorization.
Definition: Amesos2_Superludist_def.hpp:387