Actual source code: superlu.c

  1: #define PETSCMAT_DLL

  3: /*  -------------------------------------------------------------------- 

  5:      This file implements a subclass of the SeqAIJ matrix class that uses
  6:      the SuperLU sparse solver. You can use this as a starting point for 
  7:      implementing your own subclass of a PETSc matrix class.

  9:      This demonstrates a way to make an implementation inheritence of a PETSc
 10:      matrix type. This means constructing a new matrix type (SuperLU) by changing some
 11:      of the methods of the previous type (SeqAIJ), adding additional data, and possibly
 12:      additional method. (See any book on object oriented programming).
 13: */

 15: /*
 16:      Defines the data structure for the base matrix type (SeqAIJ)
 17: */
 18:  #include ../src/mat/impls/aij/seq/aij.h

 20: /*
 21:      SuperLU include files
 22: */
 24: #if defined(PETSC_USE_COMPLEX)
 25: #include "slu_zdefs.h"
 26: #else
 27: #include "slu_ddefs.h"
 28: #endif  
 29: #include "slu_util.h"

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

 52:   /* Flag to clean up (non-global) SuperLU objects during Destroy */
 53:   PetscTruth CleanUpSuperLU;
 54: } Mat_SuperLU;


 66: /*
 67:     Utility function
 68: */
 71: PetscErrorCode MatFactorInfo_SuperLU(Mat A,PetscViewer viewer)
 72: {
 73:   Mat_SuperLU       *lu= (Mat_SuperLU*)A->spptr;
 74:   PetscErrorCode    ierr;
 75:   superlu_options_t options;

 78:   /* check if matrix is superlu_dist type */
 79:   if (A->ops->solve != MatSolve_SuperLU) return(0);

 81:   options = lu->options;
 82:   PetscViewerASCIIPrintf(viewer,"SuperLU run parameters:\n");
 83:   PetscViewerASCIIPrintf(viewer,"  Equil: %s\n",(options.Equil != NO) ? "YES": "NO");
 84:   PetscViewerASCIIPrintf(viewer,"  ColPerm: %D\n",options.ColPerm);
 85:   PetscViewerASCIIPrintf(viewer,"  IterRefine: %D\n",options.IterRefine);
 86:   PetscViewerASCIIPrintf(viewer,"  SymmetricMode: %s\n",(options.SymmetricMode != NO) ? "YES": "NO");
 87:   PetscViewerASCIIPrintf(viewer,"  DiagPivotThresh: %g\n",options.DiagPivotThresh);
 88:   PetscViewerASCIIPrintf(viewer,"  PivotGrowth: %s\n",(options.PivotGrowth != NO) ? "YES": "NO");
 89:   PetscViewerASCIIPrintf(viewer,"  ConditionNumber: %s\n",(options.ConditionNumber != NO) ? "YES": "NO");
 90:   PetscViewerASCIIPrintf(viewer,"  RowPerm: %D\n",options.RowPerm);
 91:   PetscViewerASCIIPrintf(viewer,"  ReplaceTinyPivot: %s\n",(options.ReplaceTinyPivot != NO) ? "YES": "NO");
 92:   PetscViewerASCIIPrintf(viewer,"  PrintStat: %s\n",(options.PrintStat != NO) ? "YES": "NO");
 93:   PetscViewerASCIIPrintf(viewer,"  lwork: %D\n",lu->lwork);
 94:   if (A->factor == MAT_FACTOR_ILU){
 95:     PetscViewerASCIIPrintf(viewer,"  ILU_DropTol: %g\n",options.ILU_DropTol);
 96:     PetscViewerASCIIPrintf(viewer,"  ILU_FillTol: %g\n",options.ILU_FillTol);
 97:     PetscViewerASCIIPrintf(viewer,"  ILU_FillFactor: %g\n",options.ILU_FillFactor);
 98:     PetscViewerASCIIPrintf(viewer,"  ILU_DropRule: %D\n",options.ILU_DropRule);
 99:     PetscViewerASCIIPrintf(viewer,"  ILU_Norm: %D\n",options.ILU_Norm);
100:     PetscViewerASCIIPrintf(viewer,"  ILU_MILU: %D\n",options.ILU_MILU);
101:   }
102:   return(0);
103: }

