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 3.0 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;

 49:   /* 
 50:      This is where the methods for the superclass (SeqAIJ) are kept while we 
 51:      reset the pointers in the function table to the new (SuperLU) versions
 52:   */
 53:   PetscErrorCode (*MatDuplicate)(Mat,MatDuplicateOption,Mat*);
 54:   PetscErrorCode (*MatView)(Mat,PetscViewer);
 55:   PetscErrorCode (*MatAssemblyEnd)(Mat,MatAssemblyType);
 56:   PetscErrorCode (*MatLUFactorSymbolic)(Mat,IS,IS,MatFactorInfo*,Mat*);
 57:   PetscErrorCode (*MatDestroy)(Mat);

 59:   /* Flag to clean up (non-global) SuperLU objects during Destroy */
 60:   PetscTruth CleanUpSuperLU;
 61: } Mat_SuperLU;


 73: /*
 74:     Takes a SuperLU matrix (that is a SeqAIJ matrix with the additional SuperLU data-structures
 75:    and methods) and converts it back to a regular SeqAIJ matrix.
 76: */
 80: PetscErrorCode  MatConvert_SuperLU_SeqAIJ(Mat A,MatType type,MatReuse reuse,Mat *newmat)
 81: {
 83:   Mat            B=*newmat;
 84:   Mat_SuperLU    *lu=(Mat_SuperLU *)A->spptr;

 87:   if (reuse == MAT_INITIAL_MATRIX) {
 88:     MatDuplicate(A,MAT_COPY_VALUES,&B);
 89:   }
 90:   /* Reset the original SeqAIJ function pointers */
 91:   B->ops->duplicate        = lu->MatDuplicate;
 92:   B->ops->view             = lu->MatView;
 93:   B->ops->assemblyend      = lu->MatAssemblyEnd;
 94:   B->ops->lufactorsymbolic = lu->MatLUFactorSymbolic;
 95:   B->ops->destroy          = lu->MatDestroy;

 97:   PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaij_superlu_C","",PETSC_NULL);
 98:   PetscObjectComposeFunction((PetscObject)B,"MatConvert_superlu_seqaij_C","",PETSC_NULL);

100:   /* change the type name back to its original value */
101:   PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJ);
102:   *newmat = B;
103:   return(0);
104: }

110: PetscErrorCode  MatConvert_SeqAIJ_SuperLU(Mat A,MatType type,MatReuse reuse,Mat *newmat)
111: {
113:   Mat            B=*newmat;
114:   Mat_SuperLU    *lu;

117:   if (reuse == MAT_INITIAL_MATRIX) {
118:     MatDuplicate(A,MAT_COPY_VALUES,&B);
119:   }

121:   PetscNew(Mat_SuperLU,&lu);
122:   /* save the original SeqAIJ methods that we are changing */
123:   lu->MatDuplicate         = A->ops->duplicate;
124:   lu->MatView              = A->ops->view;
125:   lu->MatAssemblyEnd       = A->ops->assemblyend;
126:   lu->MatLUFactorSymbolic  = A->ops->lufactorsymbolic;
127:   lu->MatDestroy           = A->ops->destroy;
128:   lu->CleanUpSuperLU       = PETSC_FALSE;

130:   /* add to the matrix the location for all the SuperLU data is to be stored */
131:   B->spptr                 = (void*)lu;

133:   /* set the methods in the function table to the SuperLU versions */
134:   B->ops->duplicate        = MatDuplicate_SuperLU;
135:   B->ops->view             = MatView_SuperLU;
136:   B->ops->assemblyend      = MatAssemblyEnd_SuperLU;
137:   B->ops->lufactorsymbolic = MatLUFactorSymbolic_SuperLU;
138:   B->ops->choleskyfactorsymbolic = 0;
139:   B->ops->destroy          = MatDestroy_SuperLU;

141:   PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_superlu_C",
142:                                            "MatConvert_SeqAIJ_SuperLU",MatConvert_SeqAIJ_SuperLU);
143:   PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_superlu_seqaij_C",
144:                                            "MatConvert_SuperLU_SeqAIJ",MatConvert_SuperLU_SeqAIJ);
145:   PetscInfo(0,"Using SuperLU for SeqAIJ LU factorization and solves.\n");
146:   PetscObjectChangeTypeName((PetscObject)B,MATSUPERLU);
147:   *newmat = B;
148:   return(0);
149: }

