MFEM  v4.5.1
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
superlu.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2022, Lawrence Livermore National Security, LLC. Produced
2 // at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3 // LICENSE and NOTICE for details. LLNL-CODE-806117.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability visit https://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the BSD-3 license. We welcome feedback and contributions, see file
10 // CONTRIBUTING.md for details.
11 
12 #include "../config/config.hpp"
13 
14 #ifdef MFEM_USE_SUPERLU
15 #ifdef MFEM_USE_MPI
16 
17 #include "superlu.hpp"
18 
19 // SuperLU headers
20 #include "superlu_defs.h"
21 #include "superlu_ddefs.h"
22 
23 #if XSDK_INDEX_SIZE == 64
24 #error "SuperLUDist has been built with 64bit integers. This is not supported"
25 #endif
26 
27 // For now, it is assumed that HYPRE_BigInt is int.
28 #if defined(HYPRE_BIGINT) || defined(HYPRE_MIXEDINT)
29 #error "SuperLUDist support requires HYPRE_BigInt == int, for now."
30 #endif
31 
32 #if SUPERLU_DIST_MAJOR_VERSION > 6 || \
33  (SUPERLU_DIST_MAJOR_VERSION == 6 && SUPERLU_DIST_MINOR_VERSION > 2)
34 #define ScalePermstruct_t dScalePermstruct_t
35 #define LUstruct_t dLUstruct_t
36 #define SOLVEstruct_t dSOLVEstruct_t
37 #define ScalePermstructFree dScalePermstructFree
38 #define Destroy_LU dDestroy_LU
39 #define LUstructFree dLUstructFree
40 #define LUstructInit dLUstructInit
41 #endif
42 
43 
44 using namespace std;
45 
46 namespace mfem
47 {
48 unsigned int superlu_internal::sqrti( const unsigned int & a )
49 {
50  unsigned int a_ = a;
51  unsigned int rem = 0;
52  unsigned int root = 0;
53  unsigned short len = sizeof(int); len <<= 2;
54  unsigned short shift = (unsigned short)((len<<1) - 2);
55 
56  for (int i=0; i<len; i++)
57  {
58  root <<= 1;
59  rem = ((rem << 2) + (a_ >> shift));
60  a_ <<= 2;
61  root ++;
62  if (root <= rem)
63  {
64  rem -= root;
65  root++;
66  }
67  else
68  {
69  root--;
70  }
71  }
72  return (root >> 1);
73 }
74 
75 SuperLURowLocMatrix::SuperLURowLocMatrix(MPI_Comm comm,
76  int num_loc_rows, int first_loc_row,
77  int glob_nrows, int glob_ncols,
78  int *I, int *J, double *data)
79  : comm_(comm),
80  rowLocPtr_(NULL)
81 {
82  // Set mfem::Operator member data
83  height = num_loc_rows;
84  width = num_loc_rows;
85 
86  // Allocate SuperLU's SuperMatrix struct
87  rowLocPtr_ = new SuperMatrix;
88  SuperMatrix * A = (SuperMatrix*)rowLocPtr_;
89 
90  A->Store = NULL;
91 
92  int m = glob_nrows;
93  int n = glob_ncols;
94  int nnz_loc = I[num_loc_rows];
95  int m_loc = num_loc_rows;
96  int fst_row = first_loc_row;
97 
98  double * nzval = NULL;
99  int * colind = NULL;
100  int * rowptr = NULL;
101 
102  if ( !(nzval = doubleMalloc_dist(nnz_loc)) )
103  {
104  ABORT("Malloc fails for nzval[].");
105  }
106  for (int i=0; i<nnz_loc; i++)
107  {
108  nzval[i] = data[i];
109  }
110 
111  if ( !(colind = intMalloc_dist(nnz_loc)) )
112  {
113  ABORT("Malloc fails for colind[].");
114  }
115  for (int i=0; i<nnz_loc; i++)
116  {
117  colind[i] = J[i];
118  }
119 
120  if ( !(rowptr = intMalloc_dist(m_loc+1)) )
121  {
122  ABORT("Malloc fails for rowptr[].");
123  }
124  for (int i=0; i<=m_loc; i++)
125  {
126  rowptr[i] = I[i];
127  }
128 
129  // Assign he matrix data to SuperLU's SuperMatrix structure
130  dCreate_CompRowLoc_Matrix_dist(A, m, n, nnz_loc, m_loc, fst_row,
131  nzval, colind, rowptr,
132  SLU_NR_loc, SLU_D, SLU_GE);
133 }
134 
136  : comm_(hypParMat.GetComm()),
137  rowLocPtr_(NULL)
138 {
139  rowLocPtr_ = new SuperMatrix;
140  SuperMatrix * A = (SuperMatrix*)rowLocPtr_;
141 
142  A->Store = NULL;
143 
144  // First cast the parameter to a hypre_ParCSRMatrix
145  hypre_ParCSRMatrix * parcsr_op =
146  (hypre_ParCSRMatrix *)const_cast<HypreParMatrix&>(hypParMat);
147 
148  MFEM_ASSERT(parcsr_op != NULL,"SuperLU: const_cast failed in SetOperator");
149 
150  // Create the SuperMatrix A by borrowing the internal data from a
151  // hypre_CSRMatrix.
152  hypParMat.HostRead();
153  hypre_CSRMatrix * csr_op = hypre_MergeDiagAndOffd(parcsr_op);
154  hypParMat.HypreRead();
155  hypre_CSRMatrixSetDataOwner(csr_op,0);
156 #if MFEM_HYPRE_VERSION >= 21600
157  // For now, this method assumes that HYPRE_BigInt is int. Also, csr_op->num_cols
158  // is of type HYPRE_Int, so if we want to check for big indices in
159  // csr_op->big_j, we'll have to check all entries and that check will only be
160  // necessary in HYPRE_MIXEDINT mode which is not supported at the moment.
161  hypre_CSRMatrixBigJtoJ(csr_op);
162 #endif
163 
164  int m = parcsr_op->global_num_rows;
165  int n = parcsr_op->global_num_cols;
166  int fst_row = parcsr_op->first_row_index;
167  int nnz_loc = csr_op->num_nonzeros;
168  int m_loc = csr_op->num_rows;
169 
170  height = m_loc;
171  width = m_loc;
172 
173  double * nzval = csr_op->data;
174  int * colind = csr_op->j;
175  int * rowptr = NULL;
176 
177  // The "i" array cannot be stolen from the hypre_CSRMatrix so we'll copy it
178  if ( !(rowptr = intMalloc_dist(m_loc+1)) )
179  {
180  ABORT("Malloc fails for rowptr[].");
181  }
182  for (int i=0; i<=m_loc; i++)
183  {
184  rowptr[i] = (csr_op->i)[i];
185  }
186 
187  // Everything has been copied or abducted so delete the structure
188  hypre_CSRMatrixDestroy(csr_op);
189 
190  // Assign he matrix data to SuperLU's SuperMatrix structure
191  dCreate_CompRowLoc_Matrix_dist(A, m, n, nnz_loc, m_loc, fst_row,
192  nzval, colind, rowptr,
193  SLU_NR_loc, SLU_D, SLU_GE);
194 
195  // Save global number of columns (width) of the matrix
196  num_global_cols = n;
197 }
198 
200 {
201  SuperMatrix * A = (SuperMatrix*)rowLocPtr_;
202 
203  // Delete the internal data
204  Destroy_CompRowLoc_Matrix_dist(A);
205 
206  // Delete the struct
207  if ( A != NULL ) { delete A; }
208 }
209 
211  : comm_(comm),
212  APtr_(NULL),
213  optionsPtr_(NULL),
214  statPtr_(NULL),
215  ScalePermstructPtr_(NULL),
216  LUstructPtr_(NULL),
217  SOLVEstructPtr_(NULL),
218  gridPtr_(NULL),
219  berr_(NULL),
220  perm_r_(NULL),
221  nrhs_(1),
222  nprow_(0),
223  npcol_(0),
224  firstSolveWithThisA_(false),
225  gridInitialized_(false),
226  LUStructInitialized_(false)
227 {
228  this->Init();
229 }
230 
232  : comm_(A.GetComm()),
233  APtr_(&A),
234  optionsPtr_(NULL),
235  statPtr_(NULL),
236  ScalePermstructPtr_(NULL),
237  LUstructPtr_(NULL),
238  SOLVEstructPtr_(NULL),
239  gridPtr_(NULL),
240  berr_(NULL),
241  perm_r_(NULL),
242  nrhs_(1),
243  nprow_(0),
244  npcol_(0),
245  firstSolveWithThisA_(true),
246  gridInitialized_(false),
247  LUStructInitialized_(false)
248 {
249  height = A.Height();
250  width = A.Width();
251 
252  this->Init();
253 }
254 
256 {
257  superlu_dist_options_t * options = (superlu_dist_options_t*)optionsPtr_;
258  SuperLUStat_t * stat = (SuperLUStat_t*)statPtr_;
259  ScalePermstruct_t * SPstruct = (ScalePermstruct_t*)ScalePermstructPtr_;
260  LUstruct_t * LUstruct = (LUstruct_t*)LUstructPtr_;
261  SOLVEstruct_t * SOLVEstruct = (SOLVEstruct_t*)SOLVEstructPtr_;
262  gridinfo_t * grid = (gridinfo_t*)gridPtr_;
263 
264  SUPERLU_FREE(berr_);
265  PStatFree(stat);
266 
267  if ( LUStructInitialized_ )
268  {
269  ScalePermstructFree(SPstruct);
270  Destroy_LU(APtr_->GetGlobalNumColumns(), grid, LUstruct);
271  LUstructFree(LUstruct);
272  }
273 
274  if ( options->SolveInitialized )
275  {
276  dSolveFinalize(options, SOLVEstruct);
277  }
278 
279  if ( options != NULL ) { delete options; }
280  if ( stat != NULL ) { delete stat; }
281  if ( SPstruct != NULL ) { delete SPstruct; }
282  if ( LUstruct != NULL ) { delete LUstruct; }
283  if ( SOLVEstruct != NULL ) { delete SOLVEstruct; }
284  if ( grid != NULL ) { delete grid; }
285  if ( perm_r_ != NULL ) { SUPERLU_FREE(perm_r_); }
286 }
287 
288 void SuperLUSolver::Init()
289 {
290  MPI_Comm_size(comm_, &numProcs_);
291  MPI_Comm_rank(comm_, &myid_);
292 
293  optionsPtr_ = new superlu_dist_options_t;
294  statPtr_ = new SuperLUStat_t;
295  ScalePermstructPtr_ = new ScalePermstruct_t;
296  LUstructPtr_ = new LUstruct_t;
297  SOLVEstructPtr_ = new SOLVEstruct_t;
298  gridPtr_ = new gridinfo_t;
299 
300  superlu_dist_options_t * options = (superlu_dist_options_t*)optionsPtr_;
301  SuperLUStat_t * stat = (SuperLUStat_t*)statPtr_;
302 
303  if ( !(berr_ = doubleMalloc_dist(nrhs_)) )
304  {
305  ABORT("Malloc fails for berr[].");
306  }
307 
308  // Set default options
309  set_default_options_dist(options);
310 
311  // Choose nprow and npcol so that the process grid is as square as possible.
312  // If the processes cannot be divided evenly, keep the row dimension smaller
313  // than the column dimension.
314 
315  nprow_ = (int)superlu_internal::sqrti((unsigned int)numProcs_);
316  while (numProcs_ % nprow_ != 0 && nprow_ > 0)
317  {
318  nprow_--;
319  }
320 
321  npcol_ = (int)(numProcs_ / nprow_);
322  MFEM_ASSERT(nprow_ * npcol_ == numProcs_, "");
323 
324  PStatInit(stat); // Initialize the statistics variables.
325 }
326 
327 void SuperLUSolver::SetPrintStatistics( bool print_stat )
328 {
329  superlu_dist_options_t * options = (superlu_dist_options_t*)optionsPtr_;
330 
331  yes_no_t opt = print_stat?YES:NO;
332 
333  options->PrintStat = opt;
334 }
335 
337 {
338  superlu_dist_options_t * options = (superlu_dist_options_t*)optionsPtr_;
339 
340  yes_no_t opt = equil?YES:NO;
341 
342  options->Equil = opt;
343 }
344 
346 {
347  superlu_dist_options_t * options = (superlu_dist_options_t*)optionsPtr_;
348 
349  colperm_t opt = (colperm_t)col_perm;
350 
351  options->ColPerm = opt;
352 }
353 
355  Array<int> * perm )
356 {
357  superlu_dist_options_t * options = (superlu_dist_options_t*)optionsPtr_;
358 
359  rowperm_t opt = (rowperm_t)row_perm;
360 
361  options->RowPerm = opt;
362 
363  if ( opt == MY_PERMR )
364  {
365  if ( perm == NULL )
366  {
367  mfem_error("SuperLUSolver::SetRowPermutation :"
368  " permutation vector not set!");
369  }
370 
371  if ( !(perm_r_ = intMalloc_dist(perm->Size())) )
372  {
373  ABORT("Malloc fails for perm_r[].");
374  }
375  for (int i=0; i<perm->Size(); i++)
376  {
377  perm_r_[i] = (*perm)[i];
378  }
379  }
380 }
381 
383 {
384  superlu_dist_options_t * options = (superlu_dist_options_t*)optionsPtr_;
385 
386  trans_t opt = (trans_t)trans;
387 
388  options->Trans = opt;
389 }
390 
392 {
393  superlu_dist_options_t * options = (superlu_dist_options_t*)optionsPtr_;
394 
395  IterRefine_t opt = (IterRefine_t)iter_ref;
396 
397  options->IterRefine = opt;
398 }
399 
401 {
402  superlu_dist_options_t * options = (superlu_dist_options_t*)optionsPtr_;
403 
404  yes_no_t opt = rtp?YES:NO;
405 
406  options->ReplaceTinyPivot = opt;
407 }
408 
409 void SuperLUSolver::SetNumLookAheads( int num_lookaheads )
410 {
411  superlu_dist_options_t * options = (superlu_dist_options_t*)optionsPtr_;
412 
413  options->num_lookaheads = num_lookaheads;
414 }
415 
417 {
418  superlu_dist_options_t * options = (superlu_dist_options_t*)optionsPtr_;
419 
420  yes_no_t opt = etree?YES:NO;
421 
422  options->lookahead_etree = opt;
423 }
424 
426 {
427  superlu_dist_options_t * options = (superlu_dist_options_t*)optionsPtr_;
428 
429  yes_no_t opt = sym?YES:NO;
430 
431  options->SymPattern = opt;
432 }
433 
435 {
436  superlu_dist_options_t * options = (superlu_dist_options_t*)optionsPtr_;
437 
438  yes_no_t opt = par?YES:NO;
439 
440  options->ParSymbFact = opt;
441 }
442 
444 {
445  gridinfo_t * grid = (gridinfo_t*)gridPtr_;
446 
447  // Make sure the values of nprow and npcol are reasonable
448  if ( ((nprow_ * npcol_) > numProcs_) || ((nprow_ * npcol_) < 1) )
449  {
450  if ( myid_ == 0 )
451  {
452  mfem::err << "Warning: User specified nprow and npcol are such that "
453  << "(nprow * npcol) > numProcs or (nprow * npcol) < 1. "
454  << "Using default values for nprow and npcol instead."
455  << endl;
456  }
457 
458  nprow_ = (int)superlu_internal::sqrti((unsigned int)numProcs_);
459  while (numProcs_ % nprow_ != 0 && nprow_ > 0)
460  {
461  nprow_--;
462  }
463 
464  npcol_ = (int)(numProcs_ / nprow_);
465  MFEM_ASSERT(nprow_ * npcol_ == numProcs_, "");
466  }
467 
468  superlu_gridinit(comm_, nprow_, npcol_, grid);
469 
470  gridInitialized_ = true;
471 }
472 
474 {
475  if ( gridInitialized_ )
476  {
477  gridinfo_t * grid = (gridinfo_t*)gridPtr_;
478 
479  superlu_gridexit(grid);
480  }
481 
482  gridInitialized_ = false;
483 }
484 
485 void SuperLUSolver::Mult( const Vector & x, Vector & y ) const
486 {
487  MFEM_ASSERT(APtr_ != NULL,
488  "SuperLU Error: The operator must be set before"
489  " the system can be solved.");
490 
491  superlu_dist_options_t * options = (superlu_dist_options_t*)optionsPtr_;
492  SuperLUStat_t * stat = (SuperLUStat_t*)statPtr_;
493  SuperMatrix * A = (SuperMatrix*)APtr_->InternalData();
494 
495  ScalePermstruct_t * SPstruct = (ScalePermstruct_t*)ScalePermstructPtr_;
496  LUstruct_t * LUstruct = (LUstruct_t*)LUstructPtr_;
497  SOLVEstruct_t * SOLVEstruct = (SOLVEstruct_t*)SOLVEstructPtr_;
498  gridinfo_t * grid = (gridinfo_t*)gridPtr_;
499 
501  {
502  options->Fact = FACTORED; // Indicate the factored form of A is supplied.
503  }
504  else // This is the first solve with this A
505  {
506  firstSolveWithThisA_ = false;
507 
508  // Make sure that the parameters have been initialized The only parameter
509  // we might have to worry about is ScalePermstruct, if the user is
510  // supplying a row or column permutation.
511 
512  // Initialize ScalePermstruct and LUstruct.
513  SPstruct->DiagScale = NOEQUIL;
514 
515  // Transfer ownership of the row permutations if available
516  if ( perm_r_ != NULL )
517  {
518  SPstruct->perm_r = perm_r_;
519  perm_r_ = NULL;
520  }
521  else
522  {
523  if ( !(SPstruct->perm_r = intMalloc_dist(A->nrow)) )
524  {
525  ABORT("Malloc fails for perm_r[].");
526  }
527  }
528  if ( !(SPstruct->perm_c = intMalloc_dist(A->ncol)) )
529  {
530  ABORT("Malloc fails for perm_c[].");
531  }
532 
533  LUstructInit(A->ncol, LUstruct);
534  LUStructInitialized_ = true;
535  }
536 
537  // SuperLU overwrites x with y, so copy x to y and pass that to the solve
538  // routine.
539 
540  const double *xPtr = x.HostRead();
541  y = xPtr;
542  double * yPtr = y.HostReadWrite();
543 
544  int info = -1, locSize = y.Size();
545 
546  // Solve the system
547  pdgssvx(options, A, SPstruct, yPtr, locSize, nrhs_, grid,
548  LUstruct, SOLVEstruct, berr_, stat, &info);
549 
550  if ( info != 0 )
551  {
552  if ( info < 0 )
553  {
554  switch (-info)
555  {
556  case 1:
557  MFEM_ABORT("SuperLU: SuperLU options are invalid.");
558  break;
559  case 2:
560  MFEM_ABORT("SuperLU: Matrix A (in Ax=b) is invalid.");
561  break;
562  case 5:
563  MFEM_ABORT("SuperLU: Vector b dimension (in Ax=b) is invalid.");
564  break;
565  case 6:
566  MFEM_ABORT("SuperLU: Number of right-hand sides is invalid.");
567  break;
568  default:
569  MFEM_ABORT("SuperLU: Parameter with index "
570  << -info << "invalid. (1-indexed)");
571  break;
572  }
573  }
574  else if ( info <= A->ncol )
575  {
576  MFEM_ABORT("SuperLU: Found a singular matrix, U("
577  << info << "," << info << ") is exactly zero.");
578  }
579  else if ( info > A->ncol )
580  {
581  MFEM_ABORT("SuperLU: Memory allocation error with "
582  << info - A->ncol << " bytes already allocated,");
583  }
584  else
585  {
586  MFEM_ABORT("Unknown SuperLU Error");
587  }
588  }
589 }
590 
592 {
593  // Verify that we have a compatible operator
594  APtr_ = dynamic_cast<const SuperLURowLocMatrix*>(&op);
595  if ( APtr_ == NULL )
596  {
597  mfem_error("SuperLUSolver::SetOperator : not SuperLURowLocMatrix!");
598  }
599 
600  // Everything is OK so finish setting the operator
601  firstSolveWithThisA_ = true;
602 
603  // Set mfem::Operator member data
604  height = op.Height();
605  width = op.Width();
606 
607  // Initialize the processor grid if necessary
608  if (!gridInitialized_)
609  {
610  this->SetupGrid();
611  }
612 }
613 
614 } // mfem namespace
615 
616 #endif // MFEM_USE_MPI
617 #endif // MFEM_USE_SUPERLU
int Size() const
Return the logical size of the array.
Definition: array.hpp:138
void * InternalData() const
Definition: superlu.hpp:73
void trans(const Vector &u, Vector &x)
Definition: ex27.cpp:412
void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Definition: superlu.cpp:485
void * ScalePermstructPtr_
Definition: superlu.hpp:148
void HypreRead() const
Update the internal hypre_ParCSRMatrix object, A, to be in hypre memory space.
Definition: hypre.hpp:845
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
Definition: operator.hpp:73
void SetRowPermutation(superlu::RowPerm row_perm, Array< int > *perm=NULL)
Definition: superlu.cpp:354
int Size() const
Returns the size of the vector.
Definition: vector.hpp:200
void SetSymmetricPattern(bool sym)
Definition: superlu.cpp:425
void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition: superlu.cpp:591
HYPRE_BigInt GetGlobalNumColumns() const
Definition: superlu.hpp:75
void SetIterativeRefine(superlu::IterRefine iter_ref)
Definition: superlu.cpp:391
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition: operator.hpp:67
void SetColumnPermutation(superlu::ColPerm col_perm)
Definition: superlu.cpp:345
void mfem_error(const char *msg)
Function called when an error is encountered. Used by the macros MFEM_ABORT, MFEM_ASSERT, MFEM_VERIFY.
Definition: error.cpp:154
void SetNumLookAheads(int num_lookaheads)
Definition: superlu.cpp:409
void SetPrintStatistics(bool print_stat)
Definition: superlu.cpp:327
virtual const double * HostRead() const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), false).
Definition: vector.hpp:453
OutStream err(std::cerr)
Global stream used by the library for standard error output. Initially it uses the same std::streambu...
Definition: globals.hpp:71
void HostRead() const
Update the internal hypre_ParCSRMatrix object, A, to be on host.
Definition: hypre.hpp:828
A class to initialize the size of a Tensor.
Definition: dtensor.hpp:54
int height
Dimension of the output / number of rows in the matrix.
Definition: operator.hpp:27
double a
Definition: lissajous.cpp:41
void * SOLVEstructPtr_
Definition: superlu.hpp:150
void SetLookAheadElimTree(bool etree)
Definition: superlu.cpp:416
unsigned int sqrti(const unsigned int &a)
Definition: superlu.cpp:48
void SetReplaceTinyPivot(bool rtp)
Definition: superlu.cpp:400
Vector data type.
Definition: vector.hpp:60
void SetTranspose(superlu::Trans trans)
Definition: superlu.cpp:382
const SuperLURowLocMatrix * APtr_
Definition: superlu.hpp:138
void SetParSymbFact(bool par)
Definition: superlu.cpp:434
SuperLURowLocMatrix(MPI_Comm comm, int num_loc_rows, int first_loc_row, int glob_nrows, int glob_ncols, int *I, int *J, double *data)
Definition: superlu.cpp:75
Abstract operator.
Definition: operator.hpp:24
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:343
SuperLUSolver(MPI_Comm comm)
Definition: superlu.cpp:210
virtual double * HostReadWrite()
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), false).
Definition: vector.hpp:469
int width
Dimension of the input / number of columns in the matrix.
Definition: operator.hpp:28
void SetEquilibriate(bool equil)
Definition: superlu.cpp:336