Actual source code: superlu.c

petsc-master 2016-08-24
Report Typos and Errors
  2: /*  --------------------------------------------------------------------

  4:      This file implements a subclass of the SeqAIJ matrix class that uses
  5:      the SuperLU sparse solver.
  6: */

  8: /*
  9:      Defines the data structure for the base matrix type (SeqAIJ)
 10: */
 11:  #include <../src/mat/impls/aij/seq/aij.h>

 13: /*
 14:      SuperLU include files
 15: */
 16: EXTERN_C_BEGIN
 17: #if defined(PETSC_USE_COMPLEX)
 18: #if defined(PETSC_USE_REAL_SINGLE)
 19: #include <slu_cdefs.h>
 20: #else
 21: #include <slu_zdefs.h>
 22: #endif
 23: #else
 24: #if defined(PETSC_USE_REAL_SINGLE)
 25: #include <slu_sdefs.h>
 26: #else
 27: #include <slu_ddefs.h>
 28: #endif
 29: #endif
 30: #include <slu_util.h>
 31: EXTERN_C_END

 33: /*
 34:      This is the data that defines the SuperLU factored matrix type
 35: */
 36: typedef struct {
 37:   SuperMatrix       A,L,U,B,X;
 38:   superlu_options_t options;
 39:   PetscInt          *perm_c; /* column permutation vector */
 40:   PetscInt          *perm_r; /* row permutations from partial pivoting */
 41:   PetscInt          *etree;
 42:   PetscReal         *R, *C;
 43:   char              equed[1];
 44:   PetscInt          lwork;
 45:   void              *work;
 46:   PetscReal         rpg, rcond;
 47:   mem_usage_t       mem_usage;
 48:   MatStructure      flg;
 49:   SuperLUStat_t     stat;
 50:   Mat               A_dup;
 51:   PetscScalar       *rhs_dup;
 52:   GlobalLU_t        Glu;

 54:   /* Flag to clean up (non-global) SuperLU objects during Destroy */
 55:   PetscBool CleanUpSuperLU;
 56: } Mat_SuperLU;

 58: /*
 59:     Utility function
 60: */
 63: static PetscErrorCode MatFactorInfo_SuperLU(Mat A,PetscViewer viewer)
 64: {
 65:   Mat_SuperLU       *lu= (Mat_SuperLU*)A->data;
 66:   PetscErrorCode    ierr;
 67:   superlu_options_t options;

 70:   options = lu->options;

 72:   PetscViewerASCIIPrintf(viewer,"SuperLU run parameters:\n");
 73:   PetscViewerASCIIPrintf(viewer,"  Equil: %s\n",(options.Equil != NO) ? "YES" : "NO");
 74:   PetscViewerASCIIPrintf(viewer,"  ColPerm: %D\n",options.ColPerm);
 75:   PetscViewerASCIIPrintf(viewer,"  IterRefine: %D\n",options.IterRefine);
 76:   PetscViewerASCIIPrintf(viewer,"  SymmetricMode: %s\n",(options.SymmetricMode != NO) ? "YES" : "NO");
 77:   PetscViewerASCIIPrintf(viewer,"  DiagPivotThresh: %g\n",options.DiagPivotThresh);
 78:   PetscViewerASCIIPrintf(viewer,"  PivotGrowth: %s\n",(options.PivotGrowth != NO) ? "YES" : "NO");
 79:   PetscViewerASCIIPrintf(viewer,"  ConditionNumber: %s\n",(options.ConditionNumber != NO) ? "YES" : "NO");
 80:   PetscViewerASCIIPrintf(viewer,"  RowPerm: %D\n",options.RowPerm);
 81:   PetscViewerASCIIPrintf(viewer,"  ReplaceTinyPivot: %s\n",(options.ReplaceTinyPivot != NO) ? "YES" : "NO");
 82:   PetscViewerASCIIPrintf(viewer,"  PrintStat: %s\n",(options.PrintStat != NO) ? "YES" : "NO");
 83:   PetscViewerASCIIPrintf(viewer,"  lwork: %D\n",lu->lwork);
 84:   if (A->factortype == MAT_FACTOR_ILU) {
 85:     PetscViewerASCIIPrintf(viewer,"  ILU_DropTol: %g\n",options.ILU_DropTol);
 86:     PetscViewerASCIIPrintf(viewer,"  ILU_FillTol: %g\n",options.ILU_FillTol);
 87:     PetscViewerASCIIPrintf(viewer,"  ILU_FillFactor: %g\n",options.ILU_FillFactor);
 88:     PetscViewerASCIIPrintf(viewer,"  ILU_DropRule: %D\n",options.ILU_DropRule);
 89:     PetscViewerASCIIPrintf(viewer,"  ILU_Norm: %D\n",options.ILU_Norm);
 90:     PetscViewerASCIIPrintf(viewer,"  ILU_MILU: %D\n",options.ILU_MILU);
 91:   }
 92:   return(0);
 93: }

 97: PetscErrorCode MatSolve_SuperLU_Private(Mat A,Vec b,Vec x)
 98: {
 99:   Mat_SuperLU       *lu = (Mat_SuperLU*)A->data;
100:   const PetscScalar *barray;
101:   PetscScalar       *xarray;
102:   PetscErrorCode    ierr;
103:   PetscInt          info,i,n;
104:   PetscReal         ferr,berr;
105:   static PetscBool  cite = PETSC_FALSE;

108:   if (lu->lwork == -1) return(0);
109:   PetscCitationsRegister("@article{superlu99,\n  author  = {James W. Demmel and Stanley C. Eisenstat and\n             John R. Gilbert and Xiaoye S. Li and Joseph W. H. Liu},\n  title = {A supernodal approach to sparse partial pivoting},\n  journal = {SIAM J. Matrix Analysis and Applications},\n  year = {1999},\n  volume  = {20},\n  number = {3},\n  pages = {720-755}\n}\n",&cite);

111:   VecGetLocalSize(x,&n);
112:   lu->B.ncol = 1;   /* Set the number of right-hand side */
113:   if (lu->options.Equil && !lu->rhs_dup) {
114:     /* superlu overwrites b when Equil is used, thus create rhs_dup to keep user's b unchanged */
115:     PetscMalloc1(n,&lu->rhs_dup);
116:   }
117:   if (lu->options.Equil) {
118:     /* Copy b into rsh_dup */
119:     VecGetArrayRead(b,&barray);
120:     PetscMemcpy(lu->rhs_dup,barray,n*sizeof(PetscScalar));
121:     VecRestoreArrayRead(b,&barray);
122:     barray = lu->rhs_dup;
123:   } else {
124:     VecGetArrayRead(b,&barray);
125:   }
126:   VecGetArray(x,&xarray);

128: #if defined(PETSC_USE_COMPLEX)
129: #if defined(PETSC_USE_REAL_SINGLE)
130:   ((DNformat*)lu->B.Store)->nzval = (singlecomplex*)barray;
131:   ((DNformat*)lu->X.Store)->nzval = (singlecomplex*)xarray;
132: #else
133:   ((DNformat*)lu->B.Store)->nzval = (doublecomplex*)barray;
134:   ((DNformat*)lu->X.Store)->nzval = (doublecomplex*)xarray;
135: #endif
136: #else
137:   ((DNformat*)lu->B.Store)->nzval = (void*)barray;
138:   ((DNformat*)lu->X.Store)->nzval = xarray;
139: #endif

141:   lu->options.Fact = FACTORED; /* Indicate the factored form of A is supplied. */
142:   if (A->factortype == MAT_FACTOR_LU) {
143: #if defined(PETSC_USE_COMPLEX)
144: #if defined(PETSC_USE_REAL_SINGLE)
145:     PetscStackCall("SuperLU:cgssvx",cgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
146:                                      &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
147:                                      &lu->Glu, &lu->mem_usage, &lu->stat, &info));
148: #else
149:     PetscStackCall("SuperLU:zgssvx",zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
150:                                      &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
151:                                      &lu->Glu, &lu->mem_usage, &lu->stat, &info));
152: #endif
153: #else
154: #if defined(PETSC_USE_REAL_SINGLE)
155:     PetscStackCall("SuperLU:sgssvx",sgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
156:                                      &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
157:                                      &lu->Glu, &lu->mem_usage, &lu->stat, &info));
158: #else
159:     PetscStackCall("SuperLU:dgssvx",dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
160:                                      &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
161:                                      &lu->Glu,&lu->mem_usage, &lu->stat, &info));
162: #endif
163: #endif
164:   } else if (A->factortype == MAT_FACTOR_ILU) {
165: #if defined(PETSC_USE_COMPLEX)
166: #if defined(PETSC_USE_REAL_SINGLE)
167:     PetscStackCall("SuperLU:cgsisx",cgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
168:                                      &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond,
169:                                      &lu->Glu, &lu->mem_usage, &lu->stat, &info));
170: #else
171:     PetscStackCall("SuperLU:zgsisx",zgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
172:                                      &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond,
173:                                      &lu->Glu, &lu->mem_usage, &lu->stat, &info));
174: #endif
175: #else
176: #if defined(PETSC_USE_REAL_SINGLE)
177:     PetscStackCall("SuperLU:sgsisx",sgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
178:                                      &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond,
179:                                      &lu->Glu, &lu->mem_usage, &lu->stat, &info));
180: #else
181:     PetscStackCall("SuperLU:dgsisx",dgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
182:                                      &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond,
183:                                      &lu->Glu, &lu->mem_usage, &lu->stat, &info));
184: #endif
185: #endif
186:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported");
187:   if (!lu->options.Equil) {
188:     VecRestoreArrayRead(b,&barray);
189:   }
190:   VecRestoreArray(x,&xarray);

