xref: /petsc/src/mat/impls/aij/seq/superlu/superlu.c (revision 8914a3f7a19025903633f2addf24f794ed8f0f25)
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;
200*8914a3f7SHong Zhang   int           ierr,info,i;
20175af56d4SHong Zhang   SuperLUStat_t stat;
202*8914a3f7SHong Zhang   double        ferr,berr;
20314b396bbSKris Buschelman 
20414b396bbSKris Buschelman   PetscFunctionBegin;
20575af56d4SHong Zhang   /* rhs vector */
20675af56d4SHong Zhang   lu->B.ncol = 1;   /* Set the number of right-hand side */
20775af56d4SHong Zhang   ierr = VecGetArray(b,&barray);CHKERRQ(ierr);
20875af56d4SHong Zhang   ((DNformat*)lu->B.Store)->nzval = barray;
20975af56d4SHong Zhang 
21075af56d4SHong Zhang   /* solution vector */
21175af56d4SHong Zhang   ierr = VecGetArray(x,&xarray);CHKERRQ(ierr);
21275af56d4SHong Zhang   ((DNformat*)lu->X.Store)->nzval = xarray;
21375af56d4SHong Zhang 
21475af56d4SHong Zhang   /* Initialize the statistics variables. */
21575af56d4SHong Zhang   StatInit(&stat);
21675af56d4SHong Zhang 
21775af56d4SHong Zhang   lu->options.Fact  = FACTORED; /* Indicate the factored form of A is supplied. */
21875af56d4SHong Zhang   lu->options.Trans = TRANS;
219*8914a3f7SHong Zhang #if defined(PETSC_USE_COMPLEX)
220*8914a3f7SHong Zhang   zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
221*8914a3f7SHong Zhang            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
222*8914a3f7SHong Zhang            &lu->mem_usage, &stat, &info);
223*8914a3f7SHong Zhang #else
22475af56d4SHong Zhang   dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
22575af56d4SHong Zhang            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
22675af56d4SHong Zhang            &lu->mem_usage, &stat, &info);
227*8914a3f7SHong Zhang #endif
22875af56d4SHong Zhang   ierr = VecRestoreArray(b,&barray);CHKERRQ(ierr);
22975af56d4SHong Zhang   ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
23075af56d4SHong Zhang 
23175af56d4SHong Zhang   if ( info == 0 || info == lu->A.ncol+1 ) {
23275af56d4SHong Zhang     if ( lu->options.IterRefine ) {
233*8914a3f7SHong Zhang       ierr = PetscPrintf(PETSC_COMM_SELF,"Iterative Refinement:\n");
234*8914a3f7SHong Zhang       ierr = PetscPrintf(PETSC_COMM_SELF,"  %8s%8s%16s%16s\n", "rhs", "Steps", "FERR", "BERR");
23575af56d4SHong Zhang       for (i = 0; i < 1; ++i)
236*8914a3f7SHong Zhang         ierr = PetscPrintf(PETSC_COMM_SELF,"  %8d%8d%16e%16e\n", i+1, stat.RefineSteps, ferr, berr);
23775af56d4SHong Zhang     }
238*8914a3f7SHong Zhang   } else if ( info > 0 ){
239*8914a3f7SHong Zhang     if ( lu->lwork == -1 ) {
240*8914a3f7SHong Zhang       ierr = PetscPrintf(PETSC_COMM_SELF,"  ** Estimated memory: %d bytes\n", info - lu->A.ncol);
241*8914a3f7SHong Zhang     } else {
242*8914a3f7SHong Zhang       ierr = PetscPrintf(PETSC_COMM_SELF,"  Warning: gssvx() returns info %d\n",info);
243*8914a3f7SHong Zhang     }
244*8914a3f7SHong Zhang   } else if (info < 0){
245*8914a3f7SHong Zhang     SETERRQ2(1, "info = %d, the %d-th argument in gssvx() had an illegal value", info,-info);
24614b396bbSKris Buschelman   }
24714b396bbSKris Buschelman 
248*8914a3f7SHong Zhang   if ( lu->options.PrintStat ) {
249*8914a3f7SHong Zhang     ierr = PetscPrintf(PETSC_COMM_SELF,"MatSolve__SuperLU():\n");
250*8914a3f7SHong Zhang     StatPrint(&stat);
251*8914a3f7SHong Zhang   }
25275af56d4SHong Zhang   StatFree(&stat);
25375af56d4SHong Zhang   PetscFunctionReturn(0);
25475af56d4SHong Zhang }
25514b396bbSKris Buschelman 
25614b396bbSKris Buschelman #undef __FUNCT__
257f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorNumeric_SuperLU"
258f0c56d0fSKris Buschelman int MatLUFactorNumeric_SuperLU(Mat A,Mat *F)
25914b396bbSKris Buschelman {
26014b396bbSKris Buschelman   Mat_SeqAIJ    *aa = (Mat_SeqAIJ*)(A)->data;
261f0c56d0fSKris Buschelman   Mat_SuperLU   *lu = (Mat_SuperLU*)(*F)->spptr;
26275af56d4SHong Zhang   int           ierr,info;
26314b396bbSKris Buschelman   PetscTruth    flag;
26475af56d4SHong Zhang   SuperLUStat_t stat;
26575af56d4SHong Zhang   double        ferr, berr;
266*8914a3f7SHong Zhang   NCformat      *Ustore;
267*8914a3f7SHong Zhang   SCformat      *Lstore;
26814b396bbSKris Buschelman 
26914b396bbSKris Buschelman   PetscFunctionBegin;
2709ce50919SHong Zhang   if (lu->flg == SAME_NONZERO_PATTERN){ /* successing numerical factorization */
2719ce50919SHong Zhang     lu->options.Fact = SamePattern;
27275af56d4SHong Zhang     /* Ref: ~SuperLU_3.0/EXAMPLE/dlinsolx2.c */
27375af56d4SHong Zhang     Destroy_SuperMatrix_Store(&lu->A);
2749ce50919SHong Zhang     if ( lu->lwork >= 0 ) {
27575af56d4SHong Zhang       Destroy_SuperNode_Matrix(&lu->L);
27675af56d4SHong Zhang       Destroy_CompCol_Matrix(&lu->U);
27775af56d4SHong Zhang       lu->options.Fact = SamePattern;
27814b396bbSKris Buschelman     }
27975af56d4SHong Zhang   }
28014b396bbSKris Buschelman 
28175af56d4SHong Zhang   /* Create the SuperMatrix for lu->A=A^T:
28275af56d4SHong Zhang        Since SuperLU likes column-oriented matrices,we pass it the transpose,
28375af56d4SHong Zhang        and then solve A^T X = B in MatSolve(). */
28475af56d4SHong Zhang #if defined(PETSC_USE_COMPLEX)
28575af56d4SHong Zhang   zCreate_CompCol_Matrix(&lu->A,A->n,A->m,aa->nz,aa->a,aa->j,aa->i,
28675af56d4SHong Zhang                            SLU_NC,SLU_Z,SLU_GE);
28775af56d4SHong Zhang #else
28875af56d4SHong Zhang   dCreate_CompCol_Matrix(&lu->A,A->n,A->m,aa->nz,aa->a,aa->j,aa->i,
28975af56d4SHong Zhang                            SLU_NC,SLU_D,SLU_GE);
29075af56d4SHong Zhang #endif
29175af56d4SHong Zhang 
29275af56d4SHong Zhang   /* Initialize the statistics variables. */
29375af56d4SHong Zhang   StatInit(&stat);
29475af56d4SHong Zhang 
2959ce50919SHong Zhang   /* Numerical factorization */
29675af56d4SHong Zhang   lu->B.ncol = 0;  /* Indicate not to solve the system */
297*8914a3f7SHong Zhang #if defined(PETSC_USE_COMPLEX)
298*8914a3f7SHong Zhang    zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
299*8914a3f7SHong Zhang            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
300*8914a3f7SHong Zhang            &lu->mem_usage, &stat, &info);
301*8914a3f7SHong Zhang #else
30275af56d4SHong Zhang   dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
30375af56d4SHong Zhang            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
30475af56d4SHong Zhang            &lu->mem_usage, &stat, &info);
305*8914a3f7SHong Zhang #endif
30675af56d4SHong Zhang   if ( info == 0 || info == lu->A.ncol+1 ) {
307*8914a3f7SHong Zhang     if ( lu->options.PivotGrowth )
308*8914a3f7SHong Zhang       ierr = PetscPrintf(PETSC_COMM_SELF,"  Recip. pivot growth = %e\n", lu->rpg);
30975af56d4SHong Zhang     if ( lu->options.ConditionNumber )
310*8914a3f7SHong Zhang       ierr = PetscPrintf(PETSC_COMM_SELF,"  Recip. condition number = %e\n", lu->rcond);
311*8914a3f7SHong Zhang   } else if ( info > 0 ){
312*8914a3f7SHong Zhang     if ( lu->lwork == -1 ) {
313*8914a3f7SHong Zhang       ierr = PetscPrintf(PETSC_COMM_SELF,"  ** Estimated memory: %d bytes\n", info - lu->A.ncol);
314*8914a3f7SHong Zhang     } else {
315*8914a3f7SHong Zhang       ierr = PetscPrintf(PETSC_COMM_SELF,"  Warning: gssvx() returns info %d\n",info);
316*8914a3f7SHong Zhang     }
317*8914a3f7SHong Zhang   } else if (info < 0){
318*8914a3f7SHong Zhang     SETERRQ2(1, "info = %d, the %d-th argument in gssvx() had an illegal value", info,-info);
31975af56d4SHong Zhang   }
32075af56d4SHong Zhang 
321*8914a3f7SHong Zhang   if ( lu->options.PrintStat ) {
322*8914a3f7SHong Zhang     ierr = PetscPrintf(PETSC_COMM_SELF,"MatLUFactorNumeric_SuperLU():\n");
323*8914a3f7SHong Zhang     StatPrint(&stat);
324*8914a3f7SHong Zhang     Lstore = (SCformat *) lu->L.Store;
325*8914a3f7SHong Zhang     Ustore = (NCformat *) lu->U.Store;
326*8914a3f7SHong Zhang     ierr = PetscPrintf(PETSC_COMM_SELF,"  No of nonzeros in factor L = %d\n", Lstore->nnz);
327*8914a3f7SHong Zhang     ierr = PetscPrintf(PETSC_COMM_SELF,"  No of nonzeros in factor U = %d\n", Ustore->nnz);
328*8914a3f7SHong Zhang     ierr = PetscPrintf(PETSC_COMM_SELF,"  No of nonzeros in L+U = %d\n", Lstore->nnz + Ustore->nnz - lu->A.ncol);
329*8914a3f7SHong Zhang     ierr = PetscPrintf(PETSC_COMM_SELF,"  L\\U MB %.3f\ttotal MB needed %.3f\texpansions %d\n",
330*8914a3f7SHong Zhang 	       lu->mem_usage.for_lu/1e6, lu->mem_usage.total_needed/1e6,
331*8914a3f7SHong Zhang 	       lu->mem_usage.expansions);
332*8914a3f7SHong Zhang   }
33375af56d4SHong Zhang   StatFree(&stat);
33475af56d4SHong Zhang 
33575af56d4SHong Zhang   lu->flg = SAME_NONZERO_PATTERN;
33675af56d4SHong Zhang   PetscFunctionReturn(0);
33714b396bbSKris Buschelman }
33814b396bbSKris Buschelman 
33914b396bbSKris Buschelman /*
34014b396bbSKris Buschelman    Note the r permutation is ignored
34114b396bbSKris Buschelman */
34214b396bbSKris Buschelman #undef __FUNCT__
343f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_SuperLU"
344f0c56d0fSKris Buschelman int MatLUFactorSymbolic_SuperLU(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F)
34514b396bbSKris Buschelman {
34614b396bbSKris Buschelman   Mat          B;
347f0c56d0fSKris Buschelman   Mat_SuperLU  *lu;
3489ce50919SHong Zhang   int          ierr,m=A->m,n=A->n,indx;
3499ce50919SHong Zhang   PetscTruth   flg;
350*8914a3f7SHong Zhang   char         *colperm[]={"NATURAL","MMD_ATA","MMD_AT_PLUS_A","COLAMD"}; /* MY_PERMC - not supported by the petsc interface yet */
3519ce50919SHong Zhang   char         *iterrefine[]={"NOREFINE", "SINGLE", "DOUBLE", "EXTRA"};
352*8914a3f7SHong Zhang   char         *rowperm[]={"NOROWPERM", "LargeDiag"}; /* MY_PERMC - not supported by the petsc interface yet */
35314b396bbSKris Buschelman 
35414b396bbSKris Buschelman   PetscFunctionBegin;
35514b396bbSKris Buschelman 
356899d7b4fSKris Buschelman   ierr = MatCreate(A->comm,A->m,A->n,PETSC_DETERMINE,PETSC_DETERMINE,&B);CHKERRQ(ierr);
357899d7b4fSKris Buschelman   ierr = MatSetType(B,MATSUPERLU);CHKERRQ(ierr);
358899d7b4fSKris Buschelman   ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
359f0c56d0fSKris Buschelman 
360f0c56d0fSKris Buschelman   B->ops->lufactornumeric = MatLUFactorNumeric_SuperLU;
361f0c56d0fSKris Buschelman   B->ops->solve           = MatSolve_SuperLU;
36214b396bbSKris Buschelman   B->factor               = FACTOR_LU;
363899d7b4fSKris Buschelman   B->assembled            = PETSC_TRUE;  /* required by -sles_view */
36414b396bbSKris Buschelman 
365f0c56d0fSKris Buschelman   lu = (Mat_SuperLU*)(B->spptr);
3669ce50919SHong Zhang 
3679ce50919SHong Zhang   /* Set SuperLU options */
3689ce50919SHong Zhang     /* the default values for options argument:
3699ce50919SHong Zhang 	options.Fact = DOFACT;
3709ce50919SHong Zhang         options.Equil = YES;
3719ce50919SHong Zhang     	options.ColPerm = COLAMD;
3729ce50919SHong Zhang 	options.DiagPivotThresh = 1.0;
3739ce50919SHong Zhang     	options.Trans = NOTRANS;
3749ce50919SHong Zhang     	options.IterRefine = NOREFINE;
3759ce50919SHong Zhang     	options.SymmetricMode = NO;
3769ce50919SHong Zhang     	options.PivotGrowth = NO;
3779ce50919SHong Zhang     	options.ConditionNumber = NO;
3789ce50919SHong Zhang     	options.PrintStat = YES;
3799ce50919SHong Zhang     */
3809ce50919SHong Zhang   set_default_options(&lu->options);
381*8914a3f7SHong Zhang   /* equilibration causes error in solve(), thus not supported here. See dgssvx.c for possible reason. */
382*8914a3f7SHong Zhang   lu->options.Equil = NO;
3839ce50919SHong Zhang   lu->options.PrintStat = NO;
384*8914a3f7SHong Zhang   lu->lwork = 0;   /* allocate space internally by system malloc */
3859ce50919SHong Zhang 
3869ce50919SHong Zhang   ierr = PetscOptionsBegin(A->comm,A->prefix,"SuperLU Options","Mat");CHKERRQ(ierr);
3879ce50919SHong Zhang   /*
388*8914a3f7SHong Zhang   ierr = PetscOptionsLogical("-mat_superlu_equil","Equil","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr);
389*8914a3f7SHong Zhang   if (flg) lu->options.Equil = YES; -- not supported by the interface !!!
3909ce50919SHong Zhang   */
391*8914a3f7SHong Zhang   ierr = PetscOptionsEList("-mat_superlu_colperm","ColPerm","None",colperm,4,colperm[3],&indx,&flg);CHKERRQ(ierr);
3929ce50919SHong Zhang   if (flg) {lu->options.ColPerm = (colperm_t)indx;}
393*8914a3f7SHong Zhang   ierr = PetscOptionsEList("-mat_superlu_iterrefine","IterRefine","None",iterrefine,4,iterrefine[0],&indx,&flg);CHKERRQ(ierr);
3949ce50919SHong Zhang   if (flg) { lu->options.IterRefine = (IterRefine_t)indx;}
395*8914a3f7SHong Zhang   ierr = PetscOptionsLogical("-mat_superlu_symmetricmode","SymmetricMode","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr);
3969ce50919SHong Zhang   if (flg) lu->options.SymmetricMode = YES;
397*8914a3f7SHong Zhang   ierr = PetscOptionsReal("-mat_superlu_diagpivotthresh","DiagPivotThresh","None",lu->options.DiagPivotThresh,&lu->options.DiagPivotThresh,PETSC_NULL);CHKERRQ(ierr);
398*8914a3f7SHong Zhang   ierr = PetscOptionsLogical("-mat_superlu_pivotgrowth","PivotGrowth","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr);
3999ce50919SHong Zhang   if (flg) lu->options.PivotGrowth = YES;
400*8914a3f7SHong Zhang   ierr = PetscOptionsLogical("-mat_superlu_conditionnumber","ConditionNumber","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr);
4019ce50919SHong Zhang   if (flg) lu->options.ConditionNumber = YES;
402*8914a3f7SHong Zhang   ierr = PetscOptionsEList("-mat_superlu_rowperm","rowperm","None",rowperm,2,rowperm[0],&indx,&flg);CHKERRQ(ierr);
4039ce50919SHong Zhang   if (flg) {lu->options.RowPerm = (rowperm_t)indx;}
404*8914a3f7SHong Zhang   ierr = PetscOptionsLogical("-mat_superlu_replacetinypivot","ReplaceTinyPivot","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr);
4059ce50919SHong Zhang   if (flg) lu->options.ReplaceTinyPivot = YES;
406*8914a3f7SHong Zhang   ierr = PetscOptionsLogical("-mat_superlu_printstat","PrintStat","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr);
4079ce50919SHong Zhang   if (flg) lu->options.PrintStat = YES;
408*8914a3f7SHong Zhang   ierr = PetscOptionsInt("-mat_superlu_lwork","size of work array in bytes used by factorization","None",lu->lwork,&lu->lwork,PETSC_NULL);CHKERRQ(ierr);
409*8914a3f7SHong Zhang   if (lu->lwork != 0 && lu->lwork != -1){
410*8914a3f7SHong Zhang     ierr = PetscPrintf(PETSC_COMM_SELF,"   Warning: lwork %d is not supported by PETSc interface. The default lwork=0 is used.\n",lu->lwork);
411*8914a3f7SHong Zhang     lu->lwork = 0;
412*8914a3f7SHong Zhang   }
4139ce50919SHong Zhang   PetscOptionsEnd();
4149ce50919SHong Zhang 
41575af56d4SHong Zhang #ifdef SUPERLU2
4162ecb5974SKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatCreateNull","MatCreateNull_SuperLU",
417f0c56d0fSKris Buschelman                                     (void(*)(void))MatCreateNull_SuperLU);CHKERRQ(ierr);
41875af56d4SHong Zhang #endif
41914b396bbSKris Buschelman 
42075af56d4SHong Zhang   /* Allocate spaces (notice sizes are for the transpose) */
4219ce50919SHong Zhang   ierr = PetscMalloc(m*sizeof(int),&lu->etree);CHKERRQ(ierr);
4229ce50919SHong Zhang   ierr = PetscMalloc(n*sizeof(int),&lu->perm_r);CHKERRQ(ierr);
4239ce50919SHong Zhang   ierr = PetscMalloc(m*sizeof(int),&lu->perm_c);CHKERRQ(ierr);
4249ce50919SHong Zhang   ierr = PetscMalloc(n*sizeof(int),&lu->R);CHKERRQ(ierr);
4259ce50919SHong Zhang   ierr = PetscMalloc(m*sizeof(int),&lu->C);CHKERRQ(ierr);
42675af56d4SHong Zhang 
4279ce50919SHong Zhang   /* create rhs and solution x without allocate space for .Store */
42875af56d4SHong Zhang   dCreate_Dense_Matrix(&lu->B, m, 1, PETSC_NULL, m, SLU_DN, SLU_D, SLU_GE);
42975af56d4SHong Zhang   dCreate_Dense_Matrix(&lu->X, m, 1, PETSC_NULL, m, SLU_DN, SLU_D, SLU_GE);
43014b396bbSKris Buschelman 
43114b396bbSKris Buschelman   lu->flg            = DIFFERENT_NONZERO_PATTERN;
432e740cb95SKris Buschelman   lu->CleanUpSuperLU = PETSC_TRUE;
433899d7b4fSKris Buschelman 
434899d7b4fSKris Buschelman   *F = B;
4359ce50919SHong Zhang   PetscLogObjectMemory(B,(A->m+A->n)*sizeof(int)+sizeof(Mat_SuperLU));
43614b396bbSKris Buschelman   PetscFunctionReturn(0);
43714b396bbSKris Buschelman }
43814b396bbSKris Buschelman 
43914b396bbSKris Buschelman /* used by -sles_view */
44014b396bbSKris Buschelman #undef __FUNCT__
441f0c56d0fSKris Buschelman #define __FUNCT__ "MatFactorInfo_SuperLU"
442f0c56d0fSKris Buschelman int MatFactorInfo_SuperLU(Mat A,PetscViewer viewer)
44314b396bbSKris Buschelman {
444f0c56d0fSKris Buschelman   Mat_SuperLU       *lu= (Mat_SuperLU*)A->spptr;
445f0c56d0fSKris Buschelman   int               ierr;
4469ce50919SHong Zhang   superlu_options_t options;
447f0c56d0fSKris Buschelman 
448f0c56d0fSKris Buschelman   PetscFunctionBegin;
4499ce50919SHong Zhang   /* check if matrix is superlu_dist type */
4509ce50919SHong Zhang   if (A->ops->solve != MatSolve_SuperLU) PetscFunctionReturn(0);
4519ce50919SHong Zhang 
4529ce50919SHong Zhang   options = lu->options;
453f0c56d0fSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"SuperLU run parameters:\n");CHKERRQ(ierr);
4549ce50919SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Equil: %s\n",(options.Equil != NO) ? "YES": "NO");CHKERRQ(ierr);
4559ce50919SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ColPerm: %d\n",options.ColPerm);CHKERRQ(ierr);
4569ce50919SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  IterRefine: %d\n",options.IterRefine);CHKERRQ(ierr);
4579ce50919SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  SymmetricMode: %s\n",(options.SymmetricMode != NO) ? "YES": "NO");CHKERRQ(ierr);
4589ce50919SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  DiagPivotThresh: %g\n",options.DiagPivotThresh);CHKERRQ(ierr);
4599ce50919SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  PivotGrowth: %s\n",(options.PivotGrowth != NO) ? "YES": "NO");CHKERRQ(ierr);
4609ce50919SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ConditionNumber: %s\n",(options.ConditionNumber != NO) ? "YES": "NO");CHKERRQ(ierr);
4619ce50919SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  RowPerm: %d\n",options.RowPerm);CHKERRQ(ierr);
4629ce50919SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ReplaceTinyPivot: %s\n",(options.ReplaceTinyPivot != NO) ? "YES": "NO");CHKERRQ(ierr);
4639ce50919SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  PrintStat: %s\n",(options.PrintStat != NO) ? "YES": "NO");CHKERRQ(ierr);
464*8914a3f7SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  lwork: %d\n",lu->lwork);CHKERRQ(ierr);
4659ce50919SHong Zhang 
466f0c56d0fSKris Buschelman   PetscFunctionReturn(0);
467f0c56d0fSKris Buschelman }
468f0c56d0fSKris Buschelman 
469f0c56d0fSKris Buschelman #undef __FUNCT__
470f0c56d0fSKris Buschelman #define __FUNCT__ "MatDuplicate_SuperLU"
471f0c56d0fSKris Buschelman int MatDuplicate_SuperLU(Mat A, MatDuplicateOption op, Mat *M) {
47214b396bbSKris Buschelman   int         ierr;
4738f340917SKris Buschelman   Mat_SuperLU *lu=(Mat_SuperLU *)A->spptr;
4748f340917SKris Buschelman 
47514b396bbSKris Buschelman   PetscFunctionBegin;
4768f340917SKris Buschelman   ierr = (*lu->MatDuplicate)(A,op,M);CHKERRQ(ierr);
477f0c56d0fSKris Buschelman   ierr = MatConvert_SeqAIJ_SuperLU(*M,MATSUPERLU,M);CHKERRQ(ierr);
4783f953163SKris Buschelman   ierr = PetscMemcpy((*M)->spptr,lu,sizeof(Mat_SuperLU));CHKERRQ(ierr);
47914b396bbSKris Buschelman   PetscFunctionReturn(0);
48014b396bbSKris Buschelman }
48114b396bbSKris Buschelman 
48214b396bbSKris Buschelman EXTERN_C_BEGIN
48314b396bbSKris Buschelman #undef __FUNCT__
484b0592601SKris Buschelman #define __FUNCT__ "MatConvert_SuperLU_SeqAIJ"
485b0592601SKris Buschelman int MatConvert_SuperLU_SeqAIJ(Mat A,MatType type,Mat *newmat) {
486fb3e25aaSKris Buschelman   /* This routine is only called to convert an unfactored PETSc-SuperLU matrix */
487fb3e25aaSKris Buschelman   /* to its base PETSc type, so we will ignore 'MatType type'. */
48814b396bbSKris Buschelman   int                  ierr;
489b0592601SKris Buschelman   Mat                  B=*newmat;
490f0c56d0fSKris Buschelman   Mat_SuperLU   *lu=(Mat_SuperLU *)A->spptr;
491b0592601SKris Buschelman 
492b0592601SKris Buschelman   PetscFunctionBegin;
493b0592601SKris Buschelman   if (B != A) {
494b0592601SKris Buschelman     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
495f0c56d0fSKris Buschelman   }
496b0592601SKris Buschelman   /* Reset the original function pointers */
497f0c56d0fSKris Buschelman   B->ops->duplicate        = lu->MatDuplicate;
498b0592601SKris Buschelman   B->ops->view             = lu->MatView;
499b0592601SKris Buschelman   B->ops->assemblyend      = lu->MatAssemblyEnd;
500b0592601SKris Buschelman   B->ops->lufactorsymbolic = lu->MatLUFactorSymbolic;
501b0592601SKris Buschelman   B->ops->destroy          = lu->MatDestroy;
502b0592601SKris Buschelman   /* lu is only a function pointer stash unless we've factored the matrix, which we haven't! */
503b0592601SKris Buschelman   ierr = PetscFree(lu);CHKERRQ(ierr);
504b0592601SKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJ);CHKERRQ(ierr);
505b0592601SKris Buschelman   *newmat = B;
506b0592601SKris Buschelman   PetscFunctionReturn(0);
507b0592601SKris Buschelman }
508b0592601SKris Buschelman EXTERN_C_END
509b0592601SKris Buschelman 
510b0592601SKris Buschelman EXTERN_C_BEGIN
511b0592601SKris Buschelman #undef __FUNCT__
512b0592601SKris Buschelman #define __FUNCT__ "MatConvert_SeqAIJ_SuperLU"
513b0592601SKris Buschelman int MatConvert_SeqAIJ_SuperLU(Mat A,MatType type,Mat *newmat) {
514fb3e25aaSKris Buschelman   /* This routine is only called to convert to MATSUPERLU */
515fb3e25aaSKris Buschelman   /* from MATSEQAIJ, so we will ignore 'MatType type'. */
516b0592601SKris Buschelman   int         ierr;
517b0592601SKris Buschelman   Mat         B=*newmat;
518f0c56d0fSKris Buschelman   Mat_SuperLU *lu;
51914b396bbSKris Buschelman 
52014b396bbSKris Buschelman   PetscFunctionBegin;
521b0592601SKris Buschelman   if (B != A) {
522b0592601SKris Buschelman     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
523b0592601SKris Buschelman   }
52414b396bbSKris Buschelman 
525f0c56d0fSKris Buschelman   ierr = PetscNew(Mat_SuperLU,&lu);CHKERRQ(ierr);
526f0c56d0fSKris Buschelman   lu->MatDuplicate         = A->ops->duplicate;
52714b396bbSKris Buschelman   lu->MatView              = A->ops->view;
52814b396bbSKris Buschelman   lu->MatAssemblyEnd       = A->ops->assemblyend;
529b0592601SKris Buschelman   lu->MatLUFactorSymbolic  = A->ops->lufactorsymbolic;
53014b396bbSKris Buschelman   lu->MatDestroy           = A->ops->destroy;
53114b396bbSKris Buschelman   lu->CleanUpSuperLU       = PETSC_FALSE;
532b0592601SKris Buschelman 
533b0592601SKris Buschelman   B->spptr                 = (void*)lu;
534f0c56d0fSKris Buschelman   B->ops->duplicate        = MatDuplicate_SuperLU;
535f0c56d0fSKris Buschelman   B->ops->view             = MatView_SuperLU;
536f0c56d0fSKris Buschelman   B->ops->assemblyend      = MatAssemblyEnd_SuperLU;
537f0c56d0fSKris Buschelman   B->ops->lufactorsymbolic = MatLUFactorSymbolic_SuperLU;
538f0c56d0fSKris Buschelman   B->ops->destroy          = MatDestroy_SuperLU;
539b0592601SKris Buschelman 
540b0592601SKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_superlu_C",
541b0592601SKris Buschelman                                            "MatConvert_SeqAIJ_SuperLU",MatConvert_SeqAIJ_SuperLU);CHKERRQ(ierr);
542b0592601SKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_superlu_seqaij_C",
543b0592601SKris Buschelman                                            "MatConvert_SuperLU_SeqAIJ",MatConvert_SuperLU_SeqAIJ);CHKERRQ(ierr);
544fb3e25aaSKris Buschelman   PetscLogInfo(0,"Using SuperLU for SeqAIJ LU factorization and solves.");
545fb3e25aaSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)B,MATSUPERLU);CHKERRQ(ierr);
546b0592601SKris Buschelman   *newmat = B;
547b0592601SKris Buschelman   PetscFunctionReturn(0);
548b0592601SKris Buschelman }
549b0592601SKris Buschelman EXTERN_C_END
550b0592601SKris Buschelman 
55124b6179bSKris Buschelman /*MC
552fafad747SKris Buschelman   MATSUPERLU - MATSUPERLU = "superlu" - A matrix type providing direct solvers (LU) for sequential matrices
55324b6179bSKris Buschelman   via the external package SuperLU.
55424b6179bSKris Buschelman 
55524b6179bSKris Buschelman   If SuperLU is installed (see the manual for
55624b6179bSKris Buschelman   instructions on how to declare the existence of external packages),
55724b6179bSKris Buschelman   a matrix type can be constructed which invokes SuperLU solvers.
55824b6179bSKris Buschelman   After calling MatCreate(...,A), simply call MatSetType(A,MATSUPERLU).
55924b6179bSKris Buschelman   This matrix type is only supported for double precision real.
56024b6179bSKris Buschelman 
56124b6179bSKris Buschelman   This matrix inherits from MATSEQAIJ.  As a result, MatSeqAIJSetPreallocation is
56228b08bd3SKris Buschelman   supported for this matrix type.  One can also call MatConvert for an inplace conversion to or from
56328b08bd3SKris Buschelman   the MATSEQAIJ type without data copy.
56424b6179bSKris Buschelman 
56524b6179bSKris Buschelman   Options Database Keys:
5660bad9183SKris Buschelman + -mat_type superlu - sets the matrix type to "superlu" during a call to MatSetFromOptions()
56724b6179bSKris Buschelman - -mat_superlu_ordering <0,1,2,3> - 0: natural ordering,
56824b6179bSKris Buschelman                                     1: MMD applied to A'*A,
56924b6179bSKris Buschelman                                     2: MMD applied to A'+A,
57024b6179bSKris Buschelman                                     3: COLAMD, approximate minimum degree column ordering
57124b6179bSKris Buschelman 
57224b6179bSKris Buschelman    Level: beginner
57324b6179bSKris Buschelman 
57424b6179bSKris Buschelman .seealso: PCLU
57524b6179bSKris Buschelman M*/
57624b6179bSKris Buschelman 
577b0592601SKris Buschelman EXTERN_C_BEGIN
578b0592601SKris Buschelman #undef __FUNCT__
579f0c56d0fSKris Buschelman #define __FUNCT__ "MatCreate_SuperLU"
580f0c56d0fSKris Buschelman int MatCreate_SuperLU(Mat A) {
581b0592601SKris Buschelman   int ierr;
582b0592601SKris Buschelman 
583b0592601SKris Buschelman   PetscFunctionBegin;
5845441df8eSKris Buschelman   /* Change type name before calling MatSetType to force proper construction of SeqAIJ and SUPERLU types */
5855441df8eSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)A,MATSUPERLU);CHKERRQ(ierr);
586b0592601SKris Buschelman   ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr);
587b0592601SKris Buschelman   ierr = MatConvert_SeqAIJ_SuperLU(A,MATSUPERLU,&A);CHKERRQ(ierr);
58814b396bbSKris Buschelman   PetscFunctionReturn(0);
58914b396bbSKris Buschelman }
59014b396bbSKris Buschelman EXTERN_C_END
591