MFEM  v4.1.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 using namespace std;
28 
29 namespace mfem
30 {
31 unsigned int superlu_internal::sqrti( const unsigned int & a )
32 {
33  unsigned int a_ = a;
34  unsigned int rem = 0;
35  unsigned int root = 0;
36  unsigned short len = sizeof(int); len <<= 2;
37  unsigned short shift = (unsigned short)((len<<1) - 2);
38 
39  for (int i=0; i<len; i++)
40  {
41  root <<= 1;
42  rem = ((rem << 2) + (a_ >> shift));
43  a_ <<= 2;
44  root ++;
45  if (root <= rem)
46  {
47  rem -= root;
48  root++;
49  }
50  else
51  {
52  root--;
53  }
54  }
55  return (root >> 1);
56 }
57 
58 SuperLURowLocMatrix::SuperLURowLocMatrix(MPI_Comm comm,
59  int num_loc_rows, int first_loc_row,
60  int glob_nrows, int glob_ncols,
61  int *I, int *J, double *data)
62  : comm_(comm),
63  rowLocPtr_(NULL)
64 {
65  // Set mfem::Operator member data
66  height = num_loc_rows;
67  width = num_loc_rows;
68 
69  // Allocate SuperLU's SuperMatrix struct
70  rowLocPtr_ = new SuperMatrix;
71  SuperMatrix * A = (SuperMatrix*)rowLocPtr_;
72 
73  A->Store = NULL;
74 
75  int m = glob_nrows;
76  int n = glob_ncols;
77  int nnz_loc = I[num_loc_rows];
78  int m_loc = num_loc_rows;
79  int fst_row = first_loc_row;
80 
81  double * nzval = NULL;
82  int * colind = NULL;
83  int * rowptr = NULL;
84 
85  if ( !(nzval = doubleMalloc_dist(nnz_loc)) )
86  {
87  ABORT("Malloc fails for nzval[].");
88  }
89  for (int i=0; i<nnz_loc; i++)
90  {
91  nzval[i] = data[i];
92  }
93 
94  if ( !(colind = intMalloc_dist(nnz_loc)) )
95  {
96  ABORT("Malloc fails for colind[].");
97  }
98  for (int i=0; i<nnz_loc; i++)
99  {
100  colind[i] = J[i];
101  }
102 
103  if ( !(rowptr = intMalloc_dist(m_loc+1)) )
104  {
105  ABORT("Malloc fails for rowptr[].");
106  }
107  for (int i=0; i<=m_loc; i++)
108  {
109  rowptr[i] = I[i];
110  }
111 
112  // Assign he matrix data to SuperLU's SuperMatrix structure
113  dCreate_CompRowLoc_Matrix_dist(A, m, n, nnz_loc, m_loc, fst_row,
114  nzval, colind, rowptr,
115  SLU_NR_loc, SLU_D, SLU_GE);
116 }
117 
119  : comm_(hypParMat.GetComm()),
120  rowLocPtr_(NULL)
121 {
122  rowLocPtr_ = new SuperMatrix;
123  SuperMatrix * A = (SuperMatrix*)rowLocPtr_;
124 
125  A->Store = NULL;
126 
127  // First cast the parameter to a hypre_ParCSRMatrix
128  hypre_ParCSRMatrix * parcsr_op =
129  (hypre_ParCSRMatrix *)const_cast<HypreParMatrix&>(hypParMat);
130 
131  MFEM_ASSERT(parcsr_op != NULL,"SuperLU: const_cast failed in SetOperator");
132 
133  // Create the SuperMatrix A by borrowing the internal data from a
134  // hypre_CSRMatrix.
135  hypre_CSRMatrix * csr_op = hypre_MergeDiagAndOffd(parcsr_op);
136  hypre_CSRMatrixSetDataOwner(csr_op,0);
137 #if MFEM_HYPRE_VERSION >= 21600
138  MFEM_VERIFY(csr_op->num_rows < INT_MAX,"SuperLU: number of local rows "
139  "is too large to store as an integer.");
140  hypre_CSRMatrixBigJtoJ(csr_op);
141 #endif
142 
143  int m = parcsr_op->global_num_rows;
144  int n = parcsr_op->global_num_cols;
145  int fst_row = parcsr_op->first_row_index;
146  int nnz_loc = csr_op->num_nonzeros;
147  int m_loc = csr_op->num_rows;
148 
149  height = m_loc;
150  width = m_loc;
151 
152  double * nzval = csr_op->data;
153  int * colind = csr_op->j;
154  int * rowptr = NULL;
155 
156  // The "i" array cannot be stolen from the hypre_CSRMatrix so we'll copy it
157  if ( !(rowptr = intMalloc_dist(m_loc+1)) )
158  {
159  ABORT("Malloc fails for rowptr[].");
160  }
161  for (int i=0; i<=m_loc; i++)
162  {
163  rowptr[i] = (csr_op->i)[i];
164  }
165 
166  // Everything has been copied or abducted so delete the structure
167  hypre_CSRMatrixDestroy(csr_op);
168 
169  // Assign he matrix data to SuperLU's SuperMatrix structure
170  dCreate_CompRowLoc_Matrix_dist(A, m, n, nnz_loc, m_loc, fst_row,
171  nzval, colind, rowptr,
172  SLU_NR_loc, SLU_D, SLU_GE);
173 }
174 
176 {
177  SuperMatrix * A = (SuperMatrix*)rowLocPtr_;
178 
179  // Delete the internal data
180  Destroy_CompRowLoc_Matrix_dist(A);
181 
182  // Delete the struct
183  if ( A != NULL ) { delete A; }
184 }
185 
187  : comm_(comm),
188  APtr_(NULL),
189  optionsPtr_(NULL),
190  statPtr_(NULL),
191  ScalePermstructPtr_(NULL),
192  LUstructPtr_(NULL),
193  SOLVEstructPtr_(NULL),
194  gridPtr_(NULL),
195  berr_(NULL),
196  perm_r_(NULL),
197  nrhs_(1),
198  nprow_(0),
199  npcol_(0),
200  firstSolveWithThisA_(false),
201  gridInitialized_(false),
202  LUStructInitialized_(false)
203 {
204  this->Init();
205 }
206 
208  : comm_(A.GetComm()),
209  APtr_(&A),
210  optionsPtr_(NULL),
211  statPtr_(NULL),
212  ScalePermstructPtr_(NULL),
213  LUstructPtr_(NULL),
214  SOLVEstructPtr_(NULL),
215  gridPtr_(NULL),
216  berr_(NULL),
217  perm_r_(NULL),
218  nrhs_(1),
219  nprow_(0),
220  npcol_(0),
221  firstSolveWithThisA_(true),
222  gridInitialized_(false),
223  LUStructInitialized_(false)
224 {
225  height = A.Height();
226  width = A.Width();
227 
228  this->Init();
229 }
230 
232 {
233  superlu_dist_options_t * options = (superlu_dist_options_t*)optionsPtr_;
234  SuperLUStat_t * stat = (SuperLUStat_t*)statPtr_;
235  ScalePermstruct_t * SPstruct = (ScalePermstruct_t*)ScalePermstructPtr_;
236  LUstruct_t * LUstruct = (LUstruct_t*)LUstructPtr_;
237  SOLVEstruct_t * SOLVEstruct = (SOLVEstruct_t*)SOLVEstructPtr_;
238  gridinfo_t * grid = (gridinfo_t*)gridPtr_;
239 
240  SUPERLU_FREE(berr_);
241  PStatFree(stat);
242 
243  if ( LUStructInitialized_ )
244  {
245  ScalePermstructFree(SPstruct);
246  Destroy_LU(width, grid, LUstruct);
247  LUstructFree(LUstruct);
248  }
249 
250  if ( options->SolveInitialized )
251  {
252  dSolveFinalize(options, SOLVEstruct);
253  }
254 
255  if ( options != NULL ) { delete options; }
256  if ( stat != NULL ) { delete stat; }
257  if ( SPstruct != NULL ) { delete SPstruct; }
258  if ( LUstruct != NULL ) { delete LUstruct; }
259  if ( SOLVEstruct != NULL ) { delete SOLVEstruct; }
260  if ( grid != NULL ) { delete grid; }
261  if ( perm_r_ != NULL ) { SUPERLU_FREE(perm_r_); }
262 }
263 
264 void SuperLUSolver::Init()
265 {
266  MPI_Comm_size(comm_, &numProcs_);
267  MPI_Comm_rank(comm_, &myid_);
268 
269  optionsPtr_ = new superlu_dist_options_t;
270  statPtr_ = new SuperLUStat_t;
271  ScalePermstructPtr_ = new ScalePermstruct_t;
272  LUstructPtr_ = new LUstruct_t;
273  SOLVEstructPtr_ = new SOLVEstruct_t;
274  gridPtr_ = new gridinfo_t;
275 
276  superlu_dist_options_t * options = (superlu_dist_options_t*)optionsPtr_;
277  SuperLUStat_t * stat = (SuperLUStat_t*)statPtr_;
278 
279  if ( !(berr_ = doubleMalloc_dist(nrhs_)) )
280  {
281  ABORT("Malloc fails for berr[].");
282  }
283 
284  // Set default options
285  set_default_options_dist(options);
286 
287  options->ParSymbFact = YES;
288  options->ColPerm = NATURAL;
289 
290  // Choose nprow and npcol so that the process grid is as square as possible.
291  // If the processes cannot be divided evenly, keep the row dimension smaller
292  // than the column dimension.
293 
294  nprow_ = (int)superlu_internal::sqrti((unsigned int)numProcs_);
295  while (numProcs_ % nprow_ != 0 && nprow_ > 0)
296  {
297  nprow_--;
298  }
299 
300  npcol_ = (int)(numProcs_ / nprow_);
301  MFEM_ASSERT(nprow_ * npcol_ == numProcs_, "");
302 
303  PStatInit(stat); // Initialize the statistics variables.
304 }
305 
306 void SuperLUSolver::SetPrintStatistics( bool print_stat )
307 {
308  superlu_dist_options_t * options = (superlu_dist_options_t*)optionsPtr_;
309 
310  yes_no_t opt = print_stat?YES:NO;
311 
312  options->PrintStat = opt;
313 }
314 
316 {
317  superlu_dist_options_t * options = (superlu_dist_options_t*)optionsPtr_;
318 
319  yes_no_t opt = equil?YES:NO;
320 
321  options->Equil = opt;
322 }
323 
325 {
326  superlu_dist_options_t * options = (superlu_dist_options_t*)optionsPtr_;
327 
328  colperm_t opt = (colperm_t)col_perm;
329 
330  options->ColPerm = opt;
331 }
332 
334  Array<int> * perm )
335 {
336  superlu_dist_options_t * options = (superlu_dist_options_t*)optionsPtr_;
337 
338  rowperm_t opt = (rowperm_t)row_perm;
339 
340  options->RowPerm = opt;
341 
342  if ( opt == MY_PERMR )
343  {
344  if ( perm == NULL )
345  {
346  mfem_error("SuperLUSolver::SetRowPermutation :"
347  " permutation vector not set!");
348  }
349 
350  if ( !(perm_r_ = intMalloc_dist(perm->Size())) )
351  {
352  ABORT("Malloc fails for perm_r[].");
353  }
354  for (int i=0; i<perm->Size(); i++)
355  {
356  perm_r_[i] = (*perm)[i];
357  }
358  }
359 }
360 
362 {
363  superlu_dist_options_t * options = (superlu_dist_options_t*)optionsPtr_;
364 
365  trans_t opt = (trans_t)trans;
366 
367  options->Trans = opt;
368 }
369 
371 {
372  superlu_dist_options_t * options = (superlu_dist_options_t*)optionsPtr_;
373 
374  IterRefine_t opt = (IterRefine_t)iter_ref;
375 
376  options->IterRefine = opt;
377 }
378 
380 {
381  superlu_dist_options_t * options = (superlu_dist_options_t*)optionsPtr_;
382 
383  yes_no_t opt = rtp?YES:NO;
384 
385  options->ReplaceTinyPivot = opt;
386 }
387 
388 void SuperLUSolver::SetNumLookAheads( int num_lookaheads )
389 {
390  superlu_dist_options_t * options = (superlu_dist_options_t*)optionsPtr_;
391 
392  options->num_lookaheads = num_lookaheads;
393 }
394 
396 {
397  superlu_dist_options_t * options = (superlu_dist_options_t*)optionsPtr_;
398 
399  yes_no_t opt = etree?YES:NO;
400 
401  options->lookahead_etree = opt;
402 }
403 
405 {
406  superlu_dist_options_t * options = (superlu_dist_options_t*)optionsPtr_;
407 
408  yes_no_t opt = sym?YES:NO;
409 
410  options->SymPattern = opt;
411 }
412 
414 {
415  gridinfo_t * grid = (gridinfo_t*)gridPtr_;
416 
417  // Make sure the values of nprow and npcol are reasonable
418  if ( ((nprow_ * npcol_) > numProcs_) || ((nprow_ * npcol_) < 1) )
419  {
420  if ( myid_ == 0 )
421  {
422  mfem::err << "Warning: User specified nprow and npcol are such that "
423  << "(nprow * npcol) > numProcs or (nprow * npcol) < 1. "
424  << "Using default values for nprow and npcol instead."
425  << endl;
426  }
427 
428  nprow_ = (int)superlu_internal::sqrti((unsigned int)numProcs_);
429  while (numProcs_ % nprow_ != 0 && nprow_ > 0)
430  {
431  nprow_--;
432  }
433 
434  npcol_ = (int)(numProcs_ / nprow_);
435  MFEM_ASSERT(nprow_ * npcol_ == numProcs_, "");
436  }
437 
438  superlu_gridinit(comm_, nprow_, npcol_, grid);
439 
440  gridInitialized_ = true;
441 }
442 
444 {
445  if ( gridInitialized_ )
446  {
447  gridinfo_t * grid = (gridinfo_t*)gridPtr_;
448 
449  superlu_gridexit(grid);
450  }
451 
452  gridInitialized_ = false;
453 }
454 
455 void SuperLUSolver::Mult( const Vector & x, Vector & y ) const
456 {
457  MFEM_ASSERT(APtr_ != NULL,
458  "SuperLU Error: The operator must be set before"
459  " the system can be solved.");
460 
461  superlu_dist_options_t * options = (superlu_dist_options_t*)optionsPtr_;
462  SuperLUStat_t * stat = (SuperLUStat_t*)statPtr_;
463  SuperMatrix * A = (SuperMatrix*)APtr_->InternalData();
464 
465  ScalePermstruct_t * SPstruct = (ScalePermstruct_t*)ScalePermstructPtr_;
466  LUstruct_t * LUstruct = (LUstruct_t*)LUstructPtr_;
467  SOLVEstruct_t * SOLVEstruct = (SOLVEstruct_t*)SOLVEstructPtr_;
468  gridinfo_t * grid = (gridinfo_t*)gridPtr_;
469 
471  {
472  options->Fact = FACTORED; // Indicate the factored form of A is supplied.
473  }
474  else // This is the first solve with this A
475  {
476  firstSolveWithThisA_ = false;
477 
478  // Make sure that the parameters have been initialized The only parameter
479  // we might have to worry about is ScalePermstruct, if the user is
480  // supplying a row or column permutation.
481 
482  // Initialize ScalePermstruct and LUstruct.
483  SPstruct->DiagScale = NOEQUIL;
484 
485  // Transfer ownership of the row permutations if available
486  if ( perm_r_ != NULL )
487  {
488  SPstruct->perm_r = perm_r_;
489  perm_r_ = NULL;
490  }
491  else
492  {
493  if ( !(SPstruct->perm_r = intMalloc_dist(A->nrow)) )
494  {
495  ABORT("Malloc fails for perm_r[].");
496  }
497  }
498  if ( !(SPstruct->perm_c = intMalloc_dist(A->ncol)) )
499  {
500  ABORT("Malloc fails for perm_c[].");
501  }
502 
503  LUstructInit(A->ncol, LUstruct);
504  LUStructInitialized_ = true;
505  }
506 
507  // SuperLU overwrites x with y, so copy x to y and pass that to the solve
508  // routine.
509 
510  y = x;
511 
512  double* yPtr = (double*)y;
513  int info = -1, locSize = y.Size();
514 
515  // Solve the system
516  pdgssvx(options, A, SPstruct, yPtr, locSize, nrhs_, grid,
517  LUstruct, SOLVEstruct, berr_, stat, &info);
518 
519  if ( info != 0 )
520  {
521  if ( info <= A->ncol )
522  {
523  MFEM_ABORT("SuperLU: Found a singular matrix, U("
524  << info << "," << info << ") is exactly zero.");
525  }
526  else if ( info > A->ncol )
527  {
528  MFEM_ABORT("SuperLU: Memory allocation error with "
529  << info - A->ncol << " bytes already allocated,");
530  }
531  else
532  {
533  MFEM_ABORT("Unknown SuperLU Error");
534  }
535  }
536 }
537 
539 {
540  // Verify that we have a compatible operator
541  APtr_ = dynamic_cast<const SuperLURowLocMatrix*>(&op);
542  if ( APtr_ == NULL )
543  {
544  mfem_error("SuperLUSolver::SetOperator : not SuperLURowLocMatrix!");
545  }
546 
547  // Everything is OK so finish setting the operator
548  firstSolveWithThisA_ = true;
549 
550  // Set mfem::Operator member data
551  height = op.Height();
552  width = op.Width();
553 
554  // Initialize the processor grid if necessary
555  if (!gridInitialized_)
556  {
557  this->SetupGrid();
558  }
559 }
560 
561 } // mfem namespace
562 
563 #endif // MFEM_USE_MPI
564 #endif // MFEM_USE_SUPERLU
int Size() const
Logical size of the array.
Definition: array.hpp:124
void * InternalData() const
Definition: superlu.hpp:69
void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Definition: superlu.cpp:455
void * ScalePermstructPtr_
Definition: superlu.hpp:140
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
Definition: operator.hpp:63
void SetRowPermutation(superlu::RowPerm row_perm, Array< int > *perm=NULL)
Definition: superlu.cpp:333
int Size() const
Returns the size of the vector.
Definition: vector.hpp:157
void SetSymmetricPattern(bool sym)
Definition: superlu.cpp:404
void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition: superlu.cpp:538
void SetIterativeRefine(superlu::IterRefine iter_ref)
Definition: superlu.cpp:370
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition: operator.hpp:57
void SetColumnPermutation(superlu::ColPerm col_perm)
Definition: superlu.cpp:324
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:388
void SetPrintStatistics(bool print_stat)
Definition: superlu.cpp:306
void trans(const Vector &x, Vector &p)
Definition: toroid.cpp:239
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:142
void SetLookAheadElimTree(bool etree)
Definition: superlu.cpp:395
unsigned int sqrti(const unsigned int &a)
Definition: superlu.cpp:31
void SetReplaceTinyPivot(bool rtp)
Definition: superlu.cpp:379
Vector data type.
Definition: vector.hpp:48
void SetTranspose(superlu::Trans trans)
Definition: superlu.cpp:361
const SuperLURowLocMatrix * APtr_
Definition: superlu.hpp:130
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:58
Abstract operator.
Definition: operator.hpp:24
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:187
SuperLUSolver(MPI_Comm comm)
Definition: superlu.cpp:186
int width
Dimension of the input / number of columns in the matrix.
Definition: operator.hpp:28
void SetEquilibriate(bool equil)
Definition: superlu.cpp:315