105: /*
106:     These are the methods provided to REPLACE the corresponding methods of the 
107:    SeqAIJ matrix class. Other methods of SeqAIJ are not replaced
108: */
111: PetscErrorCode MatLUFactorNumeric_SuperLU(Mat F,Mat A,const MatFactorInfo *info)
112: {
113:   Mat_SuperLU    *lu = (Mat_SuperLU*)(F)->spptr;
114:   Mat_SeqAIJ     *aa;
116:   PetscInt       sinfo;
117:   PetscReal      ferr, berr;
118:   NCformat       *Ustore;
119:   SCformat       *Lstore;
120: 
122:   if (lu->flg == SAME_NONZERO_PATTERN){ /* successing numerical factorization */
123:     lu->options.Fact = SamePattern;
124:     /* Ref: ~SuperLU_3.0/EXAMPLE/dlinsolx2.c */
125:     Destroy_SuperMatrix_Store(&lu->A);
126:     if (lu->options.Equil){
127:       MatCopy_SeqAIJ(A,lu->A_dup,SAME_NONZERO_PATTERN);
128:     }
129:     if ( lu->lwork >= 0 ) {
130:       Destroy_SuperNode_Matrix(&lu->L);
131:       Destroy_CompCol_Matrix(&lu->U);
132:       lu->options.Fact = SamePattern;
133:     }
134:   }

136:   /* Create the SuperMatrix for lu->A=A^T:
137:        Since SuperLU likes column-oriented matrices,we pass it the transpose,
138:        and then solve A^T X = B in MatSolve(). */
139:   if (lu->options.Equil){
140:     aa = (Mat_SeqAIJ*)(lu->A_dup)->data;
141:   } else {
142:     aa = (Mat_SeqAIJ*)(A)->data;
143:   }
144: #if defined(PETSC_USE_COMPLEX)
145:   zCreate_CompCol_Matrix(&lu->A,A->cmap->n,A->rmap->n,aa->nz,(doublecomplex*)aa->a,aa->j,aa->i,
146:                            SLU_NC,SLU_Z,SLU_GE);
147: #else
148:   dCreate_CompCol_Matrix(&lu->A,A->cmap->n,A->rmap->n,aa->nz,aa->a,aa->j,aa->i,
149:                            SLU_NC,SLU_D,SLU_GE);
150: #endif

152:   /* Numerical factorization */
153:   lu->B.ncol = 0;  /* Indicate not to solve the system */
154:   if (F->factor == MAT_FACTOR_LU){
155: #if defined(PETSC_USE_COMPLEX)
156:     zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
157:            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
158:            &lu->mem_usage, &lu->stat, &sinfo);
159: #else
160:     dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
161:            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
162:            &lu->mem_usage, &lu->stat, &sinfo);
163: #endif
164:   } else if (F->factor == MAT_FACTOR_ILU){
165:     /* Compute the incomplete factorization, condition number and pivot growth */
166: #if defined(PETSC_USE_COMPLEX)
167:     zgsisx(&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->mem_usage, &lu->stat, &sinfo);
170: #else
171:     dgsisx(&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->mem_usage, &lu->stat, &sinfo);
174: #endif
175:   } else {
176:     SETERRQ(PETSC_ERR_SUP,"Factor type not supported");
177:   }
178:   if ( !sinfo || sinfo == lu->A.ncol+1 ) {
179:     if ( lu->options.PivotGrowth )
180:       PetscPrintf(PETSC_COMM_SELF,"  Recip. pivot growth = %e\n", lu->rpg);
181:     if ( lu->options.ConditionNumber )
182:       PetscPrintf(PETSC_COMM_SELF,"  Recip. condition number = %e\n", lu->rcond);
183:   } else if ( sinfo > 0 ){
184:     if ( lu->lwork == -1 ) {
185:       PetscPrintf(PETSC_COMM_SELF,"  ** Estimated memory: %D bytes\n", sinfo - lu->A.ncol);
186:     } else {
187:       PetscPrintf(PETSC_COMM_SELF,"  Warning: gssvx() returns info %D\n",sinfo);
188:     }
189:   } else { /* sinfo < 0 */
190:     SETERRQ2(PETSC_ERR_LIB, "info = %D, the %D-th argument in gssvx() had an illegal value", sinfo,-sinfo);
191:   }