152: /*
153:     Utility function
154: */
157: PetscErrorCode MatFactorInfo_SuperLU(Mat A,PetscViewer viewer)
158: {
159:   Mat_SuperLU       *lu= (Mat_SuperLU*)A->spptr;
160:   PetscErrorCode    ierr;
161:   superlu_options_t options;

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

167:   options = lu->options;
168:   PetscViewerASCIIPrintf(viewer,"SuperLU run parameters:\n");
169:   PetscViewerASCIIPrintf(viewer,"  Equil: %s\n",(options.Equil != NO) ? "YES": "NO");
170:   PetscViewerASCIIPrintf(viewer,"  ColPerm: %D\n",options.ColPerm);
171:   PetscViewerASCIIPrintf(viewer,"  IterRefine: %D\n",options.IterRefine);
172:   PetscViewerASCIIPrintf(viewer,"  SymmetricMode: %s\n",(options.SymmetricMode != NO) ? "YES": "NO");
173:   PetscViewerASCIIPrintf(viewer,"  DiagPivotThresh: %g\n",options.DiagPivotThresh);
174:   PetscViewerASCIIPrintf(viewer,"  PivotGrowth: %s\n",(options.PivotGrowth != NO) ? "YES": "NO");
175:   PetscViewerASCIIPrintf(viewer,"  ConditionNumber: %s\n",(options.ConditionNumber != NO) ? "YES": "NO");
176:   PetscViewerASCIIPrintf(viewer,"  RowPerm: %D\n",options.RowPerm);
177:   PetscViewerASCIIPrintf(viewer,"  ReplaceTinyPivot: %s\n",(options.ReplaceTinyPivot != NO) ? "YES": "NO");
178:   PetscViewerASCIIPrintf(viewer,"  PrintStat: %s\n",(options.PrintStat != NO) ? "YES": "NO");
179:   PetscViewerASCIIPrintf(viewer,"  lwork: %D\n",lu->lwork);

181:   return(0);
182: }

184: /*
185:     These are the methods provided to REPLACE the corresponding methods of the 
186:    SeqAIJ matrix class. Other methods of SeqAIJ are not replaced
187: */
190: PetscErrorCode MatLUFactorNumeric_SuperLU(Mat A,MatFactorInfo *info,Mat *F)
191: {
192:   Mat_SeqAIJ     *aa = (Mat_SeqAIJ*)(A)->data;
193:   Mat_SuperLU    *lu = (Mat_SuperLU*)(*F)->spptr;
195:   PetscInt       sinfo;
196:   SuperLUStat_t  stat;
197:   PetscReal      ferr, berr;
198:   NCformat       *Ustore;
199:   SCformat       *Lstore;
200: 
202:   if (lu->flg == SAME_NONZERO_PATTERN){ /* successing numerical factorization */
203:     lu->options.Fact = SamePattern;
204:     /* Ref: ~SuperLU_3.0/EXAMPLE/dlinsolx2.c */
205:     Destroy_SuperMatrix_Store(&lu->A);
206:     if ( lu->lwork >= 0 ) {
207:       Destroy_SuperNode_Matrix(&lu->L);
208:       Destroy_CompCol_Matrix(&lu->U);
209:       lu->options.Fact = SamePattern;
210:     }
211:   }

213:   /* Create the SuperMatrix for lu->A=A^T:
214:        Since SuperLU likes column-oriented matrices,we pass it the transpose,
215:        and then solve A^T X = B in MatSolve(). */
216: #if defined(PETSC_USE_COMPLEX)
217:   zCreate_CompCol_Matrix(&lu->A,A->cmap.n,A->rmap.n,aa->nz,(doublecomplex*)aa->a,aa->j,aa->i,
218:                            SLU_NC,SLU_Z,SLU_GE);
219: #else
220:   dCreate_CompCol_Matrix(&lu->A,A->cmap.n,A->rmap.n,aa->nz,aa->a,aa->j,aa->i,
221:                            SLU_NC,SLU_D,SLU_GE);
222: #endif
223: 
224:   /* Initialize the statistics variables. */
225:   StatInit(&stat);

227:   /* Numerical factorization */
228:   lu->B.ncol = 0;  /* Indicate not to solve the system */
229: #if defined(PETSC_USE_COMPLEX)
230:    zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
231:            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
232:            &lu->mem_usage, &stat, &sinfo);
233: #else
234:   dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
235:            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
236:            &lu->mem_usage, &stat, &sinfo);
237: #endif
238:   if ( !sinfo || sinfo == lu->A.ncol+1 ) {
239:     if ( lu->options.PivotGrowth )
240:       PetscPrintf(PETSC_COMM_SELF,"  Recip. pivot growth = %e\n", lu->rpg);
241:     if ( lu->options.ConditionNumber )
242:       PetscPrintf(PETSC_COMM_SELF,"  Recip. condition number = %e\n", lu->rcond);
243:   } else if ( sinfo > 0 ){
244:     if ( lu->lwork == -1 ) {
245:       PetscPrintf(PETSC_COMM_SELF,"  ** Estimated memory: %D bytes\n", sinfo - lu->A.ncol);
246:     } else {
247:       PetscPrintf(PETSC_COMM_SELF,"  Warning: gssvx() returns info %D\n",sinfo);
248:     }
249:   } else { /* sinfo < 0 */
250:     SETERRQ2(PETSC_ERR_LIB, "info = %D, the %D-th argument in gssvx() had an illegal value", sinfo,-sinfo);
251:   }

