Amesos2 - Direct Sparse Solver Interfaces  Version of the Day
Amesos2_Superludist_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 
52 #ifndef AMESOS2_SUPERLUDIST_DEF_HPP
53 #define AMESOS2_SUPERLUDIST_DEF_HPP
54 
55 #include <Teuchos_Tuple.hpp>
56 #include <Teuchos_StandardParameterEntryValidators.hpp>
57 #include <Teuchos_DefaultMpiComm.hpp>
58 
61 #include "Amesos2_Util.hpp"
62 
63 
64 namespace Amesos2 {
65 
66 
67  template <class Matrix, class Vector>
68  Superludist<Matrix,Vector>::Superludist(Teuchos::RCP<const Matrix> A,
69  Teuchos::RCP<Vector> X,
70  Teuchos::RCP<const Vector> B)
71  : SolverCore<Amesos2::Superludist,Matrix,Vector>(A, X, B)
72  , nzvals_() // initialization to empty arrays
73  , colind_()
74  , rowptr_()
75  , bvals_()
76  , xvals_()
77  , in_grid_(false)
78  , is_contiguous_(true)
79  {
80  using Teuchos::Comm;
81  // It's OK to depend on MpiComm explicitly here, because
82  // SuperLU_DIST requires MPI anyway.
83  using Teuchos::MpiComm;
84  using Teuchos::outArg;
85  using Teuchos::ParameterList;
86  using Teuchos::parameterList;
87  using Teuchos::RCP;
88  using Teuchos::rcp;
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;
94 
96  // Set up the SuperLU_DIST processor grid //
98 
99  RCP<const Comm<int> > comm = this->getComm ();
100  const int myRank = comm->getRank ();
101  const int numProcs = comm->getSize ();
102 
103  SLUD::int_t nprow, npcol;
104  get_default_grid_size (numProcs, nprow, npcol);
105 
106  {
107  // FIXME (mfh 16 Dec 2014) getComm() just returns
108  // matrixA_->getComm(), so it's not clear why we need to ask for
109  // the matrix's communicator separately here.
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);
116  // FIXME (mfh 16 Dec 2014) If the matrix's communicator is a
117  // SerialComm, we probably could just use MPI_COMM_SELF here.
118  // I'm not sure if SuperLU_DIST is smart enough to handle that
119  // case, though.
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 "
130  "MPI_COMM_NULL.");
131  MPI_Comm rawMpiComm = (* (matMpiComm->getRawMpiComm ())) ();
132  data_.mat_comm = rawMpiComm;
133  // This looks a bit like ScaLAPACK's grid initialization (which
134  // technically takes place in the BLACS, not in ScaLAPACK
135  // proper). See http://netlib.org/scalapack/slug/node34.html.
136  // The main difference is that SuperLU_DIST depends explicitly
137  // on MPI, while the BLACS hides its communication protocol.
138  SLUD::superlu_gridinit(data_.mat_comm, nprow, npcol, &(data_.grid));
139  }
140 
142  // Set some default parameters. //
143  // //
144  // Must do this after grid has been created in //
145  // case user specifies the nprow and npcol parameters //
147  SLUD::set_default_options_dist(&data_.options);
148 
149  RCP<ParameterList> default_params =
150  parameterList (* (this->getValidParameters ()));
151  this->setParameters (default_params);
152 
153  // Set some internal options
154  data_.options.Fact = SLUD::DOFACT;
155  data_.equed = SLUD::NOEQUIL; // No equilibration has yet been performed
156  data_.options.SolveInitialized = SLUD::NO;
157  data_.options.RefineInitialized = SLUD::NO;
158  data_.rowequ = false;
159  data_.colequ = false;
160  data_.perm_r.resize(this->globalNumRows_);
161  data_.perm_c.resize(this->globalNumCols_);
162 
164  // Set up a communicator for the parallel column ordering and //
165  // parallel symbolic factorization. //
167  data_.symb_comm = MPI_COMM_NULL;
168 
169  // domains is the next power of 2 less than nprow*npcol. This
170  // value will be used for creating an MPI communicator for the
171  // pre-ordering and symbolic factorization methods.
172  data_.domains = (int) ( pow(2.0, floor(log10((double)nprow*npcol)/log10(2.0))) );
173 
174  const int color = (myRank < data_.domains) ? 0 : MPI_UNDEFINED;
175  MPI_Comm_split (data_.mat_comm, color, myRank, &(data_.symb_comm));
176 
178  // Set up a row Map that only includes processes that are in the
179  // SuperLU process grid. This will be used for redistributing A.
181 
182  // mfh 16 Dec 2014: We could use createWeightedContigMapWithNode
183  // with myProcParticipates as the weight, but that costs an extra
184  // all-reduce.
185 
186  // Set to 1 if I am in the grid, and I get some of the matrix rows.
187  int myProcParticipates = 0;
188  if (myRank < nprow * npcol) {
189  in_grid_ = true;
190  myProcParticipates = 1;
191  }
192 
193  // Compute how many processes in the communicator belong to the
194  // process grid.
195  int numParticipatingProcs = 0;
196  reduceAll<int, int> (*comm, REDUCE_SUM, myProcParticipates,
197  outArg (numParticipatingProcs));
198  TEUCHOS_TEST_FOR_EXCEPTION(
199  this->globalNumRows_ != 0 && numParticipatingProcs == 0,
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 << ".");
204 
205  // Divide up the rows among the participating processes.
206  size_t myNumRows = 0;
207  {
208  const GO GNR = static_cast<GO> (this->globalNumRows_);
209  const GO quotient = (numParticipatingProcs == 0) ? static_cast<GO> (0) :
210  GNR / static_cast<GO> (numParticipatingProcs);
211  const GO remainder =
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);
216  }
217 
218  // TODO: might only need to initialize if parallel symbolic factorization is requested.
219  const GO indexBase = this->matrixA_->getRowMap ()->getIndexBase ();
221  rcp (new map_type (this->globalNumRows_, myNumRows, indexBase, comm));
222 
224  // Do some other initialization //
226 
227  data_.A.Store = NULL;
228  function_map::LUstructInit(this->globalNumRows_, this->globalNumCols_, &(data_.lu));
229  SLUD::PStatInit(&(data_.stat));
230  // We do not use ScalePermstructInit because we will use our own
231  // arrays for storing perm_r and perm_c
232  data_.scale_perm.perm_r = data_.perm_r.getRawPtr();
233  data_.scale_perm.perm_c = data_.perm_c.getRawPtr();
234  }
235 
236 
237  template <class Matrix, class Vector>
239  {
240  /* Free SuperLU_DIST data_types
241  * - Matrices
242  * - Vectors
243  * - Stat object
244  * - ScalePerm, LUstruct, grid, and solve objects
245  *
246  * Note: the function definitions are the same regardless whether
247  * complex or real, so we arbitrarily use the D namespace
248  */
249  if ( this->status_.getNumPreOrder() > 0 ){
250  free( data_.sizes );
251  free( data_.fstVtxSep );
252  }
253 
254  // Cleanup old matrix store memory if it's non-NULL. Our
255  // Teuchos::Array's will destroy rowind, colptr, and nzval for us
256  if( data_.A.Store != NULL ){
257  SLUD::Destroy_SuperMatrix_Store_dist( &(data_.A) );
258  }
259 
260  // LU data is initialized in numericFactorization_impl()
261  if ( this->status_.getNumNumericFact() > 0 ){
262  function_map::Destroy_LU(this->globalNumRows_, &(data_.grid), &(data_.lu));
263  }
264  function_map::LUstructFree(&(data_.lu));
265 
266  // If a symbolic factorization is ever performed without a
267  // follow-up numericfactorization, there are some arrays in the
268  // Pslu_freeable struct which will never be free'd by
269  // SuperLU_DIST.
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 );
275  }
276  if ( data_.pslu_freeable.xusub != NULL ){
277  free( data_.pslu_freeable.xusub );
278  free( data_.pslu_freeable.usub );
279  }
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 );
284  }
285  free( data_.pslu_freeable.globToLoc );
286  }
287 
288  SLUD::PStatFree( &(data_.stat) ) ;
289 
290  // Teuchos::Arrays will free R, C, perm_r, and perm_c
291  // SLUD::D::ScalePermstructFree(&(data_.scale_perm));
292 
293  if ( data_.options.SolveInitialized == SLUD::YES )
294  function_map::SolveFinalize(&(data_.options), &(data_.solve_struct));
295 
296  SLUD::superlu_gridexit(&(data_.grid)); // TODO: are there any
297  // cases where grid
298  // wouldn't be initialized?
299 
300  if ( data_.symb_comm != MPI_COMM_NULL ) MPI_Comm_free(&(data_.symb_comm));
301  }
302 
303  template<class Matrix, class Vector>
304  int
306  {
307  // We will always use the NATURAL row ordering to avoid the
308  // sequential bottleneck present when doing any other row
309  // ordering scheme from SuperLU_DIST
310  //
311  // Set perm_r to be the natural ordering
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;
314 
315  // loadA_impl(); // Refresh matrix values
316 
317  if( in_grid_ ){
318  // If this function has been called at least once, then the
319  // sizes, and fstVtxSep arrays were allocated in
320  // get_perm_c_parmetis. Delete them before calling that
321  // function again. These arrays will also be dealloc'd in the
322  // deconstructor.
323  if( this->status_.getNumPreOrder() > 0 ){
324  free( data_.sizes );
325  free( data_.fstVtxSep );
326  }
327 #ifdef HAVE_AMESOS2_TIMERS
328  Teuchos::TimeMonitor preOrderTime( this->timers_.preOrderTime_ );
329 #endif
330 
331  float info = 0.0;
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) );
337 
338  TEUCHOS_TEST_FOR_EXCEPTION( info > 0.0,
339  std::runtime_error,
340  "SuperLU_DIST pre-ordering ran out of memory after allocating "
341  << info << " bytes of memory" );
342  }
343 
344  // Ordering will be applied directly before numeric factorization,
345  // after we have a chance to get updated coefficients from the
346  // matrix
347 
348  return EXIT_SUCCESS;
349  }
350 
351 
352 
353  template <class Matrix, class Vector>
354  int
356  {
357  // loadA_impl(); // Refresh matrix values
358 
359  if( in_grid_ ){
360 
361 #ifdef HAVE_AMESOS2_TIMERS
362  Teuchos::TimeMonitor symFactTime( this->timers_.symFactTime_ );
363 #endif
364 
365  float info = 0.0;
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),
371  &(data_.mem_usage));
372 
373  TEUCHOS_TEST_FOR_EXCEPTION( info > 0.0,
374  std::runtime_error,
375  "SuperLU_DIST symbolic factorization ran out of memory after"
376  " allocating " << info << " bytes of memory" );
377  }
378  same_symbolic_ = false;
379  same_solve_struct_ = false;
380 
381  return EXIT_SUCCESS;
382  }
383 
384 
385  template <class Matrix, class Vector>
386  int
388  using Teuchos::as;
389 
390  // loadA_impl(); // Refresh the matrix values
391 
392  // if( data_.options.Equil == SLUD::YES ){
393  // // Apply the scalings computed in preOrdering
394  // function_map::laqgs(&(data_.A), data_.R.getRawPtr(),
395  // data_.C.getRawPtr(), data_.rowcnd, data_.colcnd,
396  // data_.amax, &(data_.equed));
397 
398  // data_.rowequ = (data_.equed == SLUD::ROW) || (data_.equed == SLUD::BOTH);
399  // data_.colequ = (data_.equed == SLUD::COL) || (data_.equed == SLUD::BOTH);
400  // }
401 
402  if( in_grid_ ){
403  // Apply the column ordering, so that AC is the column-permuted A, and compute etree
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]];
406 
407  // Distribute data from the symbolic factorization
408  if( same_symbolic_ ){
409  // Note: with the SamePattern_SameRowPerm options, it does not
410  // matter that the glu_freeable member has never been
411  // initialized, because it is never accessed. It is a
412  // placeholder arg. The real work is done in data_.lu
413  function_map::pdistribute(SLUD::SamePattern_SameRowPerm,
414  as<SLUD::int_t>(this->globalNumRows_), // aka "n"
415  &(data_.A), &(data_.scale_perm),
416  &(data_.glu_freeable), &(data_.lu),
417  &(data_.grid));
418  } else {
419  function_map::dist_psymbtonum(SLUD::DOFACT,
420  as<SLUD::int_t>(this->globalNumRows_), // aka "n"
421  &(data_.A), &(data_.scale_perm),
422  &(data_.pslu_freeable), &(data_.lu),
423  &(data_.grid));
424  }
425 
426  // Retrieve the normI of A (required by gstrf).
427  double anorm = function_map::plangs((char *)"I", &(data_.A), &(data_.grid));
428 
429  int info = 0;
430  {
431 #ifdef HAVE_AMESOS2_TIMERS
432  Teuchos::TimeMonitor numFactTimer(this->timers_.numFactTime_);
433 #endif
434 
435  function_map::gstrf(&(data_.options), this->globalNumRows_,
436  this->globalNumCols_, anorm, &(data_.lu),
437  &(data_.grid), &(data_.stat), &info);
438  }
439 
440  // Check output
441  TEUCHOS_TEST_FOR_EXCEPTION( info > 0,
442  std::runtime_error,
443  "L and U factors have been computed but U("
444  << info << "," << info << ") is exactly zero "
445  "(i.e. U is singular)");
446  }
447 
448  // The other option, that info_st < 0, denotes invalid parameters
449  // to the function, but we'll assume for now that that won't
450  // happen.
451 
452  data_.options.Fact = SLUD::FACTORED;
453  same_symbolic_ = true;
454 
455  return EXIT_SUCCESS;
456  }
457 
458 
459  template <class Matrix, class Vector>
460  int
462  const Teuchos::Ptr<const MultiVecAdapter<Vector> > B) const
463  {
464  using Teuchos::as;
465 
466  // local_len_rhs is how many of the multivector rows belong to
467  // this processor in the SuperLU_DIST processor grid.
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();
471 
472  // make sure our multivector storage is sized appropriately
473  bvals_.resize(nrhs * local_len_rhs);
474  xvals_.resize(nrhs * local_len_rhs);
475 
476  // We assume the global length of the two vectors have already been
477  // checked for compatibility
478 
479  { // get the values from B
480 #ifdef HAVE_AMESOS2_TIMERS
481  Teuchos::TimeMonitor convTimer(this->timers_.vecConvTime_);
482 #endif
483 
484  {
485  // The input dense matrix for B should be distributed in the
486  // same manner as the superlu_dist matrix. That is, if a
487  // processor has m_loc rows of A, then it should also have
488  // m_loc rows of B (and the same rows). We accomplish this by
489  // distributing the multivector rows with the same Map that
490  // the matrix A's rows are distributed.
491 #ifdef HAVE_AMESOS2_TIMERS
492  Teuchos::TimeMonitor redistTimer(this->timers_.vecRedistTime_);
493 #endif
494 
495  // get grid-distributed mv data. The multivector data will be
496  // distributed across the processes in the SuperLU_DIST grid.
497  typedef Util::get_1d_copy_helper<MultiVecAdapter<Vector>,slu_type> copy_helper;
498  copy_helper::do_get(B,
499  bvals_(),
500  local_len_rhs,
501  Teuchos::ptrInArg(*superlu_rowmap_));
502  }
503  } // end block for conversion time
504 
505  if( in_grid_ ){
506  // if( data_.options.trans == SLUD::NOTRANS ){
507  // if( data_.rowequ ){ // row equilibration has been done on AC
508  // // scale bxvals_ by diag(R)
509  // Util::scale(bxvals_(), as<size_t>(len_rhs), ldbx_, data_.R(),
510  // SLUD::slu_mt_mult<slu_type,magnitude_type>());
511  // }
512  // } else if( data_.colequ ){ // column equilibration has been done on AC
513  // // scale bxvals_ by diag(C)
514  // Util::scale(bxvals_(), as<size_t>(len_rhs), ldbx_, data_.C(),
515  // SLUD::slu_mt_mult<slu_type,magnitude_type>());
516  // }
517 
518  // Initialize the SOLVEstruct_t.
519  //
520  // We are able to reuse the solve struct if we have not changed
521  // the sparsity pattern of L and U since the last solve
522  if( !same_solve_struct_ ){
523  if( data_.options.SolveInitialized == SLUD::YES ){
524  function_map::SolveFinalize(&(data_.options), &(data_.solve_struct));
525  }
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));
529  // Flag that we can reuse this solve_struct unless another
530  // symbolicFactorization is called between here and the next
531  // solve.
532  same_solve_struct_ = true;
533  }
534 
535  int ierr = 0; // returned error code
536  {
537 #ifdef HAVE_AMESOS2_TIMERS
538  Teuchos::TimeMonitor solveTimer(this->timers_.solveTime_);
539 #endif
540 
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);
546  } // end block for solve time
547 
548  TEUCHOS_TEST_FOR_EXCEPTION( ierr < 0,
549  std::runtime_error,
550  "Argument " << -ierr << " to gstrs had an illegal value" );
551 
552  // "Un-scale" the solution so that it is a solution of the original system
553  // if( data_.options.trans == SLUD::NOTRANS ){
554  // if( data_.colequ ){ // column equilibration has been done on AC
555  // // scale bxvals_ by diag(C)
556  // Util::scale(bxvals_(), as<size_t>(len_rhs), ldbx_, data_.C(),
557  // SLUD::slu_mt_mult<slu_type,magnitude_type>());
558  // }
559  // } else if( data_.rowequ ){ // row equilibration has been done on AC
560  // // scale bxvals_ by diag(R)
561  // Util::scale(bxvals_(), as<size_t>(len_rhs), ldbx_, data_.R(),
562  // SLUD::slu_mt_mult<slu_type,magnitude_type>());
563  // }
564  { // permute B to a solution of the original system
565 #ifdef HAVE_AMESOS2_TIMERS
566  Teuchos::TimeMonitor redistTimer(this->timers_.vecRedistTime_);
567 #endif
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,
575  as<int>(nrhs),
576  &(data_.grid));
577  }
578  }
579 
580  /* Update X's global values */
581  {
582 #ifdef HAVE_AMESOS2_TIMERS
583  Teuchos::TimeMonitor redistTimer(this->timers_.vecRedistTime_);
584 #endif
585 
586  typedef Util::put_1d_data_helper<MultiVecAdapter<Vector>,slu_type> put_helper;
587  put_helper::do_put(X,
588  xvals_(),
589  local_len_rhs,
590  Teuchos::ptrInArg(*superlu_rowmap_));
591  }
592 
593  return EXIT_SUCCESS;
594  }
595 
596 
597  template <class Matrix, class Vector>
598  bool
600  {
601  // SuperLU_DIST requires square matrices
602  return( this->globalNumRows_ == this->globalNumCols_ );
603  }
604 
605 
606  template <class Matrix, class Vector>
607  void
608  Superludist<Matrix,Vector>::setParameters_impl(const Teuchos::RCP<Teuchos::ParameterList> & parameterList )
609  {
610  using Teuchos::as;
611  using Teuchos::RCP;
612  using Teuchos::getIntegralValue;
613  using Teuchos::ParameterEntryValidator;
614 
615  RCP<const Teuchos::ParameterList> valid_params = getValidParameters_impl();
616 
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" );
622 
623  SLUD::int_t nprow = parameterList->template get<SLUD::int_t>("nprow");
624  SLUD::int_t npcol = parameterList->template get<SLUD::int_t>("npcol");
625 
626  TEUCHOS_TEST_FOR_EXCEPTION( nprow * npcol > this->getComm()->getSize(),
627  std::invalid_argument,
628  "nprow and npcol combination invalid" );
629 
630  if( (npcol != data_.grid.npcol) || (nprow != data_.grid.nprow) ){
631  // De-allocate the default grid that was initialized in the constructor
632  SLUD::superlu_gridexit(&(data_.grid));
633  // Create a new grid
634  SLUD::superlu_gridinit(data_.mat_comm, nprow, npcol, &(data_.grid));
635  } // else our grid has not changed size since the last initialization
636  }
637 
638  TEUCHOS_TEST_FOR_EXCEPTION( this->control_.useTranspose_,
639  std::invalid_argument,
640  "SuperLU_DIST does not support solving the tranpose system" );
641 
642  data_.options.Trans = SLUD::NOTRANS; // should always be set this way;
643 
644  // TODO: Uncomment when supported
645  // bool equil = parameterList->get<bool>("Equil", true);
646  // data_.options.Equil = equil ? SLUD::YES : SLUD::NO;
647  data_.options.Equil = SLUD::NO;
648 
649  if( parameterList->isParameter("ColPerm") ){
650  RCP<const ParameterEntryValidator> colperm_validator = valid_params->getEntry("ColPerm").validator();
651  parameterList->getEntry("ColPerm").setValidator(colperm_validator);
652 
653  data_.options.ColPerm = getIntegralValue<SLUD::colperm_t>(*parameterList, "ColPerm");
654  }
655 
656  // Always use the "NOROWPERM" option to avoid a serial bottleneck
657  // with the weighted bipartite matching algorithm used for the
658  // "LargeDiag" RowPerm. Note the inconsistency with the SuperLU
659  // User guide (which states that the value should be "NATURAL").
660  data_.options.RowPerm = SLUD::NOROWPERM;
661 
662  // TODO: Uncomment when supported
663  // if( parameterList->isParameter("IterRefine") ){
664  // RCP<const ParameterEntryValidator> iter_refine_validator = valid_params->getEntry("IterRefine").validator();
665  // parameterList->getEntry("IterRefine").setValidator(iter_refine_validator);
666 
667  // data_.options.IterRefine = getIntegralValue<SLUD::IterRefine_t>(*parameterList, "IterRefine");
668  // }
669  data_.options.IterRefine = SLUD::NOREFINE;
670 
671  bool replace_tiny = parameterList->get<bool>("ReplaceTinyPivot", true);
672  data_.options.ReplaceTinyPivot = replace_tiny ? SLUD::YES : SLUD::NO;
673 
674  if( parameterList->isParameter("IsContiguous") ){
675  is_contiguous_ = parameterList->get<bool>("IsContiguous");
676  }
677  }
678 
679 
680  template <class Matrix, class Vector>
681  Teuchos::RCP<const Teuchos::ParameterList>
683  {
684  using std::string;
685  using Teuchos::tuple;
686  using Teuchos::ParameterList;
687  using Teuchos::EnhancedNumberValidator;
688  using Teuchos::setStringToIntegralParameter;
689  using Teuchos::stringToIntegralParameterEntryValidator;
690 
691  static Teuchos::RCP<const Teuchos::ParameterList> valid_params;
692 
693  if( is_null(valid_params) ){
694  Teuchos::RCP<Teuchos::ParameterList> pl = Teuchos::parameterList();
695 
696  Teuchos::RCP<EnhancedNumberValidator<SLUD::int_t> > col_row_validator
697  = Teuchos::rcp( new EnhancedNumberValidator<SLUD::int_t>() );
698  col_row_validator->setMin(1);
699 
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);
706 
707  // validator will catch any value besides NOTRANS
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),
713  pl.getRawPtr());
714 
715  // TODO: uncomment when supported
716  // pl->set("Equil", false, "Whether to equilibrate the system before solve");
717 
718  // TODO: uncomment when supported
719  // setStringToIntegralParameter<SLUD::IterRefine_t>("IterRefine", "NOREFINE",
720  // "Type of iterative refinement to use",
721  // tuple<string>("NOREFINE", "DOUBLE"),
722  // tuple<string>("Do not use iterative refinement",
723  // "Do double iterative refinement"),
724  // tuple<SLUD::IterRefine_t>(SLUD::NOREFINE,
725  // SLUD::DOUBLE),
726  // pl.getRawPtr());
727 
728  pl->set("ReplaceTinyPivot", true,
729  "Specifies whether to replace tiny diagonals during LU factorization");
730 
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,
738  SLUD::PARMETIS),
739  pl.getRawPtr());
740 
741  pl->set("IsContiguous", true, "Whether GIDs contiguous");
742 
743  valid_params = pl;
744  }
745 
746  return valid_params;
747  }
748 
749 
750  template <class Matrix, class Vector>
751  void
753  SLUD::int_t& nprow,
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++;
760  nprow = npcol = --r; // fall back to square grid
761  c = nprocs / r;
762  while( (r--)*c != nprocs ){
763  c = nprocs / r; // note integer division
764  }
765  ++r;
766  // prefer the square grid over a single row (which will only happen
767  // in the case of a prime nprocs
768  if( r > 1 || nprocs < 9){ // nprocs < 9 is a heuristic for the small cases
769  nprow = r;
770  npcol = c;
771  }
772  }
773 
774 
775  template <class Matrix, class Vector>
776  bool
778  // Extract the necessary information from mat and call SLU function
779  using Teuchos::Array;
780  using Teuchos::ArrayView;
781  using Teuchos::ptrInArg;
782  using Teuchos::as;
783 
784  using SLUD::int_t;
785 
786 #ifdef HAVE_AMESOS2_TIMERS
787  Teuchos::TimeMonitor convTimer(this->timers_.mtxConvTime_);
788 #endif
789 
790  // Cleanup old store memory if it's non-NULL
791  if( data_.A.Store != NULL ){
792  SLUD::Destroy_SuperMatrix_Store_dist( &(data_.A) );
793  data_.A.Store = NULL;
794  }
795 
796  Teuchos::RCP<const MatrixAdapter<Matrix> > redist_mat
797  = this->matrixA_->get(ptrInArg(*superlu_rowmap_));
798 
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());
803  g_cols = g_rows; // we deal with square matrices
804  fst_global_row = as<int_t>(superlu_rowmap_->getMinGlobalIndex());
805 
806  nzvals_.resize(l_nnz);
807  colind_.resize(l_nnz);
808  rowptr_.resize(l_rows + 1);
809 
810  int_t nnz_ret = 0;
811  {
812 #ifdef HAVE_AMESOS2_TIMERS
813  Teuchos::TimeMonitor mtxRedistTimer( this->timers_.mtxRedistTime_ );
814 #endif
815 
818  slu_type, int_t, int_t >::do_get(redist_mat.ptr(),
819  nzvals_(), colind_(), rowptr_(),
820  nnz_ret,
821  ptrInArg(*superlu_rowmap_),
822  ROOTED,
823  ARBITRARY);
824  }
825 
826  TEUCHOS_TEST_FOR_EXCEPTION( nnz_ret != l_nnz,
827  std::runtime_error,
828  "Did not get the expected number of non-zero vals");
829 
830  // Get the SLU data type for this type of matrix
831  SLUD::Dtype_t dtype = type_map::dtype;
832 
833  if( in_grid_ ){
834  function_map::create_CompRowLoc_Matrix(&(data_.A),
835  g_rows, g_cols,
836  l_nnz, l_rows, fst_global_row,
837  nzvals_.getRawPtr(),
838  colind_.getRawPtr(),
839  rowptr_.getRawPtr(),
840  SLUD::SLU_NR_loc,
841  dtype, SLUD::SLU_GE);
842  }
843 
844  return true;
845 }
846 
847 
848  template<class Matrix, class Vector>
849  const char* Superludist<Matrix,Vector>::name = "SuperLU_DIST";
850 
851 
852 } // end namespace Amesos2
853 
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 > &parameterList)
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-&gt;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 > &parameterList)
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&#39;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