MFEM  v4.6.0
Finite element discretization library
superlu.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2023, 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 header
20 #include "superlu_ddefs.h"
21 
22 #if XSDK_INDEX_SIZE == 64 && !(defined(HYPRE_BIGINT) || defined(HYPRE_MIXEDINT))
23 #error "Mismatch between HYPRE (32bit) and SuperLU (64bit) integer types"
24 #endif
25 #if XSDK_INDEX_SIZE == 32 && (defined(HYPRE_BIGINT) || defined(HYPRE_MIXEDINT))
26 #error "Mismatch between HYPRE (64bit) and SuperLU (32bit) integer types"
27 #endif
28 
29 #if SUPERLU_DIST_MAJOR_VERSION > 6 || \
30  (SUPERLU_DIST_MAJOR_VERSION == 6 && SUPERLU_DIST_MINOR_VERSION >= 3)
31 #define ScalePermstruct_t dScalePermstruct_t
32 #define LUstruct_t dLUstruct_t
33 #define SOLVEstruct_t dSOLVEstruct_t
34 #define ZeroLblocks dZeroLblocks
35 #define ZeroUblocks dZeroUblocks
36 #define Destroy_LU dDestroy_LU
37 #define SolveFinalize dSolveFinalize
38 #define ScalePermstructInit dScalePermstructInit
39 #define ScalePermstructFree dScalePermstructFree
40 #define LUstructFree dLUstructFree
41 #define LUstructInit dLUstructInit
42 #endif
43 
44 #if SUPERLU_DIST_MAJOR_VERSION > 7 || \
45  (SUPERLU_DIST_MAJOR_VERSION == 7 && SUPERLU_DIST_MINOR_VERSION >= 2)
46 #define DeAllocLlu_3d dDeAllocLlu_3d
47 #define DeAllocGlu_3d dDeAllocGlu_3d
48 #define Destroy_A3d_gathered_on_2d dDestroy_A3d_gathered_on_2d
49 #endif
50 
51 unsigned int sqrti(unsigned int a)
52 {
53  unsigned int rem = 0;
54  unsigned int root = 0;
55  unsigned short len = sizeof(int); len <<= 2;
56  unsigned short shift = (unsigned short)((len << 1) - 2);
57 
58  for (int i = 0; i < len; i++)
59  {
60  root <<= 1;
61  rem = ((rem << 2) + (a >> shift));
62  a <<= 2;
63  root ++;
64  if (root <= rem)
65  {
66  rem -= root;
67  root++;
68  }
69  else
70  {
71  root--;
72  }
73  }
74  return (root >> 1);
75 }
76 
77 int GetGridRows(MPI_Comm comm, int npdep)
78 {
79  int np;
80  MPI_Comm_size(comm, &np);
81  MFEM_VERIFY(npdep > 0 && np % npdep == 0 && !(npdep & (npdep - 1)),
82  "SuperLUSolver: 3D partition depth must be a power of two "
83  "and evenly divide the number of processors!");
84  int nr = (int)sqrti((unsigned int)(np / npdep));
85  while (np % nr != 0 && nr > 0)
86  {
87  nr--;
88  }
89  MFEM_VERIFY(nr > 0,
90  "SuperLUSolver: Unable to determine processor grid for np = " << np);
91  return nr;
92 }
93 
94 int GetGridCols(MPI_Comm comm, int npdep, int nr)
95 {
96  int np;
97  MPI_Comm_size(comm, &np);
98  int nc = np / (nr * npdep);
99  MFEM_VERIFY(nr * nc * npdep == np,
100  "SuperLUSolver: Impossible processor partition!");
101  return nc;
102 }
103 
104 namespace mfem
105 {
106 
108  int num_loc_rows,
109  HYPRE_BigInt first_loc_row,
110  HYPRE_BigInt glob_nrows,
111  HYPRE_BigInt glob_ncols,
112  int *I, HYPRE_BigInt *J,
113  double *data)
114  : comm_(comm)
115 {
116  // Set mfem::Operator member data
117  height = num_loc_rows;
118  width = num_loc_rows;
119 
120  // Allocate SuperLU's SuperMatrix struct
121  rowLocPtr_ = new SuperMatrix;
122  SuperMatrix *A = (SuperMatrix *)rowLocPtr_;
123  A->Store = NULL;
124 
125  int_t m = glob_nrows;
126  int_t n = glob_ncols;
127  int_t nnz_loc = I[num_loc_rows];
128  int_t m_loc = num_loc_rows;
129  int_t fst_row = first_loc_row;
130 
131  double *nzval = NULL;
132  int_t *colind = NULL;
133  int_t *rowptr = NULL;
134 
135  if (!(nzval = doubleMalloc_dist(nnz_loc)))
136  {
137  MFEM_ABORT("SuperLURowLocMatrix: Malloc failed for nzval!");
138  }
139  for (int_t i = 0; i < nnz_loc; i++)
140  {
141  nzval[i] = data[i];
142  }
143 
144  if (!(colind = intMalloc_dist(nnz_loc)))
145  {
146  MFEM_ABORT("SuperLURowLocMatrix: Malloc failed for colind!")
147  }
148  for (int_t i = 0; i < nnz_loc; i++)
149  {
150  colind[i] = J[i];
151  }
152 
153  if (!(rowptr = intMalloc_dist(m_loc+1)))
154  {
155  MFEM_ABORT("SuperLURowLocMatrix: Malloc failed for rowptr!")
156  }
157  for (int_t i = 0; i <= m_loc; i++)
158  {
159  rowptr[i] = I[i];
160  }
161 
162  // Assign the matrix data to SuperLU's SuperMatrix structure
163  dCreate_CompRowLoc_Matrix_dist(A, m, n, nnz_loc, m_loc, fst_row,
164  nzval, colind, rowptr,
165  SLU_NR_loc, SLU_D, SLU_GE);
166 
167  // Save global number of rows and columns of the matrix
168  num_global_rows_ = m;
169  num_global_cols_ = n;
170 }
171 
173 {
174  const HypreParMatrix *APtr = dynamic_cast<const HypreParMatrix *>(&op);
175  MFEM_VERIFY(APtr, "Not a compatible matrix type");
176  comm_ = APtr->GetComm();
177 
178  // Set mfem::Operator member data
179  height = op.Height();
180  width = op.Width();
181 
182  // Allocate SuperLU's SuperMatrix struct
183  rowLocPtr_ = new SuperMatrix;
184  SuperMatrix *A = (SuperMatrix *)rowLocPtr_;
185  A->Store = NULL;
186 
187  // First cast the parameter to a hypre_ParCSRMatrix
188  hypre_ParCSRMatrix *parcsr_op =
189  (hypre_ParCSRMatrix *)const_cast<HypreParMatrix &>(*APtr);
190 
191  // Create the SuperMatrix A by taking the internal data from a
192  // hypre_CSRMatrix
193  APtr->HostRead();
194  hypre_CSRMatrix *csr_op = hypre_MergeDiagAndOffd(parcsr_op);
195  APtr->HypreRead();
196  HYPRE_Int *Iptr = csr_op->i;
197 #if MFEM_HYPRE_VERSION >= 21600
198  HYPRE_BigInt *Jptr = csr_op->big_j;
199 #else
200  HYPRE_Int *Jptr = csr_op->j;
201 #endif
202  int_t m = parcsr_op->global_num_rows;
203  int_t n = parcsr_op->global_num_cols;
204  int_t fst_row = parcsr_op->first_row_index;
205  int_t nnz_loc = csr_op->num_nonzeros;
206  int_t m_loc = csr_op->num_rows;
207 
208  // We copy the data from the hypre_CSRMatrix because SuperLU_DIST will
209  // free the memory assuming it has been allocated with its *Malloc_dist
210  // wrappers
211  double *nzval = NULL;
212  int_t *colind = NULL;
213  int_t *rowptr = NULL;
214 
215  if (!(nzval = doubleMalloc_dist(nnz_loc)))
216  {
217  MFEM_ABORT("SuperLURowLocMatrix: Malloc failed for nzval!");
218  }
219  for (int_t i = 0; i < nnz_loc; i++)
220  {
221  nzval[i] = csr_op->data[i];
222  }
223 
224  if (!(colind = intMalloc_dist(nnz_loc)))
225  {
226  MFEM_ABORT("SuperLURowLocMatrix: Malloc failed for colind!")
227  }
228  for (int_t i = 0; i < nnz_loc; i++)
229  {
230  colind[i] = Jptr[i];
231  }
232 
233  if (!(rowptr = intMalloc_dist(m_loc+1)))
234  {
235  MFEM_ABORT("SuperLURowLocMatrix: Malloc failed for rowptr!")
236  }
237  for (int_t i = 0; i <= m_loc; i++)
238  {
239  rowptr[i] = Iptr[i];
240  }
241 
242  // Assign the matrix data to SuperLU's SuperMatrix structure
243  dCreate_CompRowLoc_Matrix_dist(A, m, n, nnz_loc, m_loc, fst_row,
244  nzval, colind, rowptr,
245  SLU_NR_loc, SLU_D, SLU_GE);
246 
247  // Everything has been copied so delete the structure
248  hypre_CSRMatrixDestroy(csr_op);
249 
250  // Save global number of rows and columns of the matrix
251  num_global_rows_ = m;
252  num_global_cols_ = n;
253 }
254 
256 {
257  SuperMatrix *A = (SuperMatrix *)rowLocPtr_;
258  Destroy_CompRowLoc_Matrix_dist(A);
259  delete A;
260 }
261 
262 SuperLUSolver::SuperLUSolver(MPI_Comm comm, int npdep)
263  : nprow_(GetGridRows(comm, npdep)),
264  npcol_(GetGridCols(comm, npdep, nprow_)),
265  npdep_(npdep),
266  APtr_(NULL),
267  nrhs_(0)
268 {
269  Init(comm);
270 }
271 
273  : SuperLUSolver(A.GetComm(), npdep)
274 {
275  SetOperator(A);
276 }
277 
279 {
280  superlu_dist_options_t *options = (superlu_dist_options_t *)optionsPtr_;
281 
282  ScalePermstruct_t *ScalePermstruct = (ScalePermstruct_t *)ScalePermstructPtr_;
283  LUstruct_t *LUstruct = (LUstruct_t *)LUstructPtr_;
284  SOLVEstruct_t *SOLVEstruct = (SOLVEstruct_t *)SOLVEstructPtr_;
285 
286 #if SUPERLU_DIST_MAJOR_VERSION > 7 || \
287  (SUPERLU_DIST_MAJOR_VERSION == 7 && SUPERLU_DIST_MINOR_VERSION >= 2)
288  if (npdep_ > 1)
289  {
290  gridinfo3d_t *grid3d = (gridinfo3d_t *)gridPtr_;
291 
292  if (APtr_)
293  {
294  if (grid3d->zscp.Iam == 0)
295  {
296  // Process layer 0
297  Destroy_LU(APtr_->GetGlobalNumColumns(), &(grid3d->grid2d),
298  LUstruct);
299  SolveFinalize(options, SOLVEstruct);
300  }
301  else
302  {
303  // Process layers not equal 0
304  DeAllocLlu_3d(APtr_->GetGlobalNumColumns(), LUstruct, grid3d);
305  DeAllocGlu_3d(LUstruct);
306  }
307  Destroy_A3d_gathered_on_2d(SOLVEstruct, grid3d);
308  ScalePermstructFree(ScalePermstruct);
309  LUstructFree(LUstruct);
310  }
311 
312  superlu_gridexit3d(grid3d);
313  delete grid3d;
314  }
315  else
316 #endif
317  {
318  gridinfo_t *grid = (gridinfo_t *)gridPtr_;
319 
320  if (APtr_)
321  {
322  Destroy_LU(APtr_->GetGlobalNumColumns(), grid, LUstruct);
323  SolveFinalize(options, SOLVEstruct);
324  ScalePermstructFree(ScalePermstruct);
325  LUstructFree(LUstruct);
326  }
327 
328  superlu_gridexit(grid);
329  delete grid;
330  }
331 
332  delete options;
333  delete ScalePermstruct;
334  delete LUstruct;
335  delete SOLVEstruct;
336 }
337 
338 void SuperLUSolver::Init(MPI_Comm comm)
339 {
340  optionsPtr_ = new superlu_dist_options_t;
341  ScalePermstructPtr_ = new ScalePermstruct_t;
342  LUstructPtr_ = new LUstruct_t;
343  SOLVEstructPtr_ = new SOLVEstruct_t;
344 
345  // Initialize process grid
346 #if SUPERLU_DIST_MAJOR_VERSION > 7 || \
347  (SUPERLU_DIST_MAJOR_VERSION == 7 && SUPERLU_DIST_MINOR_VERSION >= 2)
348  if (npdep_ > 1)
349  {
350  gridPtr_ = new gridinfo3d_t;
351  superlu_gridinit3d(comm, nprow_, npcol_, npdep_, (gridinfo3d_t *)gridPtr_);
352  }
353  else
354 #endif
355  {
356  gridPtr_ = new gridinfo_t;
357  MFEM_VERIFY(npdep_ == 1,
358  "SuperLUSolver: 3D partitioning is only available for "
359  "SuperLU_DIST version >= 7.2.0!");
360  superlu_gridinit(comm, nprow_, npcol_, (gridinfo_t *)gridPtr_);
361  }
362 
363  // Set default options:
364  // options.Fact = DOFACT;
365  // options.Equil = YES;
366  // options.ColPerm = METIS_AT_PLUS_A;
367  // options.RowPerm = LargeDiag_MC64;
368  // options.ReplaceTinyPivot = NO;
369  // options.Trans = NOTRANS;
370  // options.IterRefine = SLU_DOUBLE;
371  // options.SolveInitialized = NO;
372  // options.RefineInitialized = NO;
373  // options.PrintStat = YES;
374  superlu_dist_options_t *options = (superlu_dist_options_t *)optionsPtr_;
375  set_default_options_dist(options);
376 #if SUPERLU_DIST_MAJOR_VERSION > 7 || \
377  (SUPERLU_DIST_MAJOR_VERSION == 7 && SUPERLU_DIST_MINOR_VERSION >= 2)
378  if (npdep_ > 1)
379  {
380  options->Algo3d = YES;
381  }
382 #endif
383 }
384 
386 {
387  superlu_dist_options_t *options = (superlu_dist_options_t *)optionsPtr_;
388  yes_no_t opt = print_stat ? YES : NO;
389  options->PrintStat = opt;
390 }
391 
393 {
394  superlu_dist_options_t *options = (superlu_dist_options_t *)optionsPtr_;
395  yes_no_t opt = equil ? YES : NO;
396  options->Equil = opt;
397 }
398 
400 {
401  superlu_dist_options_t *options = (superlu_dist_options_t *)optionsPtr_;
402  colperm_t opt = (colperm_t)col_perm;
403  if (opt == MY_PERMC)
404  {
405  MFEM_ABORT("SuperLUSolver::SetColumnPermutation does not yet support "
406  "MY_PERMC!");
407  }
408  else if (opt == PARMETIS)
409  {
410  options->ParSymbFact = YES;
411  }
412  options->ColPerm = opt;
413 }
414 
416 {
417  superlu_dist_options_t *options = (superlu_dist_options_t *)optionsPtr_;
418  rowperm_t opt = (rowperm_t)row_perm;
419  if (opt == MY_PERMR)
420  {
421  MFEM_ABORT("SuperLUSolver::SetRowPermutation does not yet support "
422  "MY_PERMR!");
423  }
424  options->RowPerm = opt;
425 }
426 
428 {
429  superlu_dist_options_t *options = (superlu_dist_options_t *)optionsPtr_;
430  IterRefine_t opt = (IterRefine_t)iter_ref;
431  options->IterRefine = opt;
432 }
433 
435 {
436  superlu_dist_options_t *options = (superlu_dist_options_t *)optionsPtr_;
437  yes_no_t opt = rtp ? YES : NO;
438  options->ReplaceTinyPivot = opt;
439 }
440 
441 void SuperLUSolver::SetNumLookAheads(int num_lookaheads)
442 {
443  superlu_dist_options_t *options = (superlu_dist_options_t *)optionsPtr_;
444  options->num_lookaheads = num_lookaheads;
445 }
446 
448 {
449  superlu_dist_options_t *options = (superlu_dist_options_t *)optionsPtr_;
450  yes_no_t opt = etree ? YES : NO;
451  options->lookahead_etree = opt;
452 }
453 
455 {
456  superlu_dist_options_t *options = (superlu_dist_options_t *)optionsPtr_;
457  yes_no_t opt = sym ? YES : NO;
458  options->SymPattern = opt;
459 }
460 
462 {
463  superlu_dist_options_t *options = (superlu_dist_options_t *)optionsPtr_;
464  yes_no_t opt = par ? YES : NO;
465  options->ParSymbFact = opt;
466 }
467 
469 {
470  superlu_dist_options_t *options = (superlu_dist_options_t *)optionsPtr_;
471  fact_t opt = (fact_t)fact;
472  options->Fact = opt;
473 }
474 
476 {
477  // Verify that we have a compatible operator
478  bool LUStructInitialized = (APtr_ != NULL);
479  APtr_ = dynamic_cast<const SuperLURowLocMatrix *>(&op);
480  MFEM_VERIFY(APtr_, "SuperLUSolver::SetOperator: Not a SuperLURowLocMatrix!");
481 
482  superlu_dist_options_t *options = (superlu_dist_options_t *)optionsPtr_;
483 
484  ScalePermstruct_t *ScalePermstruct = (ScalePermstruct_t *)ScalePermstructPtr_;
485  LUstruct_t *LUstruct = (LUstruct_t *)LUstructPtr_;
486 
487  gridinfo_t *grid;
488 #if SUPERLU_DIST_MAJOR_VERSION > 7 || \
489  (SUPERLU_DIST_MAJOR_VERSION == 7 && SUPERLU_DIST_MINOR_VERSION >= 2)
490  gridinfo3d_t *grid3d = NULL;
491  if (npdep_ > 1)
492  {
493  grid3d = (gridinfo3d_t *)gridPtr_;
494  grid = NULL;
495  }
496  else
497 #endif
498  {
499  grid = (gridinfo_t *)gridPtr_;
500  }
501 
502  // Set mfem::Operator member data
503  MFEM_VERIFY(!LUStructInitialized ||
504  (height == op.Height() && width == op.Width()),
505  "SuperLUSolver::SetOperator: Inconsistent new matrix size!");
506  height = op.Height();
507  width = op.Width();
508 
509  if (!LUStructInitialized)
510  {
511  // Initialize ScalePermstruct and LUstruct once for all operators (must
512  // have same dimensions)
513  ScalePermstructInit(APtr_->GetGlobalNumRows(),
514  APtr_->GetGlobalNumColumns(), ScalePermstruct);
515  LUstructInit(APtr_->GetGlobalNumColumns(), LUstruct);
516  options->Fact = DOFACT;
517  }
518  else
519  {
520  // A previous matrix has already been set and factored
521  switch (options->Fact)
522  {
523  case DOFACT:
524  MFEM_ABORT("SuperLUSolver::SetOperator: Previous matrix was never used!");
525  break;
527  {
528  // Just zero the LU factors
529 #if SUPERLU_DIST_MAJOR_VERSION > 7 || \
530 (SUPERLU_DIST_MAJOR_VERSION == 7 && SUPERLU_DIST_MINOR_VERSION >= 2)
531  if (npdep_ > 1)
532  {
533  if (grid3d->zscp.Iam == 0)
534  {
535  ZeroLblocks(grid3d->iam, APtr_->GetGlobalNumColumns(),
536  &(grid3d->grid2d), LUstruct);
537  ZeroUblocks(grid3d->iam, APtr_->GetGlobalNumColumns(),
538  &(grid3d->grid2d), LUstruct);
539  }
540  }
541  else
542 #endif
543  {
544  ZeroLblocks(grid->iam, APtr_->GetGlobalNumColumns(),
545  grid, LUstruct);
546  ZeroUblocks(grid->iam, APtr_->GetGlobalNumColumns(),
547  grid, LUstruct);
548  }
549  }
550  break;
551  case SamePattern:
552  case FACTORED:
553  {
554  // Delete factors from the prior factorization
555 #if SUPERLU_DIST_MAJOR_VERSION > 7 || \
556 (SUPERLU_DIST_MAJOR_VERSION == 7 && SUPERLU_DIST_MINOR_VERSION >= 2)
557  if (npdep_ > 1)
558  {
559  if (grid3d->zscp.Iam == 0)
560  {
561  Destroy_LU(APtr_->GetGlobalNumColumns(), &(grid3d->grid2d),
562  LUstruct);
563  }
564  else
565  {
566  DeAllocLlu_3d(APtr_->GetGlobalNumColumns(), LUstruct,
567  grid3d);
568  DeAllocGlu_3d(LUstruct);
569  }
570  }
571  else
572 #endif
573  {
574  Destroy_LU(APtr_->GetGlobalNumColumns(), grid, LUstruct);
575  }
576  }
577  break;
578  default:
579  MFEM_ABORT("SuperLUSolver::SetOperator: Unexpected value for "
580  "options->Fact!");
581  break;
582  }
583  if (options->Fact == FACTORED) { options->Fact = DOFACT; }
584  }
585 }
586 
587 void SuperLUSolver::Mult(const Vector &x, Vector &y) const
588 {
590  Array<Vector *> Y(1);
591  X[0] = &x;
592  Y[0] = &y;
593  ArrayMult(X, Y);
594 }
595 
597  Array<Vector *> &Y) const
598 {
599  MFEM_ASSERT(APtr_ != NULL,
600  "SuperLU Error: The operator must be set before"
601  " the system can be solved.");
602  SuperMatrix *A = (SuperMatrix *)APtr_->InternalData();
603  superlu_dist_options_t *options = (superlu_dist_options_t *)optionsPtr_;
604 
605  ScalePermstruct_t *ScalePermstruct = (ScalePermstruct_t *)ScalePermstructPtr_;
606  LUstruct_t *LUstruct = (LUstruct_t *)LUstructPtr_;
607  SOLVEstruct_t *SOLVEstruct = (SOLVEstruct_t *)SOLVEstructPtr_;
608 
609  gridinfo_t *grid;
610 #if SUPERLU_DIST_MAJOR_VERSION > 7 || \
611  (SUPERLU_DIST_MAJOR_VERSION == 7 && SUPERLU_DIST_MINOR_VERSION >= 2)
612  gridinfo3d_t *grid3d = NULL;
613  if (npdep_ > 1)
614  {
615  grid3d = (gridinfo3d_t *)gridPtr_;
616  grid = NULL;
617  }
618  else
619 #endif
620  {
621  grid = (gridinfo_t *)gridPtr_;
622  }
623 
624  // SuperLU overwrites x with y, so copy x to y and pass that to the solve
625  // routine. Due to issues with repeated solves and changes in the number
626  // of RHS vectors, this is not supported.
627  MFEM_ASSERT(X.Size() == Y.Size(),
628  "Number of columns mismatch in SuperLUSolver::Mult!");
629  MFEM_VERIFY(nrhs_ < 1 || nrhs_ == X.Size(),
630  "SuperLUSolver does not support multiple solves with different "
631  "numbers of RHS vectors!");
632  int ldx = Height();
633  if (X.Size() == 1)
634  {
635  MFEM_ASSERT(X[0] && Y[0], "Missing Vector in SuperLUSolver::Mult!");
636  sol_.MakeRef(*Y[0], 0, Y[0]->Size());
637  sol_ = *X[0];
638  nrhs_ = 1;
639  }
640  else
641  {
642  if (nrhs_ < 1)
643  {
644  MFEM_ASSERT(X[0], "Missing Vector in SuperLUSolver::Mult!");
645  sol_.SetSize(X.Size() * ldx, *X[0]);
646  nrhs_ = X.Size();
647  }
648  for (int i = 0; i < nrhs_; i++)
649  {
650  MFEM_ASSERT(X[i], "Missing Vector in SuperLUSolver::Mult!");
651  Vector s(sol_, i * ldx, ldx);
652  s = *X[i];
653  }
654  }
655 
656  // Solve the system
657  double *B = sol_.HostReadWrite(), *berr;
658  if (!(berr = doubleMalloc_dist(nrhs_)))
659  {
660  MFEM_ABORT("SuperLUSolver::Mult: Malloc failed for berr!");
661  }
662  SuperLUStat_t stat;
663  PStatInit(&stat);
664  int info = -1;
665 #if SUPERLU_DIST_MAJOR_VERSION > 7 || \
666  (SUPERLU_DIST_MAJOR_VERSION == 7 && SUPERLU_DIST_MINOR_VERSION >= 2)
667  if (npdep_ > 1)
668  {
669  pdgssvx3d(options, A, ScalePermstruct, B, ldx, nrhs_,
670  grid3d, LUstruct, SOLVEstruct, berr, &stat, &info);
671  }
672  else
673 #endif
674  {
675  pdgssvx(options, A, ScalePermstruct, B, ldx, nrhs_,
676  grid, LUstruct, SOLVEstruct, berr, &stat, &info);
677  }
678  HandleError(info);
679  SUPERLU_FREE(berr);
680  PStatFree(&stat);
681  options->Fact = FACTORED;
682 
683  // Copy solution into output (no need to do anything for single RHS since
684  // solution is written directly into output Vector)
685  if (nrhs_ == 1)
686  {
687  sol_.SyncAliasMemory(*Y[0]);
688  }
689  else
690  {
691  for (int i = 0; i < nrhs_; i++)
692  {
693  MFEM_ASSERT(Y[i], "Missing Vector in SuperLUSolver::Mult!");
694  Vector s(sol_, i * ldx, ldx);
695  *Y[i] = s;
696  }
697  }
698 }
699 
701 {
702  // Set flag for transpose solve
703  superlu_dist_options_t *options = (superlu_dist_options_t *)optionsPtr_;
704  options->Trans = TRANS;
705  Mult(x, y);
706 
707  // Reset the flag
708  options->Trans = NOTRANS;
709 }
710 
712  Array<Vector *> &Y) const
713 {
714  // Set flag for transpose solve
715  superlu_dist_options_t *options = (superlu_dist_options_t *)optionsPtr_;
716  options->Trans = TRANS;
717  ArrayMult(X, Y);
718 
719  // Reset the flag
720  options->Trans = NOTRANS;
721 }
722 
723 void SuperLUSolver::HandleError(int info) const
724 {
725  if (info != 0)
726  {
727  SuperMatrix *A = (SuperMatrix *)APtr_->InternalData();
728  if (info < 0)
729  {
730  switch (-info)
731  {
732  case 1:
733  MFEM_ABORT("SuperLUSolver: SuperLU options are invalid!");
734  break;
735  case 2:
736  MFEM_ABORT("SuperLUSolver: Matrix A (in Ax=b) is invalid!");
737  break;
738  case 5:
739  MFEM_ABORT("SuperLUSolver: Vector b dimension (in Ax=b) is "
740  "invalid!");
741  break;
742  case 6:
743  MFEM_ABORT("SuperLUSolver: Number of right-hand sides is "
744  "invalid!");
745  break;
746  default:
747  MFEM_ABORT("SuperLUSolver: Parameter with index "
748  << -info << "invalid (1-indexed)!");
749  break;
750  }
751  }
752  else if (info <= A->ncol)
753  {
754  MFEM_ABORT("SuperLUSolver: Found a singular matrix, U("
755  << info << "," << info << ") is exactly zero!");
756  }
757  else if (info > A->ncol)
758  {
759  MFEM_ABORT("SuperLUSolver: Memory allocation error with "
760  << info - A->ncol << " bytes already allocated!");
761  }
762  else
763  {
764  MFEM_ABORT("Unknown SuperLU error: info = " << info << "!");
765  }
766  }
767 }
768 
769 } // namespace mfem
770 
771 #endif // MFEM_USE_MPI
772 #endif // MFEM_USE_SUPERLU
void SetRowPermutation(superlu::RowPerm row_perm)
Definition: superlu.cpp:415
void SetFact(superlu::Fact fact)
Definition: superlu.cpp:468
int GetGridRows(MPI_Comm comm, int npdep)
Definition: superlu.cpp:77
void ArrayMultTranspose(const Array< const Vector *> &X, Array< Vector *> &Y) const
Action of the transpose operator on a matrix: Y=A^t(X).
Definition: superlu.cpp:711
void * ScalePermstructPtr_
Definition: superlu.hpp:151
const int nprow_
Definition: superlu.hpp:130
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:517
HYPRE_BigInt GetGlobalNumRows() const
Definition: superlu.hpp:72
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
Definition: operator.hpp:72
void HostRead() const
Update the internal hypre_ParCSRMatrix object, A, to be on host.
Definition: hypre.hpp:837
void SetSymmetricPattern(bool sym)
Definition: superlu.cpp:454
HYPRE_BigInt GetGlobalNumColumns() const
Definition: superlu.hpp:74
void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Definition: superlu.cpp:587
void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition: superlu.cpp:475
void * InternalData() const
Definition: superlu.hpp:68
void MultTranspose(const Vector &x, Vector &y) const
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
Definition: superlu.cpp:700
unsigned int sqrti(unsigned int a)
Definition: superlu.cpp:51
void SetIterativeRefine(superlu::IterRefine iter_ref)
Definition: superlu.cpp:427
void SetColumnPermutation(superlu::ColPerm col_perm)
Definition: superlu.cpp:399
const int npdep_
Definition: superlu.hpp:130
void SetNumLookAheads(int num_lookaheads)
Definition: superlu.cpp:441
void SetPrintStatistics(bool print_stat)
Definition: superlu.cpp:385
MPI_Comm GetComm() const
MPI communicator.
Definition: hypre.hpp:534
HYPRE_Int HYPRE_BigInt
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition: operator.hpp:66
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:153
void SetLookAheadElimTree(bool etree)
Definition: superlu.cpp:447
void SyncAliasMemory(const Vector &v) const
Update the alias memory location of the vector to match v.
Definition: vector.hpp:238
void ArrayMult(const Array< const Vector *> &X, Array< Vector *> &Y) const
Operator application on a matrix: Y=A(X).
Definition: superlu.cpp:596
void HypreRead() const
Update the internal hypre_ParCSRMatrix object, A, to be in hypre memory space.
Definition: hypre.hpp:854
void SetReplaceTinyPivot(bool rtp)
Definition: superlu.cpp:434
int Size() const
Return the logical size of the array.
Definition: array.hpp:141
Vector data type.
Definition: vector.hpp:58
SuperLURowLocMatrix(MPI_Comm comm, int num_loc_rows, HYPRE_BigInt first_loc_row, HYPRE_BigInt glob_nrows, HYPRE_BigInt glob_ncols, int *I, HYPRE_BigInt *J, double *data)
Definition: superlu.cpp:107
void MakeRef(Vector &base, int offset, int size)
Reset the Vector to be a reference to a sub-vector of base.
Definition: vector.hpp:581
const SuperLURowLocMatrix * APtr_
Definition: superlu.hpp:140
void SetParSymbFact(bool par)
Definition: superlu.cpp:461
RefCoord s[3]
int GetGridCols(MPI_Comm comm, int npdep, int nr)
Definition: superlu.cpp:94
Abstract operator.
Definition: operator.hpp:24
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:343
const int npcol_
Definition: superlu.hpp:130
virtual double * HostReadWrite()
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), false).
Definition: vector.hpp:473
int width
Dimension of the input / number of columns in the matrix.
Definition: operator.hpp:28
SuperLUSolver(MPI_Comm comm, int npdep=1)
Definition: superlu.cpp:262
void SetEquilibriate(bool equil)
Definition: superlu.cpp:392