253:   if ( lu->options.PrintStat ) {
254:     PetscPrintf(PETSC_COMM_SELF,"MatLUFactorNumeric_SuperLU():\n");
255:     StatPrint(&stat);
256:     Lstore = (SCformat *) lu->L.Store;
257:     Ustore = (NCformat *) lu->U.Store;
258:     PetscPrintf(PETSC_COMM_SELF,"  No of nonzeros in factor L = %D\n", Lstore->nnz);
259:     PetscPrintf(PETSC_COMM_SELF,"  No of nonzeros in factor U = %D\n", Ustore->nnz);
260:     PetscPrintf(PETSC_COMM_SELF,"  No of nonzeros in L+U = %D\n", Lstore->nnz + Ustore->nnz - lu->A.ncol);
261:     PetscPrintf(PETSC_COMM_SELF,"  L\\U MB %.3f\ttotal MB needed %.3f\texpansions %D\n",
262:                lu->mem_usage.for_lu/1e6, lu->mem_usage.total_needed/1e6,
263:                lu->mem_usage.expansions);
264:   }
265:   StatFree(&stat);

267:   lu->flg = SAME_NONZERO_PATTERN;
268:   return(0);
269: }

273: PetscErrorCode MatDestroy_SuperLU(Mat A)
274: {
276:   Mat_SuperLU    *lu=(Mat_SuperLU*)A->spptr;

279:   if (lu->CleanUpSuperLU) { /* Free the SuperLU datastructures */
280:     Destroy_SuperMatrix_Store(&lu->A);
281:     Destroy_SuperMatrix_Store(&lu->B);
282:     Destroy_SuperMatrix_Store(&lu->X);

284:     PetscFree(lu->etree);
285:     PetscFree(lu->perm_r);
286:     PetscFree(lu->perm_c);
287:     PetscFree(lu->R);
288:     PetscFree(lu->C);
289:     if ( lu->lwork >= 0 ) {
290:       Destroy_SuperNode_Matrix(&lu->L);
291:       Destroy_CompCol_Matrix(&lu->U);
292:     }
293:   }
294:   MatConvert_SuperLU_SeqAIJ(A,MATSEQAIJ,MAT_REUSE_MATRIX,&A);
295:   (*A->ops->destroy)(A);
296:   return(0);
297: }

301: PetscErrorCode MatView_SuperLU(Mat A,PetscViewer viewer)
302: {
303:   PetscErrorCode    ierr;
304:   PetscTruth        iascii;
305:   PetscViewerFormat format;
306:   Mat_SuperLU       *lu=(Mat_SuperLU*)(A->spptr);

309:   (*lu->MatView)(A,viewer);

311:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
312:   if (iascii) {
313:     PetscViewerGetFormat(viewer,&format);
314:     if (format == PETSC_VIEWER_ASCII_INFO) {
315:       MatFactorInfo_SuperLU(A,viewer);
316:     }
317:   }
318:   return(0);
319: }