192:   if (!info || info == lu->A.ncol+1) {
193:     if (lu->options.IterRefine) {
194:       PetscPrintf(PETSC_COMM_SELF,"Iterative Refinement:\n");
195:       PetscPrintf(PETSC_COMM_SELF,"  %8s%8s%16s%16s\n", "rhs", "Steps", "FERR", "BERR");
196:       for (i = 0; i < 1; ++i) {
197:         PetscPrintf(PETSC_COMM_SELF,"  %8d%8d%16e%16e\n", i+1, lu->stat.RefineSteps, ferr, berr);
198:       }
199:     }
200:   } else if (info > 0) {
201:     if (lu->lwork == -1) {
202:       PetscPrintf(PETSC_COMM_SELF,"  ** Estimated memory: %D bytes\n", info - lu->A.ncol);
203:     } else {
204:       PetscPrintf(PETSC_COMM_SELF,"  Warning: gssvx() returns info %D\n",info);
205:     }
206:   } else if (info < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_LIB, "info = %D, the %D-th argument in gssvx() had an illegal value", info,-info);

208:   if (lu->options.PrintStat) {
209:     PetscPrintf(PETSC_COMM_SELF,"MatSolve__SuperLU():\n");
210:     PetscStackCall("SuperLU:StatPrint",StatPrint(&lu->stat));
211:   }
212:   return(0);
213: }

