xref: /petsc/src/mat/impls/aij/seq/superlu/superlu.c (revision 5fe6150d55853ab08dea4f53f899a46ee82b79b7)
114b396bbSKris Buschelman /*$Id: superlu.c,v 1.10 2001/08/15 15:56:50 bsmith Exp $*/
214b396bbSKris Buschelman 
314b396bbSKris Buschelman /*
475af56d4SHong Zhang         Provides an interface to the SuperLU 3.0 sparse solver
514b396bbSKris Buschelman */
614b396bbSKris Buschelman 
714b396bbSKris Buschelman #include "src/mat/impls/aij/seq/aij.h"
814b396bbSKris Buschelman 
914b396bbSKris Buschelman EXTERN_C_BEGIN
1014b396bbSKris Buschelman #if defined(PETSC_USE_COMPLEX)
1114b396bbSKris Buschelman #include "zsp_defs.h"
1214b396bbSKris Buschelman #else
1314b396bbSKris Buschelman #include "dsp_defs.h"
1414b396bbSKris Buschelman #endif
1514b396bbSKris Buschelman #include "util.h"
1614b396bbSKris Buschelman EXTERN_C_END
1714b396bbSKris Buschelman 
1814b396bbSKris Buschelman typedef struct {
1975af56d4SHong Zhang   SuperMatrix       A,L,U,B,X;
2075af56d4SHong Zhang   superlu_options_t options;
2175af56d4SHong Zhang   int               *perm_c; /* column permutation vector */
2275af56d4SHong Zhang   int               *perm_r; /* row permutations from partial pivoting */
2375af56d4SHong Zhang   int               *etree;
2475af56d4SHong Zhang   double            *R, *C;
2575af56d4SHong Zhang   char              equed[1];
2675af56d4SHong Zhang   int               lwork;
2775af56d4SHong Zhang   void              *work;
2875af56d4SHong Zhang   double            rpg, rcond;
2975af56d4SHong Zhang   mem_usage_t       mem_usage;
3075af56d4SHong Zhang   MatStructure      flg;
3114b396bbSKris Buschelman 
3214b396bbSKris Buschelman   /* A few function pointers for inheritance */
33f0c56d0fSKris Buschelman   int (*MatDuplicate)(Mat,MatDuplicateOption,Mat*);
3414b396bbSKris Buschelman   int (*MatView)(Mat,PetscViewer);
3514b396bbSKris Buschelman   int (*MatAssemblyEnd)(Mat,MatAssemblyType);
36b0592601SKris Buschelman   int (*MatLUFactorSymbolic)(Mat,IS,IS,MatFactorInfo*,Mat*);
3714b396bbSKris Buschelman   int (*MatDestroy)(Mat);
3814b396bbSKris Buschelman 
3914b396bbSKris Buschelman   /* Flag to clean up (non-global) SuperLU objects during Destroy */
4014b396bbSKris Buschelman   PetscTruth CleanUpSuperLU;
41f0c56d0fSKris Buschelman } Mat_SuperLU;
4214b396bbSKris Buschelman 
4314b396bbSKris Buschelman 
44f0c56d0fSKris Buschelman EXTERN int MatFactorInfo_SuperLU(Mat,PetscViewer);
45f0c56d0fSKris Buschelman EXTERN int MatLUFactorSymbolic_SuperLU(Mat,IS,IS,MatFactorInfo*,Mat*);
46b0592601SKris Buschelman 
47b0592601SKris Buschelman EXTERN_C_BEGIN
48b0592601SKris Buschelman EXTERN int MatConvert_SuperLU_SeqAIJ(Mat,MatType,Mat*);
49b0592601SKris Buschelman EXTERN int MatConvert_SeqAIJ_SuperLU(Mat,MatType,Mat*);
50b0592601SKris Buschelman EXTERN_C_END
5114b396bbSKris Buschelman 
5214b396bbSKris Buschelman #undef __FUNCT__
53f0c56d0fSKris Buschelman #define __FUNCT__ "MatDestroy_SuperLU"
54f0c56d0fSKris Buschelman int MatDestroy_SuperLU(Mat A)
5514b396bbSKris Buschelman {
56460c4903SKris Buschelman   int         ierr;
57f0c56d0fSKris Buschelman   Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr;
5814b396bbSKris Buschelman 
5914b396bbSKris Buschelman   PetscFunctionBegin;
6075af56d4SHong Zhang   if (lu->CleanUpSuperLU) { /* Free the SuperLU datastructures */
6175af56d4SHong Zhang     Destroy_SuperMatrix_Store(&lu->A);
6275af56d4SHong Zhang     Destroy_SuperMatrix_Store(&lu->B);
6375af56d4SHong Zhang     Destroy_SuperMatrix_Store(&lu->X);
649ce50919SHong Zhang 
659ce50919SHong Zhang     ierr = PetscFree(lu->etree);CHKERRQ(ierr);
669ce50919SHong Zhang     ierr = PetscFree(lu->perm_r);CHKERRQ(ierr);
679ce50919SHong Zhang     ierr = PetscFree(lu->perm_c);CHKERRQ(ierr);
689ce50919SHong Zhang     ierr = PetscFree(lu->R);CHKERRQ(ierr);
699ce50919SHong Zhang     ierr = PetscFree(lu->C);CHKERRQ(ierr);
7075af56d4SHong Zhang     if ( lu->lwork >= 0 ) {
71fb3e25aaSKris Buschelman       Destroy_SuperNode_Matrix(&lu->L);
72fb3e25aaSKris Buschelman       Destroy_CompCol_Matrix(&lu->U);
7375af56d4SHong Zhang     }
74fb3e25aaSKris Buschelman   }
75b0592601SKris Buschelman   ierr = MatConvert_SuperLU_SeqAIJ(A,MATSEQAIJ,&A);CHKERRQ(ierr);
76fb3e25aaSKris Buschelman   ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
7714b396bbSKris Buschelman   PetscFunctionReturn(0);
7814b396bbSKris Buschelman }
7914b396bbSKris Buschelman 
8014b396bbSKris Buschelman #undef __FUNCT__
81f0c56d0fSKris Buschelman #define __FUNCT__ "MatView_SuperLU"
82f0c56d0fSKris Buschelman int MatView_SuperLU(Mat A,PetscViewer viewer)
8314b396bbSKris Buschelman {
8414b396bbSKris Buschelman   int               ierr;
8514b396bbSKris Buschelman   PetscTruth        isascii;
8614b396bbSKris Buschelman   PetscViewerFormat format;
87f0c56d0fSKris Buschelman   Mat_SuperLU       *lu=(Mat_SuperLU*)(A->spptr);
8814b396bbSKris Buschelman 
8914b396bbSKris Buschelman   PetscFunctionBegin;
9014b396bbSKris Buschelman   ierr = (*lu->MatView)(A,viewer);CHKERRQ(ierr);
9114b396bbSKris Buschelman 
9214b396bbSKris Buschelman   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr);
9314b396bbSKris Buschelman   if (isascii) {
9414b396bbSKris Buschelman     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
9514b396bbSKris Buschelman     if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) {
96f0c56d0fSKris Buschelman       ierr = MatFactorInfo_SuperLU(A,viewer);CHKERRQ(ierr);
9714b396bbSKris Buschelman     }
9814b396bbSKris Buschelman   }
9914b396bbSKris Buschelman   PetscFunctionReturn(0);
10014b396bbSKris Buschelman }
10114b396bbSKris Buschelman 
10214b396bbSKris Buschelman #undef __FUNCT__
103f0c56d0fSKris Buschelman #define __FUNCT__ "MatAssemblyEnd_SuperLU"
104f0c56d0fSKris Buschelman int MatAssemblyEnd_SuperLU(Mat A,MatAssemblyType mode) {
10514b396bbSKris Buschelman   int         ierr;
106f0c56d0fSKris Buschelman   Mat_SuperLU *lu=(Mat_SuperLU*)(A->spptr);
10714b396bbSKris Buschelman 
10814b396bbSKris Buschelman   PetscFunctionBegin;
10914b396bbSKris Buschelman   ierr = (*lu->MatAssemblyEnd)(A,mode);CHKERRQ(ierr);
110b0592601SKris Buschelman 
111b0592601SKris Buschelman   lu->MatLUFactorSymbolic  = A->ops->lufactorsymbolic;
112f0c56d0fSKris Buschelman   A->ops->lufactorsymbolic = MatLUFactorSymbolic_SuperLU;
11314b396bbSKris Buschelman   PetscFunctionReturn(0);
11414b396bbSKris Buschelman }
11514b396bbSKris Buschelman 
11675af56d4SHong Zhang /* This function was written for SuperLU 2.0 by Matthew Knepley. Not tested for SuperLU 3.0! */
11775af56d4SHong Zhang #ifdef SuperLU2
11814b396bbSKris Buschelman #include "src/mat/impls/dense/seq/dense.h"
11914b396bbSKris Buschelman #undef __FUNCT__
120f0c56d0fSKris Buschelman #define __FUNCT__ "MatCreateNull_SuperLU"
121f0c56d0fSKris Buschelman int MatCreateNull_SuperLU(Mat A,Mat *nullMat)
12214b396bbSKris Buschelman {
123f0c56d0fSKris Buschelman   Mat_SuperLU   *lu = (Mat_SuperLU*)A->spptr;
12414b396bbSKris Buschelman   int           numRows = A->m,numCols = A->n;
12514b396bbSKris Buschelman   SCformat      *Lstore;
12614b396bbSKris Buschelman   int           numNullCols,size;
12775af56d4SHong Zhang   SuperLUStat_t stat;
12814b396bbSKris Buschelman #if defined(PETSC_USE_COMPLEX)
12914b396bbSKris Buschelman   doublecomplex *nullVals,*workVals;
13014b396bbSKris Buschelman #else
13114b396bbSKris Buschelman   PetscScalar   *nullVals,*workVals;
13214b396bbSKris Buschelman #endif
13314b396bbSKris Buschelman   int           row,newRow,col,newCol,block,b,ierr;
13414b396bbSKris Buschelman 
13514b396bbSKris Buschelman   PetscFunctionBegin;
13614b396bbSKris Buschelman   PetscValidHeaderSpecific(A,MAT_COOKIE);
13714b396bbSKris Buschelman   PetscValidPointer(nullMat);
13814b396bbSKris Buschelman   if (!A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Unfactored matrix");
13914b396bbSKris Buschelman   numNullCols = numCols - numRows;
14014b396bbSKris Buschelman   if (numNullCols < 0) SETERRQ(PETSC_ERR_ARG_WRONG,"Function only applies to underdetermined problems");
14114b396bbSKris Buschelman   /* Create the null matrix */
14214b396bbSKris Buschelman   ierr = MatCreateSeqDense(A->comm,numRows,numNullCols,PETSC_NULL,nullMat);CHKERRQ(ierr);
14314b396bbSKris Buschelman   if (numNullCols == 0) {
14414b396bbSKris Buschelman     ierr = MatAssemblyBegin(*nullMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
14514b396bbSKris Buschelman     ierr = MatAssemblyEnd(*nullMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
14614b396bbSKris Buschelman     PetscFunctionReturn(0);
14714b396bbSKris Buschelman   }
14814b396bbSKris Buschelman #if defined(PETSC_USE_COMPLEX)
14914b396bbSKris Buschelman   nullVals = (doublecomplex*)((Mat_SeqDense*)(*nullMat)->data)->v;
15014b396bbSKris Buschelman #else
15114b396bbSKris Buschelman   nullVals = ((Mat_SeqDense*)(*nullMat)->data)->v;
15214b396bbSKris Buschelman #endif
15314b396bbSKris Buschelman   /* Copy in the columns */
15414b396bbSKris Buschelman   Lstore = (SCformat*)lu->L.Store;
15514b396bbSKris Buschelman   for(block = 0; block <= Lstore->nsuper; block++) {
15614b396bbSKris Buschelman     newRow = Lstore->sup_to_col[block];
15714b396bbSKris Buschelman     size   = Lstore->sup_to_col[block+1] - Lstore->sup_to_col[block];
15814b396bbSKris Buschelman     for(col = Lstore->rowind_colptr[newRow]; col < Lstore->rowind_colptr[newRow+1]; col++) {
15914b396bbSKris Buschelman       newCol = Lstore->rowind[col];
16014b396bbSKris Buschelman       if (newCol >= numRows) {
16114b396bbSKris Buschelman         for(b = 0; b < size; b++)
16214b396bbSKris Buschelman #if defined(PETSC_USE_COMPLEX)
16314b396bbSKris Buschelman           nullVals[(newCol-numRows)*numRows+newRow+b] = ((doublecomplex*)Lstore->nzval)[Lstore->nzval_colptr[newRow+b]+col];
16414b396bbSKris Buschelman #else
16514b396bbSKris Buschelman           nullVals[(newCol-numRows)*numRows+newRow+b] = ((double*)Lstore->nzval)[Lstore->nzval_colptr[newRow+b]+col];
16614b396bbSKris Buschelman #endif
16714b396bbSKris Buschelman       }
16814b396bbSKris Buschelman     }
16914b396bbSKris Buschelman   }
17014b396bbSKris Buschelman   /* Permute rhs to form P^T_c B */
17114b396bbSKris Buschelman   ierr = PetscMalloc(numRows*sizeof(double),&workVals);CHKERRQ(ierr);
17214b396bbSKris Buschelman   for(b = 0; b < numNullCols; b++) {
17314b396bbSKris Buschelman     for(row = 0; row < numRows; row++) workVals[lu->perm_c[row]] = nullVals[b*numRows+row];
17414b396bbSKris Buschelman     for(row = 0; row < numRows; row++) nullVals[b*numRows+row]   = workVals[row];
17514b396bbSKris Buschelman   }
17614b396bbSKris Buschelman   /* Backward solve the upper triangle A x = b */
17714b396bbSKris Buschelman   for(b = 0; b < numNullCols; b++) {
17814b396bbSKris Buschelman #if defined(PETSC_USE_COMPLEX)
17975af56d4SHong Zhang     sp_ztrsv("L","T","U",&lu->L,&lu->U,&nullVals[b*numRows],&stat,&ierr);
18014b396bbSKris Buschelman #else
18175af56d4SHong Zhang     sp_dtrsv("L","T","U",&lu->L,&lu->U,&nullVals[b*numRows],&stat,&ierr);
18214b396bbSKris Buschelman #endif
18314b396bbSKris Buschelman     if (ierr < 0)
18414b396bbSKris Buschelman       SETERRQ1(PETSC_ERR_ARG_WRONG,"The argument %d was invalid",-ierr);
18514b396bbSKris Buschelman   }
18614b396bbSKris Buschelman   ierr = PetscFree(workVals);CHKERRQ(ierr);
18714b396bbSKris Buschelman 
18814b396bbSKris Buschelman   ierr = MatAssemblyBegin(*nullMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
18914b396bbSKris Buschelman   ierr = MatAssemblyEnd(*nullMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
19014b396bbSKris Buschelman   PetscFunctionReturn(0);
19114b396bbSKris Buschelman }
19275af56d4SHong Zhang #endif
19314b396bbSKris Buschelman 
19414b396bbSKris Buschelman #undef __FUNCT__
195f0c56d0fSKris Buschelman #define __FUNCT__ "MatSolve_SuperLU"
196f0c56d0fSKris Buschelman int MatSolve_SuperLU(Mat A,Vec b,Vec x)
19714b396bbSKris Buschelman {
198f0c56d0fSKris Buschelman   Mat_SuperLU   *lu = (Mat_SuperLU*)A->spptr;
19975af56d4SHong Zhang   PetscScalar   *barray,*xarray;
2008914a3f7SHong Zhang   int           ierr,info,i;
20175af56d4SHong Zhang   SuperLUStat_t stat;
2028914a3f7SHong Zhang   double        ferr,berr;
20314b396bbSKris Buschelman 
20414b396bbSKris Buschelman   PetscFunctionBegin;
20575af56d4SHong Zhang   lu->B.ncol = 1;   /* Set the number of right-hand side */
20675af56d4SHong Zhang   ierr = VecGetArray(b,&barray);CHKERRQ(ierr);
20775af56d4SHong Zhang   ierr = VecGetArray(x,&xarray);CHKERRQ(ierr);
208*5fe6150dSHong Zhang 
209*5fe6150dSHong Zhang #if defined(PETSC_USE_COMPLEX)
210*5fe6150dSHong Zhang   ((DNformat*)lu->B.Store)->nzval = (doublecomplex*)barray;
211*5fe6150dSHong Zhang   ((DNformat*)lu->X.Store)->nzval = (doublecomplex*)xarray;
212*5fe6150dSHong Zhang #else
213*5fe6150dSHong Zhang   ((DNformat*)lu->B.Store)->nzval = barray;
21475af56d4SHong Zhang   ((DNformat*)lu->X.Store)->nzval = xarray;
215*5fe6150dSHong Zhang #endif
21675af56d4SHong Zhang 
21775af56d4SHong Zhang   /* Initialize the statistics variables. */
21875af56d4SHong Zhang   StatInit(&stat);
21975af56d4SHong Zhang 
22075af56d4SHong Zhang   lu->options.Fact  = FACTORED; /* Indicate the factored form of A is supplied. */
22175af56d4SHong Zhang   lu->options.Trans = TRANS;
2228914a3f7SHong Zhang #if defined(PETSC_USE_COMPLEX)
2238914a3f7SHong Zhang   zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
2248914a3f7SHong Zhang            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
2258914a3f7SHong Zhang            &lu->mem_usage, &stat, &info);
2268914a3f7SHong Zhang #else
22775af56d4SHong Zhang   dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
22875af56d4SHong Zhang            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
22975af56d4SHong Zhang            &lu->mem_usage, &stat, &info);
2308914a3f7SHong Zhang #endif
23175af56d4SHong Zhang   ierr = VecRestoreArray(b,&barray);CHKERRQ(ierr);
23275af56d4SHong Zhang   ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
23375af56d4SHong Zhang 
23475af56d4SHong Zhang   if ( info == 0 || info == lu->A.ncol+1 ) {
23575af56d4SHong Zhang     if ( lu->options.IterRefine ) {
2368914a3f7SHong Zhang       ierr = PetscPrintf(PETSC_COMM_SELF,"Iterative Refinement:\n");
2378914a3f7SHong Zhang       ierr = PetscPrintf(PETSC_COMM_SELF,"  %8s%8s%16s%16s\n", "rhs", "Steps", "FERR", "BERR");
23875af56d4SHong Zhang       for (i = 0; i < 1; ++i)
2398914a3f7SHong Zhang         ierr = PetscPrintf(PETSC_COMM_SELF,"  %8d%8d%16e%16e\n", i+1, stat.RefineSteps, ferr, berr);
24075af56d4SHong Zhang     }
2418914a3f7SHong Zhang   } else if ( info > 0 ){
2428914a3f7SHong Zhang     if ( lu->lwork == -1 ) {
2438914a3f7SHong Zhang       ierr = PetscPrintf(PETSC_COMM_SELF,"  ** Estimated memory: %d bytes\n", info - lu->A.ncol);
2448914a3f7SHong Zhang     } else {
2458914a3f7SHong Zhang       ierr = PetscPrintf(PETSC_COMM_SELF,"  Warning: gssvx() returns info %d\n",info);
2468914a3f7SHong Zhang     }
2478914a3f7SHong Zhang   } else if (info < 0){
2488914a3f7SHong Zhang     SETERRQ2(1, "info = %d, the %d-th argument in gssvx() had an illegal value", info,-info);
24914b396bbSKris Buschelman   }
25014b396bbSKris Buschelman 
2518914a3f7SHong Zhang   if ( lu->options.PrintStat ) {
2528914a3f7SHong Zhang     ierr = PetscPrintf(PETSC_COMM_SELF,"MatSolve__SuperLU():\n");
2538914a3f7SHong Zhang     StatPrint(&stat);
2548914a3f7SHong Zhang   }
25575af56d4SHong Zhang   StatFree(&stat);
25675af56d4SHong Zhang   PetscFunctionReturn(0);
25775af56d4SHong Zhang }
25814b396bbSKris Buschelman 
25914b396bbSKris Buschelman #undef __FUNCT__
260f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorNumeric_SuperLU"
261f0c56d0fSKris Buschelman int MatLUFactorNumeric_SuperLU(Mat A,Mat *F)
26214b396bbSKris Buschelman {
26314b396bbSKris Buschelman   Mat_SeqAIJ    *aa = (Mat_SeqAIJ*)(A)->data;
264f0c56d0fSKris Buschelman   Mat_SuperLU   *lu = (Mat_SuperLU*)(*F)->spptr;
26575af56d4SHong Zhang   int           ierr,info;
26614b396bbSKris Buschelman   PetscTruth    flag;
26775af56d4SHong Zhang   SuperLUStat_t stat;
26875af56d4SHong Zhang   double        ferr, berr;
2698914a3f7SHong Zhang   NCformat      *Ustore;
2708914a3f7SHong Zhang   SCformat      *Lstore;
27114b396bbSKris Buschelman 
27214b396bbSKris Buschelman   PetscFunctionBegin;
2739ce50919SHong Zhang   if (lu->flg == SAME_NONZERO_PATTERN){ /* successing numerical factorization */
2749ce50919SHong Zhang     lu->options.Fact = SamePattern;
27575af56d4SHong Zhang     /* Ref: ~SuperLU_3.0/EXAMPLE/dlinsolx2.c */
27675af56d4SHong Zhang     Destroy_SuperMatrix_Store(&lu->A);
2779ce50919SHong Zhang     if ( lu->lwork >= 0 ) {
27875af56d4SHong Zhang       Destroy_SuperNode_Matrix(&lu->L);
27975af56d4SHong Zhang       Destroy_CompCol_Matrix(&lu->U);
28075af56d4SHong Zhang       lu->options.Fact = SamePattern;
28114b396bbSKris Buschelman     }
28275af56d4SHong Zhang   }
28314b396bbSKris Buschelman 
28475af56d4SHong Zhang   /* Create the SuperMatrix for lu->A=A^T:
28575af56d4SHong Zhang        Since SuperLU likes column-oriented matrices,we pass it the transpose,
28675af56d4SHong Zhang        and then solve A^T X = B in MatSolve(). */
28775af56d4SHong Zhang #if defined(PETSC_USE_COMPLEX)
288*5fe6150dSHong Zhang   zCreate_CompCol_Matrix(&lu->A,A->n,A->m,aa->nz,(doublecomplex*)aa->a,aa->j,aa->i,
28975af56d4SHong Zhang                            SLU_NC,SLU_Z,SLU_GE);
29075af56d4SHong Zhang #else
29175af56d4SHong Zhang   dCreate_CompCol_Matrix(&lu->A,A->n,A->m,aa->nz,aa->a,aa->j,aa->i,
29275af56d4SHong Zhang                            SLU_NC,SLU_D,SLU_GE);
29375af56d4SHong Zhang #endif
29475af56d4SHong Zhang 
29575af56d4SHong Zhang   /* Initialize the statistics variables. */
29675af56d4SHong Zhang   StatInit(&stat);
29775af56d4SHong Zhang 
2989ce50919SHong Zhang   /* Numerical factorization */
29975af56d4SHong Zhang   lu->B.ncol = 0;  /* Indicate not to solve the system */
3008914a3f7SHong Zhang #if defined(PETSC_USE_COMPLEX)
3018914a3f7SHong Zhang    zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
3028914a3f7SHong Zhang            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
3038914a3f7SHong Zhang            &lu->mem_usage, &stat, &info);
3048914a3f7SHong Zhang #else
30575af56d4SHong Zhang   dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
30675af56d4SHong Zhang            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
30775af56d4SHong Zhang            &lu->mem_usage, &stat, &info);
3088914a3f7SHong Zhang #endif
30975af56d4SHong Zhang   if ( info == 0 || info == lu->A.ncol+1 ) {
3108914a3f7SHong Zhang     if ( lu->options.PivotGrowth )
3118914a3f7SHong Zhang       ierr = PetscPrintf(PETSC_COMM_SELF,"  Recip. pivot growth = %e\n", lu->rpg);
31275af56d4SHong Zhang     if ( lu->options.ConditionNumber )
3138914a3f7SHong Zhang       ierr = PetscPrintf(PETSC_COMM_SELF,"  Recip. condition number = %e\n", lu->rcond);
3148914a3f7SHong Zhang   } else if ( info > 0 ){
3158914a3f7SHong Zhang     if ( lu->lwork == -1 ) {
3168914a3f7SHong Zhang       ierr = PetscPrintf(PETSC_COMM_SELF,"  ** Estimated memory: %d bytes\n", info - lu->A.ncol);
3178914a3f7SHong Zhang     } else {
3188914a3f7SHong Zhang       ierr = PetscPrintf(PETSC_COMM_SELF,"  Warning: gssvx() returns info %d\n",info);
3198914a3f7SHong Zhang     }
3208914a3f7SHong Zhang   } else if (info < 0){
3218914a3f7SHong Zhang     SETERRQ2(1, "info = %d, the %d-th argument in gssvx() had an illegal value", info,-info);
32275af56d4SHong Zhang   }
32375af56d4SHong Zhang 
3248914a3f7SHong Zhang   if ( lu->options.PrintStat ) {
3258914a3f7SHong Zhang     ierr = PetscPrintf(PETSC_COMM_SELF,"MatLUFactorNumeric_SuperLU():\n");
3268914a3f7SHong Zhang     StatPrint(&stat);
3278914a3f7SHong Zhang     Lstore = (SCformat *) lu->L.Store;
3288914a3f7SHong Zhang     Ustore = (NCformat *) lu->U.Store;
3298914a3f7SHong Zhang     ierr = PetscPrintf(PETSC_COMM_SELF,"  No of nonzeros in factor L = %d\n", Lstore->nnz);
3308914a3f7SHong Zhang     ierr = PetscPrintf(PETSC_COMM_SELF,"  No of nonzeros in factor U = %d\n", Ustore->nnz);
3318914a3f7SHong Zhang     ierr = PetscPrintf(PETSC_COMM_SELF,"  No of nonzeros in L+U = %d\n", Lstore->nnz + Ustore->nnz - lu->A.ncol);
3328914a3f7SHong Zhang     ierr = PetscPrintf(PETSC_COMM_SELF,"  L\\U MB %.3f\ttotal MB needed %.3f\texpansions %d\n",
3338914a3f7SHong Zhang 	       lu->mem_usage.for_lu/1e6, lu->mem_usage.total_needed/1e6,
3348914a3f7SHong Zhang 	       lu->mem_usage.expansions);
3358914a3f7SHong Zhang   }
33675af56d4SHong Zhang   StatFree(&stat);
33775af56d4SHong Zhang 
33875af56d4SHong Zhang   lu->flg = SAME_NONZERO_PATTERN;
33975af56d4SHong Zhang   PetscFunctionReturn(0);
34014b396bbSKris Buschelman }
34114b396bbSKris Buschelman 
34214b396bbSKris Buschelman /*
34314b396bbSKris Buschelman    Note the r permutation is ignored
34414b396bbSKris Buschelman */
34514b396bbSKris Buschelman #undef __FUNCT__
346f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_SuperLU"
347f0c56d0fSKris Buschelman int MatLUFactorSymbolic_SuperLU(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F)
34814b396bbSKris Buschelman {
34914b396bbSKris Buschelman   Mat          B;
350f0c56d0fSKris Buschelman   Mat_SuperLU  *lu;
3519ce50919SHong Zhang   int          ierr,m=A->m,n=A->n,indx;
3529ce50919SHong Zhang   PetscTruth   flg;
3538914a3f7SHong Zhang   char         *colperm[]={"NATURAL","MMD_ATA","MMD_AT_PLUS_A","COLAMD"}; /* MY_PERMC - not supported by the petsc interface yet */
3549ce50919SHong Zhang   char         *iterrefine[]={"NOREFINE", "SINGLE", "DOUBLE", "EXTRA"};
3558914a3f7SHong Zhang   char         *rowperm[]={"NOROWPERM", "LargeDiag"}; /* MY_PERMC - not supported by the petsc interface yet */
35614b396bbSKris Buschelman 
35714b396bbSKris Buschelman   PetscFunctionBegin;
35814b396bbSKris Buschelman 
359899d7b4fSKris Buschelman   ierr = MatCreate(A->comm,A->m,A->n,PETSC_DETERMINE,PETSC_DETERMINE,&B);CHKERRQ(ierr);
360899d7b4fSKris Buschelman   ierr = MatSetType(B,MATSUPERLU);CHKERRQ(ierr);
361899d7b4fSKris Buschelman   ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
362f0c56d0fSKris Buschelman 
363f0c56d0fSKris Buschelman   B->ops->lufactornumeric = MatLUFactorNumeric_SuperLU;
364f0c56d0fSKris Buschelman   B->ops->solve           = MatSolve_SuperLU;
36514b396bbSKris Buschelman   B->factor               = FACTOR_LU;
36694b7f48cSBarry Smith   B->assembled            = PETSC_TRUE;  /* required by -ksp_view */
36714b396bbSKris Buschelman 
368f0c56d0fSKris Buschelman   lu = (Mat_SuperLU*)(B->spptr);
3699ce50919SHong Zhang 
3709ce50919SHong Zhang   /* Set SuperLU options */
3719ce50919SHong Zhang     /* the default values for options argument:
3729ce50919SHong Zhang 	options.Fact = DOFACT;
3739ce50919SHong Zhang         options.Equil = YES;
3749ce50919SHong Zhang     	options.ColPerm = COLAMD;
3759ce50919SHong Zhang 	options.DiagPivotThresh = 1.0;
3769ce50919SHong Zhang     	options.Trans = NOTRANS;
3779ce50919SHong Zhang     	options.IterRefine = NOREFINE;
3789ce50919SHong Zhang     	options.SymmetricMode = NO;
3799ce50919SHong Zhang     	options.PivotGrowth = NO;
3809ce50919SHong Zhang     	options.ConditionNumber = NO;
3819ce50919SHong Zhang     	options.PrintStat = YES;
3829ce50919SHong Zhang     */
3839ce50919SHong Zhang   set_default_options(&lu->options);
3848914a3f7SHong Zhang   /* equilibration causes error in solve(), thus not supported here. See dgssvx.c for possible reason. */
3858914a3f7SHong Zhang   lu->options.Equil = NO;
3869ce50919SHong Zhang   lu->options.PrintStat = NO;
3878914a3f7SHong Zhang   lu->lwork = 0;   /* allocate space internally by system malloc */
3889ce50919SHong Zhang 
3899ce50919SHong Zhang   ierr = PetscOptionsBegin(A->comm,A->prefix,"SuperLU Options","Mat");CHKERRQ(ierr);
3909ce50919SHong Zhang   /*
3918914a3f7SHong Zhang   ierr = PetscOptionsLogical("-mat_superlu_equil","Equil","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr);
3928914a3f7SHong Zhang   if (flg) lu->options.Equil = YES; -- not supported by the interface !!!
3939ce50919SHong Zhang   */
3948914a3f7SHong Zhang   ierr = PetscOptionsEList("-mat_superlu_colperm","ColPerm","None",colperm,4,colperm[3],&indx,&flg);CHKERRQ(ierr);
3959ce50919SHong Zhang   if (flg) {lu->options.ColPerm = (colperm_t)indx;}
3968914a3f7SHong Zhang   ierr = PetscOptionsEList("-mat_superlu_iterrefine","IterRefine","None",iterrefine,4,iterrefine[0],&indx,&flg);CHKERRQ(ierr);
3979ce50919SHong Zhang   if (flg) { lu->options.IterRefine = (IterRefine_t)indx;}
3988914a3f7SHong Zhang   ierr = PetscOptionsLogical("-mat_superlu_symmetricmode","SymmetricMode","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr);
3999ce50919SHong Zhang   if (flg) lu->options.SymmetricMode = YES;
4008914a3f7SHong Zhang   ierr = PetscOptionsReal("-mat_superlu_diagpivotthresh","DiagPivotThresh","None",lu->options.DiagPivotThresh,&lu->options.DiagPivotThresh,PETSC_NULL);CHKERRQ(ierr);
4018914a3f7SHong Zhang   ierr = PetscOptionsLogical("-mat_superlu_pivotgrowth","PivotGrowth","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr);
4029ce50919SHong Zhang   if (flg) lu->options.PivotGrowth = YES;
4038914a3f7SHong Zhang   ierr = PetscOptionsLogical("-mat_superlu_conditionnumber","ConditionNumber","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr);
4049ce50919SHong Zhang   if (flg) lu->options.ConditionNumber = YES;
4058914a3f7SHong Zhang   ierr = PetscOptionsEList("-mat_superlu_rowperm","rowperm","None",rowperm,2,rowperm[0],&indx,&flg);CHKERRQ(ierr);
4069ce50919SHong Zhang   if (flg) {lu->options.RowPerm = (rowperm_t)indx;}
4078914a3f7SHong Zhang   ierr = PetscOptionsLogical("-mat_superlu_replacetinypivot","ReplaceTinyPivot","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr);
4089ce50919SHong Zhang   if (flg) lu->options.ReplaceTinyPivot = YES;
4098914a3f7SHong Zhang   ierr = PetscOptionsLogical("-mat_superlu_printstat","PrintStat","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr);
4109ce50919SHong Zhang   if (flg) lu->options.PrintStat = YES;
4118914a3f7SHong Zhang   ierr = PetscOptionsInt("-mat_superlu_lwork","size of work array in bytes used by factorization","None",lu->lwork,&lu->lwork,PETSC_NULL);CHKERRQ(ierr);
412*5fe6150dSHong Zhang   if (lu->lwork > 0 ){
413*5fe6150dSHong Zhang     ierr = PetscMalloc(lu->lwork,&lu->work);CHKERRQ(ierr);
414*5fe6150dSHong Zhang   } else if (lu->lwork != 0 && lu->lwork != -1){
4158914a3f7SHong Zhang     ierr = PetscPrintf(PETSC_COMM_SELF,"   Warning: lwork %d is not supported by PETSc interface. The default lwork=0 is used.\n",lu->lwork);
4168914a3f7SHong Zhang     lu->lwork = 0;
4178914a3f7SHong Zhang   }
4189ce50919SHong Zhang   PetscOptionsEnd();
4199ce50919SHong Zhang 
42075af56d4SHong Zhang #ifdef SUPERLU2
4212ecb5974SKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatCreateNull","MatCreateNull_SuperLU",
422f0c56d0fSKris Buschelman                                     (void(*)(void))MatCreateNull_SuperLU);CHKERRQ(ierr);
42375af56d4SHong Zhang #endif
42414b396bbSKris Buschelman 
42575af56d4SHong Zhang   /* Allocate spaces (notice sizes are for the transpose) */
4269ce50919SHong Zhang   ierr = PetscMalloc(m*sizeof(int),&lu->etree);CHKERRQ(ierr);
4279ce50919SHong Zhang   ierr = PetscMalloc(n*sizeof(int),&lu->perm_r);CHKERRQ(ierr);
4289ce50919SHong Zhang   ierr = PetscMalloc(m*sizeof(int),&lu->perm_c);CHKERRQ(ierr);
4299ce50919SHong Zhang   ierr = PetscMalloc(n*sizeof(int),&lu->R);CHKERRQ(ierr);
4309ce50919SHong Zhang   ierr = PetscMalloc(m*sizeof(int),&lu->C);CHKERRQ(ierr);
43175af56d4SHong Zhang 
4329ce50919SHong Zhang   /* create rhs and solution x without allocate space for .Store */
433*5fe6150dSHong Zhang #if defined(PETSC_USE_COMPLEX)
434*5fe6150dSHong Zhang   zCreate_Dense_Matrix(&lu->B, m, 1, PETSC_NULL, m, SLU_DN, SLU_D, SLU_GE);
435*5fe6150dSHong Zhang   zCreate_Dense_Matrix(&lu->X, m, 1, PETSC_NULL, m, SLU_DN, SLU_D, SLU_GE);
436*5fe6150dSHong Zhang #else
43775af56d4SHong Zhang   dCreate_Dense_Matrix(&lu->B, m, 1, PETSC_NULL, m, SLU_DN, SLU_D, SLU_GE);
43875af56d4SHong Zhang   dCreate_Dense_Matrix(&lu->X, m, 1, PETSC_NULL, m, SLU_DN, SLU_D, SLU_GE);
439*5fe6150dSHong Zhang #endif
44014b396bbSKris Buschelman 
44114b396bbSKris Buschelman   lu->flg            = DIFFERENT_NONZERO_PATTERN;
442e740cb95SKris Buschelman   lu->CleanUpSuperLU = PETSC_TRUE;
443899d7b4fSKris Buschelman 
444899d7b4fSKris Buschelman   *F = B;
4459ce50919SHong Zhang   PetscLogObjectMemory(B,(A->m+A->n)*sizeof(int)+sizeof(Mat_SuperLU));
44614b396bbSKris Buschelman   PetscFunctionReturn(0);
44714b396bbSKris Buschelman }
44814b396bbSKris Buschelman 
44994b7f48cSBarry Smith /* used by -ksp_view */
45014b396bbSKris Buschelman #undef __FUNCT__
451f0c56d0fSKris Buschelman #define __FUNCT__ "MatFactorInfo_SuperLU"
452f0c56d0fSKris Buschelman int MatFactorInfo_SuperLU(Mat A,PetscViewer viewer)
45314b396bbSKris Buschelman {
454f0c56d0fSKris Buschelman   Mat_SuperLU       *lu= (Mat_SuperLU*)A->spptr;
455f0c56d0fSKris Buschelman   int               ierr;
4569ce50919SHong Zhang   superlu_options_t options;
457f0c56d0fSKris Buschelman 
458f0c56d0fSKris Buschelman   PetscFunctionBegin;
4599ce50919SHong Zhang   /* check if matrix is superlu_dist type */
4609ce50919SHong Zhang   if (A->ops->solve != MatSolve_SuperLU) PetscFunctionReturn(0);
4619ce50919SHong Zhang 
4629ce50919SHong Zhang   options = lu->options;
463f0c56d0fSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"SuperLU run parameters:\n");CHKERRQ(ierr);
4649ce50919SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Equil: %s\n",(options.Equil != NO) ? "YES": "NO");CHKERRQ(ierr);
4659ce50919SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ColPerm: %d\n",options.ColPerm);CHKERRQ(ierr);
4669ce50919SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  IterRefine: %d\n",options.IterRefine);CHKERRQ(ierr);
4679ce50919SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  SymmetricMode: %s\n",(options.SymmetricMode != NO) ? "YES": "NO");CHKERRQ(ierr);
4689ce50919SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  DiagPivotThresh: %g\n",options.DiagPivotThresh);CHKERRQ(ierr);
4699ce50919SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  PivotGrowth: %s\n",(options.PivotGrowth != NO) ? "YES": "NO");CHKERRQ(ierr);
4709ce50919SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ConditionNumber: %s\n",(options.ConditionNumber != NO) ? "YES": "NO");CHKERRQ(ierr);
4719ce50919SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  RowPerm: %d\n",options.RowPerm);CHKERRQ(ierr);
4729ce50919SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ReplaceTinyPivot: %s\n",(options.ReplaceTinyPivot != NO) ? "YES": "NO");CHKERRQ(ierr);
4739ce50919SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  PrintStat: %s\n",(options.PrintStat != NO) ? "YES": "NO");CHKERRQ(ierr);
4748914a3f7SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  lwork: %d\n",lu->lwork);CHKERRQ(ierr);
4759ce50919SHong Zhang 
476f0c56d0fSKris Buschelman   PetscFunctionReturn(0);
477f0c56d0fSKris Buschelman }
478f0c56d0fSKris Buschelman 
479f0c56d0fSKris Buschelman #undef __FUNCT__
480f0c56d0fSKris Buschelman #define __FUNCT__ "MatDuplicate_SuperLU"
481f0c56d0fSKris Buschelman int MatDuplicate_SuperLU(Mat A, MatDuplicateOption op, Mat *M) {
48214b396bbSKris Buschelman   int         ierr;
4838f340917SKris Buschelman   Mat_SuperLU *lu=(Mat_SuperLU *)A->spptr;
4848f340917SKris Buschelman 
48514b396bbSKris Buschelman   PetscFunctionBegin;
4868f340917SKris Buschelman   ierr = (*lu->MatDuplicate)(A,op,M);CHKERRQ(ierr);
487f0c56d0fSKris Buschelman   ierr = MatConvert_SeqAIJ_SuperLU(*M,MATSUPERLU,M);CHKERRQ(ierr);
4883f953163SKris Buschelman   ierr = PetscMemcpy((*M)->spptr,lu,sizeof(Mat_SuperLU));CHKERRQ(ierr);
48914b396bbSKris Buschelman   PetscFunctionReturn(0);
49014b396bbSKris Buschelman }
49114b396bbSKris Buschelman 
49214b396bbSKris Buschelman EXTERN_C_BEGIN
49314b396bbSKris Buschelman #undef __FUNCT__
494b0592601SKris Buschelman #define __FUNCT__ "MatConvert_SuperLU_SeqAIJ"
495b0592601SKris Buschelman int MatConvert_SuperLU_SeqAIJ(Mat A,MatType type,Mat *newmat) {
496fb3e25aaSKris Buschelman   /* This routine is only called to convert an unfactored PETSc-SuperLU matrix */
497fb3e25aaSKris Buschelman   /* to its base PETSc type, so we will ignore 'MatType type'. */
49814b396bbSKris Buschelman   int                  ierr;
499b0592601SKris Buschelman   Mat                  B=*newmat;
500f0c56d0fSKris Buschelman   Mat_SuperLU   *lu=(Mat_SuperLU *)A->spptr;
501b0592601SKris Buschelman 
502b0592601SKris Buschelman   PetscFunctionBegin;
503b0592601SKris Buschelman   if (B != A) {
504b0592601SKris Buschelman     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
505f0c56d0fSKris Buschelman   }
506b0592601SKris Buschelman   /* Reset the original function pointers */
507f0c56d0fSKris Buschelman   B->ops->duplicate        = lu->MatDuplicate;
508b0592601SKris Buschelman   B->ops->view             = lu->MatView;
509b0592601SKris Buschelman   B->ops->assemblyend      = lu->MatAssemblyEnd;
510b0592601SKris Buschelman   B->ops->lufactorsymbolic = lu->MatLUFactorSymbolic;
511b0592601SKris Buschelman   B->ops->destroy          = lu->MatDestroy;
512b0592601SKris Buschelman   /* lu is only a function pointer stash unless we've factored the matrix, which we haven't! */
513b0592601SKris Buschelman   ierr = PetscFree(lu);CHKERRQ(ierr);
514b0592601SKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJ);CHKERRQ(ierr);
515b0592601SKris Buschelman   *newmat = B;
516b0592601SKris Buschelman   PetscFunctionReturn(0);
517b0592601SKris Buschelman }
518b0592601SKris Buschelman EXTERN_C_END
519b0592601SKris Buschelman 
520b0592601SKris Buschelman EXTERN_C_BEGIN
521b0592601SKris Buschelman #undef __FUNCT__
522b0592601SKris Buschelman #define __FUNCT__ "MatConvert_SeqAIJ_SuperLU"
523b0592601SKris Buschelman int MatConvert_SeqAIJ_SuperLU(Mat A,MatType type,Mat *newmat) {
524fb3e25aaSKris Buschelman   /* This routine is only called to convert to MATSUPERLU */
525fb3e25aaSKris Buschelman   /* from MATSEQAIJ, so we will ignore 'MatType type'. */
526b0592601SKris Buschelman   int         ierr;
527b0592601SKris Buschelman   Mat         B=*newmat;
528f0c56d0fSKris Buschelman   Mat_SuperLU *lu;
52914b396bbSKris Buschelman 
53014b396bbSKris Buschelman   PetscFunctionBegin;
531b0592601SKris Buschelman   if (B != A) {
532b0592601SKris Buschelman     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
533b0592601SKris Buschelman   }
53414b396bbSKris Buschelman 
535f0c56d0fSKris Buschelman   ierr = PetscNew(Mat_SuperLU,&lu);CHKERRQ(ierr);
536f0c56d0fSKris Buschelman   lu->MatDuplicate         = A->ops->duplicate;
53714b396bbSKris Buschelman   lu->MatView              = A->ops->view;
53814b396bbSKris Buschelman   lu->MatAssemblyEnd       = A->ops->assemblyend;
539b0592601SKris Buschelman   lu->MatLUFactorSymbolic  = A->ops->lufactorsymbolic;
54014b396bbSKris Buschelman   lu->MatDestroy           = A->ops->destroy;
54114b396bbSKris Buschelman   lu->CleanUpSuperLU       = PETSC_FALSE;
542b0592601SKris Buschelman 
543b0592601SKris Buschelman   B->spptr                 = (void*)lu;
544f0c56d0fSKris Buschelman   B->ops->duplicate        = MatDuplicate_SuperLU;
545f0c56d0fSKris Buschelman   B->ops->view             = MatView_SuperLU;
546f0c56d0fSKris Buschelman   B->ops->assemblyend      = MatAssemblyEnd_SuperLU;
547f0c56d0fSKris Buschelman   B->ops->lufactorsymbolic = MatLUFactorSymbolic_SuperLU;
548f0c56d0fSKris Buschelman   B->ops->destroy          = MatDestroy_SuperLU;
549b0592601SKris Buschelman 
550b0592601SKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_superlu_C",
551b0592601SKris Buschelman                                            "MatConvert_SeqAIJ_SuperLU",MatConvert_SeqAIJ_SuperLU);CHKERRQ(ierr);
552b0592601SKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_superlu_seqaij_C",
553b0592601SKris Buschelman                                            "MatConvert_SuperLU_SeqAIJ",MatConvert_SuperLU_SeqAIJ);CHKERRQ(ierr);
554fb3e25aaSKris Buschelman   PetscLogInfo(0,"Using SuperLU for SeqAIJ LU factorization and solves.");
555fb3e25aaSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)B,MATSUPERLU);CHKERRQ(ierr);
556b0592601SKris Buschelman   *newmat = B;
557b0592601SKris Buschelman   PetscFunctionReturn(0);
558b0592601SKris Buschelman }
559b0592601SKris Buschelman EXTERN_C_END
560b0592601SKris Buschelman 
56124b6179bSKris Buschelman /*MC
562fafad747SKris Buschelman   MATSUPERLU - MATSUPERLU = "superlu" - A matrix type providing direct solvers (LU) for sequential matrices
56324b6179bSKris Buschelman   via the external package SuperLU.
56424b6179bSKris Buschelman 
56524b6179bSKris Buschelman   If SuperLU is installed (see the manual for
56624b6179bSKris Buschelman   instructions on how to declare the existence of external packages),
56724b6179bSKris Buschelman   a matrix type can be constructed which invokes SuperLU solvers.
56824b6179bSKris Buschelman   After calling MatCreate(...,A), simply call MatSetType(A,MATSUPERLU).
56924b6179bSKris Buschelman   This matrix type is only supported for double precision real.
57024b6179bSKris Buschelman 
57124b6179bSKris Buschelman   This matrix inherits from MATSEQAIJ.  As a result, MatSeqAIJSetPreallocation is
57228b08bd3SKris Buschelman   supported for this matrix type.  One can also call MatConvert for an inplace conversion to or from
57328b08bd3SKris Buschelman   the MATSEQAIJ type without data copy.
57424b6179bSKris Buschelman 
57524b6179bSKris Buschelman   Options Database Keys:
5760bad9183SKris Buschelman + -mat_type superlu - sets the matrix type to "superlu" during a call to MatSetFromOptions()
57724b6179bSKris Buschelman - -mat_superlu_ordering <0,1,2,3> - 0: natural ordering,
57824b6179bSKris Buschelman                                     1: MMD applied to A'*A,
57924b6179bSKris Buschelman                                     2: MMD applied to A'+A,
58024b6179bSKris Buschelman                                     3: COLAMD, approximate minimum degree column ordering
58124b6179bSKris Buschelman 
58224b6179bSKris Buschelman    Level: beginner
58324b6179bSKris Buschelman 
58424b6179bSKris Buschelman .seealso: PCLU
58524b6179bSKris Buschelman M*/
58624b6179bSKris Buschelman 
587b0592601SKris Buschelman EXTERN_C_BEGIN
588b0592601SKris Buschelman #undef __FUNCT__
589f0c56d0fSKris Buschelman #define __FUNCT__ "MatCreate_SuperLU"
590f0c56d0fSKris Buschelman int MatCreate_SuperLU(Mat A) {
591b0592601SKris Buschelman   int ierr;
592b0592601SKris Buschelman 
593b0592601SKris Buschelman   PetscFunctionBegin;
5945441df8eSKris Buschelman   /* Change type name before calling MatSetType to force proper construction of SeqAIJ and SUPERLU types */
5955441df8eSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)A,MATSUPERLU);CHKERRQ(ierr);
596b0592601SKris Buschelman   ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr);
597b0592601SKris Buschelman   ierr = MatConvert_SeqAIJ_SuperLU(A,MATSUPERLU,&A);CHKERRQ(ierr);
59814b396bbSKris Buschelman   PetscFunctionReturn(0);
59914b396bbSKris Buschelman }
60014b396bbSKris Buschelman EXTERN_C_END
601