193:   if ( lu->options.PrintStat ) {
194:     PetscPrintf(PETSC_COMM_SELF,"MatLUFactorNumeric_SuperLU():\n");
195:     StatPrint(&lu->stat);
196:     Lstore = (SCformat *) lu->L.Store;
197:     Ustore = (NCformat *) lu->U.Store;
198:     PetscPrintf(PETSC_COMM_SELF,"  No of nonzeros in factor L = %D\n", Lstore->nnz);
199:     PetscPrintf(PETSC_COMM_SELF,"  No of nonzeros in factor U = %D\n", Ustore->nnz);
200:     PetscPrintf(PETSC_COMM_SELF,"  No of nonzeros in L+U = %D\n", Lstore->nnz + Ustore->nnz - lu->A.ncol);
201:     PetscPrintf(PETSC_COMM_SELF,"  L\\U MB %.3f\ttotal MB needed %.3f\n",
202:                lu->mem_usage.for_lu/1e6, lu->mem_usage.total_needed/1e6);
203:   }

205:   lu->flg = SAME_NONZERO_PATTERN;
206:   (F)->ops->solve          = MatSolve_SuperLU;
207:   (F)->ops->solvetranspose = MatSolveTranspose_SuperLU;
208:   return(0);
209: }

213: PetscErrorCode MatDestroy_SuperLU(Mat A)
214: {
216:   Mat_SuperLU    *lu=(Mat_SuperLU*)A->spptr;

219:   if (lu->CleanUpSuperLU) { /* Free the SuperLU datastructures */
220:     Destroy_SuperMatrix_Store(&lu->A);
221:     Destroy_SuperMatrix_Store(&lu->B);
222:     Destroy_SuperMatrix_Store(&lu->X);
223:     StatFree(&lu->stat);

225:     PetscFree(lu->etree);
226:     PetscFree(lu->perm_r);
227:     PetscFree(lu->perm_c);
228:     PetscFree(lu->R);
229:     PetscFree(lu->C);
230:     if ( lu->lwork >= 0 ) {
231:       Destroy_SuperNode_Matrix(&lu->L);
232:       Destroy_CompCol_Matrix(&lu->U);
233:     }
234:   }
235:   MatDestroy_SeqAIJ(A);
236:   if (lu->A_dup){MatDestroy(lu->A_dup);}
237:   if (lu->rhs_dup){PetscFree(lu->rhs_dup);}
238:   return(0);
239: }

243: PetscErrorCode MatView_SuperLU(Mat A,PetscViewer viewer)
244: {
245:   PetscErrorCode    ierr;
246:   PetscTruth        iascii;
247:   PetscViewerFormat format;

250:   MatView_SeqAIJ(A,viewer);

252:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
253:   if (iascii) {
254:     PetscViewerGetFormat(viewer,&format);
255:     if (format == PETSC_VIEWER_ASCII_INFO) {
256:       MatFactorInfo_SuperLU(A,viewer);
257:     }
258:   }
259:   return(0);
260: }