217: PetscErrorCode MatSolve_SuperLU(Mat A,Vec b,Vec x)
218: {
219:   Mat_SuperLU    *lu = (Mat_SuperLU*)A->data;

223:   if (A->errortype) {
224:     PetscInfo(A,"MatSolve is called with singular matrix factor, skip\n");
225:     VecSetInf(x);
226:     return(0);
227:   }

229:   lu->options.Trans = TRANS;
230:   MatSolve_SuperLU_Private(A,b,x);
231:   return(0);
232: }

236: PetscErrorCode MatSolveTranspose_SuperLU(Mat A,Vec b,Vec x)
237: {
238:   Mat_SuperLU    *lu = (Mat_SuperLU*)A->data;

242:   if (A->errortype) {
243:     PetscInfo(A,"MatSolve is called with singular matrix factor, skip\n");
244:     VecSetInf(x);
245:     return(0);
246:   }

248:   lu->options.Trans = NOTRANS;
249:   MatSolve_SuperLU_Private(A,b,x);
250:   return(0);
251: }

255: static PetscErrorCode MatLUFactorNumeric_SuperLU(Mat F,Mat A,const MatFactorInfo *info)
256: {
257:   Mat_SuperLU    *lu = (Mat_SuperLU*)F->data;
258:   Mat_SeqAIJ     *aa;
260:   PetscInt       sinfo;
261:   PetscReal      ferr, berr;
262:   NCformat       *Ustore;
263:   SCformat       *Lstore;

266:   if (lu->flg == SAME_NONZERO_PATTERN) { /* successing numerical factorization */
267:     lu->options.Fact = SamePattern;
268:     /* Ref: ~SuperLU_3.0/EXAMPLE/dlinsolx2.c */
269:     Destroy_SuperMatrix_Store(&lu->A);
270:     if (lu->options.Equil) {
271:       MatCopy_SeqAIJ(A,lu->A_dup,SAME_NONZERO_PATTERN);
272:     }
273:     if (lu->lwork >= 0) {
274:       PetscStackCall("SuperLU:Destroy_SuperNode_Matrix",Destroy_SuperNode_Matrix(&lu->L));
275:       PetscStackCall("SuperLU:Destroy_CompCol_Matrix",Destroy_CompCol_Matrix(&lu->U));
276:       lu->options.Fact = SamePattern;
277:     }
278:   }

280:   /* Create the SuperMatrix for lu->A=A^T:
281:        Since SuperLU likes column-oriented matrices,we pass it the transpose,
282:        and then solve A^T X = B in MatSolve(). */
283:   if (lu->options.Equil) {
284:     aa = (Mat_SeqAIJ*)(lu->A_dup)->data;
285:   } else {
286:     aa = (Mat_SeqAIJ*)(A)->data;
287:   }
288: #if defined(PETSC_USE_COMPLEX)
289: #if defined(PETSC_USE_REAL_SINGLE)
290:   PetscStackCall("SuperLU:cCreate_CompCol_Matrix",cCreate_CompCol_Matrix(&lu->A,A->cmap->n,A->rmap->n,aa->nz,(singlecomplex*)aa->a,aa->j,aa->i,SLU_NC,SLU_C,SLU_GE));
291: #else
292:   PetscStackCall("SuperLU:zCreate_CompCol_Matrix",zCreate_CompCol_Matrix(&lu->A,A->cmap->n,A->rmap->n,aa->nz,(doublecomplex*)aa->a,aa->j,aa->i,SLU_NC,SLU_Z,SLU_GE));
293: #endif
294: #else
295: #if defined(PETSC_USE_REAL_SINGLE)
296:   PetscStackCall("SuperLU:sCreate_CompCol_Matrix",sCreate_CompCol_Matrix(&lu->A,A->cmap->n,A->rmap->n,aa->nz,aa->a,aa->j,aa->i,SLU_NC,SLU_S,SLU_GE));
297: #else
298:   PetscStackCall("SuperLU:dCreate_CompCol_Matrix",dCreate_CompCol_Matrix(&lu->A,A->cmap->n,A->rmap->n,aa->nz,aa->a,aa->j,aa->i,SLU_NC,SLU_D,SLU_GE));
299: #endif
300: #endif

302:   /* Numerical factorization */
303:   lu->B.ncol = 0;  /* Indicate not to solve the system */
304:   if (F->factortype == MAT_FACTOR_LU) {
305: #if defined(PETSC_USE_COMPLEX)
306: #if defined(PETSC_USE_REAL_SINGLE)
307:     PetscStackCall("SuperLU:cgssvx",cgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
308:                                      &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
309:                                      &lu->Glu, &lu->mem_usage, &lu->stat, &sinfo));
310: #else
311:     PetscStackCall("SuperLU:zgssvx",zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
312:                                      &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
313:                                      &lu->Glu, &lu->mem_usage, &lu->stat, &sinfo));
314: #endif
315: #else
316: #if defined(PETSC_USE_REAL_SINGLE)
317:     PetscStackCall("SuperLU:sgssvx",sgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
318:                                      &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
319:                                      &lu->Glu, &lu->mem_usage, &lu->stat, &sinfo));
320: #else
321:     PetscStackCall("SuperLU:dgssvx",dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
322:                                      &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
323:                                      &lu->Glu,&lu->mem_usage, &lu->stat, &sinfo));
324: #endif
325: #endif
326:   } else if (F->factortype == MAT_FACTOR_ILU) {
327:     /* Compute the incomplete factorization, condition number and pivot growth */
328: #if defined(PETSC_USE_COMPLEX)
329: #if defined(PETSC_USE_REAL_SINGLE)
330:     PetscStackCall("SuperLU:cgsisx",cgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r,lu->etree, lu->equed, lu->R, lu->C,
331:                                      &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond,
332:                                      &lu->Glu, &lu->mem_usage, &lu->stat, &sinfo));
333: #else
334:     PetscStackCall("SuperLU:zgsisx",zgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r,lu->etree, lu->equed, lu->R, lu->C,
335:                                      &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond,
336:                                      &lu->Glu, &lu->mem_usage, &lu->stat, &sinfo));
337: #endif
338: #else
339: #if defined(PETSC_USE_REAL_SINGLE)
340:     PetscStackCall("SuperLU:sgsisx",sgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
341:                                      &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond,
342:                                      &lu->Glu, &lu->mem_usage, &lu->stat, &sinfo));
343: #else
344:     PetscStackCall("SuperLU:dgsisx",dgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
345:                                      &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond,
346:                                      &lu->Glu, &lu->mem_usage, &lu->stat, &sinfo));
347: #endif
348: #endif
349:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported");
350:   if (!sinfo || sinfo == lu->A.ncol+1) {
351:     if (lu->options.PivotGrowth) {
352:       PetscPrintf(PETSC_COMM_SELF,"  Recip. pivot growth = %e\n", lu->rpg);
353:     }
354:     if (lu->options.ConditionNumber) {
355:       PetscPrintf(PETSC_COMM_SELF,"  Recip. condition number = %e\n", lu->rcond);
356:     }
357:   } else if (sinfo > 0) {
358:     if (A->erroriffailure) {
359:       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot in row %D",sinfo);
360:     } else {
361:       if (sinfo <= lu->A.ncol) {
362:         if (lu->options.ILU_FillTol == 0.0) {
363:           F->errortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
364:         }
365:         PetscInfo2(F,"Number of zero pivots %D, ILU_FillTol %g\n",sinfo,lu->options.ILU_FillTol);
366:       } else if (sinfo == lu->A.ncol + 1) {
367:         /* 
368:          U is nonsingular, but RCOND is less than machine
369:                        precision, meaning that the matrix is singular to
370:                        working precision. Nevertheless, the solution and
371:                        error bounds are computed because there are a number
372:                        of situations where the computed solution can be more
373:                        accurate than the value of RCOND would suggest.
374:          */
375:         PetscInfo1(F,"Matrix factor U is nonsingular, but is singular to working precision. The solution is computed. info %D",sinfo);
376:       } else { /* sinfo > lu->A.ncol + 1 */
377:         F->errortype = MAT_FACTOR_OUTMEMORY;
378:         PetscInfo1(F,"Number of bytes allocated when memory allocation fails %D\n",sinfo);
379:       }
380:     }
381:   } else SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_LIB, "info = %D, the %D-th argument in gssvx() had an illegal value", sinfo,-sinfo);

