xref: /petsc/src/mat/impls/aij/seq/superlu/superlu.c (revision 5a46d3fad41c2283670f018990932d4d60b33ea0)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
214b396bbSKris Buschelman 
3*5a46d3faSBarry Smith /*  --------------------------------------------------------------------
4*5a46d3faSBarry Smith 
5*5a46d3faSBarry Smith      This file implements a subclass of the SeqAIJ matrix class that uses
6*5a46d3faSBarry Smith      the SuperLU 3.0 sparse solver. You can use this as a starting point for
7*5a46d3faSBarry Smith      implementing your own subclass of a PETSc matrix class.
8*5a46d3faSBarry Smith 
9*5a46d3faSBarry Smith      This demonstrates a way to make an implementation inheritence of a PETSc
10*5a46d3faSBarry Smith      matrix type. This means constructing a new matrix type (SuperLU) by changing some
11*5a46d3faSBarry Smith      of the methods of the previous type (SeqAIJ), adding additional data, and possibly
12*5a46d3faSBarry Smith      additional method. (See any book on object oriented programming).
1314b396bbSKris Buschelman */
1414b396bbSKris Buschelman 
15*5a46d3faSBarry Smith /*
16*5a46d3faSBarry Smith      Defines the data structure for the base matrix type (SeqAIJ)
17*5a46d3faSBarry Smith */
1814b396bbSKris Buschelman #include "src/mat/impls/aij/seq/aij.h"
1914b396bbSKris Buschelman 
20*5a46d3faSBarry Smith /*
21*5a46d3faSBarry Smith      SuperLU include files
22*5a46d3faSBarry Smith */
2314b396bbSKris Buschelman EXTERN_C_BEGIN
2414b396bbSKris Buschelman #if defined(PETSC_USE_COMPLEX)
2566f57492SHong Zhang #include "slu_zdefs.h"
2614b396bbSKris Buschelman #else
2766f57492SHong Zhang #include "slu_ddefs.h"
2814b396bbSKris Buschelman #endif
2966f57492SHong Zhang #include "slu_util.h"
3014b396bbSKris Buschelman EXTERN_C_END
3114b396bbSKris Buschelman 
32*5a46d3faSBarry Smith /*
33*5a46d3faSBarry Smith      This is the data we are "ADDING" to the SeqAIJ matrix type to get the SuperLU matrix type
34*5a46d3faSBarry Smith */
3514b396bbSKris Buschelman typedef struct {
3675af56d4SHong Zhang   SuperMatrix       A,L,U,B,X;
3775af56d4SHong Zhang   superlu_options_t options;
38da7a1d00SHong Zhang   PetscInt          *perm_c; /* column permutation vector */
39da7a1d00SHong Zhang   PetscInt          *perm_r; /* row permutations from partial pivoting */
40da7a1d00SHong Zhang   PetscInt          *etree;
41da7a1d00SHong Zhang   PetscReal         *R, *C;
4275af56d4SHong Zhang   char              equed[1];
43da7a1d00SHong Zhang   PetscInt          lwork;
4475af56d4SHong Zhang   void              *work;
45da7a1d00SHong Zhang   PetscReal         rpg, rcond;
4675af56d4SHong Zhang   mem_usage_t       mem_usage;
4775af56d4SHong Zhang   MatStructure      flg;
4814b396bbSKris Buschelman 
49*5a46d3faSBarry Smith   /*
50*5a46d3faSBarry Smith      This is where the methods for the superclass (SeqAIJ) are kept while we
51*5a46d3faSBarry Smith      reset the pointers in the function table to the new (SuperLU) versions
52*5a46d3faSBarry Smith   */
536849ba73SBarry Smith   PetscErrorCode (*MatDuplicate)(Mat,MatDuplicateOption,Mat*);
546849ba73SBarry Smith   PetscErrorCode (*MatView)(Mat,PetscViewer);
556849ba73SBarry Smith   PetscErrorCode (*MatAssemblyEnd)(Mat,MatAssemblyType);
566849ba73SBarry Smith   PetscErrorCode (*MatLUFactorSymbolic)(Mat,IS,IS,MatFactorInfo*,Mat*);
576849ba73SBarry Smith   PetscErrorCode (*MatDestroy)(Mat);
5814b396bbSKris Buschelman 
5914b396bbSKris Buschelman   /* Flag to clean up (non-global) SuperLU objects during Destroy */
6014b396bbSKris Buschelman   PetscTruth CleanUpSuperLU;
61f0c56d0fSKris Buschelman } Mat_SuperLU;
6214b396bbSKris Buschelman 
63*5a46d3faSBarry Smith /*
64*5a46d3faSBarry Smith     Takes a SuperLU matrix (that is a SeqAIJ matrix with the additional SuperLU data-structures
65*5a46d3faSBarry Smith    and methods) and converts it back to a regular SeqAIJ matrix.
66*5a46d3faSBarry Smith */
67*5a46d3faSBarry Smith EXTERN_C_BEGIN
68*5a46d3faSBarry Smith #undef __FUNCT__
69*5a46d3faSBarry Smith #define __FUNCT__ "MatConvert_SuperLU_SeqAIJ"
70*5a46d3faSBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SuperLU_SeqAIJ(Mat A,MatType type,MatReuse reuse,Mat *newmat)
71*5a46d3faSBarry Smith {
72*5a46d3faSBarry Smith   PetscErrorCode ierr;
73*5a46d3faSBarry Smith   Mat            B=*newmat;
74*5a46d3faSBarry Smith   Mat_SuperLU    *lu=(Mat_SuperLU *)A->spptr;
7514b396bbSKris Buschelman 
76*5a46d3faSBarry Smith   PetscFunctionBegin;
77*5a46d3faSBarry Smith   if (reuse == MAT_INITIAL_MATRIX) {
78*5a46d3faSBarry Smith     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
79*5a46d3faSBarry Smith   }
80*5a46d3faSBarry Smith   /* Reset the original SeqAIJ function pointers */
81*5a46d3faSBarry Smith   B->ops->duplicate        = lu->MatDuplicate;
82*5a46d3faSBarry Smith   B->ops->view             = lu->MatView;
83*5a46d3faSBarry Smith   B->ops->assemblyend      = lu->MatAssemblyEnd;
84*5a46d3faSBarry Smith   B->ops->lufactorsymbolic = lu->MatLUFactorSymbolic;
85*5a46d3faSBarry Smith   B->ops->destroy          = lu->MatDestroy;
86*5a46d3faSBarry Smith 
87*5a46d3faSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaij_superlu_C","",PETSC_NULL);CHKERRQ(ierr);
88*5a46d3faSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_superlu_seqaij_C","",PETSC_NULL);CHKERRQ(ierr);
89*5a46d3faSBarry Smith 
90*5a46d3faSBarry Smith   /* change the type name back to its original value */
91*5a46d3faSBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJ);CHKERRQ(ierr);
92*5a46d3faSBarry Smith   *newmat = B;
93*5a46d3faSBarry Smith   PetscFunctionReturn(0);
94*5a46d3faSBarry Smith }
95*5a46d3faSBarry Smith EXTERN_C_END
96b0592601SKris Buschelman 
97b0592601SKris Buschelman EXTERN_C_BEGIN
98*5a46d3faSBarry Smith #undef __FUNCT__
99*5a46d3faSBarry Smith #define __FUNCT__ "MatConvert_SeqAIJ_SuperLU"
100*5a46d3faSBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqAIJ_SuperLU(Mat A,MatType type,MatReuse reuse,Mat *newmat)
101*5a46d3faSBarry Smith {
102*5a46d3faSBarry Smith   PetscErrorCode ierr;
103*5a46d3faSBarry Smith   Mat            B=*newmat;
104*5a46d3faSBarry Smith   Mat_SuperLU    *lu;
105*5a46d3faSBarry Smith 
106*5a46d3faSBarry Smith   PetscFunctionBegin;
107*5a46d3faSBarry Smith   if (reuse == MAT_INITIAL_MATRIX) {
108*5a46d3faSBarry Smith     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
109*5a46d3faSBarry Smith   }
110*5a46d3faSBarry Smith 
111*5a46d3faSBarry Smith   ierr = PetscNew(Mat_SuperLU,&lu);CHKERRQ(ierr);
112*5a46d3faSBarry Smith   /* save the original SeqAIJ methods that we are changing */
113*5a46d3faSBarry Smith   lu->MatDuplicate         = A->ops->duplicate;
114*5a46d3faSBarry Smith   lu->MatView              = A->ops->view;
115*5a46d3faSBarry Smith   lu->MatAssemblyEnd       = A->ops->assemblyend;
116*5a46d3faSBarry Smith   lu->MatLUFactorSymbolic  = A->ops->lufactorsymbolic;
117*5a46d3faSBarry Smith   lu->MatDestroy           = A->ops->destroy;
118*5a46d3faSBarry Smith   lu->CleanUpSuperLU       = PETSC_FALSE;
119*5a46d3faSBarry Smith 
120*5a46d3faSBarry Smith   /* add to the matrix the location for all the SuperLU data is to be stored */
121*5a46d3faSBarry Smith   B->spptr                 = (void*)lu;
122*5a46d3faSBarry Smith 
123*5a46d3faSBarry Smith   /* set the methods in the function table to the SuperLU versions */
124*5a46d3faSBarry Smith   B->ops->duplicate        = MatDuplicate_SuperLU;
125*5a46d3faSBarry Smith   B->ops->view             = MatView_SuperLU;
126*5a46d3faSBarry Smith   B->ops->assemblyend      = MatAssemblyEnd_SuperLU;
127*5a46d3faSBarry Smith   B->ops->lufactorsymbolic = MatLUFactorSymbolic_SuperLU;
128*5a46d3faSBarry Smith   B->ops->choleskyfactorsymbolic = 0;
129*5a46d3faSBarry Smith   B->ops->destroy          = MatDestroy_SuperLU;
130*5a46d3faSBarry Smith 
131*5a46d3faSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_superlu_C",
132*5a46d3faSBarry Smith                                            "MatConvert_SeqAIJ_SuperLU",MatConvert_SeqAIJ_SuperLU);CHKERRQ(ierr);
133*5a46d3faSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_superlu_seqaij_C",
134*5a46d3faSBarry Smith                                            "MatConvert_SuperLU_SeqAIJ",MatConvert_SuperLU_SeqAIJ);CHKERRQ(ierr);
135*5a46d3faSBarry Smith   ierr = PetscInfo(0,"Using SuperLU for SeqAIJ LU factorization and solves.\n");CHKERRQ(ierr);
136*5a46d3faSBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)B,MATSUPERLU);CHKERRQ(ierr);
137*5a46d3faSBarry Smith   *newmat = B;
138*5a46d3faSBarry Smith   PetscFunctionReturn(0);
139*5a46d3faSBarry Smith }
140b0592601SKris Buschelman EXTERN_C_END
14114b396bbSKris Buschelman 
142*5a46d3faSBarry Smith /*
143*5a46d3faSBarry Smith     Utility function
144*5a46d3faSBarry Smith */
145*5a46d3faSBarry Smith #undef __FUNCT__
146*5a46d3faSBarry Smith #define __FUNCT__ "MatFactorInfo_SuperLU"
147*5a46d3faSBarry Smith PetscErrorCode MatFactorInfo_SuperLU(Mat A,PetscViewer viewer)
148*5a46d3faSBarry Smith {
149*5a46d3faSBarry Smith   Mat_SuperLU       *lu= (Mat_SuperLU*)A->spptr;
150*5a46d3faSBarry Smith   PetscErrorCode    ierr;
151*5a46d3faSBarry Smith   superlu_options_t options;
152*5a46d3faSBarry Smith 
153*5a46d3faSBarry Smith   PetscFunctionBegin;
154*5a46d3faSBarry Smith   /* check if matrix is superlu_dist type */
155*5a46d3faSBarry Smith   if (A->ops->solve != MatSolve_SuperLU) PetscFunctionReturn(0);
156*5a46d3faSBarry Smith 
157*5a46d3faSBarry Smith   options = lu->options;
158*5a46d3faSBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"SuperLU run parameters:\n");CHKERRQ(ierr);
159*5a46d3faSBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"  Equil: %s\n",(options.Equil != NO) ? "YES": "NO");CHKERRQ(ierr);
160*5a46d3faSBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"  ColPerm: %D\n",options.ColPerm);CHKERRQ(ierr);
161*5a46d3faSBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"  IterRefine: %D\n",options.IterRefine);CHKERRQ(ierr);
162*5a46d3faSBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"  SymmetricMode: %s\n",(options.SymmetricMode != NO) ? "YES": "NO");CHKERRQ(ierr);
163*5a46d3faSBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"  DiagPivotThresh: %g\n",options.DiagPivotThresh);CHKERRQ(ierr);
164*5a46d3faSBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"  PivotGrowth: %s\n",(options.PivotGrowth != NO) ? "YES": "NO");CHKERRQ(ierr);
165*5a46d3faSBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"  ConditionNumber: %s\n",(options.ConditionNumber != NO) ? "YES": "NO");CHKERRQ(ierr);
166*5a46d3faSBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"  RowPerm: %D\n",options.RowPerm);CHKERRQ(ierr);
167*5a46d3faSBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"  ReplaceTinyPivot: %s\n",(options.ReplaceTinyPivot != NO) ? "YES": "NO");CHKERRQ(ierr);
168*5a46d3faSBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"  PrintStat: %s\n",(options.PrintStat != NO) ? "YES": "NO");CHKERRQ(ierr);
169*5a46d3faSBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"  lwork: %D\n",lu->lwork);CHKERRQ(ierr);
170*5a46d3faSBarry Smith 
171*5a46d3faSBarry Smith   PetscFunctionReturn(0);
172*5a46d3faSBarry Smith }
173*5a46d3faSBarry Smith 
174*5a46d3faSBarry Smith /*
175*5a46d3faSBarry Smith     These are the methods provided to REPLACE the corresponding methods of the
176*5a46d3faSBarry Smith    SeqAIJ matrix class. Other methods of SeqAIJ are not replaced
177*5a46d3faSBarry Smith */
178*5a46d3faSBarry Smith #undef __FUNCT__
179*5a46d3faSBarry Smith #define __FUNCT__ "MatLUFactorNumeric_SuperLU"
180*5a46d3faSBarry Smith PetscErrorCode MatLUFactorNumeric_SuperLU(Mat A,MatFactorInfo *info,Mat *F)
181*5a46d3faSBarry Smith {
182*5a46d3faSBarry Smith   Mat_SeqAIJ     *aa = (Mat_SeqAIJ*)(A)->data;
183*5a46d3faSBarry Smith   Mat_SuperLU    *lu = (Mat_SuperLU*)(*F)->spptr;
184*5a46d3faSBarry Smith   PetscErrorCode ierr;
185*5a46d3faSBarry Smith   PetscInt       sinfo;
186*5a46d3faSBarry Smith   SuperLUStat_t  stat;
187*5a46d3faSBarry Smith   PetscReal      ferr, berr;
188*5a46d3faSBarry Smith   NCformat       *Ustore;
189*5a46d3faSBarry Smith   SCformat       *Lstore;
190*5a46d3faSBarry Smith 
191*5a46d3faSBarry Smith   PetscFunctionBegin;
192*5a46d3faSBarry Smith   if (lu->flg == SAME_NONZERO_PATTERN){ /* successing numerical factorization */
193*5a46d3faSBarry Smith     lu->options.Fact = SamePattern;
194*5a46d3faSBarry Smith     /* Ref: ~SuperLU_3.0/EXAMPLE/dlinsolx2.c */
195*5a46d3faSBarry Smith     Destroy_SuperMatrix_Store(&lu->A);
196*5a46d3faSBarry Smith     if ( lu->lwork >= 0 ) {
197*5a46d3faSBarry Smith       Destroy_SuperNode_Matrix(&lu->L);
198*5a46d3faSBarry Smith       Destroy_CompCol_Matrix(&lu->U);
199*5a46d3faSBarry Smith       lu->options.Fact = SamePattern;
200*5a46d3faSBarry Smith     }
201*5a46d3faSBarry Smith   }
202*5a46d3faSBarry Smith 
203*5a46d3faSBarry Smith   /* Create the SuperMatrix for lu->A=A^T:
204*5a46d3faSBarry Smith        Since SuperLU likes column-oriented matrices,we pass it the transpose,
205*5a46d3faSBarry Smith        and then solve A^T X = B in MatSolve(). */
206*5a46d3faSBarry Smith #if defined(PETSC_USE_COMPLEX)
207*5a46d3faSBarry Smith   zCreate_CompCol_Matrix(&lu->A,A->cmap.n,A->rmap.n,aa->nz,(doublecomplex*)aa->a,aa->j,aa->i,
208*5a46d3faSBarry Smith                            SLU_NC,SLU_Z,SLU_GE);
209*5a46d3faSBarry Smith #else
210*5a46d3faSBarry Smith   dCreate_CompCol_Matrix(&lu->A,A->cmap.n,A->rmap.n,aa->nz,aa->a,aa->j,aa->i,
211*5a46d3faSBarry Smith                            SLU_NC,SLU_D,SLU_GE);
212*5a46d3faSBarry Smith #endif
213*5a46d3faSBarry Smith 
214*5a46d3faSBarry Smith   /* Initialize the statistics variables. */
215*5a46d3faSBarry Smith   StatInit(&stat);
216*5a46d3faSBarry Smith 
217*5a46d3faSBarry Smith   /* Numerical factorization */
218*5a46d3faSBarry Smith   lu->B.ncol = 0;  /* Indicate not to solve the system */
219*5a46d3faSBarry Smith #if defined(PETSC_USE_COMPLEX)
220*5a46d3faSBarry Smith    zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
221*5a46d3faSBarry Smith            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
222*5a46d3faSBarry Smith            &lu->mem_usage, &stat, &sinfo);
223*5a46d3faSBarry Smith #else
224*5a46d3faSBarry Smith   dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
225*5a46d3faSBarry Smith            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
226*5a46d3faSBarry Smith            &lu->mem_usage, &stat, &sinfo);
227*5a46d3faSBarry Smith #endif
228*5a46d3faSBarry Smith   if ( !sinfo || sinfo == lu->A.ncol+1 ) {
229*5a46d3faSBarry Smith     if ( lu->options.PivotGrowth )
230*5a46d3faSBarry Smith       ierr = PetscPrintf(PETSC_COMM_SELF,"  Recip. pivot growth = %e\n", lu->rpg);
231*5a46d3faSBarry Smith     if ( lu->options.ConditionNumber )
232*5a46d3faSBarry Smith       ierr = PetscPrintf(PETSC_COMM_SELF,"  Recip. condition number = %e\n", lu->rcond);
233*5a46d3faSBarry Smith   } else if ( sinfo > 0 ){
234*5a46d3faSBarry Smith     if ( lu->lwork == -1 ) {
235*5a46d3faSBarry Smith       ierr = PetscPrintf(PETSC_COMM_SELF,"  ** Estimated memory: %D bytes\n", sinfo - lu->A.ncol);
236*5a46d3faSBarry Smith     } else {
237*5a46d3faSBarry Smith       ierr = PetscPrintf(PETSC_COMM_SELF,"  Warning: gssvx() returns info %D\n",sinfo);
238*5a46d3faSBarry Smith     }
239*5a46d3faSBarry Smith   } else { /* sinfo < 0 */
240*5a46d3faSBarry Smith     SETERRQ2(PETSC_ERR_LIB, "info = %D, the %D-th argument in gssvx() had an illegal value", sinfo,-sinfo);
241*5a46d3faSBarry Smith   }
242*5a46d3faSBarry Smith 
243*5a46d3faSBarry Smith   if ( lu->options.PrintStat ) {
244*5a46d3faSBarry Smith     ierr = PetscPrintf(PETSC_COMM_SELF,"MatLUFactorNumeric_SuperLU():\n");
245*5a46d3faSBarry Smith     StatPrint(&stat);
246*5a46d3faSBarry Smith     Lstore = (SCformat *) lu->L.Store;
247*5a46d3faSBarry Smith     Ustore = (NCformat *) lu->U.Store;
248*5a46d3faSBarry Smith     ierr = PetscPrintf(PETSC_COMM_SELF,"  No of nonzeros in factor L = %D\n", Lstore->nnz);
249*5a46d3faSBarry Smith     ierr = PetscPrintf(PETSC_COMM_SELF,"  No of nonzeros in factor U = %D\n", Ustore->nnz);
250*5a46d3faSBarry Smith     ierr = PetscPrintf(PETSC_COMM_SELF,"  No of nonzeros in L+U = %D\n", Lstore->nnz + Ustore->nnz - lu->A.ncol);
251*5a46d3faSBarry Smith     ierr = PetscPrintf(PETSC_COMM_SELF,"  L\\U MB %.3f\ttotal MB needed %.3f\texpansions %D\n",
252*5a46d3faSBarry Smith 	       lu->mem_usage.for_lu/1e6, lu->mem_usage.total_needed/1e6,
253*5a46d3faSBarry Smith 	       lu->mem_usage.expansions);
254*5a46d3faSBarry Smith   }
255*5a46d3faSBarry Smith   StatFree(&stat);
256*5a46d3faSBarry Smith 
257*5a46d3faSBarry Smith   lu->flg = SAME_NONZERO_PATTERN;
258*5a46d3faSBarry Smith   PetscFunctionReturn(0);
259*5a46d3faSBarry Smith }
260*5a46d3faSBarry Smith 
26114b396bbSKris Buschelman #undef __FUNCT__
262f0c56d0fSKris Buschelman #define __FUNCT__ "MatDestroy_SuperLU"
263dfbe8321SBarry Smith PetscErrorCode MatDestroy_SuperLU(Mat A)
26414b396bbSKris Buschelman {
265dfbe8321SBarry Smith   PetscErrorCode ierr;
266f0c56d0fSKris Buschelman   Mat_SuperLU    *lu=(Mat_SuperLU*)A->spptr;
26714b396bbSKris Buschelman 
26814b396bbSKris Buschelman   PetscFunctionBegin;
26975af56d4SHong Zhang   if (lu->CleanUpSuperLU) { /* Free the SuperLU datastructures */
27075af56d4SHong Zhang     Destroy_SuperMatrix_Store(&lu->A);
27175af56d4SHong Zhang     Destroy_SuperMatrix_Store(&lu->B);
27275af56d4SHong Zhang     Destroy_SuperMatrix_Store(&lu->X);
2739ce50919SHong Zhang 
2749ce50919SHong Zhang     ierr = PetscFree(lu->etree);CHKERRQ(ierr);
2759ce50919SHong Zhang     ierr = PetscFree(lu->perm_r);CHKERRQ(ierr);
2769ce50919SHong Zhang     ierr = PetscFree(lu->perm_c);CHKERRQ(ierr);
2779ce50919SHong Zhang     ierr = PetscFree(lu->R);CHKERRQ(ierr);
2789ce50919SHong Zhang     ierr = PetscFree(lu->C);CHKERRQ(ierr);
27975af56d4SHong Zhang     if ( lu->lwork >= 0 ) {
280fb3e25aaSKris Buschelman       Destroy_SuperNode_Matrix(&lu->L);
281fb3e25aaSKris Buschelman       Destroy_CompCol_Matrix(&lu->U);
28275af56d4SHong Zhang     }
283fb3e25aaSKris Buschelman   }
284ceb03754SKris Buschelman   ierr = MatConvert_SuperLU_SeqAIJ(A,MATSEQAIJ,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr);
285fb3e25aaSKris Buschelman   ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
28614b396bbSKris Buschelman   PetscFunctionReturn(0);
28714b396bbSKris Buschelman }
28814b396bbSKris Buschelman 
28914b396bbSKris Buschelman #undef __FUNCT__
290f0c56d0fSKris Buschelman #define __FUNCT__ "MatView_SuperLU"
291dfbe8321SBarry Smith PetscErrorCode MatView_SuperLU(Mat A,PetscViewer viewer)
29214b396bbSKris Buschelman {
293dfbe8321SBarry Smith   PetscErrorCode    ierr;
29432077d6dSBarry Smith   PetscTruth        iascii;
29514b396bbSKris Buschelman   PetscViewerFormat format;
296f0c56d0fSKris Buschelman   Mat_SuperLU       *lu=(Mat_SuperLU*)(A->spptr);
29714b396bbSKris Buschelman 
29814b396bbSKris Buschelman   PetscFunctionBegin;
29914b396bbSKris Buschelman   ierr = (*lu->MatView)(A,viewer);CHKERRQ(ierr);
30014b396bbSKris Buschelman 
30132077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
30232077d6dSBarry Smith   if (iascii) {
30314b396bbSKris Buschelman     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
3042f59403fSHong Zhang     if (format == PETSC_VIEWER_ASCII_INFO) {
305f0c56d0fSKris Buschelman       ierr = MatFactorInfo_SuperLU(A,viewer);CHKERRQ(ierr);
30614b396bbSKris Buschelman     }
30714b396bbSKris Buschelman   }
30814b396bbSKris Buschelman   PetscFunctionReturn(0);
30914b396bbSKris Buschelman }
31014b396bbSKris Buschelman 
31114b396bbSKris Buschelman #undef __FUNCT__
312f0c56d0fSKris Buschelman #define __FUNCT__ "MatAssemblyEnd_SuperLU"
313dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_SuperLU(Mat A,MatAssemblyType mode) {
314dfbe8321SBarry Smith   PetscErrorCode ierr;
315f0c56d0fSKris Buschelman   Mat_SuperLU    *lu=(Mat_SuperLU*)(A->spptr);
31614b396bbSKris Buschelman 
31714b396bbSKris Buschelman   PetscFunctionBegin;
31814b396bbSKris Buschelman   ierr = (*lu->MatAssemblyEnd)(A,mode);CHKERRQ(ierr);
319b0592601SKris Buschelman   lu->MatLUFactorSymbolic  = A->ops->lufactorsymbolic;
320f0c56d0fSKris Buschelman   A->ops->lufactorsymbolic = MatLUFactorSymbolic_SuperLU;
32114b396bbSKris Buschelman   PetscFunctionReturn(0);
32214b396bbSKris Buschelman }
32314b396bbSKris Buschelman 
32414b396bbSKris Buschelman 
32514b396bbSKris Buschelman #undef __FUNCT__
326c29e884eSHong Zhang #define __FUNCT__ "MatSolve_SuperLU_Private"
327c29e884eSHong Zhang PetscErrorCode MatSolve_SuperLU_Private(Mat A,Vec b,Vec x)
32814b396bbSKris Buschelman {
329f0c56d0fSKris Buschelman   Mat_SuperLU    *lu = (Mat_SuperLU*)A->spptr;
33075af56d4SHong Zhang   PetscScalar    *barray,*xarray;
331dfbe8321SBarry Smith   PetscErrorCode ierr;
332da7a1d00SHong Zhang   PetscInt       info,i;
33375af56d4SHong Zhang   SuperLUStat_t  stat;
334da7a1d00SHong Zhang   PetscReal      ferr,berr;
33514b396bbSKris Buschelman 
33614b396bbSKris Buschelman   PetscFunctionBegin;
337937a6b0eSHong Zhang   if ( lu->lwork == -1 ) {
338937a6b0eSHong Zhang     PetscFunctionReturn(0);
339937a6b0eSHong Zhang   }
34075af56d4SHong Zhang   lu->B.ncol = 1;   /* Set the number of right-hand side */
34175af56d4SHong Zhang   ierr = VecGetArray(b,&barray);CHKERRQ(ierr);
34275af56d4SHong Zhang   ierr = VecGetArray(x,&xarray);CHKERRQ(ierr);
3435fe6150dSHong Zhang 
3445fe6150dSHong Zhang #if defined(PETSC_USE_COMPLEX)
3455fe6150dSHong Zhang   ((DNformat*)lu->B.Store)->nzval = (doublecomplex*)barray;
3465fe6150dSHong Zhang   ((DNformat*)lu->X.Store)->nzval = (doublecomplex*)xarray;
3475fe6150dSHong Zhang #else
3485fe6150dSHong Zhang   ((DNformat*)lu->B.Store)->nzval = barray;
34975af56d4SHong Zhang   ((DNformat*)lu->X.Store)->nzval = xarray;
3505fe6150dSHong Zhang #endif
35175af56d4SHong Zhang 
35275af56d4SHong Zhang   /* Initialize the statistics variables. */
35375af56d4SHong Zhang   StatInit(&stat);
35475af56d4SHong Zhang 
35575af56d4SHong Zhang   lu->options.Fact  = FACTORED; /* Indicate the factored form of A is supplied. */
3568914a3f7SHong Zhang #if defined(PETSC_USE_COMPLEX)
3578914a3f7SHong Zhang   zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
3588914a3f7SHong Zhang            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
3598914a3f7SHong Zhang            &lu->mem_usage, &stat, &info);
3608914a3f7SHong Zhang #else
36175af56d4SHong Zhang   dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
36275af56d4SHong Zhang            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
36375af56d4SHong Zhang            &lu->mem_usage, &stat, &info);
3648914a3f7SHong Zhang #endif
36575af56d4SHong Zhang   ierr = VecRestoreArray(b,&barray);CHKERRQ(ierr);
36675af56d4SHong Zhang   ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
36775af56d4SHong Zhang 
368958c9bccSBarry Smith   if ( !info || info == lu->A.ncol+1 ) {
36975af56d4SHong Zhang     if ( lu->options.IterRefine ) {
3708914a3f7SHong Zhang       ierr = PetscPrintf(PETSC_COMM_SELF,"Iterative Refinement:\n");
3718914a3f7SHong Zhang       ierr = PetscPrintf(PETSC_COMM_SELF,"  %8s%8s%16s%16s\n", "rhs", "Steps", "FERR", "BERR");
37275af56d4SHong Zhang       for (i = 0; i < 1; ++i)
3738914a3f7SHong Zhang         ierr = PetscPrintf(PETSC_COMM_SELF,"  %8d%8d%16e%16e\n", i+1, stat.RefineSteps, ferr, berr);
37475af56d4SHong Zhang     }
3758914a3f7SHong Zhang   } else if ( info > 0 ){
3768914a3f7SHong Zhang     if ( lu->lwork == -1 ) {
37777431f27SBarry Smith       ierr = PetscPrintf(PETSC_COMM_SELF,"  ** Estimated memory: %D bytes\n", info - lu->A.ncol);
3788914a3f7SHong Zhang     } else {
37977431f27SBarry Smith       ierr = PetscPrintf(PETSC_COMM_SELF,"  Warning: gssvx() returns info %D\n",info);
3808914a3f7SHong Zhang     }
3818914a3f7SHong Zhang   } else if (info < 0){
38277431f27SBarry Smith     SETERRQ2(PETSC_ERR_LIB, "info = %D, the %D-th argument in gssvx() had an illegal value", info,-info);
38314b396bbSKris Buschelman   }
38414b396bbSKris Buschelman 
3858914a3f7SHong Zhang   if ( lu->options.PrintStat ) {
3868914a3f7SHong Zhang     ierr = PetscPrintf(PETSC_COMM_SELF,"MatSolve__SuperLU():\n");
3878914a3f7SHong Zhang     StatPrint(&stat);
3888914a3f7SHong Zhang   }
38975af56d4SHong Zhang   StatFree(&stat);
39075af56d4SHong Zhang   PetscFunctionReturn(0);
39175af56d4SHong Zhang }
39214b396bbSKris Buschelman 
39314b396bbSKris Buschelman #undef __FUNCT__
394c29e884eSHong Zhang #define __FUNCT__ "MatSolve_SuperLU"
395c29e884eSHong Zhang PetscErrorCode MatSolve_SuperLU(Mat A,Vec b,Vec x)
396c29e884eSHong Zhang {
397c29e884eSHong Zhang   Mat_SuperLU    *lu = (Mat_SuperLU*)A->spptr;
398c29e884eSHong Zhang   PetscErrorCode ierr;
399c29e884eSHong Zhang 
400c29e884eSHong Zhang   PetscFunctionBegin;
401c29e884eSHong Zhang   lu->options.Trans = TRANS;
402c29e884eSHong Zhang   ierr = MatSolve_SuperLU_Private(A,b,x);CHKERRQ(ierr);
403c29e884eSHong Zhang   PetscFunctionReturn(0);
404c29e884eSHong Zhang }
405c29e884eSHong Zhang 
406c29e884eSHong Zhang #undef __FUNCT__
407c7c1fe80SHong Zhang #define __FUNCT__ "MatSolveTranspose_SuperLU"
408c7c1fe80SHong Zhang PetscErrorCode MatSolveTranspose_SuperLU(Mat A,Vec b,Vec x)
409c7c1fe80SHong Zhang {
410c7c1fe80SHong Zhang   Mat_SuperLU    *lu = (Mat_SuperLU*)A->spptr;
411c7c1fe80SHong Zhang   PetscErrorCode ierr;
412c7c1fe80SHong Zhang 
413c7c1fe80SHong Zhang   PetscFunctionBegin;
414c7c1fe80SHong Zhang   lu->options.Trans = NOTRANS;
415c29e884eSHong Zhang   ierr = MatSolve_SuperLU_Private(A,b,x);CHKERRQ(ierr);
416c7c1fe80SHong Zhang   PetscFunctionReturn(0);
417c7c1fe80SHong Zhang }
418c7c1fe80SHong Zhang 
41914b396bbSKris Buschelman 
42014b396bbSKris Buschelman /*
42114b396bbSKris Buschelman    Note the r permutation is ignored
42214b396bbSKris Buschelman */
42314b396bbSKris Buschelman #undef __FUNCT__
424f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_SuperLU"
425dfbe8321SBarry Smith PetscErrorCode MatLUFactorSymbolic_SuperLU(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F)
42614b396bbSKris Buschelman {
42714b396bbSKris Buschelman   Mat            B;
428f0c56d0fSKris Buschelman   Mat_SuperLU    *lu;
4296849ba73SBarry Smith   PetscErrorCode ierr;
4302a4c71feSBarry Smith   PetscInt       m=A->rmap.n,n=A->cmap.n,indx;
4319ce50919SHong Zhang   PetscTruth     flg;
4325ca28756SHong Zhang   const char   *colperm[]={"NATURAL","MMD_ATA","MMD_AT_PLUS_A","COLAMD"}; /* MY_PERMC - not supported by the petsc interface yet */
4335ca28756SHong Zhang   const char   *iterrefine[]={"NOREFINE", "SINGLE", "DOUBLE", "EXTRA"};
4345ca28756SHong Zhang   const char   *rowperm[]={"NOROWPERM", "LargeDiag"}; /* MY_PERMC - not supported by the petsc interface yet */
43514b396bbSKris Buschelman 
43614b396bbSKris Buschelman   PetscFunctionBegin;
437f69a0ea3SMatthew Knepley   ierr = MatCreate(A->comm,&B);CHKERRQ(ierr);
4382a4c71feSBarry Smith   ierr = MatSetSizes(B,A->rmap.n,A->cmap.n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
439be5d1d56SKris Buschelman   ierr = MatSetType(B,A->type_name);CHKERRQ(ierr);
440899d7b4fSKris Buschelman   ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
441f0c56d0fSKris Buschelman 
442f0c56d0fSKris Buschelman   B->ops->lufactornumeric = MatLUFactorNumeric_SuperLU;
443f0c56d0fSKris Buschelman   B->ops->solve           = MatSolve_SuperLU;
444c7c1fe80SHong Zhang   B->ops->solvetranspose  = MatSolveTranspose_SuperLU;
44514b396bbSKris Buschelman   B->factor               = FACTOR_LU;
44694b7f48cSBarry Smith   B->assembled            = PETSC_TRUE;  /* required by -ksp_view */
44714b396bbSKris Buschelman 
448f0c56d0fSKris Buschelman   lu = (Mat_SuperLU*)(B->spptr);
4499ce50919SHong Zhang 
4509ce50919SHong Zhang   /* Set SuperLU options */
4519ce50919SHong Zhang     /* the default values for options argument:
4529ce50919SHong Zhang 	options.Fact = DOFACT;
4539ce50919SHong Zhang         options.Equil = YES;
4549ce50919SHong Zhang     	options.ColPerm = COLAMD;
4559ce50919SHong Zhang 	options.DiagPivotThresh = 1.0;
4569ce50919SHong Zhang     	options.Trans = NOTRANS;
4579ce50919SHong Zhang     	options.IterRefine = NOREFINE;
4589ce50919SHong Zhang     	options.SymmetricMode = NO;
4599ce50919SHong Zhang     	options.PivotGrowth = NO;
4609ce50919SHong Zhang     	options.ConditionNumber = NO;
4619ce50919SHong Zhang     	options.PrintStat = YES;
4629ce50919SHong Zhang     */
4639ce50919SHong Zhang   set_default_options(&lu->options);
4648914a3f7SHong Zhang   /* equilibration causes error in solve(), thus not supported here. See dgssvx.c for possible reason. */
4658914a3f7SHong Zhang   lu->options.Equil = NO;
4669ce50919SHong Zhang   lu->options.PrintStat = NO;
4678914a3f7SHong Zhang   lu->lwork = 0;   /* allocate space internally by system malloc */
4689ce50919SHong Zhang 
4699ce50919SHong Zhang   ierr = PetscOptionsBegin(A->comm,A->prefix,"SuperLU Options","Mat");CHKERRQ(ierr);
4709ce50919SHong Zhang   /*
4714dc4c822SBarry Smith   ierr = PetscOptionsTruth("-mat_superlu_equil","Equil","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr);
4728914a3f7SHong Zhang   if (flg) lu->options.Equil = YES; -- not supported by the interface !!!
4739ce50919SHong Zhang   */
4748914a3f7SHong Zhang   ierr = PetscOptionsEList("-mat_superlu_colperm","ColPerm","None",colperm,4,colperm[3],&indx,&flg);CHKERRQ(ierr);
4759ce50919SHong Zhang   if (flg) {lu->options.ColPerm = (colperm_t)indx;}
4768914a3f7SHong Zhang   ierr = PetscOptionsEList("-mat_superlu_iterrefine","IterRefine","None",iterrefine,4,iterrefine[0],&indx,&flg);CHKERRQ(ierr);
4779ce50919SHong Zhang   if (flg) { lu->options.IterRefine = (IterRefine_t)indx;}
4784dc4c822SBarry Smith   ierr = PetscOptionsTruth("-mat_superlu_symmetricmode","SymmetricMode","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr);
4799ce50919SHong Zhang   if (flg) lu->options.SymmetricMode = YES;
4808914a3f7SHong Zhang   ierr = PetscOptionsReal("-mat_superlu_diagpivotthresh","DiagPivotThresh","None",lu->options.DiagPivotThresh,&lu->options.DiagPivotThresh,PETSC_NULL);CHKERRQ(ierr);
4814dc4c822SBarry Smith   ierr = PetscOptionsTruth("-mat_superlu_pivotgrowth","PivotGrowth","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr);
4829ce50919SHong Zhang   if (flg) lu->options.PivotGrowth = YES;
4834dc4c822SBarry Smith   ierr = PetscOptionsTruth("-mat_superlu_conditionnumber","ConditionNumber","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr);
4849ce50919SHong Zhang   if (flg) lu->options.ConditionNumber = YES;
4858914a3f7SHong Zhang   ierr = PetscOptionsEList("-mat_superlu_rowperm","rowperm","None",rowperm,2,rowperm[0],&indx,&flg);CHKERRQ(ierr);
4869ce50919SHong Zhang   if (flg) {lu->options.RowPerm = (rowperm_t)indx;}
4874dc4c822SBarry Smith   ierr = PetscOptionsTruth("-mat_superlu_replacetinypivot","ReplaceTinyPivot","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr);
4889ce50919SHong Zhang   if (flg) lu->options.ReplaceTinyPivot = YES;
4894dc4c822SBarry Smith   ierr = PetscOptionsTruth("-mat_superlu_printstat","PrintStat","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr);
4909ce50919SHong Zhang   if (flg) lu->options.PrintStat = YES;
4918914a3f7SHong Zhang   ierr = PetscOptionsInt("-mat_superlu_lwork","size of work array in bytes used by factorization","None",lu->lwork,&lu->lwork,PETSC_NULL);CHKERRQ(ierr);
4925fe6150dSHong Zhang   if (lu->lwork > 0 ){
4935fe6150dSHong Zhang     ierr = PetscMalloc(lu->lwork,&lu->work);CHKERRQ(ierr);
4945fe6150dSHong Zhang   } else if (lu->lwork != 0 && lu->lwork != -1){
49577431f27SBarry Smith     ierr = PetscPrintf(PETSC_COMM_SELF,"   Warning: lwork %D is not supported by SUPERLU. The default lwork=0 is used.\n",lu->lwork);
4968914a3f7SHong Zhang     lu->lwork = 0;
4978914a3f7SHong Zhang   }
4989ce50919SHong Zhang   PetscOptionsEnd();
4999ce50919SHong Zhang 
50075af56d4SHong Zhang #ifdef SUPERLU2
5012ecb5974SKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatCreateNull","MatCreateNull_SuperLU",
502f0c56d0fSKris Buschelman                                     (void(*)(void))MatCreateNull_SuperLU);CHKERRQ(ierr);
50375af56d4SHong Zhang #endif
50414b396bbSKris Buschelman 
50575af56d4SHong Zhang   /* Allocate spaces (notice sizes are for the transpose) */
506da7a1d00SHong Zhang   ierr = PetscMalloc(m*sizeof(PetscInt),&lu->etree);CHKERRQ(ierr);
507da7a1d00SHong Zhang   ierr = PetscMalloc(n*sizeof(PetscInt),&lu->perm_r);CHKERRQ(ierr);
508da7a1d00SHong Zhang   ierr = PetscMalloc(m*sizeof(PetscInt),&lu->perm_c);CHKERRQ(ierr);
509da7a1d00SHong Zhang   ierr = PetscMalloc(n*sizeof(PetscInt),&lu->R);CHKERRQ(ierr);
510da7a1d00SHong Zhang   ierr = PetscMalloc(m*sizeof(PetscInt),&lu->C);CHKERRQ(ierr);
51175af56d4SHong Zhang 
5129ce50919SHong Zhang   /* create rhs and solution x without allocate space for .Store */
5135fe6150dSHong Zhang #if defined(PETSC_USE_COMPLEX)
514937a6b0eSHong Zhang   zCreate_Dense_Matrix(&lu->B, m, 1, PETSC_NULL, m, SLU_DN, SLU_Z, SLU_GE);
515937a6b0eSHong Zhang   zCreate_Dense_Matrix(&lu->X, m, 1, PETSC_NULL, m, SLU_DN, SLU_Z, SLU_GE);
5165fe6150dSHong Zhang #else
51775af56d4SHong Zhang   dCreate_Dense_Matrix(&lu->B, m, 1, PETSC_NULL, m, SLU_DN, SLU_D, SLU_GE);
51875af56d4SHong Zhang   dCreate_Dense_Matrix(&lu->X, m, 1, PETSC_NULL, m, SLU_DN, SLU_D, SLU_GE);
5195fe6150dSHong Zhang #endif
52014b396bbSKris Buschelman 
52114b396bbSKris Buschelman   lu->flg            = DIFFERENT_NONZERO_PATTERN;
522e740cb95SKris Buschelman   lu->CleanUpSuperLU = PETSC_TRUE;
523899d7b4fSKris Buschelman 
524899d7b4fSKris Buschelman   *F = B;
5252a4c71feSBarry Smith   ierr = PetscLogObjectMemory(B,(A->rmap.n+A->cmap.n)*sizeof(PetscInt)+sizeof(Mat_SuperLU));CHKERRQ(ierr);
52614b396bbSKris Buschelman   PetscFunctionReturn(0);
52714b396bbSKris Buschelman }
52814b396bbSKris Buschelman 
529f0c56d0fSKris Buschelman 
530f0c56d0fSKris Buschelman #undef __FUNCT__
531f0c56d0fSKris Buschelman #define __FUNCT__ "MatDuplicate_SuperLU"
532dfbe8321SBarry Smith PetscErrorCode MatDuplicate_SuperLU(Mat A, MatDuplicateOption op, Mat *M) {
533da7a1d00SHong Zhang   PetscErrorCode ierr;
5348f340917SKris Buschelman   Mat_SuperLU    *lu=(Mat_SuperLU *)A->spptr;
5358f340917SKris Buschelman 
53614b396bbSKris Buschelman   PetscFunctionBegin;
5378f340917SKris Buschelman   ierr = (*lu->MatDuplicate)(A,op,M);CHKERRQ(ierr);
5383f953163SKris Buschelman   ierr = PetscMemcpy((*M)->spptr,lu,sizeof(Mat_SuperLU));CHKERRQ(ierr);
53914b396bbSKris Buschelman   PetscFunctionReturn(0);
54014b396bbSKris Buschelman }
54114b396bbSKris Buschelman 
542b0592601SKris Buschelman 
54324b6179bSKris Buschelman /*MC
544fafad747SKris Buschelman   MATSUPERLU - MATSUPERLU = "superlu" - A matrix type providing direct solvers (LU) for sequential matrices
54524b6179bSKris Buschelman   via the external package SuperLU.
54624b6179bSKris Buschelman 
54724b6179bSKris Buschelman   If SuperLU is installed (see the manual for
54824b6179bSKris Buschelman   instructions on how to declare the existence of external packages),
54924b6179bSKris Buschelman   a matrix type can be constructed which invokes SuperLU solvers.
55024b6179bSKris Buschelman   After calling MatCreate(...,A), simply call MatSetType(A,MATSUPERLU).
55124b6179bSKris Buschelman 
55224b6179bSKris Buschelman   This matrix inherits from MATSEQAIJ.  As a result, MatSeqAIJSetPreallocation is
55328b08bd3SKris Buschelman   supported for this matrix type.  One can also call MatConvert for an inplace conversion to or from
55428b08bd3SKris Buschelman   the MATSEQAIJ type without data copy.
55524b6179bSKris Buschelman 
55624b6179bSKris Buschelman   Options Database Keys:
5570bad9183SKris Buschelman + -mat_type superlu - sets the matrix type to "superlu" during a call to MatSetFromOptions()
55851678082SBarry Smith . -mat_superlu_ordering <0,1,2,3> - 0: natural ordering,
55924b6179bSKris Buschelman                                     1: MMD applied to A'*A,
56024b6179bSKris Buschelman                                     2: MMD applied to A'+A,
56124b6179bSKris Buschelman                                     3: COLAMD, approximate minimum degree column ordering
56251678082SBarry Smith . -mat_superlu_iterrefine - have SuperLU do iterative refinement after the triangular solve
56351678082SBarry Smith                           choices: NOREFINE, SINGLE, DOUBLE, EXTRA; default is NOREFINE
56451678082SBarry Smith - -mat_superlu_printstat - print SuperLU statistics about the factorization
56524b6179bSKris Buschelman 
56624b6179bSKris Buschelman    Level: beginner
56724b6179bSKris Buschelman 
56824b6179bSKris Buschelman .seealso: PCLU
56924b6179bSKris Buschelman M*/
57024b6179bSKris Buschelman 
571*5a46d3faSBarry Smith /*
572*5a46d3faSBarry Smith     Constructor for the new derived matrix class. It simply creates the base
573*5a46d3faSBarry Smith    matrix class and then adds the additional information/methods needed by SuperLU.
574*5a46d3faSBarry Smith */
575b0592601SKris Buschelman EXTERN_C_BEGIN
576b0592601SKris Buschelman #undef __FUNCT__
577f0c56d0fSKris Buschelman #define __FUNCT__ "MatCreate_SuperLU"
578be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_SuperLU(Mat A)
579dfbe8321SBarry Smith {
580dfbe8321SBarry Smith   PetscErrorCode ierr;
581b0592601SKris Buschelman 
582b0592601SKris Buschelman   PetscFunctionBegin;
583b0592601SKris Buschelman   ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr);
584ceb03754SKris Buschelman   ierr = MatConvert_SeqAIJ_SuperLU(A,MATSUPERLU,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr);
58514b396bbSKris Buschelman   PetscFunctionReturn(0);
58614b396bbSKris Buschelman }
58714b396bbSKris Buschelman EXTERN_C_END
588