265: PetscErrorCode MatSolve_SuperLU_Private(Mat A,Vec b,Vec x)
266: {
267:   Mat_SuperLU    *lu = (Mat_SuperLU*)A->spptr;
268:   PetscScalar    *barray,*xarray;
270:   PetscInt       info,i,n=x->map->n;
271:   PetscReal      ferr,berr;
272: 
274:   if ( lu->lwork == -1 ) {
275:     return(0);
276:   }

278:   lu->B.ncol = 1;   /* Set the number of right-hand side */
279:   if (lu->options.Equil && !lu->rhs_dup){
280:     /* superlu overwrites b when Equil is used, thus create rhs_dup to keep user's b unchanged */
281:     PetscMalloc(n*sizeof(PetscScalar),&lu->rhs_dup);
282:   }
283:   if (lu->options.Equil){
284:     /* Copy b into rsh_dup */
285:     VecGetArray(b,&barray);
286:     PetscMemcpy(lu->rhs_dup,barray,n*sizeof(PetscScalar));
287:     VecRestoreArray(b,&barray);
288:     barray = lu->rhs_dup;
289:   } else {
290:     VecGetArray(b,&barray);
291:   }
292:   VecGetArray(x,&xarray);

294: #if defined(PETSC_USE_COMPLEX)
295:   ((DNformat*)lu->B.Store)->nzval = (doublecomplex*)barray;
296:   ((DNformat*)lu->X.Store)->nzval = (doublecomplex*)xarray;
297: #else
298:   ((DNformat*)lu->B.Store)->nzval = barray;
299:   ((DNformat*)lu->X.Store)->nzval = xarray;
300: #endif

302:   lu->options.Fact = FACTORED; /* Indicate the factored form of A is supplied. */
303:   if (A->factor == MAT_FACTOR_LU){
304: #if defined(PETSC_USE_COMPLEX)
305:     zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
306:            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
307:            &lu->mem_usage, &lu->stat, &info);
308: #else
309:     dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
310:            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
311:            &lu->mem_usage, &lu->stat, &info);
312: #endif
313:   } else if (A->factor == MAT_FACTOR_ILU){
314: #if defined(PETSC_USE_COMPLEX)
315:     zgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
316:            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond,
317:            &lu->mem_usage, &lu->stat, &info);
318: #else
319:     dgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
320:            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond,
321:            &lu->mem_usage, &lu->stat, &info);
322: #endif
323:   } else {
324:     SETERRQ(PETSC_ERR_SUP,"Factor type not supported");
325:   }
326:   if (!lu->options.Equil){
327:     VecRestoreArray(b,&barray);
328:   }
329:   VecRestoreArray(x,&xarray);

331:   if ( !info || info == lu->A.ncol+1 ) {
332:     if ( lu->options.IterRefine ) {
333:       PetscPrintf(PETSC_COMM_SELF,"Iterative Refinement:\n");
334:       PetscPrintf(PETSC_COMM_SELF,"  %8s%8s%16s%16s\n", "rhs", "Steps", "FERR", "BERR");
335:       for (i = 0; i < 1; ++i)
336:         PetscPrintf(PETSC_COMM_SELF,"  %8d%8d%16e%16e\n", i+1, lu->stat.RefineSteps, ferr, berr);
337:     }
338:   } else if ( info > 0 ){
339:     if ( lu->lwork == -1 ) {
340:       PetscPrintf(PETSC_COMM_SELF,"  ** Estimated memory: %D bytes\n", info - lu->A.ncol);
341:     } else {
342:       PetscPrintf(PETSC_COMM_SELF,"  Warning: gssvx() returns info %D\n",info);
343:     }
344:   } else if (info < 0){
345:     SETERRQ2(PETSC_ERR_LIB, "info = %D, the %D-th argument in gssvx() had an illegal value", info,-info);
346:   }

348:   if ( lu->options.PrintStat ) {
349:     PetscPrintf(PETSC_COMM_SELF,"MatSolve__SuperLU():\n");
350:     StatPrint(&lu->stat);
351:   }
352:   return(0);
353: }

357: PetscErrorCode MatSolve_SuperLU(Mat A,Vec b,Vec x)
358: {
359:   Mat_SuperLU    *lu = (Mat_SuperLU*)A->spptr;

363:   lu->options.Trans = TRANS;
364:   MatSolve_SuperLU_Private(A,b,x);
365:   return(0);
366: }

370: PetscErrorCode MatSolveTranspose_SuperLU(Mat A,Vec b,Vec x)
371: {
372:   Mat_SuperLU    *lu = (Mat_SuperLU*)A->spptr;

376:   lu->options.Trans = NOTRANS;
377:   MatSolve_SuperLU_Private(A,b,x);
378:   return(0);
379: }