383:   if (lu->options.PrintStat) {
384:     PetscPrintf(PETSC_COMM_SELF,"MatLUFactorNumeric_SuperLU():\n");
385:     PetscStackCall("SuperLU:StatPrint",StatPrint(&lu->stat));
386:     Lstore = (SCformat*) lu->L.Store;
387:     Ustore = (NCformat*) lu->U.Store;
388:     PetscPrintf(PETSC_COMM_SELF,"  No of nonzeros in factor L = %D\n", Lstore->nnz);
389:     PetscPrintf(PETSC_COMM_SELF,"  No of nonzeros in factor U = %D\n", Ustore->nnz);
390:     PetscPrintf(PETSC_COMM_SELF,"  No of nonzeros in L+U = %D\n", Lstore->nnz + Ustore->nnz - lu->A.ncol);
391:     PetscPrintf(PETSC_COMM_SELF,"  L\\U MB %.3f\ttotal MB needed %.3f\n",
392:                          lu->mem_usage.for_lu/1e6, lu->mem_usage.total_needed/1e6);
393:   }

395:   lu->flg                = SAME_NONZERO_PATTERN;
396:   F->ops->solve          = MatSolve_SuperLU;
397:   F->ops->solvetranspose = MatSolveTranspose_SuperLU;
398:   F->ops->matsolve       = NULL;
399:   return(0);
400: }

