MFEM  v4.2.0
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-2020, 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_Int is int.
28 #ifdef HYPRE_BIGINT
29 #error "SuperLUDist support requires HYPRE_Int == 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  hypre_CSRMatrix * csr_op = hypre_MergeDiagAndOffd(parcsr_op);
153  hypre_CSRMatrixSetDataOwner(csr_op,0);
154 #if MFEM_HYPRE_VERSION >= 21600
155  // For now, this method assumes that HYPRE_Int is int. Also, csr_op->num_cols
156  // is of type HYPRE_Int, so if we want to check for big indices in
157  // csr_op->big_j, we'll have to check all entries and that check will only be
158  // necessary in HYPRE_MIXEDINT mode which is not supported at the moment.
159  hypre_CSRMatrixBigJtoJ(csr_op);
160 #endif
161 
162  int m = parcsr_op->global_num_rows;
163  int n = parcsr_op->global_num_cols;
164  int fst_row = parcsr_op->first_row_index;
165  int nnz_loc = csr_op->num_nonzeros;
166  int m_loc = csr_op->num_rows;
167 
168  height = m_loc;
169  width = m_loc;
170 
171  double * nzval = csr_op->data;
172  int * colind = csr_op->j;
173  int * rowptr = NULL;
174 
175  // The "i" array cannot be stolen from the hypre_CSRMatrix so we'll copy it
176  if ( !(rowptr = intMalloc_dist(m_loc+1)) )
177  {
178  ABORT("Malloc fails for rowptr[].");
179  }
180  for (int i=0; i<=m_loc; i++)
181  {
182  rowptr[i] = (csr_op->i)[i];
183  }
184 
185  // Everything has been copied or abducted so delete the structure
186  hypre_CSRMatrixDestroy(csr_op);
187 
188  // Assign he matrix data to SuperLU's SuperMatrix structure
189  dCreate_CompRowLoc_Matrix_dist(A, m, n, nnz_loc, m_loc, fst_row,
190  nzval, colind, rowptr,
191  SLU_NR_loc, SLU_D, SLU_GE);
192 }
193 
195 {
196  SuperMatrix * A = (SuperMatrix*)rowLocPtr_;
197 
198  // Delete the internal data
199  Destroy_CompRowLoc_Matrix_dist(A);
200 
201  // Delete the struct
202  if ( A != NULL ) { delete A; }
203 }
204 
206  : comm_(comm),
207  APtr_(NULL),
208  optionsPtr_(NULL),
209  statPtr_(NULL),
210  ScalePermstructPtr_(NULL),
211  LUstructPtr_(NULL),
212  SOLVEstructPtr_(NULL),
213  gridPtr_(NULL),
214  berr_(NULL),
215  perm_r_(NULL),
216  nrhs_(1),
217  nprow_(0),
218  npcol_(0),
219  firstSolveWithThisA_(false),
220  gridInitialized_(false),
221  LUStructInitialized_(false)
222 {
223  this->Init();
224 }
225 
227  : comm_(A.GetComm()),
228  APtr_(&A),
229  optionsPtr_(NULL),
230  statPtr_(NULL),
231  ScalePermstructPtr_(NULL),
232  LUstructPtr_(NULL),
233  SOLVEstructPtr_(NULL),
234  gridPtr_(NULL),
235  berr_(NULL),
236  perm_r_(NULL),
237  nrhs_(1),
238  nprow_(0),
239  npcol_(0),
240  firstSolveWithThisA_(true),
241  gridInitialized_(false),
242  LUStructInitialized_(false)
243 {
244  height = A.Height();
245  width = A.Width();
246 
247  this->Init();
248 }
249 
251 {
252  superlu_dist_options_t * options = (superlu_dist_options_t*)optionsPtr_;
253  SuperLUStat_t * stat = (SuperLUStat_t*)statPtr_;
254  ScalePermstruct_t * SPstruct = (ScalePermstruct_t*)ScalePermstructPtr_;
255  LUstruct_t * LUstruct = (LUstruct_t*)LUstructPtr_;
256  SOLVEstruct_t * SOLVEstruct = (SOLVEstruct_t*)SOLVEstructPtr_;
257  gridinfo_t * grid = (gridinfo_t*)gridPtr_;
258 
259  SUPERLU_FREE(berr_);
260  PStatFree(stat);
261 
262  if ( LUStructInitialized_ )
263  {
264  ScalePermstructFree(SPstruct);
265  Destroy_LU(width, grid, LUstruct);
266  LUstructFree(LUstruct);
267  }
268 
269  if ( options->SolveInitialized )
270  {
271  dSolveFinalize(options, SOLVEstruct);
272  }
273 
274  if ( options != NULL ) { delete options; }
275  if ( stat != NULL ) { delete stat; }
276  if ( SPstruct != NULL ) { delete SPstruct; }
277  if ( LUstruct != NULL ) { delete LUstruct; }
278  if ( SOLVEstruct != NULL ) { delete SOLVEstruct; }
279  if ( grid != NULL ) { delete grid; }
280  if ( perm_r_ != NULL ) { SUPERLU_FREE(perm_r_); }
281 }
282 
283 void SuperLUSolver::Init()
284 {
285  MPI_Comm_size(comm_, &numProcs_);
286  MPI_Comm_rank(comm_, &myid_);
287 
288  optionsPtr_ = new superlu_dist_options_t;
289  statPtr_ = new SuperLUStat_t;
290  ScalePermstructPtr_ = new ScalePermstruct_t;
291  LUstructPtr_ = new LUstruct_t;
292  SOLVEstructPtr_ = new SOLVEstruct_t;
293  gridPtr_ = new gridinfo_t;
294 
295  superlu_dist_options_t * options = (superlu_dist_options_t*)optionsPtr_;
296  SuperLUStat_t * stat = (SuperLUStat_t*)statPtr_;
297 
298  if ( !(berr_ = doubleMalloc_dist(nrhs_)) )
299  {
300  ABORT("Malloc fails for berr[].");
301  }
302 
303  // Set default options
304  set_default_options_dist(options);
305 
306  // Choose nprow and npcol so that the process grid is as square as possible.
307  // If the processes cannot be divided evenly, keep the row dimension smaller
308  // than the column dimension.
309 
310  nprow_ = (int)superlu_internal::sqrti((unsigned int)numProcs_);
311  while (numProcs_ % nprow_ != 0 && nprow_ > 0)
312  {
313  nprow_--;
314  }
315 
316  npcol_ = (int)(numProcs_ / nprow_);
317  MFEM_ASSERT(nprow_ * npcol_ == numProcs_, "");
318 
319  PStatInit(stat); // Initialize the statistics variables.
320 }
321 
322 void SuperLUSolver::SetPrintStatistics( bool print_stat )
323 {
324  superlu_dist_options_t * options = (superlu_dist_options_t*)optionsPtr_;
325 
326  yes_no_t opt = print_stat?YES:NO;
327 
328  options->PrintStat = opt;
329 }
330 
332 {
333  superlu_dist_options_t * options = (superlu_dist_options_t*)optionsPtr_;
334 
335  yes_no_t opt = equil?YES:NO;
336 
337  options->Equil = opt;
338 }
339 
341 {
342  superlu_dist_options_t * options = (superlu_dist_options_t*)optionsPtr_;
343 
344  colperm_t opt = (colperm_t)col_perm;
345 
346  options->ColPerm = opt;
347 }
348 
350  Array<int> * perm )
351 {
352  superlu_dist_options_t * options = (superlu_dist_options_t*)optionsPtr_;
353 
354  rowperm_t opt = (rowperm_t)row_perm;
355 
356  options->RowPerm = opt;
357 
358  if ( opt == MY_PERMR )
359  {
360  if ( perm == NULL )
361  {
362  mfem_error("SuperLUSolver::SetRowPermutation :"
363  " permutation vector not set!");
364  }
365 
366  if ( !(perm_r_ = intMalloc_dist(perm->Size())) )
367  {
368  ABORT("Malloc fails for perm_r[].");
369  }
370  for (int i=0; i<perm->Size(); i++)
371  {
372  perm_r_[i] = (*perm)[i];
373  }
374  }
375 }
376 
378 {
379  superlu_dist_options_t * options = (superlu_dist_options_t*)optionsPtr_;
380 
381  trans_t opt = (trans_t)trans;
382 
383  options->Trans = opt;
384 }
385 
387 {
388  superlu_dist_options_t * options = (superlu_dist_options_t*)optionsPtr_;
389 
390  IterRefine_t opt = (IterRefine_t)iter_ref;
391 
392  options->IterRefine = opt;
393 }
394 
396 {
397  superlu_dist_options_t * options = (superlu_dist_options_t*)optionsPtr_;
398 
399  yes_no_t opt = rtp?YES:NO;
400 
401  options->ReplaceTinyPivot = opt;
402 }
403 
404 void SuperLUSolver::SetNumLookAheads( int num_lookaheads )
405 {
406  superlu_dist_options_t * options = (superlu_dist_options_t*)optionsPtr_;
407 
408  options->num_lookaheads = num_lookaheads;
409 }
410 
412 {
413  superlu_dist_options_t * options = (superlu_dist_options_t*)optionsPtr_;
414 
415  yes_no_t opt = etree?YES:NO;
416 
417  options->lookahead_etree = opt;
418 }
419 
421 {
422  superlu_dist_options_t * options = (superlu_dist_options_t*)optionsPtr_;
423 
424  yes_no_t opt = sym?YES:NO;
425 
426  options->SymPattern = opt;
427 }
428 
430 {
431  gridinfo_t * grid = (gridinfo_t*)gridPtr_;
432 
433  // Make sure the values of nprow and npcol are reasonable
434  if ( ((nprow_ * npcol_) > numProcs_) || ((nprow_ * npcol_) < 1) )
435  {
436  if ( myid_ == 0 )
437  {
438  mfem::err << "Warning: User specified nprow and npcol are such that "
439  << "(nprow * npcol) > numProcs or (nprow * npcol) < 1. "
440  << "Using default values for nprow and npcol instead."
441  << endl;
442  }
443 
444  nprow_ = (int)superlu_internal::sqrti((unsigned int)numProcs_);
445  while (numProcs_ % nprow_ != 0 && nprow_ > 0)
446  {
447  nprow_--;
448  }
449 
450  npcol_ = (int)(numProcs_ / nprow_);
451  MFEM_ASSERT(nprow_ * npcol_ == numProcs_, "");
452  }
453 
454  superlu_gridinit(comm_, nprow_, npcol_, grid);
455 
456  gridInitialized_ = true;
457 }
458 
460 {
461  if ( gridInitialized_ )
462  {
463  gridinfo_t * grid = (gridinfo_t*)gridPtr_;
464 
465  superlu_gridexit(grid);
466  }
467 
468  gridInitialized_ = false;
469 }
470 
471 void SuperLUSolver::Mult( const Vector & x, Vector & y ) const
472 {
473  MFEM_ASSERT(APtr_ != NULL,
474  "SuperLU Error: The operator must be set before"
475  " the system can be solved.");
476 
477  superlu_dist_options_t * options = (superlu_dist_options_t*)optionsPtr_;
478  SuperLUStat_t * stat = (SuperLUStat_t*)statPtr_;
479  SuperMatrix * A = (SuperMatrix*)APtr_->InternalData();
480 
481  ScalePermstruct_t * SPstruct = (ScalePermstruct_t*)ScalePermstructPtr_;
482  LUstruct_t * LUstruct = (LUstruct_t*)LUstructPtr_;
483  SOLVEstruct_t * SOLVEstruct = (SOLVEstruct_t*)SOLVEstructPtr_;
484  gridinfo_t * grid = (gridinfo_t*)gridPtr_;
485 
487  {
488  options->Fact = FACTORED; // Indicate the factored form of A is supplied.
489  }
490  else // This is the first solve with this A
491  {
492  firstSolveWithThisA_ = false;
493 
494  // Make sure that the parameters have been initialized The only parameter
495  // we might have to worry about is ScalePermstruct, if the user is
496  // supplying a row or column permutation.
497 
498  // Initialize ScalePermstruct and LUstruct.
499  SPstruct->DiagScale = NOEQUIL;
500 
501  // Transfer ownership of the row permutations if available
502  if ( perm_r_ != NULL )
503  {
504  SPstruct->perm_r = perm_r_;
505  perm_r_ = NULL;
506  }
507  else
508  {
509  if ( !(SPstruct->perm_r = intMalloc_dist(A->nrow)) )
510  {
511  ABORT("Malloc fails for perm_r[].");
512  }
513  }
514  if ( !(SPstruct->perm_c = intMalloc_dist(A->ncol)) )
515  {
516  ABORT("Malloc fails for perm_c[].");
517  }
518 
519  LUstructInit(A->ncol, LUstruct);
520  LUStructInitialized_ = true;
521  }
522 
523  // SuperLU overwrites x with y, so copy x to y and pass that to the solve
524  // routine.
525 
526  y = x;
527 
528  double* yPtr = (double*)y;
529  int info = -1, locSize = y.Size();
530 
531  // Solve the system
532  pdgssvx(options, A, SPstruct, yPtr, locSize, nrhs_, grid,
533  LUstruct, SOLVEstruct, berr_, stat, &info);
534 
535  if ( info != 0 )
536  {
537  if ( info < 0 )
538  {
539  switch (-info)
540  {
541  case 1:
542  MFEM_ABORT("SuperLU: SuperLU options are invalid.");
543  break;
544  case 2:
545  MFEM_ABORT("SuperLU: Matrix A (in Ax=b) is invalid.");
546  break;
547  case 5:
548  MFEM_ABORT("SuperLU: Vector b dimension (in Ax=b) is invalid.");
549  break;
550  case 6:
551  MFEM_ABORT("SuperLU: Number of right-hand sides is invalid.");
552  break;
553  default:
554  MFEM_ABORT("SuperLU: Parameter with index "
555  << -info << "invalid. (1-indexed)");
556  break;
557  }
558  }
559  else if ( info <= A->ncol )
560  {
561  MFEM_ABORT("SuperLU: Found a singular matrix, U("
562  << info << "," << info << ") is exactly zero.");
563  }
564  else if ( info > A->ncol )
565  {
566  MFEM_ABORT("SuperLU: Memory allocation error with "
567  << info - A->ncol << " bytes already allocated,");
568  }
569  else
570  {
571  MFEM_ABORT("Unknown SuperLU Error");
572  }
573  }
574 }
575 
577 {
578  // Verify that we have a compatible operator
579  APtr_ = dynamic_cast<const SuperLURowLocMatrix*>(&op);
580  if ( APtr_ == NULL )
581  {
582  mfem_error("SuperLUSolver::SetOperator : not SuperLURowLocMatrix!");
583  }
584 
585  // Everything is OK so finish setting the operator
586  firstSolveWithThisA_ = true;
587 
588  // Set mfem::Operator member data
589  height = op.Height();
590  width = op.Width();
591 
592  // Initialize the processor grid if necessary
593  if (!gridInitialized_)
594  {
595  this->SetupGrid();
596  }
597 }
598 
599 } // mfem namespace
600 
601 #endif // MFEM_USE_MPI
602 #endif // MFEM_USE_SUPERLU
int Size() const
Return the logical size of the array.
Definition: array.hpp:124
void * InternalData() const
Definition: superlu.hpp:73
void trans(const Vector &u, Vector &x)
Definition: ex27.cpp:421
void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Definition: superlu.cpp:471
void * ScalePermstructPtr_
Definition: superlu.hpp:144
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
Definition: operator.hpp:71
void SetRowPermutation(superlu::RowPerm row_perm, Array< int > *perm=NULL)
Definition: superlu.cpp:349
int Size() const
Returns the size of the vector.
Definition: vector.hpp:160
void SetSymmetricPattern(bool sym)
Definition: superlu.cpp:420
void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition: superlu.cpp:576
void SetIterativeRefine(superlu::IterRefine iter_ref)
Definition: superlu.cpp:386
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition: operator.hpp:65
void SetColumnPermutation(superlu::ColPerm col_perm)
Definition: superlu.cpp:340
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:153
void SetNumLookAheads(int num_lookaheads)
Definition: superlu.cpp:404
void SetPrintStatistics(bool print_stat)
Definition: superlu.cpp:322
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
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:146
void SetLookAheadElimTree(bool etree)
Definition: superlu.cpp:411
unsigned int sqrti(const unsigned int &a)
Definition: superlu.cpp:48
void SetReplaceTinyPivot(bool rtp)
Definition: superlu.cpp:395
Vector data type.
Definition: vector.hpp:51
void SetTranspose(superlu::Trans trans)
Definition: superlu.cpp:377
const SuperLURowLocMatrix * APtr_
Definition: superlu.hpp:134
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:181
SuperLUSolver(MPI_Comm comm)
Definition: superlu.cpp:205
int width
Dimension of the input / number of columns in the matrix.
Definition: operator.hpp:28
void SetEquilibriate(bool equil)
Definition: superlu.cpp:331