xref: /petsc/src/mat/impls/aij/seq/superlu/superlu.c (revision be5d1d56a128fdbca06f8d9818f1d611ccde2ba2)
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
485ca28756SHong Zhang EXTERN int MatConvert_SuperLU_SeqAIJ(Mat,const MatType,Mat*);
495ca28756SHong Zhang EXTERN int MatConvert_SeqAIJ_SuperLU(Mat,const 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");
141*be5d1d56SKris Buschelman   /* Create the null matrix using MATSEQDENSE explicitly */
142878740d9SKris Buschelman   ierr = MatCreate(A->comm,numRows,numNullCols,numRows,numNullCols,nullMat);CHKERRQ(ierr);
143878740d9SKris Buschelman   ierr = MatSetType(*nullMat,MATSEQDENSE);CHKERRQ(ierr);
144878740d9SKris Buschelman   ierr = MatSeqDenseSetPreallocation(*nullMat,PETSC_NULL);CHKERRQ(ierr);
14514b396bbSKris Buschelman   if (numNullCols == 0) {
14614b396bbSKris Buschelman     ierr = MatAssemblyBegin(*nullMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
14714b396bbSKris Buschelman     ierr = MatAssemblyEnd(*nullMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
14814b396bbSKris Buschelman     PetscFunctionReturn(0);
14914b396bbSKris Buschelman   }
15014b396bbSKris Buschelman #if defined(PETSC_USE_COMPLEX)
15114b396bbSKris Buschelman   nullVals = (doublecomplex*)((Mat_SeqDense*)(*nullMat)->data)->v;
15214b396bbSKris Buschelman #else
15314b396bbSKris Buschelman   nullVals = ((Mat_SeqDense*)(*nullMat)->data)->v;
15414b396bbSKris Buschelman #endif
15514b396bbSKris Buschelman   /* Copy in the columns */
15614b396bbSKris Buschelman   Lstore = (SCformat*)lu->L.Store;
15714b396bbSKris Buschelman   for(block = 0; block <= Lstore->nsuper; block++) {
15814b396bbSKris Buschelman     newRow = Lstore->sup_to_col[block];
15914b396bbSKris Buschelman     size   = Lstore->sup_to_col[block+1] - Lstore->sup_to_col[block];
16014b396bbSKris Buschelman     for(col = Lstore->rowind_colptr[newRow]; col < Lstore->rowind_colptr[newRow+1]; col++) {
16114b396bbSKris Buschelman       newCol = Lstore->rowind[col];
16214b396bbSKris Buschelman       if (newCol >= numRows) {
16314b396bbSKris Buschelman         for(b = 0; b < size; b++)
16414b396bbSKris Buschelman #if defined(PETSC_USE_COMPLEX)
16514b396bbSKris Buschelman           nullVals[(newCol-numRows)*numRows+newRow+b] = ((doublecomplex*)Lstore->nzval)[Lstore->nzval_colptr[newRow+b]+col];
16614b396bbSKris Buschelman #else
16714b396bbSKris Buschelman           nullVals[(newCol-numRows)*numRows+newRow+b] = ((double*)Lstore->nzval)[Lstore->nzval_colptr[newRow+b]+col];
16814b396bbSKris Buschelman #endif
16914b396bbSKris Buschelman       }
17014b396bbSKris Buschelman     }
17114b396bbSKris Buschelman   }
17214b396bbSKris Buschelman   /* Permute rhs to form P^T_c B */
17314b396bbSKris Buschelman   ierr = PetscMalloc(numRows*sizeof(double),&workVals);CHKERRQ(ierr);
17414b396bbSKris Buschelman   for(b = 0; b < numNullCols; b++) {
17514b396bbSKris Buschelman     for(row = 0; row < numRows; row++) workVals[lu->perm_c[row]] = nullVals[b*numRows+row];
17614b396bbSKris Buschelman     for(row = 0; row < numRows; row++) nullVals[b*numRows+row]   = workVals[row];
17714b396bbSKris Buschelman   }
17814b396bbSKris Buschelman   /* Backward solve the upper triangle A x = b */
17914b396bbSKris Buschelman   for(b = 0; b < numNullCols; b++) {
18014b396bbSKris Buschelman #if defined(PETSC_USE_COMPLEX)
18175af56d4SHong Zhang     sp_ztrsv("L","T","U",&lu->L,&lu->U,&nullVals[b*numRows],&stat,&ierr);
18214b396bbSKris Buschelman #else
18375af56d4SHong Zhang     sp_dtrsv("L","T","U",&lu->L,&lu->U,&nullVals[b*numRows],&stat,&ierr);
18414b396bbSKris Buschelman #endif
18514b396bbSKris Buschelman     if (ierr < 0)
18614b396bbSKris Buschelman       SETERRQ1(PETSC_ERR_ARG_WRONG,"The argument %d was invalid",-ierr);
18714b396bbSKris Buschelman   }
18814b396bbSKris Buschelman   ierr = PetscFree(workVals);CHKERRQ(ierr);
18914b396bbSKris Buschelman 
19014b396bbSKris Buschelman   ierr = MatAssemblyBegin(*nullMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
19114b396bbSKris Buschelman   ierr = MatAssemblyEnd(*nullMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
19214b396bbSKris Buschelman   PetscFunctionReturn(0);
19314b396bbSKris Buschelman }
19475af56d4SHong Zhang #endif
19514b396bbSKris Buschelman 
19614b396bbSKris Buschelman #undef __FUNCT__
197f0c56d0fSKris Buschelman #define __FUNCT__ "MatSolve_SuperLU"
198f0c56d0fSKris Buschelman int MatSolve_SuperLU(Mat A,Vec b,Vec x)
19914b396bbSKris Buschelman {
200f0c56d0fSKris Buschelman   Mat_SuperLU   *lu = (Mat_SuperLU*)A->spptr;
20175af56d4SHong Zhang   PetscScalar   *barray,*xarray;
2028914a3f7SHong Zhang   int           ierr,info,i;
20375af56d4SHong Zhang   SuperLUStat_t stat;
2048914a3f7SHong Zhang   double        ferr,berr;
20514b396bbSKris Buschelman 
20614b396bbSKris Buschelman   PetscFunctionBegin;
207937a6b0eSHong Zhang   if ( lu->lwork == -1 ) {
208937a6b0eSHong Zhang     PetscFunctionReturn(0);
209937a6b0eSHong Zhang   }
21075af56d4SHong Zhang   lu->B.ncol = 1;   /* Set the number of right-hand side */
21175af56d4SHong Zhang   ierr = VecGetArray(b,&barray);CHKERRQ(ierr);
21275af56d4SHong Zhang   ierr = VecGetArray(x,&xarray);CHKERRQ(ierr);
2135fe6150dSHong Zhang 
2145fe6150dSHong Zhang #if defined(PETSC_USE_COMPLEX)
2155fe6150dSHong Zhang   ((DNformat*)lu->B.Store)->nzval = (doublecomplex*)barray;
2165fe6150dSHong Zhang   ((DNformat*)lu->X.Store)->nzval = (doublecomplex*)xarray;
2175fe6150dSHong Zhang #else
2185fe6150dSHong Zhang   ((DNformat*)lu->B.Store)->nzval = barray;
21975af56d4SHong Zhang   ((DNformat*)lu->X.Store)->nzval = xarray;
2205fe6150dSHong Zhang #endif
22175af56d4SHong Zhang 
22275af56d4SHong Zhang   /* Initialize the statistics variables. */
22375af56d4SHong Zhang   StatInit(&stat);
22475af56d4SHong Zhang 
22575af56d4SHong Zhang   lu->options.Fact  = FACTORED; /* Indicate the factored form of A is supplied. */
22675af56d4SHong Zhang   lu->options.Trans = TRANS;
2278914a3f7SHong Zhang #if defined(PETSC_USE_COMPLEX)
2288914a3f7SHong Zhang   zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
2298914a3f7SHong Zhang            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
2308914a3f7SHong Zhang            &lu->mem_usage, &stat, &info);
2318914a3f7SHong Zhang #else
23275af56d4SHong Zhang   dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
23375af56d4SHong Zhang            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
23475af56d4SHong Zhang            &lu->mem_usage, &stat, &info);
2358914a3f7SHong Zhang #endif
23675af56d4SHong Zhang   ierr = VecRestoreArray(b,&barray);CHKERRQ(ierr);
23775af56d4SHong Zhang   ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
23875af56d4SHong Zhang 
23975af56d4SHong Zhang   if ( info == 0 || info == lu->A.ncol+1 ) {
24075af56d4SHong Zhang     if ( lu->options.IterRefine ) {
2418914a3f7SHong Zhang       ierr = PetscPrintf(PETSC_COMM_SELF,"Iterative Refinement:\n");
2428914a3f7SHong Zhang       ierr = PetscPrintf(PETSC_COMM_SELF,"  %8s%8s%16s%16s\n", "rhs", "Steps", "FERR", "BERR");
24375af56d4SHong Zhang       for (i = 0; i < 1; ++i)
2448914a3f7SHong Zhang         ierr = PetscPrintf(PETSC_COMM_SELF,"  %8d%8d%16e%16e\n", i+1, stat.RefineSteps, ferr, berr);
24575af56d4SHong Zhang     }
2468914a3f7SHong Zhang   } else if ( info > 0 ){
2478914a3f7SHong Zhang     if ( lu->lwork == -1 ) {
2488914a3f7SHong Zhang       ierr = PetscPrintf(PETSC_COMM_SELF,"  ** Estimated memory: %d bytes\n", info - lu->A.ncol);
2498914a3f7SHong Zhang     } else {
2508914a3f7SHong Zhang       ierr = PetscPrintf(PETSC_COMM_SELF,"  Warning: gssvx() returns info %d\n",info);
2518914a3f7SHong Zhang     }
2528914a3f7SHong Zhang   } else if (info < 0){
2538914a3f7SHong Zhang     SETERRQ2(1, "info = %d, the %d-th argument in gssvx() had an illegal value", info,-info);
25414b396bbSKris Buschelman   }
25514b396bbSKris Buschelman 
2568914a3f7SHong Zhang   if ( lu->options.PrintStat ) {
2578914a3f7SHong Zhang     ierr = PetscPrintf(PETSC_COMM_SELF,"MatSolve__SuperLU():\n");
2588914a3f7SHong Zhang     StatPrint(&stat);
2598914a3f7SHong Zhang   }
26075af56d4SHong Zhang   StatFree(&stat);
26175af56d4SHong Zhang   PetscFunctionReturn(0);
26275af56d4SHong Zhang }
26314b396bbSKris Buschelman 
26414b396bbSKris Buschelman #undef __FUNCT__
265f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorNumeric_SuperLU"
266f0c56d0fSKris Buschelman int MatLUFactorNumeric_SuperLU(Mat A,Mat *F)
26714b396bbSKris Buschelman {
26814b396bbSKris Buschelman   Mat_SeqAIJ    *aa = (Mat_SeqAIJ*)(A)->data;
269f0c56d0fSKris Buschelman   Mat_SuperLU   *lu = (Mat_SuperLU*)(*F)->spptr;
27075af56d4SHong Zhang   int           ierr,info;
27114b396bbSKris Buschelman   PetscTruth    flag;
27275af56d4SHong Zhang   SuperLUStat_t stat;
27375af56d4SHong Zhang   double        ferr, berr;
2748914a3f7SHong Zhang   NCformat      *Ustore;
2758914a3f7SHong Zhang   SCformat      *Lstore;
27614b396bbSKris Buschelman 
27714b396bbSKris Buschelman   PetscFunctionBegin;
2789ce50919SHong Zhang   if (lu->flg == SAME_NONZERO_PATTERN){ /* successing numerical factorization */
2799ce50919SHong Zhang     lu->options.Fact = SamePattern;
28075af56d4SHong Zhang     /* Ref: ~SuperLU_3.0/EXAMPLE/dlinsolx2.c */
28175af56d4SHong Zhang     Destroy_SuperMatrix_Store(&lu->A);
2829ce50919SHong Zhang     if ( lu->lwork >= 0 ) {
28375af56d4SHong Zhang       Destroy_SuperNode_Matrix(&lu->L);
28475af56d4SHong Zhang       Destroy_CompCol_Matrix(&lu->U);
28575af56d4SHong Zhang       lu->options.Fact = SamePattern;
28614b396bbSKris Buschelman     }
28775af56d4SHong Zhang   }
28814b396bbSKris Buschelman 
28975af56d4SHong Zhang   /* Create the SuperMatrix for lu->A=A^T:
29075af56d4SHong Zhang        Since SuperLU likes column-oriented matrices,we pass it the transpose,
29175af56d4SHong Zhang        and then solve A^T X = B in MatSolve(). */
29275af56d4SHong Zhang #if defined(PETSC_USE_COMPLEX)
2935fe6150dSHong Zhang   zCreate_CompCol_Matrix(&lu->A,A->n,A->m,aa->nz,(doublecomplex*)aa->a,aa->j,aa->i,
29475af56d4SHong Zhang                            SLU_NC,SLU_Z,SLU_GE);
29575af56d4SHong Zhang #else
29675af56d4SHong Zhang   dCreate_CompCol_Matrix(&lu->A,A->n,A->m,aa->nz,aa->a,aa->j,aa->i,
29775af56d4SHong Zhang                            SLU_NC,SLU_D,SLU_GE);
29875af56d4SHong Zhang #endif
29975af56d4SHong Zhang 
30075af56d4SHong Zhang   /* Initialize the statistics variables. */
30175af56d4SHong Zhang   StatInit(&stat);
30275af56d4SHong Zhang 
3039ce50919SHong Zhang   /* Numerical factorization */
30475af56d4SHong Zhang   lu->B.ncol = 0;  /* Indicate not to solve the system */
3058914a3f7SHong Zhang #if defined(PETSC_USE_COMPLEX)
3068914a3f7SHong Zhang    zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
3078914a3f7SHong Zhang            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
3088914a3f7SHong Zhang            &lu->mem_usage, &stat, &info);
3098914a3f7SHong Zhang #else
31075af56d4SHong Zhang   dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
31175af56d4SHong Zhang            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
31275af56d4SHong Zhang            &lu->mem_usage, &stat, &info);
3138914a3f7SHong Zhang #endif
31475af56d4SHong Zhang   if ( info == 0 || info == lu->A.ncol+1 ) {
3158914a3f7SHong Zhang     if ( lu->options.PivotGrowth )
3168914a3f7SHong Zhang       ierr = PetscPrintf(PETSC_COMM_SELF,"  Recip. pivot growth = %e\n", lu->rpg);
31775af56d4SHong Zhang     if ( lu->options.ConditionNumber )
3188914a3f7SHong Zhang       ierr = PetscPrintf(PETSC_COMM_SELF,"  Recip. condition number = %e\n", lu->rcond);
3198914a3f7SHong Zhang   } else if ( info > 0 ){
3208914a3f7SHong Zhang     if ( lu->lwork == -1 ) {
3218914a3f7SHong Zhang       ierr = PetscPrintf(PETSC_COMM_SELF,"  ** Estimated memory: %d bytes\n", info - lu->A.ncol);
3228914a3f7SHong Zhang     } else {
3238914a3f7SHong Zhang       ierr = PetscPrintf(PETSC_COMM_SELF,"  Warning: gssvx() returns info %d\n",info);
3248914a3f7SHong Zhang     }
325937a6b0eSHong Zhang   } else { /* info < 0 */
3268914a3f7SHong Zhang     SETERRQ2(1, "info = %d, the %d-th argument in gssvx() had an illegal value", info,-info);
32775af56d4SHong Zhang   }
32875af56d4SHong Zhang 
3298914a3f7SHong Zhang   if ( lu->options.PrintStat ) {
3308914a3f7SHong Zhang     ierr = PetscPrintf(PETSC_COMM_SELF,"MatLUFactorNumeric_SuperLU():\n");
3318914a3f7SHong Zhang     StatPrint(&stat);
3328914a3f7SHong Zhang     Lstore = (SCformat *) lu->L.Store;
3338914a3f7SHong Zhang     Ustore = (NCformat *) lu->U.Store;
3348914a3f7SHong Zhang     ierr = PetscPrintf(PETSC_COMM_SELF,"  No of nonzeros in factor L = %d\n", Lstore->nnz);
3358914a3f7SHong Zhang     ierr = PetscPrintf(PETSC_COMM_SELF,"  No of nonzeros in factor U = %d\n", Ustore->nnz);
3368914a3f7SHong Zhang     ierr = PetscPrintf(PETSC_COMM_SELF,"  No of nonzeros in L+U = %d\n", Lstore->nnz + Ustore->nnz - lu->A.ncol);
3378914a3f7SHong Zhang     ierr = PetscPrintf(PETSC_COMM_SELF,"  L\\U MB %.3f\ttotal MB needed %.3f\texpansions %d\n",
3388914a3f7SHong Zhang 	       lu->mem_usage.for_lu/1e6, lu->mem_usage.total_needed/1e6,
3398914a3f7SHong Zhang 	       lu->mem_usage.expansions);
3408914a3f7SHong Zhang   }
34175af56d4SHong Zhang   StatFree(&stat);
34275af56d4SHong Zhang 
34375af56d4SHong Zhang   lu->flg = SAME_NONZERO_PATTERN;
34475af56d4SHong Zhang   PetscFunctionReturn(0);
34514b396bbSKris Buschelman }
34614b396bbSKris Buschelman 
34714b396bbSKris Buschelman /*
34814b396bbSKris Buschelman    Note the r permutation is ignored
34914b396bbSKris Buschelman */
35014b396bbSKris Buschelman #undef __FUNCT__
351f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_SuperLU"
352f0c56d0fSKris Buschelman int MatLUFactorSymbolic_SuperLU(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F)
35314b396bbSKris Buschelman {
35414b396bbSKris Buschelman   Mat          B;
355f0c56d0fSKris Buschelman   Mat_SuperLU  *lu;
3569ce50919SHong Zhang   int          ierr,m=A->m,n=A->n,indx;
3579ce50919SHong Zhang   PetscTruth   flg;
3585ca28756SHong Zhang   const char   *colperm[]={"NATURAL","MMD_ATA","MMD_AT_PLUS_A","COLAMD"}; /* MY_PERMC - not supported by the petsc interface yet */
3595ca28756SHong Zhang   const char   *iterrefine[]={"NOREFINE", "SINGLE", "DOUBLE", "EXTRA"};
3605ca28756SHong Zhang   const char   *rowperm[]={"NOROWPERM", "LargeDiag"}; /* MY_PERMC - not supported by the petsc interface yet */
36114b396bbSKris Buschelman 
36214b396bbSKris Buschelman   PetscFunctionBegin;
36314b396bbSKris Buschelman 
364899d7b4fSKris Buschelman   ierr = MatCreate(A->comm,A->m,A->n,PETSC_DETERMINE,PETSC_DETERMINE,&B);CHKERRQ(ierr);
365*be5d1d56SKris Buschelman   ierr = MatSetType(B,A->type_name);CHKERRQ(ierr);
366899d7b4fSKris Buschelman   ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
367f0c56d0fSKris Buschelman 
368f0c56d0fSKris Buschelman   B->ops->lufactornumeric = MatLUFactorNumeric_SuperLU;
369f0c56d0fSKris Buschelman   B->ops->solve           = MatSolve_SuperLU;
37014b396bbSKris Buschelman   B->factor               = FACTOR_LU;
37194b7f48cSBarry Smith   B->assembled            = PETSC_TRUE;  /* required by -ksp_view */
37214b396bbSKris Buschelman 
373f0c56d0fSKris Buschelman   lu = (Mat_SuperLU*)(B->spptr);
3749ce50919SHong Zhang 
3759ce50919SHong Zhang   /* Set SuperLU options */
3769ce50919SHong Zhang     /* the default values for options argument:
3779ce50919SHong Zhang 	options.Fact = DOFACT;
3789ce50919SHong Zhang         options.Equil = YES;
3799ce50919SHong Zhang     	options.ColPerm = COLAMD;
3809ce50919SHong Zhang 	options.DiagPivotThresh = 1.0;
3819ce50919SHong Zhang     	options.Trans = NOTRANS;
3829ce50919SHong Zhang     	options.IterRefine = NOREFINE;
3839ce50919SHong Zhang     	options.SymmetricMode = NO;
3849ce50919SHong Zhang     	options.PivotGrowth = NO;
3859ce50919SHong Zhang     	options.ConditionNumber = NO;
3869ce50919SHong Zhang     	options.PrintStat = YES;
3879ce50919SHong Zhang     */
3889ce50919SHong Zhang   set_default_options(&lu->options);
3898914a3f7SHong Zhang   /* equilibration causes error in solve(), thus not supported here. See dgssvx.c for possible reason. */
3908914a3f7SHong Zhang   lu->options.Equil = NO;
3919ce50919SHong Zhang   lu->options.PrintStat = NO;
3928914a3f7SHong Zhang   lu->lwork = 0;   /* allocate space internally by system malloc */
3939ce50919SHong Zhang 
3949ce50919SHong Zhang   ierr = PetscOptionsBegin(A->comm,A->prefix,"SuperLU Options","Mat");CHKERRQ(ierr);
3959ce50919SHong Zhang   /*
3968914a3f7SHong Zhang   ierr = PetscOptionsLogical("-mat_superlu_equil","Equil","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr);
3978914a3f7SHong Zhang   if (flg) lu->options.Equil = YES; -- not supported by the interface !!!
3989ce50919SHong Zhang   */
3998914a3f7SHong Zhang   ierr = PetscOptionsEList("-mat_superlu_colperm","ColPerm","None",colperm,4,colperm[3],&indx,&flg);CHKERRQ(ierr);
4009ce50919SHong Zhang   if (flg) {lu->options.ColPerm = (colperm_t)indx;}
4018914a3f7SHong Zhang   ierr = PetscOptionsEList("-mat_superlu_iterrefine","IterRefine","None",iterrefine,4,iterrefine[0],&indx,&flg);CHKERRQ(ierr);
4029ce50919SHong Zhang   if (flg) { lu->options.IterRefine = (IterRefine_t)indx;}
4038914a3f7SHong Zhang   ierr = PetscOptionsLogical("-mat_superlu_symmetricmode","SymmetricMode","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr);
4049ce50919SHong Zhang   if (flg) lu->options.SymmetricMode = YES;
4058914a3f7SHong Zhang   ierr = PetscOptionsReal("-mat_superlu_diagpivotthresh","DiagPivotThresh","None",lu->options.DiagPivotThresh,&lu->options.DiagPivotThresh,PETSC_NULL);CHKERRQ(ierr);
4068914a3f7SHong Zhang   ierr = PetscOptionsLogical("-mat_superlu_pivotgrowth","PivotGrowth","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr);
4079ce50919SHong Zhang   if (flg) lu->options.PivotGrowth = YES;
4088914a3f7SHong Zhang   ierr = PetscOptionsLogical("-mat_superlu_conditionnumber","ConditionNumber","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr);
4099ce50919SHong Zhang   if (flg) lu->options.ConditionNumber = YES;
4108914a3f7SHong Zhang   ierr = PetscOptionsEList("-mat_superlu_rowperm","rowperm","None",rowperm,2,rowperm[0],&indx,&flg);CHKERRQ(ierr);
4119ce50919SHong Zhang   if (flg) {lu->options.RowPerm = (rowperm_t)indx;}
4128914a3f7SHong Zhang   ierr = PetscOptionsLogical("-mat_superlu_replacetinypivot","ReplaceTinyPivot","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr);
4139ce50919SHong Zhang   if (flg) lu->options.ReplaceTinyPivot = YES;
4148914a3f7SHong Zhang   ierr = PetscOptionsLogical("-mat_superlu_printstat","PrintStat","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr);
4159ce50919SHong Zhang   if (flg) lu->options.PrintStat = YES;
4168914a3f7SHong Zhang   ierr = PetscOptionsInt("-mat_superlu_lwork","size of work array in bytes used by factorization","None",lu->lwork,&lu->lwork,PETSC_NULL);CHKERRQ(ierr);
4175fe6150dSHong Zhang   if (lu->lwork > 0 ){
4185fe6150dSHong Zhang     ierr = PetscMalloc(lu->lwork,&lu->work);CHKERRQ(ierr);
4195fe6150dSHong Zhang   } else if (lu->lwork != 0 && lu->lwork != -1){
420937a6b0eSHong Zhang     ierr = PetscPrintf(PETSC_COMM_SELF,"   Warning: lwork %d is not supported by SUPERLU. The default lwork=0 is used.\n",lu->lwork);
4218914a3f7SHong Zhang     lu->lwork = 0;
4228914a3f7SHong Zhang   }
4239ce50919SHong Zhang   PetscOptionsEnd();
4249ce50919SHong Zhang 
42575af56d4SHong Zhang #ifdef SUPERLU2
4262ecb5974SKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatCreateNull","MatCreateNull_SuperLU",
427f0c56d0fSKris Buschelman                                     (void(*)(void))MatCreateNull_SuperLU);CHKERRQ(ierr);
42875af56d4SHong Zhang #endif
42914b396bbSKris Buschelman 
43075af56d4SHong Zhang   /* Allocate spaces (notice sizes are for the transpose) */
4319ce50919SHong Zhang   ierr = PetscMalloc(m*sizeof(int),&lu->etree);CHKERRQ(ierr);
4329ce50919SHong Zhang   ierr = PetscMalloc(n*sizeof(int),&lu->perm_r);CHKERRQ(ierr);
4339ce50919SHong Zhang   ierr = PetscMalloc(m*sizeof(int),&lu->perm_c);CHKERRQ(ierr);
4349ce50919SHong Zhang   ierr = PetscMalloc(n*sizeof(int),&lu->R);CHKERRQ(ierr);
4359ce50919SHong Zhang   ierr = PetscMalloc(m*sizeof(int),&lu->C);CHKERRQ(ierr);
43675af56d4SHong Zhang 
4379ce50919SHong Zhang   /* create rhs and solution x without allocate space for .Store */
4385fe6150dSHong Zhang #if defined(PETSC_USE_COMPLEX)
439937a6b0eSHong Zhang   zCreate_Dense_Matrix(&lu->B, m, 1, PETSC_NULL, m, SLU_DN, SLU_Z, SLU_GE);
440937a6b0eSHong Zhang   zCreate_Dense_Matrix(&lu->X, m, 1, PETSC_NULL, m, SLU_DN, SLU_Z, SLU_GE);
4415fe6150dSHong Zhang #else
44275af56d4SHong Zhang   dCreate_Dense_Matrix(&lu->B, m, 1, PETSC_NULL, m, SLU_DN, SLU_D, SLU_GE);
44375af56d4SHong Zhang   dCreate_Dense_Matrix(&lu->X, m, 1, PETSC_NULL, m, SLU_DN, SLU_D, SLU_GE);
4445fe6150dSHong Zhang #endif
44514b396bbSKris Buschelman 
44614b396bbSKris Buschelman   lu->flg            = DIFFERENT_NONZERO_PATTERN;
447e740cb95SKris Buschelman   lu->CleanUpSuperLU = PETSC_TRUE;
448899d7b4fSKris Buschelman 
449899d7b4fSKris Buschelman   *F = B;
4509ce50919SHong Zhang   PetscLogObjectMemory(B,(A->m+A->n)*sizeof(int)+sizeof(Mat_SuperLU));
45114b396bbSKris Buschelman   PetscFunctionReturn(0);
45214b396bbSKris Buschelman }
45314b396bbSKris Buschelman 
45494b7f48cSBarry Smith /* used by -ksp_view */
45514b396bbSKris Buschelman #undef __FUNCT__
456f0c56d0fSKris Buschelman #define __FUNCT__ "MatFactorInfo_SuperLU"
457f0c56d0fSKris Buschelman int MatFactorInfo_SuperLU(Mat A,PetscViewer viewer)
45814b396bbSKris Buschelman {
459f0c56d0fSKris Buschelman   Mat_SuperLU       *lu= (Mat_SuperLU*)A->spptr;
460f0c56d0fSKris Buschelman   int               ierr;
4619ce50919SHong Zhang   superlu_options_t options;
462f0c56d0fSKris Buschelman 
463f0c56d0fSKris Buschelman   PetscFunctionBegin;
4649ce50919SHong Zhang   /* check if matrix is superlu_dist type */
4659ce50919SHong Zhang   if (A->ops->solve != MatSolve_SuperLU) PetscFunctionReturn(0);
4669ce50919SHong Zhang 
4679ce50919SHong Zhang   options = lu->options;
468f0c56d0fSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"SuperLU run parameters:\n");CHKERRQ(ierr);
4699ce50919SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  Equil: %s\n",(options.Equil != NO) ? "YES": "NO");CHKERRQ(ierr);
4709ce50919SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ColPerm: %d\n",options.ColPerm);CHKERRQ(ierr);
4719ce50919SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  IterRefine: %d\n",options.IterRefine);CHKERRQ(ierr);
4729ce50919SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  SymmetricMode: %s\n",(options.SymmetricMode != NO) ? "YES": "NO");CHKERRQ(ierr);
4739ce50919SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  DiagPivotThresh: %g\n",options.DiagPivotThresh);CHKERRQ(ierr);
4749ce50919SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  PivotGrowth: %s\n",(options.PivotGrowth != NO) ? "YES": "NO");CHKERRQ(ierr);
4759ce50919SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ConditionNumber: %s\n",(options.ConditionNumber != NO) ? "YES": "NO");CHKERRQ(ierr);
4769ce50919SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  RowPerm: %d\n",options.RowPerm);CHKERRQ(ierr);
4779ce50919SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ReplaceTinyPivot: %s\n",(options.ReplaceTinyPivot != NO) ? "YES": "NO");CHKERRQ(ierr);
4789ce50919SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  PrintStat: %s\n",(options.PrintStat != NO) ? "YES": "NO");CHKERRQ(ierr);
4798914a3f7SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  lwork: %d\n",lu->lwork);CHKERRQ(ierr);
4809ce50919SHong Zhang 
481f0c56d0fSKris Buschelman   PetscFunctionReturn(0);
482f0c56d0fSKris Buschelman }
483f0c56d0fSKris Buschelman 
484f0c56d0fSKris Buschelman #undef __FUNCT__
485f0c56d0fSKris Buschelman #define __FUNCT__ "MatDuplicate_SuperLU"
486f0c56d0fSKris Buschelman int MatDuplicate_SuperLU(Mat A, MatDuplicateOption op, Mat *M) {
48714b396bbSKris Buschelman   int         ierr;
4888f340917SKris Buschelman   Mat_SuperLU *lu=(Mat_SuperLU *)A->spptr;
4898f340917SKris Buschelman 
49014b396bbSKris Buschelman   PetscFunctionBegin;
4918f340917SKris Buschelman   ierr = (*lu->MatDuplicate)(A,op,M);CHKERRQ(ierr);
492f0c56d0fSKris Buschelman   ierr = MatConvert_SeqAIJ_SuperLU(*M,MATSUPERLU,M);CHKERRQ(ierr);
4933f953163SKris Buschelman   ierr = PetscMemcpy((*M)->spptr,lu,sizeof(Mat_SuperLU));CHKERRQ(ierr);
49414b396bbSKris Buschelman   PetscFunctionReturn(0);
49514b396bbSKris Buschelman }
49614b396bbSKris Buschelman 
49714b396bbSKris Buschelman EXTERN_C_BEGIN
49814b396bbSKris Buschelman #undef __FUNCT__
499b0592601SKris Buschelman #define __FUNCT__ "MatConvert_SuperLU_SeqAIJ"
5008e9aea5cSBarry Smith int MatConvert_SuperLU_SeqAIJ(Mat A,const MatType type,Mat *newmat) {
501fb3e25aaSKris Buschelman   /* This routine is only called to convert an unfactored PETSc-SuperLU matrix */
502fb3e25aaSKris Buschelman   /* to its base PETSc type, so we will ignore 'MatType type'. */
50314b396bbSKris Buschelman   int                  ierr;
504b0592601SKris Buschelman   Mat                  B=*newmat;
505f0c56d0fSKris Buschelman   Mat_SuperLU   *lu=(Mat_SuperLU *)A->spptr;
506b0592601SKris Buschelman 
507b0592601SKris Buschelman   PetscFunctionBegin;
508b0592601SKris Buschelman   if (B != A) {
509b0592601SKris Buschelman     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
510f0c56d0fSKris Buschelman   }
511b0592601SKris Buschelman   /* Reset the original function pointers */
512f0c56d0fSKris Buschelman   B->ops->duplicate        = lu->MatDuplicate;
513b0592601SKris Buschelman   B->ops->view             = lu->MatView;
514b0592601SKris Buschelman   B->ops->assemblyend      = lu->MatAssemblyEnd;
515b0592601SKris Buschelman   B->ops->lufactorsymbolic = lu->MatLUFactorSymbolic;
516b0592601SKris Buschelman   B->ops->destroy          = lu->MatDestroy;
517b0592601SKris Buschelman   /* lu is only a function pointer stash unless we've factored the matrix, which we haven't! */
518b0592601SKris Buschelman   ierr = PetscFree(lu);CHKERRQ(ierr);
519b0592601SKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJ);CHKERRQ(ierr);
520b0592601SKris Buschelman   *newmat = B;
521b0592601SKris Buschelman   PetscFunctionReturn(0);
522b0592601SKris Buschelman }
523b0592601SKris Buschelman EXTERN_C_END
524b0592601SKris Buschelman 
525b0592601SKris Buschelman EXTERN_C_BEGIN
526b0592601SKris Buschelman #undef __FUNCT__
527b0592601SKris Buschelman #define __FUNCT__ "MatConvert_SeqAIJ_SuperLU"
5288e9aea5cSBarry Smith int MatConvert_SeqAIJ_SuperLU(Mat A,const MatType type,Mat *newmat) {
529fb3e25aaSKris Buschelman   /* This routine is only called to convert to MATSUPERLU */
530fb3e25aaSKris Buschelman   /* from MATSEQAIJ, so we will ignore 'MatType type'. */
531b0592601SKris Buschelman   int         ierr;
532b0592601SKris Buschelman   Mat         B=*newmat;
533f0c56d0fSKris Buschelman   Mat_SuperLU *lu;
53414b396bbSKris Buschelman 
53514b396bbSKris Buschelman   PetscFunctionBegin;
536b0592601SKris Buschelman   if (B != A) {
537b0592601SKris Buschelman     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
538b0592601SKris Buschelman   }
53914b396bbSKris Buschelman 
540f0c56d0fSKris Buschelman   ierr = PetscNew(Mat_SuperLU,&lu);CHKERRQ(ierr);
541f0c56d0fSKris Buschelman   lu->MatDuplicate         = A->ops->duplicate;
54214b396bbSKris Buschelman   lu->MatView              = A->ops->view;
54314b396bbSKris Buschelman   lu->MatAssemblyEnd       = A->ops->assemblyend;
544b0592601SKris Buschelman   lu->MatLUFactorSymbolic  = A->ops->lufactorsymbolic;
54514b396bbSKris Buschelman   lu->MatDestroy           = A->ops->destroy;
54614b396bbSKris Buschelman   lu->CleanUpSuperLU       = PETSC_FALSE;
547b0592601SKris Buschelman 
548b0592601SKris Buschelman   B->spptr                 = (void*)lu;
549f0c56d0fSKris Buschelman   B->ops->duplicate        = MatDuplicate_SuperLU;
550f0c56d0fSKris Buschelman   B->ops->view             = MatView_SuperLU;
551f0c56d0fSKris Buschelman   B->ops->assemblyend      = MatAssemblyEnd_SuperLU;
552f0c56d0fSKris Buschelman   B->ops->lufactorsymbolic = MatLUFactorSymbolic_SuperLU;
553f0c56d0fSKris Buschelman   B->ops->destroy          = MatDestroy_SuperLU;
554b0592601SKris Buschelman 
555b0592601SKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_superlu_C",
556b0592601SKris Buschelman                                            "MatConvert_SeqAIJ_SuperLU",MatConvert_SeqAIJ_SuperLU);CHKERRQ(ierr);
557b0592601SKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_superlu_seqaij_C",
558b0592601SKris Buschelman                                            "MatConvert_SuperLU_SeqAIJ",MatConvert_SuperLU_SeqAIJ);CHKERRQ(ierr);
559fb3e25aaSKris Buschelman   PetscLogInfo(0,"Using SuperLU for SeqAIJ LU factorization and solves.");
560fb3e25aaSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)B,MATSUPERLU);CHKERRQ(ierr);
561b0592601SKris Buschelman   *newmat = B;
562b0592601SKris Buschelman   PetscFunctionReturn(0);
563b0592601SKris Buschelman }
564b0592601SKris Buschelman EXTERN_C_END
565b0592601SKris Buschelman 
56624b6179bSKris Buschelman /*MC
567fafad747SKris Buschelman   MATSUPERLU - MATSUPERLU = "superlu" - A matrix type providing direct solvers (LU) for sequential matrices
56824b6179bSKris Buschelman   via the external package SuperLU.
56924b6179bSKris Buschelman 
57024b6179bSKris Buschelman   If SuperLU is installed (see the manual for
57124b6179bSKris Buschelman   instructions on how to declare the existence of external packages),
57224b6179bSKris Buschelman   a matrix type can be constructed which invokes SuperLU solvers.
57324b6179bSKris Buschelman   After calling MatCreate(...,A), simply call MatSetType(A,MATSUPERLU).
57424b6179bSKris Buschelman   This matrix type is only supported for double precision real.
57524b6179bSKris Buschelman 
57624b6179bSKris Buschelman   This matrix inherits from MATSEQAIJ.  As a result, MatSeqAIJSetPreallocation is
57728b08bd3SKris Buschelman   supported for this matrix type.  One can also call MatConvert for an inplace conversion to or from
57828b08bd3SKris Buschelman   the MATSEQAIJ type without data copy.
57924b6179bSKris Buschelman 
58024b6179bSKris Buschelman   Options Database Keys:
5810bad9183SKris Buschelman + -mat_type superlu - sets the matrix type to "superlu" during a call to MatSetFromOptions()
58224b6179bSKris Buschelman - -mat_superlu_ordering <0,1,2,3> - 0: natural ordering,
58324b6179bSKris Buschelman                                     1: MMD applied to A'*A,
58424b6179bSKris Buschelman                                     2: MMD applied to A'+A,
58524b6179bSKris Buschelman                                     3: COLAMD, approximate minimum degree column ordering
58624b6179bSKris Buschelman 
58724b6179bSKris Buschelman    Level: beginner
58824b6179bSKris Buschelman 
58924b6179bSKris Buschelman .seealso: PCLU
59024b6179bSKris Buschelman M*/
59124b6179bSKris Buschelman 
592b0592601SKris Buschelman EXTERN_C_BEGIN
593b0592601SKris Buschelman #undef __FUNCT__
594f0c56d0fSKris Buschelman #define __FUNCT__ "MatCreate_SuperLU"
595f0c56d0fSKris Buschelman int MatCreate_SuperLU(Mat A) {
596b0592601SKris Buschelman   int ierr;
597b0592601SKris Buschelman 
598b0592601SKris Buschelman   PetscFunctionBegin;
5995441df8eSKris Buschelman   /* Change type name before calling MatSetType to force proper construction of SeqAIJ and SUPERLU types */
6005441df8eSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)A,MATSUPERLU);CHKERRQ(ierr);
601b0592601SKris Buschelman   ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr);
602b0592601SKris Buschelman   ierr = MatConvert_SeqAIJ_SuperLU(A,MATSUPERLU,&A);CHKERRQ(ierr);
60314b396bbSKris Buschelman   PetscFunctionReturn(0);
60414b396bbSKris Buschelman }
60514b396bbSKris Buschelman EXTERN_C_END
606