404: static PetscErrorCode MatDestroy_SuperLU(Mat A)
405: {
407:   Mat_SuperLU    *lu=(Mat_SuperLU*)A->data;

410:   if (lu->CleanUpSuperLU) { /* Free the SuperLU datastructures */
411:     PetscStackCall("SuperLU:Destroy_SuperMatrix_Store",Destroy_SuperMatrix_Store(&lu->A));
412:     PetscStackCall("SuperLU:Destroy_SuperMatrix_Store",Destroy_SuperMatrix_Store(&lu->B));
413:     PetscStackCall("SuperLU:Destroy_SuperMatrix_Store",Destroy_SuperMatrix_Store(&lu->X));
414:     PetscStackCall("SuperLU:StatFree",StatFree(&lu->stat));
415:     if (lu->lwork >= 0) {
416:       PetscStackCall("SuperLU:Destroy_SuperNode_Matrix",Destroy_SuperNode_Matrix(&lu->L));
417:       PetscStackCall("SuperLU:Destroy_CompCol_Matrix",Destroy_CompCol_Matrix(&lu->U));
418:     }
419:   }
420:   PetscFree(lu->etree);
421:   PetscFree(lu->perm_r);
422:   PetscFree(lu->perm_c);
423:   PetscFree(lu->R);
424:   PetscFree(lu->C);
425:   PetscFree(lu->rhs_dup);
426:   MatDestroy(&lu->A_dup);
427:   PetscFree(A->data);

429:   /* clear composed functions */
430:   PetscObjectComposeFunction((PetscObject)A,"MatFactorGetSolverPackage_C",NULL);
431:   PetscObjectComposeFunction((PetscObject)A,"MatSuperluSetILUDropTol_C",NULL);
432:   return(0);
433: }

437: static PetscErrorCode MatView_SuperLU(Mat A,PetscViewer viewer)
438: {
439:   PetscErrorCode    ierr;
440:   PetscBool         iascii;
441:   PetscViewerFormat format;

444:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
445:   if (iascii) {
446:     PetscViewerGetFormat(viewer,&format);
447:     if (format == PETSC_VIEWER_ASCII_INFO) {
448:       MatFactorInfo_SuperLU(A,viewer);
449:     }
450:   }
451:   return(0);
452: }

456: PetscErrorCode MatMatSolve_SuperLU(Mat A,Mat B,Mat X)
457: {
458:   Mat_SuperLU    *lu = (Mat_SuperLU*)A->data;
459:   PetscBool      flg;

463:   PetscObjectTypeCompareAny((PetscObject)B,&flg,MATSEQDENSE,MATMPIDENSE,NULL);
464:   if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix");
465:   PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,NULL);
466:   if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix");
467:   lu->options.Trans = TRANS;
468:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatMatSolve_SuperLU() is not implemented yet");
469:   return(0);
470: }