323: PetscErrorCode MatAssemblyEnd_SuperLU(Mat A,MatAssemblyType mode) {
325:   Mat_SuperLU    *lu=(Mat_SuperLU*)(A->spptr);

328:   (*lu->MatAssemblyEnd)(A,mode);
329:   lu->MatLUFactorSymbolic  = A->ops->lufactorsymbolic;
330:   A->ops->lufactorsymbolic = MatLUFactorSymbolic_SuperLU;
331:   return(0);
332: }


337: PetscErrorCode MatSolve_SuperLU_Private(Mat A,Vec b,Vec x)
338: {
339:   Mat_SuperLU    *lu = (Mat_SuperLU*)A->spptr;
340:   PetscScalar    *barray,*xarray;
342:   PetscInt       info,i;
343:   SuperLUStat_t  stat;
344:   PetscReal      ferr,berr;

347:   if ( lu->lwork == -1 ) {
348:     return(0);
349:   }
350:   lu->B.ncol = 1;   /* Set the number of right-hand side */
351:   VecGetArray(b,&barray);
352:   VecGetArray(x,&xarray);

354: #if defined(PETSC_USE_COMPLEX)
355:   ((DNformat*)lu->B.Store)->nzval = (doublecomplex*)barray;
356:   ((DNformat*)lu->X.Store)->nzval = (doublecomplex*)xarray;
357: #else
358:   ((DNformat*)lu->B.Store)->nzval = barray;
359:   ((DNformat*)lu->X.Store)->nzval = xarray;
360: #endif

362:   /* Initialize the statistics variables. */
363:   StatInit(&stat);

365:   lu->options.Fact  = FACTORED; /* Indicate the factored form of A is supplied. */
366: #if defined(PETSC_USE_COMPLEX)
367:   zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
368:            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
369:            &lu->mem_usage, &stat, &info);
370: #else
371:   dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
372:            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
373:            &lu->mem_usage, &stat, &info);
374: #endif   
375:   VecRestoreArray(b,&barray);
376:   VecRestoreArray(x,&xarray);

378:   if ( !info || info == lu->A.ncol+1 ) {
379:     if ( lu->options.IterRefine ) {
380:       PetscPrintf(PETSC_COMM_SELF,"Iterative Refinement:\n");
381:       PetscPrintf(PETSC_COMM_SELF,"  %8s%8s%16s%16s\n", "rhs", "Steps", "FERR", "BERR");
382:       for (i = 0; i < 1; ++i)
383:         PetscPrintf(PETSC_COMM_SELF,"  %8d%8d%16e%16e\n", i+1, stat.RefineSteps, ferr, berr);
384:     }
385:   } else if ( info > 0 ){
386:     if ( lu->lwork == -1 ) {
387:       PetscPrintf(PETSC_COMM_SELF,"  ** Estimated memory: %D bytes\n", info - lu->A.ncol);
388:     } else {
389:       PetscPrintf(PETSC_COMM_SELF,"  Warning: gssvx() returns info %D\n",info);
390:     }
391:   } else if (info < 0){
392:     SETERRQ2(PETSC_ERR_LIB, "info = %D, the %D-th argument in gssvx() had an illegal value", info,-info);
393:   }

395:   if ( lu->options.PrintStat ) {
396:     PetscPrintf(PETSC_COMM_SELF,"MatSolve__SuperLU():\n");
397:     StatPrint(&stat);
398:   }
399:   StatFree(&stat);
400:   return(0);
401: }

405: PetscErrorCode MatSolve_SuperLU(Mat A,Vec b,Vec x)
406: {
407:   Mat_SuperLU    *lu = (Mat_SuperLU*)A->spptr;

411:   lu->options.Trans = TRANS;
412:   MatSolve_SuperLU_Private(A,b,x);
413:   return(0);
414: }

418: PetscErrorCode MatSolveTranspose_SuperLU(Mat A,Vec b,Vec x)
419: {
420:   Mat_SuperLU    *lu = (Mat_SuperLU*)A->spptr;

424:   lu->options.Trans = NOTRANS;
425:   MatSolve_SuperLU_Private(A,b,x);
426:   return(0);
427: }


