xref: /petsc/src/mat/impls/aij/seq/superlu/superlu.c (revision 460c49030a5280f0e471e5c8a7d2f060f5488c87)
114b396bbSKris Buschelman /*$Id: superlu.c,v 1.10 2001/08/15 15:56:50 bsmith Exp $*/
214b396bbSKris Buschelman 
314b396bbSKris Buschelman /*
414b396bbSKris Buschelman         Provides an interface to the SuperLU sparse solver
514b396bbSKris Buschelman           Modified for SuperLU 2.0 by Matthew Knepley
614b396bbSKris Buschelman */
714b396bbSKris Buschelman 
814b396bbSKris Buschelman #include "src/mat/impls/aij/seq/aij.h"
914b396bbSKris Buschelman 
1014b396bbSKris Buschelman EXTERN_C_BEGIN
1114b396bbSKris Buschelman #if defined(PETSC_USE_COMPLEX)
1214b396bbSKris Buschelman #include "zsp_defs.h"
1314b396bbSKris Buschelman #else
1414b396bbSKris Buschelman #include "dsp_defs.h"
1514b396bbSKris Buschelman #endif
1614b396bbSKris Buschelman #include "util.h"
1714b396bbSKris Buschelman EXTERN_C_END
1814b396bbSKris Buschelman 
1914b396bbSKris Buschelman typedef struct {
2014b396bbSKris Buschelman   SuperMatrix  A,B,AC,L,U;
2114b396bbSKris Buschelman   int          *perm_r,*perm_c,ispec,relax,panel_size;
2214b396bbSKris Buschelman   double       pivot_threshold;
2314b396bbSKris Buschelman   NCformat     *store;
2414b396bbSKris Buschelman   MatStructure flg;
2514b396bbSKris Buschelman   PetscTruth   SuperluMatOdering;
2614b396bbSKris Buschelman 
2714b396bbSKris Buschelman   /* A few function pointers for inheritance */
2814b396bbSKris Buschelman   int (*MatView)(Mat,PetscViewer);
2914b396bbSKris Buschelman   int (*MatAssemblyEnd)(Mat,MatAssemblyType);
30b0592601SKris Buschelman   int (*MatLUFactorSymbolic)(Mat,IS,IS,MatFactorInfo*,Mat*);
3114b396bbSKris Buschelman   int (*MatDestroy)(Mat);
3214b396bbSKris Buschelman 
3314b396bbSKris Buschelman   /* Flag to clean up (non-global) SuperLU objects during Destroy */
3414b396bbSKris Buschelman   PetscTruth CleanUpSuperLU;
3514b396bbSKris Buschelman } Mat_SeqAIJ_SuperLU;
3614b396bbSKris Buschelman 
3714b396bbSKris Buschelman 
380e3434eeSKris Buschelman EXTERN int MatSeqAIJFactorInfo_SuperLU(Mat,PetscViewer);
39b0592601SKris Buschelman EXTERN int MatLUFactorSymbolic_SeqAIJ_SuperLU(Mat,IS,IS,MatFactorInfo*,Mat*);
40b0592601SKris Buschelman 
41b0592601SKris Buschelman EXTERN_C_BEGIN
42b0592601SKris Buschelman EXTERN int MatConvert_SuperLU_SeqAIJ(Mat,MatType,Mat*);
43b0592601SKris Buschelman EXTERN int MatConvert_SeqAIJ_SuperLU(Mat,MatType,Mat*);
44b0592601SKris Buschelman EXTERN_C_END
4514b396bbSKris Buschelman 
4614b396bbSKris Buschelman #undef __FUNCT__
4714b396bbSKris Buschelman #define __FUNCT__ "MatDestroy_SeqAIJ_SuperLU"
4814b396bbSKris Buschelman int MatDestroy_SeqAIJ_SuperLU(Mat A)
4914b396bbSKris Buschelman {
50*460c4903SKris Buschelman   int                ierr;
5114b396bbSKris Buschelman   Mat_SeqAIJ_SuperLU *lu = (Mat_SeqAIJ_SuperLU*)A->spptr;
52*460c4903SKris Buschelman   int                (*destroy)(Mat)=lu->MatDestroy;
5314b396bbSKris Buschelman 
5414b396bbSKris Buschelman   PetscFunctionBegin;
55b0592601SKris Buschelman   ierr = MatConvert_SuperLU_SeqAIJ(A,MATSEQAIJ,&A);CHKERRQ(ierr);
5614b396bbSKris Buschelman   ierr = (*destroy)(A);CHKERRQ(ierr);
5714b396bbSKris Buschelman   PetscFunctionReturn(0);
5814b396bbSKris Buschelman }
5914b396bbSKris Buschelman 
6014b396bbSKris Buschelman #undef __FUNCT__
6114b396bbSKris Buschelman #define __FUNCT__ "MatView_SeqAIJ_Spooles"
6214b396bbSKris Buschelman int MatView_SeqAIJ_SuperLU(Mat A,PetscViewer viewer)
6314b396bbSKris Buschelman {
6414b396bbSKris Buschelman   int                   ierr;
6514b396bbSKris Buschelman   PetscTruth            isascii;
6614b396bbSKris Buschelman   PetscViewerFormat     format;
6714b396bbSKris Buschelman   Mat_SeqAIJ_SuperLU   *lu=(Mat_SeqAIJ_SuperLU*)(A->spptr);
6814b396bbSKris Buschelman 
6914b396bbSKris Buschelman   PetscFunctionBegin;
7014b396bbSKris Buschelman   ierr = (*lu->MatView)(A,viewer);CHKERRQ(ierr);
7114b396bbSKris Buschelman 
7214b396bbSKris Buschelman   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr);
7314b396bbSKris Buschelman   if (isascii) {
7414b396bbSKris Buschelman     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
7514b396bbSKris Buschelman     if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) {
7614b396bbSKris Buschelman       ierr = MatSeqAIJFactorInfo_SuperLU(A,viewer);CHKERRQ(ierr);
7714b396bbSKris Buschelman     }
7814b396bbSKris Buschelman   }
7914b396bbSKris Buschelman   PetscFunctionReturn(0);
8014b396bbSKris Buschelman }
8114b396bbSKris Buschelman 
8214b396bbSKris Buschelman #undef __FUNCT__
8314b396bbSKris Buschelman #define __FUNCT__ "MatAssemblyEnd_SeqAIJ_SuperLU"
8414b396bbSKris Buschelman int MatAssemblyEnd_SeqAIJ_SuperLU(Mat A,MatAssemblyType mode) {
8514b396bbSKris Buschelman   int                ierr;
8614b396bbSKris Buschelman   Mat_SeqAIJ_SuperLU *lu=(Mat_SeqAIJ_SuperLU*)(A->spptr);
8714b396bbSKris Buschelman 
8814b396bbSKris Buschelman   PetscFunctionBegin;
8914b396bbSKris Buschelman   ierr = (*lu->MatAssemblyEnd)(A,mode);CHKERRQ(ierr);
90b0592601SKris Buschelman 
91b0592601SKris Buschelman   lu->MatLUFactorSymbolic  = A->ops->lufactorsymbolic;
92b0592601SKris Buschelman   A->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqAIJ_SuperLU;
9314b396bbSKris Buschelman   PetscFunctionReturn(0);
9414b396bbSKris Buschelman }
9514b396bbSKris Buschelman 
9614b396bbSKris Buschelman #include "src/mat/impls/dense/seq/dense.h"
9714b396bbSKris Buschelman #undef __FUNCT__
9814b396bbSKris Buschelman #define __FUNCT__ "MatCreateNull_SeqAIJ_SuperLU"
9914b396bbSKris Buschelman int MatCreateNull_SeqAIJ_SuperLU(Mat A,Mat *nullMat)
10014b396bbSKris Buschelman {
10114b396bbSKris Buschelman   Mat_SeqAIJ_SuperLU  *lu = (Mat_SeqAIJ_SuperLU*)A->spptr;
10214b396bbSKris Buschelman   int                 numRows = A->m,numCols = A->n;
10314b396bbSKris Buschelman   SCformat            *Lstore;
10414b396bbSKris Buschelman   int                 numNullCols,size;
10514b396bbSKris Buschelman #if defined(PETSC_USE_COMPLEX)
10614b396bbSKris Buschelman   doublecomplex       *nullVals,*workVals;
10714b396bbSKris Buschelman #else
10814b396bbSKris Buschelman   PetscScalar         *nullVals,*workVals;
10914b396bbSKris Buschelman #endif
11014b396bbSKris Buschelman   int                 row,newRow,col,newCol,block,b,ierr;
11114b396bbSKris Buschelman 
11214b396bbSKris Buschelman   PetscFunctionBegin;
11314b396bbSKris Buschelman   PetscValidHeaderSpecific(A,MAT_COOKIE);
11414b396bbSKris Buschelman   PetscValidPointer(nullMat);
11514b396bbSKris Buschelman   if (!A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Unfactored matrix");
11614b396bbSKris Buschelman   numNullCols = numCols - numRows;
11714b396bbSKris Buschelman   if (numNullCols < 0) SETERRQ(PETSC_ERR_ARG_WRONG,"Function only applies to underdetermined problems");
11814b396bbSKris Buschelman   /* Create the null matrix */
11914b396bbSKris Buschelman   ierr = MatCreateSeqDense(A->comm,numRows,numNullCols,PETSC_NULL,nullMat);CHKERRQ(ierr);
12014b396bbSKris Buschelman   if (numNullCols == 0) {
12114b396bbSKris Buschelman     ierr = MatAssemblyBegin(*nullMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
12214b396bbSKris Buschelman     ierr = MatAssemblyEnd(*nullMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
12314b396bbSKris Buschelman     PetscFunctionReturn(0);
12414b396bbSKris Buschelman   }
12514b396bbSKris Buschelman #if defined(PETSC_USE_COMPLEX)
12614b396bbSKris Buschelman   nullVals = (doublecomplex*)((Mat_SeqDense*)(*nullMat)->data)->v;
12714b396bbSKris Buschelman #else
12814b396bbSKris Buschelman   nullVals = ((Mat_SeqDense*)(*nullMat)->data)->v;
12914b396bbSKris Buschelman #endif
13014b396bbSKris Buschelman   /* Copy in the columns */
13114b396bbSKris Buschelman   Lstore = (SCformat*)lu->L.Store;
13214b396bbSKris Buschelman   for(block = 0; block <= Lstore->nsuper; block++) {
13314b396bbSKris Buschelman     newRow = Lstore->sup_to_col[block];
13414b396bbSKris Buschelman     size   = Lstore->sup_to_col[block+1] - Lstore->sup_to_col[block];
13514b396bbSKris Buschelman     for(col = Lstore->rowind_colptr[newRow]; col < Lstore->rowind_colptr[newRow+1]; col++) {
13614b396bbSKris Buschelman       newCol = Lstore->rowind[col];
13714b396bbSKris Buschelman       if (newCol >= numRows) {
13814b396bbSKris Buschelman         for(b = 0; b < size; b++)
13914b396bbSKris Buschelman #if defined(PETSC_USE_COMPLEX)
14014b396bbSKris Buschelman           nullVals[(newCol-numRows)*numRows+newRow+b] = ((doublecomplex*)Lstore->nzval)[Lstore->nzval_colptr[newRow+b]+col];
14114b396bbSKris Buschelman #else
14214b396bbSKris Buschelman           nullVals[(newCol-numRows)*numRows+newRow+b] = ((double*)Lstore->nzval)[Lstore->nzval_colptr[newRow+b]+col];
14314b396bbSKris Buschelman #endif
14414b396bbSKris Buschelman       }
14514b396bbSKris Buschelman     }
14614b396bbSKris Buschelman   }
14714b396bbSKris Buschelman   /* Permute rhs to form P^T_c B */
14814b396bbSKris Buschelman   ierr = PetscMalloc(numRows*sizeof(double),&workVals);CHKERRQ(ierr);
14914b396bbSKris Buschelman   for(b = 0; b < numNullCols; b++) {
15014b396bbSKris Buschelman     for(row = 0; row < numRows; row++) workVals[lu->perm_c[row]] = nullVals[b*numRows+row];
15114b396bbSKris Buschelman     for(row = 0; row < numRows; row++) nullVals[b*numRows+row]   = workVals[row];
15214b396bbSKris Buschelman   }
15314b396bbSKris Buschelman   /* Backward solve the upper triangle A x = b */
15414b396bbSKris Buschelman   for(b = 0; b < numNullCols; b++) {
15514b396bbSKris Buschelman #if defined(PETSC_USE_COMPLEX)
15614b396bbSKris Buschelman     sp_ztrsv("L","T","U",&lu->L,&lu->U,&nullVals[b*numRows],&ierr);
15714b396bbSKris Buschelman #else
15814b396bbSKris Buschelman     sp_dtrsv("L","T","U",&lu->L,&lu->U,&nullVals[b*numRows],&ierr);
15914b396bbSKris Buschelman #endif
16014b396bbSKris Buschelman     if (ierr < 0)
16114b396bbSKris Buschelman       SETERRQ1(PETSC_ERR_ARG_WRONG,"The argument %d was invalid",-ierr);
16214b396bbSKris Buschelman   }
16314b396bbSKris Buschelman   ierr = PetscFree(workVals);CHKERRQ(ierr);
16414b396bbSKris Buschelman 
16514b396bbSKris Buschelman   ierr = MatAssemblyBegin(*nullMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
16614b396bbSKris Buschelman   ierr = MatAssemblyEnd(*nullMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
16714b396bbSKris Buschelman   PetscFunctionReturn(0);
16814b396bbSKris Buschelman }
16914b396bbSKris Buschelman 
17014b396bbSKris Buschelman #undef __FUNCT__
17114b396bbSKris Buschelman #define __FUNCT__ "MatSolve_SeqAIJ_SuperLU"
17214b396bbSKris Buschelman int MatSolve_SeqAIJ_SuperLU(Mat A,Vec b,Vec x)
17314b396bbSKris Buschelman {
17414b396bbSKris Buschelman   Mat_SeqAIJ_SuperLU *lu = (Mat_SeqAIJ_SuperLU*)A->spptr;
17514b396bbSKris Buschelman   PetscScalar        *array;
17614b396bbSKris Buschelman   int                m,ierr;
17714b396bbSKris Buschelman 
17814b396bbSKris Buschelman   PetscFunctionBegin;
17914b396bbSKris Buschelman   ierr = VecGetLocalSize(b,&m);CHKERRQ(ierr);
18014b396bbSKris Buschelman   ierr = VecCopy(b,x);CHKERRQ(ierr);
18114b396bbSKris Buschelman   ierr = VecGetArray(x,&array);CHKERRQ(ierr);
18214b396bbSKris Buschelman   /* Create the Rhs */
18314b396bbSKris Buschelman   lu->B.Stype        = SLU_DN;
18414b396bbSKris Buschelman   lu->B.Mtype        = SLU_GE;
18514b396bbSKris Buschelman   lu->B.nrow         = m;
18614b396bbSKris Buschelman   lu->B.ncol         = 1;
18714b396bbSKris Buschelman   ((DNformat*)lu->B.Store)->lda   = m;
18814b396bbSKris Buschelman   ((DNformat*)lu->B.Store)->nzval = array;
18914b396bbSKris Buschelman #if defined(PETSC_USE_COMPLEX)
19014b396bbSKris Buschelman   lu->B.Dtype        = SLU_Z;
19114b396bbSKris Buschelman   zgstrs("T",&lu->L,&lu->U,lu->perm_r,lu->perm_c,&lu->B,&ierr);
19214b396bbSKris Buschelman #else
19314b396bbSKris Buschelman   lu->B.Dtype        = SLU_D;
19414b396bbSKris Buschelman   dgstrs("T",&lu->L,&lu->U,lu->perm_r,lu->perm_c,&lu->B,&ierr);
19514b396bbSKris Buschelman #endif
19614b396bbSKris Buschelman   if (ierr < 0) SETERRQ1(PETSC_ERR_ARG_WRONG,"The diagonal element of row %d was invalid",-ierr);
19714b396bbSKris Buschelman   ierr = VecRestoreArray(x,&array);CHKERRQ(ierr);
19814b396bbSKris Buschelman   PetscFunctionReturn(0);
19914b396bbSKris Buschelman }
20014b396bbSKris Buschelman 
20114b396bbSKris Buschelman static int StatInitCalled = 0;
20214b396bbSKris Buschelman 
20314b396bbSKris Buschelman #undef __FUNCT__
20414b396bbSKris Buschelman #define __FUNCT__ "MatLUFactorNumeric_SeqAIJ_SuperLU"
20514b396bbSKris Buschelman int MatLUFactorNumeric_SeqAIJ_SuperLU(Mat A,Mat *F)
20614b396bbSKris Buschelman {
20714b396bbSKris Buschelman   Mat_SeqAIJ         *aa = (Mat_SeqAIJ*)(A)->data;
20814b396bbSKris Buschelman   Mat_SeqAIJ_SuperLU *lu = (Mat_SeqAIJ_SuperLU*)(*F)->spptr;
209d54de34fSKris Buschelman   int                *etree,ierr;
21014b396bbSKris Buschelman   PetscTruth         flag;
21114b396bbSKris Buschelman 
21214b396bbSKris Buschelman   PetscFunctionBegin;
21314b396bbSKris Buschelman   /* Create the SuperMatrix for A^T:
21414b396bbSKris Buschelman        Since SuperLU only likes column-oriented matrices,we pass it the transpose,
21514b396bbSKris Buschelman        and then solve A^T X = B in MatSolve().
21614b396bbSKris Buschelman   */
21714b396bbSKris Buschelman 
21814b396bbSKris Buschelman   if ( lu->flg == DIFFERENT_NONZERO_PATTERN){ /* first numerical factorization */
21914b396bbSKris Buschelman     lu->A.Stype   = SLU_NC;
22014b396bbSKris Buschelman #if defined(PETSC_USE_COMPLEX)
22114b396bbSKris Buschelman     lu->A.Dtype   = SLU_Z;
22214b396bbSKris Buschelman #else
22314b396bbSKris Buschelman     lu->A.Dtype   = SLU_D;
22414b396bbSKris Buschelman #endif
22514b396bbSKris Buschelman     lu->A.Mtype   = SLU_GE;
22614b396bbSKris Buschelman     lu->A.nrow    = A->n;
22714b396bbSKris Buschelman     lu->A.ncol    = A->m;
22814b396bbSKris Buschelman 
22914b396bbSKris Buschelman     ierr = PetscMalloc(sizeof(NCformat),&lu->store);CHKERRQ(ierr);
23014b396bbSKris Buschelman     ierr = PetscMalloc(sizeof(DNformat),&lu->B.Store);CHKERRQ(ierr);
23114b396bbSKris Buschelman   }
23214b396bbSKris Buschelman   lu->store->nnz    = aa->nz;
23314b396bbSKris Buschelman   lu->store->colptr = aa->i;
23414b396bbSKris Buschelman   lu->store->rowind = aa->j;
23514b396bbSKris Buschelman   lu->store->nzval  = aa->a;
23614b396bbSKris Buschelman   lu->A.Store       = lu->store;
23714b396bbSKris Buschelman 
23814b396bbSKris Buschelman   /* Set SuperLU options */
23914b396bbSKris Buschelman   lu->relax      = sp_ienv(2);
24014b396bbSKris Buschelman   lu->panel_size = sp_ienv(1);
24114b396bbSKris Buschelman   /* We have to initialize global data or SuperLU crashes (sucky design) */
24214b396bbSKris Buschelman   if (!StatInitCalled) {
24314b396bbSKris Buschelman     StatInit(lu->panel_size,lu->relax);
24414b396bbSKris Buschelman   }
24514b396bbSKris Buschelman   StatInitCalled++;
24614b396bbSKris Buschelman 
24714b396bbSKris Buschelman   ierr = PetscOptionsBegin(A->comm,A->prefix,"SuperLU Options","Mat");CHKERRQ(ierr);
24814b396bbSKris Buschelman   /* use SuperLU mat ordeing */
24914b396bbSKris Buschelman   ierr = PetscOptionsInt("-mat_superlu_ordering","SuperLU ordering type (one of 0, 1, 2, 3)\n   0: natural ordering;\n   1: MMD applied to A'*A;\n   2: MMD applied to A'+A;\n   3: COLAMD, approximate minimum degree column ordering","None",lu->ispec,&lu->ispec,&flag);CHKERRQ(ierr);
25014b396bbSKris Buschelman   if (flag) {
25114b396bbSKris Buschelman     get_perm_c(lu->ispec, &lu->A, lu->perm_c);
25214b396bbSKris Buschelman     lu->SuperluMatOdering = PETSC_TRUE;
25314b396bbSKris Buschelman   }
25414b396bbSKris Buschelman   PetscOptionsEnd();
25514b396bbSKris Buschelman 
25614b396bbSKris Buschelman   /* Create the elimination tree */
25714b396bbSKris Buschelman   ierr = PetscMalloc(A->n*sizeof(int),&etree);CHKERRQ(ierr);
25814b396bbSKris Buschelman   sp_preorder("N",&lu->A,lu->perm_c,etree,&lu->AC);
25914b396bbSKris Buschelman   /* Factor the matrix */
26014b396bbSKris Buschelman #if defined(PETSC_USE_COMPLEX)
26114b396bbSKris Buschelman   zgstrf("N",&lu->AC,lu->pivot_threshold,0.0,lu->relax,lu->panel_size,etree,PETSC_NULL,0,lu->perm_r,lu->perm_c,&lu->L,&lu->U,&ierr);
26214b396bbSKris Buschelman #else
26314b396bbSKris Buschelman   dgstrf("N",&lu->AC,lu->pivot_threshold,0.0,lu->relax,lu->panel_size,etree,PETSC_NULL,0,lu->perm_r,lu->perm_c,&lu->L,&lu->U,&ierr);
26414b396bbSKris Buschelman #endif
26514b396bbSKris Buschelman   if (ierr < 0) {
26614b396bbSKris Buschelman     SETERRQ1(PETSC_ERR_ARG_WRONG,"The diagonal element of row %d was invalid",-ierr);
26714b396bbSKris Buschelman   } else if (ierr > 0) {
26814b396bbSKris Buschelman     if (ierr <= A->m) {
26914b396bbSKris Buschelman       SETERRQ1(PETSC_ERR_ARG_WRONG,"The diagonal element %d of U is exactly zero",ierr);
27014b396bbSKris Buschelman     } else {
27114b396bbSKris Buschelman       SETERRQ1(PETSC_ERR_ARG_WRONG,"Memory allocation failure after %d bytes were allocated",ierr-A->m);
27214b396bbSKris Buschelman     }
27314b396bbSKris Buschelman   }
27414b396bbSKris Buschelman 
27514b396bbSKris Buschelman   /* Cleanup */
27614b396bbSKris Buschelman   ierr = PetscFree(etree);CHKERRQ(ierr);
27714b396bbSKris Buschelman 
27814b396bbSKris Buschelman   lu->flg = SAME_NONZERO_PATTERN;
27914b396bbSKris Buschelman   PetscFunctionReturn(0);
28014b396bbSKris Buschelman }
28114b396bbSKris Buschelman 
28214b396bbSKris Buschelman /*
28314b396bbSKris Buschelman    Note the r permutation is ignored
28414b396bbSKris Buschelman */
28514b396bbSKris Buschelman #undef __FUNCT__
28614b396bbSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_SeqAIJ_SuperLU"
28714b396bbSKris Buschelman int MatLUFactorSymbolic_SeqAIJ_SuperLU(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F)
28814b396bbSKris Buschelman {
28914b396bbSKris Buschelman   Mat                 B;
29014b396bbSKris Buschelman   Mat_SeqAIJ_SuperLU  *lu;
29114b396bbSKris Buschelman   int                 ierr,*ca;
29214b396bbSKris Buschelman 
29314b396bbSKris Buschelman   PetscFunctionBegin;
29414b396bbSKris Buschelman 
295899d7b4fSKris Buschelman   ierr = MatCreate(A->comm,A->m,A->n,PETSC_DETERMINE,PETSC_DETERMINE,&B);CHKERRQ(ierr);
296899d7b4fSKris Buschelman   ierr = MatSetType(B,MATSUPERLU);CHKERRQ(ierr);
297899d7b4fSKris Buschelman   ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
29814b396bbSKris Buschelman   B->ops->lufactornumeric = MatLUFactorNumeric_SeqAIJ_SuperLU;
29914b396bbSKris Buschelman   B->ops->solve           = MatSolve_SeqAIJ_SuperLU;
30014b396bbSKris Buschelman   B->factor               = FACTOR_LU;
301899d7b4fSKris Buschelman   B->assembled            = PETSC_TRUE;  /* required by -sles_view */
30214b396bbSKris Buschelman 
303899d7b4fSKris Buschelman   lu = (Mat_SeqAIJ_SuperLU*)(B->spptr);
30414b396bbSKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)B,"MatCreateNull","MatCreateNull_SeqAIJ_SuperLU",
30514b396bbSKris Buschelman                                     (void(*)(void))MatCreateNull_SeqAIJ_SuperLU);CHKERRQ(ierr);
30614b396bbSKris Buschelman 
30714b396bbSKris Buschelman   /* Allocate the work arrays required by SuperLU (notice sizes are for the transpose) */
30814b396bbSKris Buschelman   ierr = PetscMalloc(A->n*sizeof(int),&lu->perm_r);CHKERRQ(ierr);
30914b396bbSKris Buschelman   ierr = PetscMalloc(A->m*sizeof(int),&lu->perm_c);CHKERRQ(ierr);
31014b396bbSKris Buschelman 
31114b396bbSKris Buschelman   /* use PETSc mat ordering */
31214b396bbSKris Buschelman   ierr = ISGetIndices(c,&ca);CHKERRQ(ierr);
31314b396bbSKris Buschelman   ierr = PetscMemcpy(lu->perm_c,ca,A->m*sizeof(int));CHKERRQ(ierr);
31414b396bbSKris Buschelman   ierr = ISRestoreIndices(c,&ca);CHKERRQ(ierr);
31514b396bbSKris Buschelman   lu->SuperluMatOdering = PETSC_FALSE;
31614b396bbSKris Buschelman 
31714b396bbSKris Buschelman   lu->pivot_threshold = info->dtcol;
31814b396bbSKris Buschelman   PetscLogObjectMemory(B,(A->m+A->n)*sizeof(int)+sizeof(Mat_SeqAIJ_SuperLU));
31914b396bbSKris Buschelman 
32014b396bbSKris Buschelman   lu->flg            = DIFFERENT_NONZERO_PATTERN;
321e740cb95SKris Buschelman   lu->CleanUpSuperLU = PETSC_TRUE;
322899d7b4fSKris Buschelman 
323899d7b4fSKris Buschelman   *F = B;
32414b396bbSKris Buschelman   PetscFunctionReturn(0);
32514b396bbSKris Buschelman }
32614b396bbSKris Buschelman 
32714b396bbSKris Buschelman /* used by -sles_view */
32814b396bbSKris Buschelman #undef __FUNCT__
32914b396bbSKris Buschelman #define __FUNCT__ "MatSeqAIJFactorInfo_SuperLU"
33014b396bbSKris Buschelman int MatSeqAIJFactorInfo_SuperLU(Mat A,PetscViewer viewer)
33114b396bbSKris Buschelman {
33214b396bbSKris Buschelman   Mat_SeqAIJ_SuperLU      *lu= (Mat_SeqAIJ_SuperLU*)A->spptr;
33314b396bbSKris Buschelman   int                     ierr;
33414b396bbSKris Buschelman   PetscFunctionBegin;
33514b396bbSKris Buschelman 
33614b396bbSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"SuperLU run parameters:\n");CHKERRQ(ierr);
33714b396bbSKris Buschelman   if(lu->SuperluMatOdering) ierr = PetscViewerASCIIPrintf(viewer,"  SuperLU mat ordering: %d\n",lu->ispec);CHKERRQ(ierr);
33814b396bbSKris Buschelman 
33914b396bbSKris Buschelman   PetscFunctionReturn(0);
34014b396bbSKris Buschelman }
34114b396bbSKris Buschelman 
34214b396bbSKris Buschelman EXTERN_C_BEGIN
34314b396bbSKris Buschelman #undef __FUNCT__
344b0592601SKris Buschelman #define __FUNCT__ "MatConvert_SuperLU_SeqAIJ"
345b0592601SKris Buschelman int MatConvert_SuperLU_SeqAIJ(Mat A,MatType type,Mat *newmat) {
34614b396bbSKris Buschelman   int                  ierr;
347b0592601SKris Buschelman   Mat                  B=*newmat;
348b0592601SKris Buschelman   Mat_SeqAIJ_SuperLU   *lu=(Mat_SeqAIJ_SuperLU *)A->spptr;
349b0592601SKris Buschelman 
350b0592601SKris Buschelman   PetscFunctionBegin;
351b0592601SKris Buschelman   if (B != A) {
352b0592601SKris Buschelman     /* This routine was inherited from SeqAIJ. */
353b0592601SKris Buschelman     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
354b0592601SKris Buschelman   } else {
355b0592601SKris Buschelman     if (lu->CleanUpSuperLU) {
356b0592601SKris Buschelman       /* We have to free the global data or SuperLU crashes (sucky design)*/
357b0592601SKris Buschelman       /* Since we don't know if more solves on other matrices may be done
358b0592601SKris Buschelman          we cannot free the yucky SuperLU global data
359b0592601SKris Buschelman          StatFree();
360b0592601SKris Buschelman       */
361b0592601SKris Buschelman 
362b0592601SKris Buschelman       /* Free the SuperLU datastructures */
363b0592601SKris Buschelman       Destroy_CompCol_Permuted(&lu->AC);
364b0592601SKris Buschelman       Destroy_SuperNode_Matrix(&lu->L);
365b0592601SKris Buschelman       Destroy_CompCol_Matrix(&lu->U);
366b0592601SKris Buschelman       ierr = PetscFree(lu->B.Store);CHKERRQ(ierr);
367b0592601SKris Buschelman       ierr = PetscFree(lu->A.Store);CHKERRQ(ierr);
368b0592601SKris Buschelman       ierr = PetscFree(lu->perm_r);CHKERRQ(ierr);
369b0592601SKris Buschelman       ierr = PetscFree(lu->perm_c);CHKERRQ(ierr);
370b0592601SKris Buschelman     }
371b0592601SKris Buschelman     /* Reset the original function pointers */
372b0592601SKris Buschelman     B->ops->view             = lu->MatView;
373b0592601SKris Buschelman     B->ops->assemblyend      = lu->MatAssemblyEnd;
374b0592601SKris Buschelman     B->ops->lufactorsymbolic = lu->MatLUFactorSymbolic;
375b0592601SKris Buschelman     B->ops->destroy          = lu->MatDestroy;
376b0592601SKris Buschelman     /* lu is only a function pointer stash unless we've factored the matrix, which we haven't! */
377b0592601SKris Buschelman     ierr = PetscFree(lu);CHKERRQ(ierr);
378b0592601SKris Buschelman     ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJ);CHKERRQ(ierr);
379b0592601SKris Buschelman   }
380b0592601SKris Buschelman   *newmat = B;
381b0592601SKris Buschelman   PetscFunctionReturn(0);
382b0592601SKris Buschelman }
383b0592601SKris Buschelman EXTERN_C_END
384b0592601SKris Buschelman 
385b0592601SKris Buschelman EXTERN_C_BEGIN
386b0592601SKris Buschelman #undef __FUNCT__
387b0592601SKris Buschelman #define __FUNCT__ "MatConvert_SeqAIJ_SuperLU"
388b0592601SKris Buschelman int MatConvert_SeqAIJ_SuperLU(Mat A,MatType type,Mat *newmat) {
389b0592601SKris Buschelman   int                ierr;
390b0592601SKris Buschelman   Mat                B=*newmat;
39114b396bbSKris Buschelman   Mat_SeqAIJ_SuperLU *lu;
39214b396bbSKris Buschelman 
39314b396bbSKris Buschelman   PetscFunctionBegin;
394b0592601SKris Buschelman   if (B != A) {
395b0592601SKris Buschelman     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
396b0592601SKris Buschelman   }
39714b396bbSKris Buschelman 
39814b396bbSKris Buschelman   ierr = PetscNew(Mat_SeqAIJ_SuperLU,&lu);CHKERRQ(ierr);
39914b396bbSKris Buschelman   lu->MatView              = A->ops->view;
40014b396bbSKris Buschelman   lu->MatAssemblyEnd       = A->ops->assemblyend;
401b0592601SKris Buschelman   lu->MatLUFactorSymbolic  = A->ops->lufactorsymbolic;
40214b396bbSKris Buschelman   lu->MatDestroy           = A->ops->destroy;
40314b396bbSKris Buschelman   lu->CleanUpSuperLU       = PETSC_FALSE;
404b0592601SKris Buschelman 
405b0592601SKris Buschelman   B->spptr                 = (void*)lu;
406b0592601SKris Buschelman   B->ops->view             = MatView_SeqAIJ_SuperLU;
407b0592601SKris Buschelman   B->ops->assemblyend      = MatAssemblyEnd_SeqAIJ_SuperLU;
408b0592601SKris Buschelman   B->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqAIJ_SuperLU;
409b0592601SKris Buschelman   B->ops->destroy          = MatDestroy_SeqAIJ_SuperLU;
410b0592601SKris Buschelman 
411b0592601SKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_superlu_C",
412b0592601SKris Buschelman                                            "MatConvert_SeqAIJ_SuperLU",MatConvert_SeqAIJ_SuperLU);CHKERRQ(ierr);
413b0592601SKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_superlu_seqaij_C",
414b0592601SKris Buschelman                                            "MatConvert_SuperLU_SeqAIJ",MatConvert_SuperLU_SeqAIJ);CHKERRQ(ierr);
415b0592601SKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)B,type);CHKERRQ(ierr);
416b0592601SKris Buschelman   *newmat = B;
417b0592601SKris Buschelman   PetscFunctionReturn(0);
418b0592601SKris Buschelman }
419b0592601SKris Buschelman EXTERN_C_END
420b0592601SKris Buschelman 
421b0592601SKris Buschelman EXTERN_C_BEGIN
422b0592601SKris Buschelman #undef __FUNCT__
423b0592601SKris Buschelman #define __FUNCT__ "MatCreate_SeqAIJ_SuperLU"
424b0592601SKris Buschelman int MatCreate_SeqAIJ_SuperLU(Mat A) {
425b0592601SKris Buschelman   int ierr;
426b0592601SKris Buschelman 
427b0592601SKris Buschelman   PetscFunctionBegin;
428b0592601SKris Buschelman   ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr);
429b0592601SKris Buschelman   ierr = MatConvert_SeqAIJ_SuperLU(A,MATSUPERLU,&A);CHKERRQ(ierr);
43014b396bbSKris Buschelman   PetscFunctionReturn(0);
43114b396bbSKris Buschelman }
43214b396bbSKris Buschelman EXTERN_C_END
433