472: /*
473:    Note the r permutation is ignored
474: */
477: static PetscErrorCode MatLUFactorSymbolic_SuperLU(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
478: {
479:   Mat_SuperLU *lu = (Mat_SuperLU*)(F->data);

482:   lu->flg                 = DIFFERENT_NONZERO_PATTERN;
483:   lu->CleanUpSuperLU      = PETSC_TRUE;
484:   F->ops->lufactornumeric = MatLUFactorNumeric_SuperLU;
485:   return(0);
486: }

490: static PetscErrorCode MatSuperluSetILUDropTol_SuperLU(Mat F,PetscReal dtol)
491: {
492:   Mat_SuperLU *lu= (Mat_SuperLU*)F->data;

495:   lu->options.ILU_DropTol = dtol;
496:   return(0);
497: }

501: /*@
502:   MatSuperluSetILUDropTol - Set SuperLU ILU drop tolerance
503:    Logically Collective on Mat

505:    Input Parameters:
506: +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-SuperLU interface
507: -  dtol - drop tolerance

509:   Options Database:
510: .   -mat_superlu_ilu_droptol <dtol>

512:    Level: beginner

514:    References:
515: .      SuperLU Users' Guide

517: .seealso: MatGetFactor()
518: @*/
519: PetscErrorCode MatSuperluSetILUDropTol(Mat F,PetscReal dtol)
520: {

526:   PetscTryMethod(F,"MatSuperluSetILUDropTol_C",(Mat,PetscReal),(F,dtol));
527:   return(0);
528: }

532: PetscErrorCode MatFactorGetSolverPackage_seqaij_superlu(Mat A,const MatSolverPackage *type)
533: {
535:   *type = MATSOLVERSUPERLU;
536:   return(0);
537: }

539: /*MC
540:   MATSOLVERSUPERLU = "superlu" - A solver package providing solvers LU and ILU for sequential matrices
541:   via the external package SuperLU.

543:   Use ./configure --download-superlu to have PETSc installed with SuperLU

545:   Use -pc_type lu -pc_factor_mat_solver_package superlu to us this direct solver

547:   Options Database Keys:
548: + -mat_superlu_equil <FALSE>            - Equil (None)
549: . -mat_superlu_colperm <COLAMD>         - (choose one of) NATURAL MMD_ATA MMD_AT_PLUS_A COLAMD
550: . -mat_superlu_iterrefine <NOREFINE>    - (choose one of) NOREFINE SINGLE DOUBLE EXTRA
551: . -mat_superlu_symmetricmode: <FALSE>   - SymmetricMode (None)
552: . -mat_superlu_diagpivotthresh <1>      - DiagPivotThresh (None)
553: . -mat_superlu_pivotgrowth <FALSE>      - PivotGrowth (None)
554: . -mat_superlu_conditionnumber <FALSE>  - ConditionNumber (None)
555: . -mat_superlu_rowperm <NOROWPERM>      - (choose one of) NOROWPERM LargeDiag
556: . -mat_superlu_replacetinypivot <FALSE> - ReplaceTinyPivot (None)
557: . -mat_superlu_printstat <FALSE>        - PrintStat (None)
558: . -mat_superlu_lwork <0>                - size of work array in bytes used by factorization (None)
559: . -mat_superlu_ilu_droptol <0>          - ILU_DropTol (None)
560: . -mat_superlu_ilu_filltol <0>          - ILU_FillTol (None)
561: . -mat_superlu_ilu_fillfactor <0>       - ILU_FillFactor (None)
562: . -mat_superlu_ilu_droprull <0>         - ILU_DropRule (None)
563: . -mat_superlu_ilu_norm <0>             - ILU_Norm (None)
564: - -mat_superlu_ilu_milu <0>             - ILU_MILU (None)

566:    Notes: Do not confuse this with MATSOLVERSUPERLU_DIST which is for parallel sparse solves

568:    Level: beginner

570: .seealso: PCLU, PCILU, MATSOLVERSUPERLU_DIST, MATSOLVERMUMPS, PCFactorSetMatSolverPackage(), MatSolverPackage
571: M*/

575: static PetscErrorCode MatGetFactor_seqaij_superlu(Mat A,MatFactorType ftype,Mat *F)
576: {
577:   Mat            B;
578:   Mat_SuperLU    *lu;
580:   PetscInt       indx,m=A->rmap->n,n=A->cmap->n;
581:   PetscBool      flg,set;
582:   PetscReal      real_input;
583:   const char     *colperm[]   ={"NATURAL","MMD_ATA","MMD_AT_PLUS_A","COLAMD"}; /* MY_PERMC - not supported by the petsc interface yet */
584:   const char     *iterrefine[]={"NOREFINE", "SINGLE", "DOUBLE", "EXTRA"};
585:   const char     *rowperm[]   ={"NOROWPERM", "LargeDiag"}; /* MY_PERMC - not supported by the petsc interface yet */

588:   MatCreate(PetscObjectComm((PetscObject)A),&B);
589:   MatSetSizes(B,A->rmap->n,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);
590:   PetscStrallocpy("superlu",&((PetscObject)B)->type_name);
591:   MatSetUp(B);
592:   if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU) {
593:     B->ops->lufactorsymbolic  = MatLUFactorSymbolic_SuperLU;
594:     B->ops->ilufactorsymbolic = MatLUFactorSymbolic_SuperLU;
595:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported");

597:   PetscFree(B->solvertype);
598:   PetscStrallocpy(MATSOLVERSUPERLU,&B->solvertype);

600:   B->ops->getinfo     = MatGetInfo_External;
601:   B->ops->destroy     = MatDestroy_SuperLU;
602:   B->ops->view        = MatView_SuperLU;
603:   B->factortype       = ftype;
604:   B->assembled        = PETSC_TRUE;           /* required by -ksp_view */
605:   B->preallocated     = PETSC_TRUE;

607:   PetscNewLog(B,&lu);

609:   if (ftype == MAT_FACTOR_LU) {
610:     set_default_options(&lu->options);
611:     /* Comments from SuperLU_4.0/SRC/dgssvx.c:
612:       "Whether or not the system will be equilibrated depends on the
613:        scaling of the matrix A, but if equilibration is used, A is
614:        overwritten by diag(R)*A*diag(C) and B by diag(R)*B
615:        (if options->Trans=NOTRANS) or diag(C)*B (if options->Trans = TRANS or CONJ)."
616:      We set 'options.Equil = NO' as default because additional space is needed for it.
617:     */
618:     lu->options.Equil = NO;
619:   } else if (ftype == MAT_FACTOR_ILU) {
620:     /* Set the default input options of ilu: */
621:     PetscStackCall("SuperLU:ilu_set_default_options",ilu_set_default_options(&lu->options));
622:   }
623:   lu->options.PrintStat = NO;

625:   /* Initialize the statistics variables. */
626:   PetscStackCall("SuperLU:StatInit",StatInit(&lu->stat));
627:   lu->lwork = 0;   /* allocate space internally by system malloc */

629:   PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"SuperLU Options","Mat");
630:   PetscOptionsBool("-mat_superlu_equil","Equil","None",(PetscBool)lu->options.Equil,(PetscBool*)&lu->options.Equil,NULL);
631:   PetscOptionsEList("-mat_superlu_colperm","ColPerm","None",colperm,4,colperm[3],&indx,&flg);
632:   if (flg) lu->options.ColPerm = (colperm_t)indx;
633:   PetscOptionsEList("-mat_superlu_iterrefine","IterRefine","None",iterrefine,4,iterrefine[0],&indx,&flg);
634:   if (flg) lu->options.IterRefine = (IterRefine_t)indx;
635:   PetscOptionsBool("-mat_superlu_symmetricmode","SymmetricMode","None",(PetscBool)lu->options.SymmetricMode,&flg,&set);
636:   if (set && flg) lu->options.SymmetricMode = YES;
637:   PetscOptionsReal("-mat_superlu_diagpivotthresh","DiagPivotThresh","None",lu->options.DiagPivotThresh,&real_input,&flg);
638:   if (flg) lu->options.DiagPivotThresh = (double) real_input;
639:   PetscOptionsBool("-mat_superlu_pivotgrowth","PivotGrowth","None",(PetscBool)lu->options.PivotGrowth,&flg,&set);
640:   if (set && flg) lu->options.PivotGrowth = YES;
641:   PetscOptionsBool("-mat_superlu_conditionnumber","ConditionNumber","None",(PetscBool)lu->options.ConditionNumber,&flg,&set);
642:   if (set && flg) lu->options.ConditionNumber = YES;
643:   PetscOptionsEList("-mat_superlu_rowperm","rowperm","None",rowperm,2,rowperm[lu->options.RowPerm],&indx,&flg);
644:   if (flg) lu->options.RowPerm = (rowperm_t)indx;
645:   PetscOptionsBool("-mat_superlu_replacetinypivot","ReplaceTinyPivot","None",(PetscBool)lu->options.ReplaceTinyPivot,&flg,&set);
646:   if (set && flg) lu->options.ReplaceTinyPivot = YES;
647:   PetscOptionsBool("-mat_superlu_printstat","PrintStat","None",(PetscBool)lu->options.PrintStat,&flg,&set);
648:   if (set && flg) lu->options.PrintStat = YES;
649:   PetscOptionsInt("-mat_superlu_lwork","size of work array in bytes used by factorization","None",lu->lwork,&lu->lwork,NULL);
650:   if (lu->lwork > 0) {
651:     /* lwork is in bytes, hence PetscMalloc() is used here, not PetscMalloc1()*/
652:     PetscMalloc(lu->lwork,&lu->work);
653:   } else if (lu->lwork != 0 && lu->lwork != -1) {
654:     PetscPrintf(PETSC_COMM_SELF,"   Warning: lwork %D is not supported by SUPERLU. The default lwork=0 is used.\n",lu->lwork);
655:     lu->lwork = 0;
656:   }
657:   /* ilu options */
658:   PetscOptionsReal("-mat_superlu_ilu_droptol","ILU_DropTol","None",lu->options.ILU_DropTol,&real_input,&flg);
659:   if (flg) lu->options.ILU_DropTol = (double) real_input;
660:   PetscOptionsReal("-mat_superlu_ilu_filltol","ILU_FillTol","None",lu->options.ILU_FillTol,&real_input,&flg);
661:   if (flg) lu->options.ILU_FillTol = (double) real_input;
662:   PetscOptionsReal("-mat_superlu_ilu_fillfactor","ILU_FillFactor","None",lu->options.ILU_FillFactor,&real_input,&flg);
663:   if (flg) lu->options.ILU_FillFactor = (double) real_input;
664:   PetscOptionsInt("-mat_superlu_ilu_droprull","ILU_DropRule","None",lu->options.ILU_DropRule,&lu->options.ILU_DropRule,NULL);
665:   PetscOptionsInt("-mat_superlu_ilu_norm","ILU_Norm","None",lu->options.ILU_Norm,&indx,&flg);
666:   if (flg) lu->options.ILU_Norm = (norm_t)indx;
667:   PetscOptionsInt("-mat_superlu_ilu_milu","ILU_MILU","None",lu->options.ILU_MILU,&indx,&flg);
668:   if (flg) lu->options.ILU_MILU = (milu_t)indx;
669:   PetscOptionsEnd();
670:   if (lu->options.Equil == YES) {
671:     /* superlu overwrites input matrix and rhs when Equil is used, thus create A_dup to keep user's A unchanged */
672:     MatDuplicate_SeqAIJ(A,MAT_COPY_VALUES,&lu->A_dup);
673:   }

675:   /* Allocate spaces (notice sizes are for the transpose) */
676:   PetscMalloc1(m,&lu->etree);
677:   PetscMalloc1(n,&lu->perm_r);
678:   PetscMalloc1(m,&lu->perm_c);
679:   PetscMalloc1(n,&lu->R);
680:   PetscMalloc1(m,&lu->C);

682:   /* create rhs and solution x without allocate space for .Store */
683: #if defined(PETSC_USE_COMPLEX)
684: #if defined(PETSC_USE_REAL_SINGLE)
685:   PetscStackCall("SuperLU:cCreate_Dense_Matrix(",cCreate_Dense_Matrix(&lu->B, m, 1, NULL, m, SLU_DN, SLU_C, SLU_GE));
686:   PetscStackCall("SuperLU:cCreate_Dense_Matrix(",cCreate_Dense_Matrix(&lu->X, m, 1, NULL, m, SLU_DN, SLU_C, SLU_GE));
687: #else
688:   PetscStackCall("SuperLU:zCreate_Dense_Matrix",zCreate_Dense_Matrix(&lu->B, m, 1, NULL, m, SLU_DN, SLU_Z, SLU_GE));
689:   PetscStackCall("SuperLU:zCreate_Dense_Matrix",zCreate_Dense_Matrix(&lu->X, m, 1, NULL, m, SLU_DN, SLU_Z, SLU_GE));
690: #endif
691: #else
692: #if defined(PETSC_USE_REAL_SINGLE)
693:   PetscStackCall("SuperLU:sCreate_Dense_Matrix",sCreate_Dense_Matrix(&lu->B, m, 1, NULL, m, SLU_DN, SLU_S, SLU_GE));
694:   PetscStackCall("SuperLU:sCreate_Dense_Matrix",sCreate_Dense_Matrix(&lu->X, m, 1, NULL, m, SLU_DN, SLU_S, SLU_GE));
695: #else
696:   PetscStackCall("SuperLU:dCreate_Dense_Matrix",dCreate_Dense_Matrix(&lu->B, m, 1, NULL, m, SLU_DN, SLU_D, SLU_GE));
697:   PetscStackCall("SuperLU:dCreate_Dense_Matrix",dCreate_Dense_Matrix(&lu->X, m, 1, NULL, m, SLU_DN, SLU_D, SLU_GE));
698: #endif
699: #endif

701:   PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_seqaij_superlu);
702:   PetscObjectComposeFunction((PetscObject)B,"MatSuperluSetILUDropTol_C",MatSuperluSetILUDropTol_SuperLU);
703:   B->data  = lu;

705:   *F       = B;
706:   return(0);
707: }

711: PETSC_EXTERN PetscErrorCode MatSolverPackageRegister_SuperLU(void)
712: {

716:   MatSolverPackageRegister(MATSOLVERSUPERLU,MATSEQAIJ,       MAT_FACTOR_LU,MatGetFactor_seqaij_superlu);
717:   MatSolverPackageRegister(MATSOLVERSUPERLU,MATSEQAIJ,       MAT_FACTOR_ILU,MatGetFactor_seqaij_superlu);
718:   return(0);
719: }