430: /*
431:    Note the r permutation is ignored
432: */
435: PetscErrorCode MatLUFactorSymbolic_SuperLU(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F)
436: {
437:   Mat            B;
438:   Mat_SuperLU    *lu;
440:   PetscInt       m=A->rmap.n,n=A->cmap.n,indx;
441:   PetscTruth     flg;
442:   const char   *colperm[]={"NATURAL","MMD_ATA","MMD_AT_PLUS_A","COLAMD"}; /* MY_PERMC - not supported by the petsc interface yet */
443:   const char   *iterrefine[]={"NOREFINE", "SINGLE", "DOUBLE", "EXTRA"};
444:   const char   *rowperm[]={"NOROWPERM", "LargeDiag"}; /* MY_PERMC - not supported by the petsc interface yet */

447:   MatCreate(A->comm,&B);
448:   MatSetSizes(B,A->rmap.n,A->cmap.n,PETSC_DETERMINE,PETSC_DETERMINE);
449:   MatSetType(B,A->type_name);
450:   MatSeqAIJSetPreallocation(B,0,PETSC_NULL);

452:   B->ops->lufactornumeric = MatLUFactorNumeric_SuperLU;
453:   B->ops->solve           = MatSolve_SuperLU;
454:   B->ops->solvetranspose  = MatSolveTranspose_SuperLU;
455:   B->factor               = FACTOR_LU;
456:   B->assembled            = PETSC_TRUE;  /* required by -ksp_view */
457: 
458:   lu = (Mat_SuperLU*)(B->spptr);

460:   /* Set SuperLU options */
461:     /* the default values for options argument:
462:         options.Fact = DOFACT;
463:         options.Equil = YES;
464:             options.ColPerm = COLAMD;
465:         options.DiagPivotThresh = 1.0;
466:             options.Trans = NOTRANS;
467:             options.IterRefine = NOREFINE;
468:             options.SymmetricMode = NO;
469:             options.PivotGrowth = NO;
470:             options.ConditionNumber = NO;
471:             options.PrintStat = YES;
472:     */
473:   set_default_options(&lu->options);
474:   /* equilibration causes error in solve(), thus not supported here. See dgssvx.c for possible reason. */
475:   lu->options.Equil = NO;
476:   lu->options.PrintStat = NO;
477:   lu->lwork = 0;   /* allocate space internally by system malloc */

479:   PetscOptionsBegin(A->comm,A->prefix,"SuperLU Options","Mat");
480:   PetscOptionsEList("-mat_superlu_colperm","ColPerm","None",colperm,4,colperm[3],&indx,&flg);
481:   if (flg) {lu->options.ColPerm = (colperm_t)indx;}
482:   PetscOptionsEList("-mat_superlu_iterrefine","IterRefine","None",iterrefine,4,iterrefine[0],&indx,&flg);
483:   if (flg) { lu->options.IterRefine = (IterRefine_t)indx;}
484:   PetscOptionsTruth("-mat_superlu_symmetricmode","SymmetricMode","None",PETSC_FALSE,&flg,0);
485:   if (flg) lu->options.SymmetricMode = YES;
486:   PetscOptionsReal("-mat_superlu_diagpivotthresh","DiagPivotThresh","None",lu->options.DiagPivotThresh,&lu->options.DiagPivotThresh,PETSC_NULL);
487:   PetscOptionsTruth("-mat_superlu_pivotgrowth","PivotGrowth","None",PETSC_FALSE,&flg,0);
488:   if (flg) lu->options.PivotGrowth = YES;
489:   PetscOptionsTruth("-mat_superlu_conditionnumber","ConditionNumber","None",PETSC_FALSE,&flg,0);
490:   if (flg) lu->options.ConditionNumber = YES;
491:   PetscOptionsEList("-mat_superlu_rowperm","rowperm","None",rowperm,2,rowperm[0],&indx,&flg);
492:   if (flg) {lu->options.RowPerm = (rowperm_t)indx;}
493:   PetscOptionsTruth("-mat_superlu_replacetinypivot","ReplaceTinyPivot","None",PETSC_FALSE,&flg,0);
494:   if (flg) lu->options.ReplaceTinyPivot = YES;
495:   PetscOptionsTruth("-mat_superlu_printstat","PrintStat","None",PETSC_FALSE,&flg,0);
496:   if (flg) lu->options.PrintStat = YES;
497:   PetscOptionsInt("-mat_superlu_lwork","size of work array in bytes used by factorization","None",lu->lwork,&lu->lwork,PETSC_NULL);
498:   if (lu->lwork > 0 ){
499:     PetscMalloc(lu->lwork,&lu->work);
500:   } else if (lu->lwork != 0 && lu->lwork != -1){
501:     PetscPrintf(PETSC_COMM_SELF,"   Warning: lwork %D is not supported by SUPERLU. The default lwork=0 is used.\n",lu->lwork);
502:     lu->lwork = 0;
503:   }
504:   PetscOptionsEnd();

506: #ifdef SUPERLU2
507:   PetscObjectComposeFunctionDynamic((PetscObject)B,"MatCreateNull","MatCreateNull_SuperLU",
508:                                     (void(*)(void))MatCreateNull_SuperLU);
509: #endif

511:   /* Allocate spaces (notice sizes are for the transpose) */
512:   PetscMalloc(m*sizeof(PetscInt),&lu->etree);
513:   PetscMalloc(n*sizeof(PetscInt),&lu->perm_r);
514:   PetscMalloc(m*sizeof(PetscInt),&lu->perm_c);
515:   PetscMalloc(n*sizeof(PetscInt),&lu->R);
516:   PetscMalloc(m*sizeof(PetscInt),&lu->C);
517: 
518:   /* create rhs and solution x without allocate space for .Store */
519: #if defined(PETSC_USE_COMPLEX)
520:   zCreate_Dense_Matrix(&lu->B, m, 1, PETSC_NULL, m, SLU_DN, SLU_Z, SLU_GE);
521:   zCreate_Dense_Matrix(&lu->X, m, 1, PETSC_NULL, m, SLU_DN, SLU_Z, SLU_GE);
522: #else
523:   dCreate_Dense_Matrix(&lu->B, m, 1, PETSC_NULL, m, SLU_DN, SLU_D, SLU_GE);
524:   dCreate_Dense_Matrix(&lu->X, m, 1, PETSC_NULL, m, SLU_DN, SLU_D, SLU_GE);
525: #endif

527:   lu->flg            = DIFFERENT_NONZERO_PATTERN;
528:   lu->CleanUpSuperLU = PETSC_TRUE;

530:   *F = B;
531:   PetscLogObjectMemory(B,(A->rmap.n+A->cmap.n)*sizeof(PetscInt)+sizeof(Mat_SuperLU));
532:   return(0);
533: }