381: /*
382:    Note the r permutation is ignored
383: */
386: PetscErrorCode MatLUFactorSymbolic_SuperLU(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
387: {
388:   Mat_SuperLU    *lu = (Mat_SuperLU*)((F)->spptr);
389: 
391:   lu->flg                   = DIFFERENT_NONZERO_PATTERN;
392:   lu->CleanUpSuperLU        = PETSC_TRUE;
393:   (F)->ops->lufactornumeric = MatLUFactorNumeric_SuperLU;
394:   return(0);
395: }

400: PetscErrorCode MatFactorGetSolverPackage_seqaij_superlu(Mat A,const MatSolverPackage *type)
401: {
403:   *type = MAT_SOLVER_SUPERLU;
404:   return(0);
405: }
407: 

409: /*MC
410:   MAT_SOLVER_SUPERLU = "superlu" - A solver package providing solvers LU and ILU for sequential matrices 
411:   via the external package SuperLU.

413:   Use config/configure.py --download-superlu to have PETSc installed with SuperLU

415:   Options Database Keys:
416: +  -mat_superlu_equil: <FALSE> Equil (None)
417: .  -mat_superlu_colperm <COLAMD> (choose one of) NATURAL MMD_ATA MMD_AT_PLUS_A COLAMD
418: .  -mat_superlu_iterrefine <NOREFINE> (choose one of) NOREFINE SINGLE DOUBLE EXTRA
419: .  -mat_superlu_symmetricmode: <FALSE> SymmetricMode (None)
420: .  -mat_superlu_diagpivotthresh <1>: DiagPivotThresh (None)
421: .  -mat_superlu_pivotgrowth: <FALSE> PivotGrowth (None)
422: .  -mat_superlu_conditionnumber: <FALSE> ConditionNumber (None)
423: .  -mat_superlu_rowperm <NOROWPERM> (choose one of) NOROWPERM LargeDiag
424: .  -mat_superlu_replacetinypivot: <FALSE> ReplaceTinyPivot (None)
425: .  -mat_superlu_printstat: <FALSE> PrintStat (None)
426: .  -mat_superlu_lwork <0>: size of work array in bytes used by factorization (None)
427: .  -mat_superlu_ilu_droptol <0>: ILU_DropTol (None)
428: .  -mat_superlu_ilu_filltol <0>: ILU_FillTol (None)
429: .  -mat_superlu_ilu_fillfactor <0>: ILU_FillFactor (None)
430: .  -mat_superlu_ilu_droprull <0>: ILU_DropRule (None)
431: .  -mat_superlu_ilu_norm <0>: ILU_Norm (None)
432: -  -mat_superlu_ilu_milu <0>: ILU_MILU (None)

434:    Notes: Do not confuse this with MAT_SOLVER_SUPERLU_DIST which is for parallel sparse solves

436:    Level: beginner

438: .seealso: PCLU, PCILU, MAT_SOLVER_SUPERLU_DIST, MAT_SOLVER_MUMPS, MAT_SOLVER_SPOOLES, PCFactorSetMatSolverPackage(), MatSolverPackage
439: M*/

444: PetscErrorCode MatGetFactor_seqaij_superlu(Mat A,MatFactorType ftype,Mat *F)
445: {
446:   Mat            B;
447:   Mat_SuperLU    *lu;
449:   PetscInt       indx,m=A->rmap->n,n=A->cmap->n;
450:   PetscTruth     flg;
451:   const char     *colperm[]={"NATURAL","MMD_ATA","MMD_AT_PLUS_A","COLAMD"}; /* MY_PERMC - not supported by the petsc interface yet */
452:   const char     *iterrefine[]={"NOREFINE", "SINGLE", "DOUBLE", "EXTRA"};
453:   const char     *rowperm[]={"NOROWPERM", "LargeDiag"}; /* MY_PERMC - not supported by the petsc interface yet */

456:   MatCreate(((PetscObject)A)->comm,&B);
457:   MatSetSizes(B,A->rmap->n,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);
458:   MatSetType(B,((PetscObject)A)->type_name);
459:   MatSeqAIJSetPreallocation(B,0,PETSC_NULL);

461:   if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU){
462:     B->ops->lufactorsymbolic  = MatLUFactorSymbolic_SuperLU;
463:     B->ops->ilufactorsymbolic = MatLUFactorSymbolic_SuperLU;
464:   } else {
465:     SETERRQ(PETSC_ERR_SUP,"Factor type not supported");
466:   }

468:   B->ops->destroy          = MatDestroy_SuperLU;
469:   B->ops->view             = MatView_SuperLU;
470:   B->factor                = ftype;
471:   B->assembled             = PETSC_TRUE;  /* required by -ksp_view */
472:   B->preallocated          = PETSC_TRUE;
473: 
474:   PetscNewLog(B,Mat_SuperLU,&lu);
475: 
476:   if (ftype == MAT_FACTOR_LU){
477:     set_default_options(&lu->options);
478:   } else if (ftype == MAT_FACTOR_ILU){
479:     /* Set the default input options of ilu:
480:         options.Fact = DOFACT;
481:         options.Equil = YES;
482:         options.ColPerm = COLAMD;
483:         options.DiagPivotThresh = 0.1; //different from complete LU
484:         options.Trans = NOTRANS;
485:         options.IterRefine = NOREFINE;
486:         options.SymmetricMode = NO;
487:         options.PivotGrowth = NO;
488:         options.ConditionNumber = NO;
489:         options.PrintStat = YES;
490:         options.RowPerm = LargeDiag;
491:         options.ILU_DropTol = 1e-4;
492:         options.ILU_FillTol = 1e-2;
493:         options.ILU_FillFactor = 10.0;
494:         options.ILU_DropRule = DROP_BASIC | DROP_AREA;
495:         options.ILU_Norm = INF_NORM;
496:         options.ILU_MILU = SMILU_2;
497:     */
498: 
499:     ilu_set_default_options(&lu->options);
500:   }
501:   /* Comments from SuperLU_4.0/SRC/dgssvx.c:
502:       "Whether or not the system will be equilibrated depends on the
503:        scaling of the matrix A, but if equilibration is used, A is
504:        overwritten by diag(R)*A*diag(C) and B by diag(R)*B
505:        (if options->Trans=NOTRANS) or diag(C)*B (if options->Trans = TRANS or CONJ)."
506:      We set 'options.Equil = NO' as default because additional space is needed for it.
507:   */
508:   lu->options.Equil     = NO;
509:   lu->options.PrintStat = NO;

511:   /* Initialize the statistics variables. */
512:   StatInit(&lu->stat);
513:   lu->lwork = 0;   /* allocate space internally by system malloc */

515:   PetscOptionsBegin(((PetscObject)A)->comm,((PetscObject)A)->prefix,"SuperLU Options","Mat");
516:     PetscOptionsTruth("-mat_superlu_equil","Equil","None",PETSC_FALSE,&flg,0);
517:     if (flg) { /* superlu overwrites input matrix and rhs when Equil is used, thus create A_dup to keep user's A unchanged */
518:       lu->options.Equil = YES;
519:       MatDuplicate_SeqAIJ(A,MAT_COPY_VALUES,&lu->A_dup);
520:     }
521:     PetscOptionsEList("-mat_superlu_colperm","ColPerm","None",colperm,4,colperm[3],&indx,&flg);
522:     if (flg) {lu->options.ColPerm = (colperm_t)indx;}
523:     PetscOptionsEList("-mat_superlu_iterrefine","IterRefine","None",iterrefine,4,iterrefine[0],&indx,&flg);
524:     if (flg) { lu->options.IterRefine = (IterRefine_t)indx;}
525:     PetscOptionsTruth("-mat_superlu_symmetricmode","SymmetricMode","None",PETSC_FALSE,&flg,0);
526:     if (flg) lu->options.SymmetricMode = YES;
527:     PetscOptionsReal("-mat_superlu_diagpivotthresh","DiagPivotThresh","None",lu->options.DiagPivotThresh,&lu->options.DiagPivotThresh,PETSC_NULL);
528:     PetscOptionsTruth("-mat_superlu_pivotgrowth","PivotGrowth","None",PETSC_FALSE,&flg,0);
529:     if (flg) lu->options.PivotGrowth = YES;
530:     PetscOptionsTruth("-mat_superlu_conditionnumber","ConditionNumber","None",PETSC_FALSE,&flg,0);
531:     if (flg) lu->options.ConditionNumber = YES;
532:     PetscOptionsEList("-mat_superlu_rowperm","rowperm","None",rowperm,2,rowperm[0],&indx,&flg);
533:     if (flg) {lu->options.RowPerm = (rowperm_t)indx;}
534:     PetscOptionsTruth("-mat_superlu_replacetinypivot","ReplaceTinyPivot","None",PETSC_FALSE,&flg,0);
535:     if (flg) lu->options.ReplaceTinyPivot = YES;
536:     PetscOptionsTruth("-mat_superlu_printstat","PrintStat","None",PETSC_FALSE,&flg,0);
537:     if (flg) lu->options.PrintStat = YES;
538:     PetscOptionsInt("-mat_superlu_lwork","size of work array in bytes used by factorization","None",lu->lwork,&lu->lwork,PETSC_NULL);
539:     if (lu->lwork > 0 ){
540:       PetscMalloc(lu->lwork,&lu->work);
541:     } else if (lu->lwork != 0 && lu->lwork != -1){
542:       PetscPrintf(PETSC_COMM_SELF,"   Warning: lwork %D is not supported by SUPERLU. The default lwork=0 is used.\n",lu->lwork);
543:       lu->lwork = 0;
544:     }
545:     /* ilu options */
546:     PetscOptionsReal("-mat_superlu_ilu_droptol","ILU_DropTol","None",lu->options.ILU_DropTol,&lu->options.ILU_DropTol,PETSC_NULL);
547:     PetscOptionsReal("-mat_superlu_ilu_filltol","ILU_FillTol","None",lu->options.ILU_FillTol,&lu->options.ILU_FillTol,PETSC_NULL);
548:     PetscOptionsReal("-mat_superlu_ilu_fillfactor","ILU_FillFactor","None",lu->options.ILU_FillFactor,&lu->options.ILU_FillFactor,PETSC_NULL);
549:     PetscOptionsInt("-mat_superlu_ilu_droprull","ILU_DropRule","None",lu->options.ILU_DropRule,&lu->options.ILU_DropRule,PETSC_NULL);
550:     PetscOptionsInt("-mat_superlu_ilu_norm","ILU_Norm","None",lu->options.ILU_Norm,&indx,&flg);
551:     if (flg){
552:       lu->options.ILU_Norm = (norm_t)indx;
553:     }
554:     PetscOptionsInt("-mat_superlu_ilu_milu","ILU_MILU","None",lu->options.ILU_MILU,&indx,&flg);
555:     if (flg){
556:       lu->options.ILU_MILU = (milu_t)indx;
557:     }
558:   PetscOptionsEnd();

560:   /* Allocate spaces (notice sizes are for the transpose) */
561:   PetscMalloc(m*sizeof(PetscInt),&lu->etree);
562:   PetscMalloc(n*sizeof(PetscInt),&lu->perm_r);
563:   PetscMalloc(m*sizeof(PetscInt),&lu->perm_c);
564:   PetscMalloc(n*sizeof(PetscScalar),&lu->R);
565:   PetscMalloc(m*sizeof(PetscScalar),&lu->C);
566: 
567:   /* create rhs and solution x without allocate space for .Store */
568: #if defined(PETSC_USE_COMPLEX)
569:   zCreate_Dense_Matrix(&lu->B, m, 1, PETSC_NULL, m, SLU_DN, SLU_Z, SLU_GE);
570:   zCreate_Dense_Matrix(&lu->X, m, 1, PETSC_NULL, m, SLU_DN, SLU_Z, SLU_GE);
571: #else
572:   dCreate_Dense_Matrix(&lu->B, m, 1, PETSC_NULL, m, SLU_DN, SLU_D, SLU_GE);
573:   dCreate_Dense_Matrix(&lu->X, m, 1, PETSC_NULL, m, SLU_DN, SLU_D, SLU_GE);
574: #endif

576: #ifdef SUPERLU2
577:   PetscObjectComposeFunctionDynamic((PetscObject)B,"MatCreateNull","MatCreateNull_SuperLU",(void(*)(void))MatCreateNull_SuperLU);
578: #endif
579:   PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_seqaij_superlu",MatFactorGetSolverPackage_seqaij_superlu);
580:   B->spptr = lu;
581:   *F = B;
582:   return(0);
583: }