538: PetscErrorCode MatDuplicate_SuperLU(Mat A, MatDuplicateOption op, Mat *M) {
540:   Mat_SuperLU    *lu=(Mat_SuperLU *)A->spptr;

543:   (*lu->MatDuplicate)(A,op,M);
544:   PetscMemcpy((*M)->spptr,lu,sizeof(Mat_SuperLU));
545:   return(0);
546: }


549: /*MC
550:   MATSUPERLU - MATSUPERLU = "superlu" - A matrix type providing direct solvers (LU) for sequential matrices 
551:   via the external package SuperLU.

553:   If SuperLU is installed (see the manual for
554:   instructions on how to declare the existence of external packages),
555:   a matrix type can be constructed which invokes SuperLU solvers.
556:   After calling MatCreate(...,A), simply call MatSetType(A,MATSUPERLU).

558:   This matrix inherits from MATSEQAIJ.  As a result, MatSeqAIJSetPreallocation is 
559:   supported for this matrix type.  One can also call MatConvert for an inplace conversion to or from 
560:   the MATSEQAIJ type without data copy.

562:   Options Database Keys:
563: + -mat_type superlu - sets the matrix type to "superlu" during a call to MatSetFromOptions()
564: . -mat_superlu_ordering <0,1,2,3> - 0: natural ordering, 
565:                                     1: MMD applied to A'*A, 
566:                                     2: MMD applied to A'+A, 
567:                                     3: COLAMD, approximate minimum degree column ordering
568: . -mat_superlu_iterrefine - have SuperLU do iterative refinement after the triangular solve
569:                           choices: NOREFINE, SINGLE, DOUBLE, EXTRA; default is NOREFINE
570: - -mat_superlu_printstat - print SuperLU statistics about the factorization

572:    Level: beginner

574: .seealso: PCLU
575: M*/

577: /*
578:     Constructor for the new derived matrix class. It simply creates the base 
579:    matrix class and then adds the additional information/methods needed by SuperLU.
580: */
584: PetscErrorCode  MatCreate_SuperLU(Mat A)
585: {

589:   MatSetType(A,MATSEQAIJ);
590:   MatConvert_SeqAIJ_SuperLU(A,MATSUPERLU,MAT_REUSE_MATRIX,&A);
591:   return(0);
592: }