xref: /petsc/src/mat/impls/aij/seq/superlu/superlu.c (revision c6db04a5321582041def2b1e244c75985478b3ef)
114b396bbSKris Buschelman 
25a46d3faSBarry Smith /*  --------------------------------------------------------------------
35a46d3faSBarry Smith 
45a46d3faSBarry Smith      This file implements a subclass of the SeqAIJ matrix class that uses
55d8b2efaSHong Zhang      the SuperLU sparse solver. You can use this as a starting point for
65a46d3faSBarry Smith      implementing your own subclass of a PETSc matrix class.
75a46d3faSBarry Smith 
85a46d3faSBarry Smith      This demonstrates a way to make an implementation inheritence of a PETSc
95a46d3faSBarry Smith      matrix type. This means constructing a new matrix type (SuperLU) by changing some
105a46d3faSBarry Smith      of the methods of the previous type (SeqAIJ), adding additional data, and possibly
115a46d3faSBarry Smith      additional method. (See any book on object oriented programming).
1214b396bbSKris Buschelman */
1314b396bbSKris Buschelman 
145a46d3faSBarry Smith /*
155a46d3faSBarry Smith      Defines the data structure for the base matrix type (SeqAIJ)
165a46d3faSBarry Smith */
17*c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h>    /*I "petscmat.h" I*/
1814b396bbSKris Buschelman 
195a46d3faSBarry Smith /*
205a46d3faSBarry Smith      SuperLU include files
215a46d3faSBarry Smith */
2214b396bbSKris Buschelman EXTERN_C_BEGIN
2314b396bbSKris Buschelman #if defined(PETSC_USE_COMPLEX)
24*c6db04a5SJed Brown #include <slu_zdefs.h>
2514b396bbSKris Buschelman #else
26*c6db04a5SJed Brown #include <slu_ddefs.h>
2714b396bbSKris Buschelman #endif
28*c6db04a5SJed Brown #include <slu_util.h>
2914b396bbSKris Buschelman EXTERN_C_END
3014b396bbSKris Buschelman 
315a46d3faSBarry Smith /*
325a46d3faSBarry Smith      This is the data we are "ADDING" to the SeqAIJ matrix type to get the SuperLU matrix type
335a46d3faSBarry Smith */
3414b396bbSKris Buschelman typedef struct {
3575af56d4SHong Zhang   SuperMatrix       A,L,U,B,X;
3675af56d4SHong Zhang   superlu_options_t options;
37da7a1d00SHong Zhang   PetscInt          *perm_c; /* column permutation vector */
38da7a1d00SHong Zhang   PetscInt          *perm_r; /* row permutations from partial pivoting */
39da7a1d00SHong Zhang   PetscInt          *etree;
40da7a1d00SHong Zhang   PetscReal         *R, *C;
4175af56d4SHong Zhang   char              equed[1];
42da7a1d00SHong Zhang   PetscInt          lwork;
4375af56d4SHong Zhang   void              *work;
44da7a1d00SHong Zhang   PetscReal         rpg, rcond;
4575af56d4SHong Zhang   mem_usage_t       mem_usage;
4675af56d4SHong Zhang   MatStructure      flg;
475d8b2efaSHong Zhang   SuperLUStat_t     stat;
48cae5a23dSHong Zhang   Mat               A_dup;
49cae5a23dSHong Zhang   PetscScalar       *rhs_dup;
5014b396bbSKris Buschelman 
5114b396bbSKris Buschelman   /* Flag to clean up (non-global) SuperLU objects during Destroy */
52ace3abfcSBarry Smith   PetscBool  CleanUpSuperLU;
53f0c56d0fSKris Buschelman } Mat_SuperLU;
5414b396bbSKris Buschelman 
55e0e586b9SHong Zhang extern PetscErrorCode MatFactorInfo_SuperLU(Mat,PetscViewer);
560481f469SBarry Smith extern PetscErrorCode MatLUFactorNumeric_SuperLU(Mat,Mat,const MatFactorInfo *);
57e0e586b9SHong Zhang extern PetscErrorCode MatDestroy_SuperLU(Mat);
58e0e586b9SHong Zhang extern PetscErrorCode MatView_SuperLU(Mat,PetscViewer);
59e0e586b9SHong Zhang extern PetscErrorCode MatAssemblyEnd_SuperLU(Mat,MatAssemblyType);
60e0e586b9SHong Zhang extern PetscErrorCode MatSolve_SuperLU(Mat,Vec,Vec);
61e0e586b9SHong Zhang extern PetscErrorCode MatSolveTranspose_SuperLU(Mat,Vec,Vec);
620481f469SBarry Smith extern PetscErrorCode MatLUFactorSymbolic_SuperLU(Mat,Mat,IS,IS,const MatFactorInfo*);
63e0e586b9SHong Zhang extern PetscErrorCode MatDuplicate_SuperLU(Mat, MatDuplicateOption, Mat *);
64e0e586b9SHong Zhang 
655a46d3faSBarry Smith /*
665a46d3faSBarry Smith     Utility function
675a46d3faSBarry Smith */
685a46d3faSBarry Smith #undef __FUNCT__
695a46d3faSBarry Smith #define __FUNCT__ "MatFactorInfo_SuperLU"
705a46d3faSBarry Smith PetscErrorCode MatFactorInfo_SuperLU(Mat A,PetscViewer viewer)
715a46d3faSBarry Smith {
725a46d3faSBarry Smith   Mat_SuperLU       *lu= (Mat_SuperLU*)A->spptr;
735a46d3faSBarry Smith   PetscErrorCode    ierr;
745a46d3faSBarry Smith   superlu_options_t options;
755a46d3faSBarry Smith 
765a46d3faSBarry Smith   PetscFunctionBegin;
775a46d3faSBarry Smith   /* check if matrix is superlu_dist type */
785a46d3faSBarry Smith   if (A->ops->solve != MatSolve_SuperLU) PetscFunctionReturn(0);
795a46d3faSBarry Smith 
805a46d3faSBarry Smith   options = lu->options;
815a46d3faSBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"SuperLU run parameters:\n");CHKERRQ(ierr);
825a46d3faSBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"  Equil: %s\n",(options.Equil != NO) ? "YES": "NO");CHKERRQ(ierr);
835a46d3faSBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"  ColPerm: %D\n",options.ColPerm);CHKERRQ(ierr);
845a46d3faSBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"  IterRefine: %D\n",options.IterRefine);CHKERRQ(ierr);
855a46d3faSBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"  SymmetricMode: %s\n",(options.SymmetricMode != NO) ? "YES": "NO");CHKERRQ(ierr);
865a46d3faSBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"  DiagPivotThresh: %g\n",options.DiagPivotThresh);CHKERRQ(ierr);
875a46d3faSBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"  PivotGrowth: %s\n",(options.PivotGrowth != NO) ? "YES": "NO");CHKERRQ(ierr);
885a46d3faSBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"  ConditionNumber: %s\n",(options.ConditionNumber != NO) ? "YES": "NO");CHKERRQ(ierr);
895a46d3faSBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"  RowPerm: %D\n",options.RowPerm);CHKERRQ(ierr);
905a46d3faSBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"  ReplaceTinyPivot: %s\n",(options.ReplaceTinyPivot != NO) ? "YES": "NO");CHKERRQ(ierr);
915a46d3faSBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"  PrintStat: %s\n",(options.PrintStat != NO) ? "YES": "NO");CHKERRQ(ierr);
925a46d3faSBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"  lwork: %D\n",lu->lwork);CHKERRQ(ierr);
93d5f3da31SBarry Smith   if (A->factortype == MAT_FACTOR_ILU){
94cffbb591SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  ILU_DropTol: %g\n",options.ILU_DropTol);CHKERRQ(ierr);
95cffbb591SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  ILU_FillTol: %g\n",options.ILU_FillTol);CHKERRQ(ierr);
96cffbb591SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  ILU_FillFactor: %g\n",options.ILU_FillFactor);CHKERRQ(ierr);
97cffbb591SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  ILU_DropRule: %D\n",options.ILU_DropRule);CHKERRQ(ierr);
98cffbb591SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  ILU_Norm: %D\n",options.ILU_Norm);CHKERRQ(ierr);
99cffbb591SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  ILU_MILU: %D\n",options.ILU_MILU);CHKERRQ(ierr);
100cffbb591SHong Zhang   }
1015a46d3faSBarry Smith   PetscFunctionReturn(0);
1025a46d3faSBarry Smith }
1035a46d3faSBarry Smith 
1045a46d3faSBarry Smith /*
1055a46d3faSBarry Smith     These are the methods provided to REPLACE the corresponding methods of the
1065a46d3faSBarry Smith    SeqAIJ matrix class. Other methods of SeqAIJ are not replaced
1075a46d3faSBarry Smith */
1085a46d3faSBarry Smith #undef __FUNCT__
1095a46d3faSBarry Smith #define __FUNCT__ "MatLUFactorNumeric_SuperLU"
1100481f469SBarry Smith PetscErrorCode MatLUFactorNumeric_SuperLU(Mat F,Mat A,const MatFactorInfo *info)
1115a46d3faSBarry Smith {
1121d5ca7f3SHong Zhang   Mat_SuperLU    *lu = (Mat_SuperLU*)F->spptr;
113cae5a23dSHong Zhang   Mat_SeqAIJ     *aa;
1145a46d3faSBarry Smith   PetscErrorCode ierr;
1155a46d3faSBarry Smith   PetscInt       sinfo;
1165a46d3faSBarry Smith   PetscReal      ferr, berr;
1175a46d3faSBarry Smith   NCformat       *Ustore;
1185a46d3faSBarry Smith   SCformat       *Lstore;
1195a46d3faSBarry Smith 
1205a46d3faSBarry Smith   PetscFunctionBegin;
1215a46d3faSBarry Smith   if (lu->flg == SAME_NONZERO_PATTERN){ /* successing numerical factorization */
1225a46d3faSBarry Smith     lu->options.Fact = SamePattern;
1235a46d3faSBarry Smith     /* Ref: ~SuperLU_3.0/EXAMPLE/dlinsolx2.c */
1245a46d3faSBarry Smith     Destroy_SuperMatrix_Store(&lu->A);
125cae5a23dSHong Zhang     if (lu->options.Equil){
126cae5a23dSHong Zhang       ierr = MatCopy_SeqAIJ(A,lu->A_dup,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
127cae5a23dSHong Zhang     }
1285a46d3faSBarry Smith     if ( lu->lwork >= 0 ) {
1295a46d3faSBarry Smith       Destroy_SuperNode_Matrix(&lu->L);
1305a46d3faSBarry Smith       Destroy_CompCol_Matrix(&lu->U);
1315a46d3faSBarry Smith       lu->options.Fact = SamePattern;
1325a46d3faSBarry Smith     }
1335a46d3faSBarry Smith   }
1345a46d3faSBarry Smith 
1355a46d3faSBarry Smith   /* Create the SuperMatrix for lu->A=A^T:
1365a46d3faSBarry Smith        Since SuperLU likes column-oriented matrices,we pass it the transpose,
1375a46d3faSBarry Smith        and then solve A^T X = B in MatSolve(). */
138cae5a23dSHong Zhang   if (lu->options.Equil){
139cae5a23dSHong Zhang     aa = (Mat_SeqAIJ*)(lu->A_dup)->data;
140cae5a23dSHong Zhang   } else {
141cae5a23dSHong Zhang     aa = (Mat_SeqAIJ*)(A)->data;
142cae5a23dSHong Zhang   }
1435a46d3faSBarry Smith #if defined(PETSC_USE_COMPLEX)
144d0f46423SBarry Smith   zCreate_CompCol_Matrix(&lu->A,A->cmap->n,A->rmap->n,aa->nz,(doublecomplex*)aa->a,aa->j,aa->i,
1455a46d3faSBarry Smith                            SLU_NC,SLU_Z,SLU_GE);
1465a46d3faSBarry Smith #else
147d0f46423SBarry Smith   dCreate_CompCol_Matrix(&lu->A,A->cmap->n,A->rmap->n,aa->nz,aa->a,aa->j,aa->i,
1485a46d3faSBarry Smith                            SLU_NC,SLU_D,SLU_GE);
1495a46d3faSBarry Smith #endif
1505a46d3faSBarry Smith 
1515a46d3faSBarry Smith   /* Numerical factorization */
1525a46d3faSBarry Smith   lu->B.ncol = 0;  /* Indicate not to solve the system */
153d5f3da31SBarry Smith   if (F->factortype == MAT_FACTOR_LU){
1545a46d3faSBarry Smith #if defined(PETSC_USE_COMPLEX)
1555a46d3faSBarry Smith     zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
1565a46d3faSBarry Smith            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
1575d8b2efaSHong Zhang            &lu->mem_usage, &lu->stat, &sinfo);
1585a46d3faSBarry Smith #else
1595a46d3faSBarry Smith     dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
1605a46d3faSBarry Smith            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
1615d8b2efaSHong Zhang            &lu->mem_usage, &lu->stat, &sinfo);
1625a46d3faSBarry Smith #endif
163d5f3da31SBarry Smith   } else if (F->factortype == MAT_FACTOR_ILU){
164cffbb591SHong Zhang     /* Compute the incomplete factorization, condition number and pivot growth */
165cffbb591SHong Zhang #if defined(PETSC_USE_COMPLEX)
166cffbb591SHong Zhang     zgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r,lu->etree, lu->equed, lu->R, lu->C,
167cffbb591SHong Zhang            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond,
1685d8b2efaSHong Zhang            &lu->mem_usage, &lu->stat, &sinfo);
169cffbb591SHong Zhang #else
170cffbb591SHong Zhang     dgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
171cffbb591SHong Zhang           &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond,
1725d8b2efaSHong Zhang           &lu->mem_usage, &lu->stat, &sinfo);
173cffbb591SHong Zhang #endif
174cffbb591SHong Zhang   } else {
175e32f2f54SBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported");
176cffbb591SHong Zhang   }
1775a46d3faSBarry Smith   if ( !sinfo || sinfo == lu->A.ncol+1 ) {
1785a46d3faSBarry Smith     if ( lu->options.PivotGrowth )
1795a46d3faSBarry Smith       ierr = PetscPrintf(PETSC_COMM_SELF,"  Recip. pivot growth = %e\n", lu->rpg);
1805a46d3faSBarry Smith     if ( lu->options.ConditionNumber )
1815a46d3faSBarry Smith       ierr = PetscPrintf(PETSC_COMM_SELF,"  Recip. condition number = %e\n", lu->rcond);
1825a46d3faSBarry Smith   } else if ( sinfo > 0 ){
1835a46d3faSBarry Smith     if ( lu->lwork == -1 ) {
1845a46d3faSBarry Smith       ierr = PetscPrintf(PETSC_COMM_SELF,"  ** Estimated memory: %D bytes\n", sinfo - lu->A.ncol);
1859308c137SBarry Smith     } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot in row %D",sinfo);
1865a46d3faSBarry Smith   } else { /* sinfo < 0 */
187e32f2f54SBarry Smith     SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_LIB, "info = %D, the %D-th argument in gssvx() had an illegal value", sinfo,-sinfo);
1885a46d3faSBarry Smith   }
1895a46d3faSBarry Smith 
1905a46d3faSBarry Smith   if ( lu->options.PrintStat ) {
1915a46d3faSBarry Smith     ierr = PetscPrintf(PETSC_COMM_SELF,"MatLUFactorNumeric_SuperLU():\n");
1925d8b2efaSHong Zhang     StatPrint(&lu->stat);
1935a46d3faSBarry Smith     Lstore = (SCformat *) lu->L.Store;
1945a46d3faSBarry Smith     Ustore = (NCformat *) lu->U.Store;
1955a46d3faSBarry Smith     ierr = PetscPrintf(PETSC_COMM_SELF,"  No of nonzeros in factor L = %D\n", Lstore->nnz);
1965a46d3faSBarry Smith     ierr = PetscPrintf(PETSC_COMM_SELF,"  No of nonzeros in factor U = %D\n", Ustore->nnz);
1975a46d3faSBarry Smith     ierr = PetscPrintf(PETSC_COMM_SELF,"  No of nonzeros in L+U = %D\n", Lstore->nnz + Ustore->nnz - lu->A.ncol);
1986da386baSHong Zhang     ierr = PetscPrintf(PETSC_COMM_SELF,"  L\\U MB %.3f\ttotal MB needed %.3f\n",
1996da386baSHong Zhang 	       lu->mem_usage.for_lu/1e6, lu->mem_usage.total_needed/1e6);
2005a46d3faSBarry Smith   }
2015a46d3faSBarry Smith 
2025a46d3faSBarry Smith   lu->flg = SAME_NONZERO_PATTERN;
2031d5ca7f3SHong Zhang   F->ops->solve          = MatSolve_SuperLU;
2041d5ca7f3SHong Zhang   F->ops->solvetranspose = MatSolveTranspose_SuperLU;
2055a46d3faSBarry Smith   PetscFunctionReturn(0);
2065a46d3faSBarry Smith }
2075a46d3faSBarry Smith 
20814b396bbSKris Buschelman #undef __FUNCT__
209f0c56d0fSKris Buschelman #define __FUNCT__ "MatDestroy_SuperLU"
210dfbe8321SBarry Smith PetscErrorCode MatDestroy_SuperLU(Mat A)
21114b396bbSKris Buschelman {
212dfbe8321SBarry Smith   PetscErrorCode ierr;
213f0c56d0fSKris Buschelman   Mat_SuperLU    *lu=(Mat_SuperLU*)A->spptr;
21414b396bbSKris Buschelman 
21514b396bbSKris Buschelman   PetscFunctionBegin;
21675af56d4SHong Zhang   if (lu->CleanUpSuperLU) { /* Free the SuperLU datastructures */
21775af56d4SHong Zhang     Destroy_SuperMatrix_Store(&lu->A);
21875af56d4SHong Zhang     Destroy_SuperMatrix_Store(&lu->B);
21975af56d4SHong Zhang     Destroy_SuperMatrix_Store(&lu->X);
2205d8b2efaSHong Zhang     StatFree(&lu->stat);
2210e742b69SHong Zhang     if ( lu->lwork >= 0 ) {
2220e742b69SHong Zhang       Destroy_SuperNode_Matrix(&lu->L);
2230e742b69SHong Zhang       Destroy_CompCol_Matrix(&lu->U);
2240e742b69SHong Zhang     }
2250e742b69SHong Zhang   }
2269ce50919SHong Zhang 
2279ce50919SHong Zhang   ierr = PetscFree(lu->etree);CHKERRQ(ierr);
2289ce50919SHong Zhang   ierr = PetscFree(lu->perm_r);CHKERRQ(ierr);
2299ce50919SHong Zhang   ierr = PetscFree(lu->perm_c);CHKERRQ(ierr);
2309ce50919SHong Zhang   ierr = PetscFree(lu->R);CHKERRQ(ierr);
2319ce50919SHong Zhang   ierr = PetscFree(lu->C);CHKERRQ(ierr);
2320e742b69SHong Zhang 
233d954db57SHong Zhang   /* clear composed functions */
234d954db57SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatFactorGetSolverPackage_C","",PETSC_NULL);CHKERRQ(ierr);
235790a96dfSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatSuperluSetILUDropTol_C","",PETSC_NULL);CHKERRQ(ierr);
236d954db57SHong Zhang 
237b24902e0SBarry Smith   ierr = MatDestroy_SeqAIJ(A);CHKERRQ(ierr);
238cae5a23dSHong Zhang   if (lu->A_dup){ierr = MatDestroy(lu->A_dup);CHKERRQ(ierr);}
239c31cb41cSBarry Smith   ierr = PetscFree(lu->rhs_dup);CHKERRQ(ierr);
24014b396bbSKris Buschelman   PetscFunctionReturn(0);
24114b396bbSKris Buschelman }
24214b396bbSKris Buschelman 
24314b396bbSKris Buschelman #undef __FUNCT__
244f0c56d0fSKris Buschelman #define __FUNCT__ "MatView_SuperLU"
245dfbe8321SBarry Smith PetscErrorCode MatView_SuperLU(Mat A,PetscViewer viewer)
24614b396bbSKris Buschelman {
247dfbe8321SBarry Smith   PetscErrorCode    ierr;
248ace3abfcSBarry Smith   PetscBool         iascii;
24914b396bbSKris Buschelman   PetscViewerFormat format;
25014b396bbSKris Buschelman 
25114b396bbSKris Buschelman   PetscFunctionBegin;
2522692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
25332077d6dSBarry Smith   if (iascii) {
25414b396bbSKris Buschelman     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
2552f59403fSHong Zhang     if (format == PETSC_VIEWER_ASCII_INFO) {
256f0c56d0fSKris Buschelman       ierr = MatFactorInfo_SuperLU(A,viewer);CHKERRQ(ierr);
25714b396bbSKris Buschelman     }
25814b396bbSKris Buschelman   }
25914b396bbSKris Buschelman   PetscFunctionReturn(0);
26014b396bbSKris Buschelman }
26114b396bbSKris Buschelman 
26214b396bbSKris Buschelman 
26314b396bbSKris Buschelman #undef __FUNCT__
264c29e884eSHong Zhang #define __FUNCT__ "MatSolve_SuperLU_Private"
265c29e884eSHong Zhang PetscErrorCode MatSolve_SuperLU_Private(Mat A,Vec b,Vec x)
26614b396bbSKris Buschelman {
267f0c56d0fSKris Buschelman   Mat_SuperLU    *lu = (Mat_SuperLU*)A->spptr;
26875af56d4SHong Zhang   PetscScalar    *barray,*xarray;
269dfbe8321SBarry Smith   PetscErrorCode ierr;
270cae5a23dSHong Zhang   PetscInt       info,i,n=x->map->n;
271da7a1d00SHong Zhang   PetscReal      ferr,berr;
27214b396bbSKris Buschelman 
27314b396bbSKris Buschelman   PetscFunctionBegin;
274937a6b0eSHong Zhang   if ( lu->lwork == -1 ) {
275937a6b0eSHong Zhang     PetscFunctionReturn(0);
276937a6b0eSHong Zhang   }
277cae5a23dSHong Zhang 
27875af56d4SHong Zhang   lu->B.ncol = 1;   /* Set the number of right-hand side */
279cae5a23dSHong Zhang   if (lu->options.Equil && !lu->rhs_dup){
280cae5a23dSHong Zhang     /* superlu overwrites b when Equil is used, thus create rhs_dup to keep user's b unchanged */
281cae5a23dSHong Zhang     ierr = PetscMalloc(n*sizeof(PetscScalar),&lu->rhs_dup);CHKERRQ(ierr);
282cae5a23dSHong Zhang   }
283cae5a23dSHong Zhang   if (lu->options.Equil){
284cae5a23dSHong Zhang     /* Copy b into rsh_dup */
28575af56d4SHong Zhang     ierr = VecGetArray(b,&barray);CHKERRQ(ierr);
286cae5a23dSHong Zhang     ierr = PetscMemcpy(lu->rhs_dup,barray,n*sizeof(PetscScalar));CHKERRQ(ierr);
287cae5a23dSHong Zhang     ierr = VecRestoreArray(b,&barray);CHKERRQ(ierr);
288cae5a23dSHong Zhang     barray = lu->rhs_dup;
289cae5a23dSHong Zhang   } else {
290cae5a23dSHong Zhang     ierr = VecGetArray(b,&barray);CHKERRQ(ierr);
291cae5a23dSHong Zhang   }
29275af56d4SHong Zhang   ierr = VecGetArray(x,&xarray);CHKERRQ(ierr);
2935fe6150dSHong Zhang 
2945fe6150dSHong Zhang #if defined(PETSC_USE_COMPLEX)
2955fe6150dSHong Zhang   ((DNformat*)lu->B.Store)->nzval = (doublecomplex*)barray;
2965fe6150dSHong Zhang   ((DNformat*)lu->X.Store)->nzval = (doublecomplex*)xarray;
2975fe6150dSHong Zhang #else
2985fe6150dSHong Zhang   ((DNformat*)lu->B.Store)->nzval = barray;
29975af56d4SHong Zhang   ((DNformat*)lu->X.Store)->nzval = xarray;
3005fe6150dSHong Zhang #endif
30175af56d4SHong Zhang 
30275af56d4SHong Zhang   lu->options.Fact = FACTORED; /* Indicate the factored form of A is supplied. */
303d5f3da31SBarry Smith   if (A->factortype == MAT_FACTOR_LU){
3048914a3f7SHong Zhang #if defined(PETSC_USE_COMPLEX)
3058914a3f7SHong Zhang     zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
3068914a3f7SHong Zhang            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
3075d8b2efaSHong Zhang            &lu->mem_usage, &lu->stat, &info);
3088914a3f7SHong Zhang #else
30975af56d4SHong Zhang     dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
31075af56d4SHong Zhang            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
3115d8b2efaSHong Zhang            &lu->mem_usage, &lu->stat, &info);
3128914a3f7SHong Zhang #endif
313d5f3da31SBarry Smith   } else if (A->factortype == MAT_FACTOR_ILU){
314cffbb591SHong Zhang #if defined(PETSC_USE_COMPLEX)
315cffbb591SHong Zhang     zgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
316cffbb591SHong Zhang            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond,
3175d8b2efaSHong Zhang            &lu->mem_usage, &lu->stat, &info);
318cffbb591SHong Zhang #else
319cffbb591SHong Zhang     dgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
320cffbb591SHong Zhang            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond,
3215d8b2efaSHong Zhang            &lu->mem_usage, &lu->stat, &info);
322cffbb591SHong Zhang #endif
323cffbb591SHong Zhang   } else {
324e32f2f54SBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported");
325cffbb591SHong Zhang   }
326cae5a23dSHong Zhang   if (!lu->options.Equil){
32775af56d4SHong Zhang     ierr = VecRestoreArray(b,&barray);CHKERRQ(ierr);
328cae5a23dSHong Zhang   }
32975af56d4SHong Zhang   ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
33075af56d4SHong Zhang 
331958c9bccSBarry Smith   if ( !info || info == lu->A.ncol+1 ) {
33275af56d4SHong Zhang     if ( lu->options.IterRefine ) {
3338914a3f7SHong Zhang       ierr = PetscPrintf(PETSC_COMM_SELF,"Iterative Refinement:\n");
3348914a3f7SHong Zhang       ierr = PetscPrintf(PETSC_COMM_SELF,"  %8s%8s%16s%16s\n", "rhs", "Steps", "FERR", "BERR");
33575af56d4SHong Zhang       for (i = 0; i < 1; ++i)
3365d8b2efaSHong Zhang         ierr = PetscPrintf(PETSC_COMM_SELF,"  %8d%8d%16e%16e\n", i+1, lu->stat.RefineSteps, ferr, berr);
33775af56d4SHong Zhang     }
3388914a3f7SHong Zhang   } else if ( info > 0 ){
3398914a3f7SHong Zhang     if ( lu->lwork == -1 ) {
34077431f27SBarry Smith       ierr = PetscPrintf(PETSC_COMM_SELF,"  ** Estimated memory: %D bytes\n", info - lu->A.ncol);
3418914a3f7SHong Zhang     } else {
34277431f27SBarry Smith       ierr = PetscPrintf(PETSC_COMM_SELF,"  Warning: gssvx() returns info %D\n",info);
3438914a3f7SHong Zhang     }
3448914a3f7SHong Zhang   } else if (info < 0){
345e32f2f54SBarry Smith     SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_LIB, "info = %D, the %D-th argument in gssvx() had an illegal value", info,-info);
34614b396bbSKris Buschelman   }
34714b396bbSKris Buschelman 
3488914a3f7SHong Zhang   if ( lu->options.PrintStat ) {
3498914a3f7SHong Zhang     ierr = PetscPrintf(PETSC_COMM_SELF,"MatSolve__SuperLU():\n");
3505d8b2efaSHong Zhang     StatPrint(&lu->stat);
3518914a3f7SHong Zhang   }
35275af56d4SHong Zhang   PetscFunctionReturn(0);
35375af56d4SHong Zhang }
35414b396bbSKris Buschelman 
35514b396bbSKris Buschelman #undef __FUNCT__
356c29e884eSHong Zhang #define __FUNCT__ "MatSolve_SuperLU"
357c29e884eSHong Zhang PetscErrorCode MatSolve_SuperLU(Mat A,Vec b,Vec x)
358c29e884eSHong Zhang {
359c29e884eSHong Zhang   Mat_SuperLU    *lu = (Mat_SuperLU*)A->spptr;
360c29e884eSHong Zhang   PetscErrorCode ierr;
361c29e884eSHong Zhang 
362c29e884eSHong Zhang   PetscFunctionBegin;
363c29e884eSHong Zhang   lu->options.Trans = TRANS;
364c29e884eSHong Zhang   ierr = MatSolve_SuperLU_Private(A,b,x);CHKERRQ(ierr);
365c29e884eSHong Zhang   PetscFunctionReturn(0);
366c29e884eSHong Zhang }
367c29e884eSHong Zhang 
368c29e884eSHong Zhang #undef __FUNCT__
369c7c1fe80SHong Zhang #define __FUNCT__ "MatSolveTranspose_SuperLU"
370c7c1fe80SHong Zhang PetscErrorCode MatSolveTranspose_SuperLU(Mat A,Vec b,Vec x)
371c7c1fe80SHong Zhang {
372c7c1fe80SHong Zhang   Mat_SuperLU    *lu = (Mat_SuperLU*)A->spptr;
373c7c1fe80SHong Zhang   PetscErrorCode ierr;
374c7c1fe80SHong Zhang 
375c7c1fe80SHong Zhang   PetscFunctionBegin;
376c7c1fe80SHong Zhang   lu->options.Trans = NOTRANS;
377c29e884eSHong Zhang   ierr = MatSolve_SuperLU_Private(A,b,x);CHKERRQ(ierr);
378c7c1fe80SHong Zhang   PetscFunctionReturn(0);
379c7c1fe80SHong Zhang }
380c7c1fe80SHong Zhang 
38114b396bbSKris Buschelman /*
38214b396bbSKris Buschelman    Note the r permutation is ignored
38314b396bbSKris Buschelman */
38414b396bbSKris Buschelman #undef __FUNCT__
385f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_SuperLU"
3860481f469SBarry Smith PetscErrorCode MatLUFactorSymbolic_SuperLU(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
38714b396bbSKris Buschelman {
3881d5ca7f3SHong Zhang   Mat_SuperLU    *lu = (Mat_SuperLU*)(F->spptr);
389b24902e0SBarry Smith 
390b24902e0SBarry Smith   PetscFunctionBegin;
391b24902e0SBarry Smith   lu->flg                 = DIFFERENT_NONZERO_PATTERN;
392b24902e0SBarry Smith   lu->CleanUpSuperLU      = PETSC_TRUE;
3931d5ca7f3SHong Zhang   F->ops->lufactornumeric = MatLUFactorNumeric_SuperLU;
394b24902e0SBarry Smith   PetscFunctionReturn(0);
395b24902e0SBarry Smith }
396b24902e0SBarry Smith 
39735bd34faSBarry Smith EXTERN_C_BEGIN
39835bd34faSBarry Smith #undef __FUNCT__
3995ccb76cbSHong Zhang #define __FUNCT__ "MatSuperluSetILUDropTol_SuperLU"
4005ccb76cbSHong Zhang PetscErrorCode MatSuperluSetILUDropTol_SuperLU(Mat F,PetscReal dtol)
4015ccb76cbSHong Zhang {
4025ccb76cbSHong Zhang   Mat_SuperLU *lu= (Mat_SuperLU*)F->spptr;
4035ccb76cbSHong Zhang 
4045ccb76cbSHong Zhang   PetscFunctionBegin;
4055ccb76cbSHong Zhang   lu->options.ILU_DropTol = dtol;
4065ccb76cbSHong Zhang   PetscFunctionReturn(0);
4075ccb76cbSHong Zhang }
4085ccb76cbSHong Zhang EXTERN_C_END
4095ccb76cbSHong Zhang 
4105ccb76cbSHong Zhang #undef __FUNCT__
4115ccb76cbSHong Zhang #define __FUNCT__ "MatSuperluSetILUDropTol"
4125ccb76cbSHong Zhang /*@
4135ccb76cbSHong Zhang   MatSuperluSetILUDropTol - Set SuperLU ILU drop tolerance
4145ccb76cbSHong Zhang    Logically Collective on Mat
4155ccb76cbSHong Zhang 
4165ccb76cbSHong Zhang    Input Parameters:
4175ccb76cbSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-SuperLU interface
4185ccb76cbSHong Zhang -  dtol - drop tolerance
4195ccb76cbSHong Zhang 
4205ccb76cbSHong Zhang   Options Database:
4215ccb76cbSHong Zhang .   -mat_superlu_ilu_droptol <dtol>
4225ccb76cbSHong Zhang 
4235ccb76cbSHong Zhang    Level: beginner
4245ccb76cbSHong Zhang 
4255ccb76cbSHong Zhang    References: SuperLU Users' Guide
4265ccb76cbSHong Zhang 
4275ccb76cbSHong Zhang .seealso: MatGetFactor()
4285ccb76cbSHong Zhang @*/
4295ccb76cbSHong Zhang PetscErrorCode MatSuperluSetILUDropTol(Mat F,PetscReal dtol)
4305ccb76cbSHong Zhang {
4315ccb76cbSHong Zhang   PetscErrorCode ierr;
4325ccb76cbSHong Zhang 
4335ccb76cbSHong Zhang   PetscFunctionBegin;
4345ccb76cbSHong Zhang   PetscValidHeaderSpecific(F,MAT_CLASSID,1);
4355ccb76cbSHong Zhang   PetscValidLogicalCollectiveInt(F,dtol,2);
4365ccb76cbSHong Zhang   ierr = PetscTryMethod(F,"MatSuperluSetILUDropTol_C",(Mat,PetscReal),(F,dtol));CHKERRQ(ierr);
4375ccb76cbSHong Zhang   PetscFunctionReturn(0);
4385ccb76cbSHong Zhang }
4395ccb76cbSHong Zhang 
4405ccb76cbSHong Zhang EXTERN_C_BEGIN
4415ccb76cbSHong Zhang #undef __FUNCT__
44235bd34faSBarry Smith #define __FUNCT__ "MatFactorGetSolverPackage_seqaij_superlu"
44335bd34faSBarry Smith PetscErrorCode MatFactorGetSolverPackage_seqaij_superlu(Mat A,const MatSolverPackage *type)
44435bd34faSBarry Smith {
44535bd34faSBarry Smith   PetscFunctionBegin;
4462692d6eeSBarry Smith   *type = MATSOLVERSUPERLU;
44735bd34faSBarry Smith   PetscFunctionReturn(0);
44835bd34faSBarry Smith }
44935bd34faSBarry Smith EXTERN_C_END
45035bd34faSBarry Smith 
451b24902e0SBarry Smith 
452b24902e0SBarry Smith /*MC
453ba992d64SSatish Balay   MATSOLVERSUPERLU = "superlu" - A solver package providing solvers LU and ILU for sequential matrices
454b24902e0SBarry Smith   via the external package SuperLU.
455b24902e0SBarry Smith 
456e2e64c6bSBarry Smith   Use ./configure --download-superlu to have PETSc installed with SuperLU
457b24902e0SBarry Smith 
458b24902e0SBarry Smith   Options Database Keys:
4592de22e2eSHong Zhang +  -mat_superlu_equil: <FALSE> Equil (None)
4602de22e2eSHong Zhang .  -mat_superlu_colperm <COLAMD> (choose one of) NATURAL MMD_ATA MMD_AT_PLUS_A COLAMD
4612de22e2eSHong Zhang .  -mat_superlu_iterrefine <NOREFINE> (choose one of) NOREFINE SINGLE DOUBLE EXTRA
4622de22e2eSHong Zhang .  -mat_superlu_symmetricmode: <FALSE> SymmetricMode (None)
4632de22e2eSHong Zhang .  -mat_superlu_diagpivotthresh <1>: DiagPivotThresh (None)
4642de22e2eSHong Zhang .  -mat_superlu_pivotgrowth: <FALSE> PivotGrowth (None)
4652de22e2eSHong Zhang .  -mat_superlu_conditionnumber: <FALSE> ConditionNumber (None)
4662de22e2eSHong Zhang .  -mat_superlu_rowperm <NOROWPERM> (choose one of) NOROWPERM LargeDiag
4672de22e2eSHong Zhang .  -mat_superlu_replacetinypivot: <FALSE> ReplaceTinyPivot (None)
4682de22e2eSHong Zhang .  -mat_superlu_printstat: <FALSE> PrintStat (None)
4692de22e2eSHong Zhang .  -mat_superlu_lwork <0>: size of work array in bytes used by factorization (None)
4702de22e2eSHong Zhang .  -mat_superlu_ilu_droptol <0>: ILU_DropTol (None)
4712de22e2eSHong Zhang .  -mat_superlu_ilu_filltol <0>: ILU_FillTol (None)
4722de22e2eSHong Zhang .  -mat_superlu_ilu_fillfactor <0>: ILU_FillFactor (None)
4732de22e2eSHong Zhang .  -mat_superlu_ilu_droprull <0>: ILU_DropRule (None)
4742de22e2eSHong Zhang .  -mat_superlu_ilu_norm <0>: ILU_Norm (None)
4752de22e2eSHong Zhang -  -mat_superlu_ilu_milu <0>: ILU_MILU (None)
476b24902e0SBarry Smith 
4772692d6eeSBarry Smith    Notes: Do not confuse this with MATSOLVERSUPERLU_DIST which is for parallel sparse solves
4785c9eb25fSBarry Smith 
479b24902e0SBarry Smith    Level: beginner
480b24902e0SBarry Smith 
481ba992d64SSatish Balay .seealso: PCLU, PCILU, MATSOLVERSUPERLU_DIST, MATSOLVERMUMPS, MATSOLVERSPOOLES, PCFactorSetMatSolverPackage(), MatSolverPackage
482b24902e0SBarry Smith M*/
483b24902e0SBarry Smith 
4845c9eb25fSBarry Smith EXTERN_C_BEGIN
485b24902e0SBarry Smith #undef __FUNCT__
486b24902e0SBarry Smith #define __FUNCT__ "MatGetFactor_seqaij_superlu"
4875c9eb25fSBarry Smith PetscErrorCode MatGetFactor_seqaij_superlu(Mat A,MatFactorType ftype,Mat *F)
488b24902e0SBarry Smith {
48914b396bbSKris Buschelman   Mat            B;
490f0c56d0fSKris Buschelman   Mat_SuperLU    *lu;
4916849ba73SBarry Smith   PetscErrorCode ierr;
4925d8b2efaSHong Zhang   PetscInt       indx,m=A->rmap->n,n=A->cmap->n;
493ace3abfcSBarry Smith   PetscBool      flg;
4945ca28756SHong Zhang   const char     *colperm[]={"NATURAL","MMD_ATA","MMD_AT_PLUS_A","COLAMD"}; /* MY_PERMC - not supported by the petsc interface yet */
4955ca28756SHong Zhang   const char     *iterrefine[]={"NOREFINE", "SINGLE", "DOUBLE", "EXTRA"};
4965ca28756SHong Zhang   const char     *rowperm[]={"NOROWPERM", "LargeDiag"}; /* MY_PERMC - not supported by the petsc interface yet */
49714b396bbSKris Buschelman 
49814b396bbSKris Buschelman   PetscFunctionBegin;
4997adad957SLisandro Dalcin   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
500d0f46423SBarry Smith   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
5017adad957SLisandro Dalcin   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
502899d7b4fSKris Buschelman   ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
503f0c56d0fSKris Buschelman 
504cffbb591SHong Zhang   if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU){
505b24902e0SBarry Smith     B->ops->lufactorsymbolic  = MatLUFactorSymbolic_SuperLU;
506cffbb591SHong Zhang     B->ops->ilufactorsymbolic = MatLUFactorSymbolic_SuperLU;
507b3a44c85SBarry Smith   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported");
508cffbb591SHong Zhang 
509b24902e0SBarry Smith   B->ops->destroy          = MatDestroy_SuperLU;
5103519fcdcSHong Zhang   B->ops->view             = MatView_SuperLU;
511d5f3da31SBarry Smith   B->factortype            = ftype;
51294b7f48cSBarry Smith   B->assembled             = PETSC_TRUE;  /* required by -ksp_view */
5135c9eb25fSBarry Smith   B->preallocated          = PETSC_TRUE;
51414b396bbSKris Buschelman 
515b24902e0SBarry Smith   ierr = PetscNewLog(B,Mat_SuperLU,&lu);CHKERRQ(ierr);
516cae5a23dSHong Zhang 
517cffbb591SHong Zhang   if (ftype == MAT_FACTOR_LU){
5189ce50919SHong Zhang     set_default_options(&lu->options);
5193d6581fbSHong Zhang     /* Comments from SuperLU_4.0/SRC/dgssvx.c:
5203d6581fbSHong Zhang       "Whether or not the system will be equilibrated depends on the
5213d6581fbSHong Zhang        scaling of the matrix A, but if equilibration is used, A is
5223d6581fbSHong Zhang        overwritten by diag(R)*A*diag(C) and B by diag(R)*B
5233d6581fbSHong Zhang        (if options->Trans=NOTRANS) or diag(C)*B (if options->Trans = TRANS or CONJ)."
5243d6581fbSHong Zhang      We set 'options.Equil = NO' as default because additional space is needed for it.
5253d6581fbSHong Zhang     */
5263d6581fbSHong Zhang     lu->options.Equil = NO;
527cffbb591SHong Zhang   } else if (ftype == MAT_FACTOR_ILU){
528cffbb591SHong Zhang     /* Set the default input options of ilu:
529cffbb591SHong Zhang 	options.Fact = DOFACT;
5303d6581fbSHong Zhang 	options.Equil = YES;           // must be YES for ilu - don't know why
531cffbb591SHong Zhang 	options.ColPerm = COLAMD;
532cffbb591SHong Zhang 	options.DiagPivotThresh = 0.1; //different from complete LU
533cffbb591SHong Zhang 	options.Trans = NOTRANS;
534cffbb591SHong Zhang 	options.IterRefine = NOREFINE;
535cffbb591SHong Zhang 	options.SymmetricMode = NO;
536cffbb591SHong Zhang 	options.PivotGrowth = NO;
537cffbb591SHong Zhang 	options.ConditionNumber = NO;
538cffbb591SHong Zhang 	options.PrintStat = YES;
539cffbb591SHong Zhang 	options.RowPerm = LargeDiag;
540cffbb591SHong Zhang 	options.ILU_DropTol = 1e-4;
541cffbb591SHong Zhang 	options.ILU_FillTol = 1e-2;
542cffbb591SHong Zhang 	options.ILU_FillFactor = 10.0;
543cffbb591SHong Zhang 	options.ILU_DropRule = DROP_BASIC | DROP_AREA;
544cffbb591SHong Zhang 	options.ILU_Norm = INF_NORM;
545cffbb591SHong Zhang 	options.ILU_MILU = SMILU_2;
546cffbb591SHong Zhang     */
547cffbb591SHong Zhang     ilu_set_default_options(&lu->options);
548cffbb591SHong Zhang   }
5499ce50919SHong Zhang   lu->options.PrintStat = NO;
5501a2e2f44SHong Zhang 
5515d8b2efaSHong Zhang   /* Initialize the statistics variables. */
5525d8b2efaSHong Zhang   StatInit(&lu->stat);
5538914a3f7SHong Zhang   lu->lwork = 0;   /* allocate space internally by system malloc */
5549ce50919SHong Zhang 
5557adad957SLisandro Dalcin   ierr = PetscOptionsBegin(((PetscObject)A)->comm,((PetscObject)A)->prefix,"SuperLU Options","Mat");CHKERRQ(ierr);
556acfcf0e5SJed Brown     ierr = PetscOptionsBool("-mat_superlu_equil","Equil","None",(PetscBool)lu->options.Equil,(PetscBool*)&lu->options.Equil,0);CHKERRQ(ierr);
5578914a3f7SHong Zhang     ierr = PetscOptionsEList("-mat_superlu_colperm","ColPerm","None",colperm,4,colperm[3],&indx,&flg);CHKERRQ(ierr);
5589ce50919SHong Zhang     if (flg) {lu->options.ColPerm = (colperm_t)indx;}
5598914a3f7SHong Zhang     ierr = PetscOptionsEList("-mat_superlu_iterrefine","IterRefine","None",iterrefine,4,iterrefine[0],&indx,&flg);CHKERRQ(ierr);
5609ce50919SHong Zhang     if (flg) { lu->options.IterRefine = (IterRefine_t)indx;}
561acfcf0e5SJed Brown     ierr = PetscOptionsBool("-mat_superlu_symmetricmode","SymmetricMode","None",(PetscBool)lu->options.SymmetricMode,&flg,0);CHKERRQ(ierr);
5629ce50919SHong Zhang     if (flg) lu->options.SymmetricMode = YES;
5638914a3f7SHong Zhang     ierr = PetscOptionsReal("-mat_superlu_diagpivotthresh","DiagPivotThresh","None",lu->options.DiagPivotThresh,&lu->options.DiagPivotThresh,PETSC_NULL);CHKERRQ(ierr);
564acfcf0e5SJed Brown     ierr = PetscOptionsBool("-mat_superlu_pivotgrowth","PivotGrowth","None",(PetscBool)lu->options.PivotGrowth,&flg,0);CHKERRQ(ierr);
5659ce50919SHong Zhang     if (flg) lu->options.PivotGrowth = YES;
566acfcf0e5SJed Brown     ierr = PetscOptionsBool("-mat_superlu_conditionnumber","ConditionNumber","None",(PetscBool)lu->options.ConditionNumber,&flg,0);CHKERRQ(ierr);
5679ce50919SHong Zhang     if (flg) lu->options.ConditionNumber = YES;
568d7ebd59bSHong Zhang     ierr = PetscOptionsEList("-mat_superlu_rowperm","rowperm","None",rowperm,2,rowperm[lu->options.RowPerm],&indx,&flg);CHKERRQ(ierr);
5699ce50919SHong Zhang     if (flg) {lu->options.RowPerm = (rowperm_t)indx;}
570acfcf0e5SJed Brown     ierr = PetscOptionsBool("-mat_superlu_replacetinypivot","ReplaceTinyPivot","None",(PetscBool)lu->options.ReplaceTinyPivot,&flg,0);CHKERRQ(ierr);
5719ce50919SHong Zhang     if (flg) lu->options.ReplaceTinyPivot = YES;
572acfcf0e5SJed Brown     ierr = PetscOptionsBool("-mat_superlu_printstat","PrintStat","None",(PetscBool)lu->options.PrintStat,&flg,0);CHKERRQ(ierr);
5739ce50919SHong Zhang     if (flg) lu->options.PrintStat = YES;
5748914a3f7SHong Zhang     ierr = PetscOptionsInt("-mat_superlu_lwork","size of work array in bytes used by factorization","None",lu->lwork,&lu->lwork,PETSC_NULL);CHKERRQ(ierr);
5755fe6150dSHong Zhang     if (lu->lwork > 0 ){
5765fe6150dSHong Zhang       ierr = PetscMalloc(lu->lwork,&lu->work);CHKERRQ(ierr);
5775fe6150dSHong Zhang     } else if (lu->lwork != 0 && lu->lwork != -1){
57877431f27SBarry Smith       ierr = PetscPrintf(PETSC_COMM_SELF,"   Warning: lwork %D is not supported by SUPERLU. The default lwork=0 is used.\n",lu->lwork);
5798914a3f7SHong Zhang       lu->lwork = 0;
5808914a3f7SHong Zhang     }
581cffbb591SHong Zhang     /* ilu options */
582cffbb591SHong Zhang     ierr = PetscOptionsReal("-mat_superlu_ilu_droptol","ILU_DropTol","None",lu->options.ILU_DropTol,&lu->options.ILU_DropTol,PETSC_NULL);CHKERRQ(ierr);
583cffbb591SHong Zhang     ierr = PetscOptionsReal("-mat_superlu_ilu_filltol","ILU_FillTol","None",lu->options.ILU_FillTol,&lu->options.ILU_FillTol,PETSC_NULL);CHKERRQ(ierr);
584cffbb591SHong Zhang     ierr = PetscOptionsReal("-mat_superlu_ilu_fillfactor","ILU_FillFactor","None",lu->options.ILU_FillFactor,&lu->options.ILU_FillFactor,PETSC_NULL);CHKERRQ(ierr);
585cffbb591SHong Zhang     ierr = PetscOptionsInt("-mat_superlu_ilu_droprull","ILU_DropRule","None",lu->options.ILU_DropRule,&lu->options.ILU_DropRule,PETSC_NULL);CHKERRQ(ierr);
586cffbb591SHong Zhang     ierr = PetscOptionsInt("-mat_superlu_ilu_norm","ILU_Norm","None",lu->options.ILU_Norm,&indx,&flg);CHKERRQ(ierr);
587cffbb591SHong Zhang     if (flg){
588cffbb591SHong Zhang       lu->options.ILU_Norm = (norm_t)indx;
589cffbb591SHong Zhang     }
590cffbb591SHong Zhang     ierr = PetscOptionsInt("-mat_superlu_ilu_milu","ILU_MILU","None",lu->options.ILU_MILU,&indx,&flg);CHKERRQ(ierr);
591cffbb591SHong Zhang     if (flg){
592cffbb591SHong Zhang       lu->options.ILU_MILU = (milu_t)indx;
593cffbb591SHong Zhang     }
5949ce50919SHong Zhang   PetscOptionsEnd();
59538abfdaeSHong Zhang   if (lu->options.Equil == YES) {
59638abfdaeSHong Zhang     /* superlu overwrites input matrix and rhs when Equil is used, thus create A_dup to keep user's A unchanged */
59738abfdaeSHong Zhang     ierr = MatDuplicate_SeqAIJ(A,MAT_COPY_VALUES,&lu->A_dup);CHKERRQ(ierr);
59838abfdaeSHong Zhang   }
5999ce50919SHong Zhang 
6005d8b2efaSHong Zhang   /* Allocate spaces (notice sizes are for the transpose) */
6015d8b2efaSHong Zhang   ierr = PetscMalloc(m*sizeof(PetscInt),&lu->etree);CHKERRQ(ierr);
6025d8b2efaSHong Zhang   ierr = PetscMalloc(n*sizeof(PetscInt),&lu->perm_r);CHKERRQ(ierr);
6035d8b2efaSHong Zhang   ierr = PetscMalloc(m*sizeof(PetscInt),&lu->perm_c);CHKERRQ(ierr);
6045d8b2efaSHong Zhang   ierr = PetscMalloc(n*sizeof(PetscScalar),&lu->R);CHKERRQ(ierr);
6055d8b2efaSHong Zhang   ierr = PetscMalloc(m*sizeof(PetscScalar),&lu->C);CHKERRQ(ierr);
6065d8b2efaSHong Zhang 
6075d8b2efaSHong Zhang   /* create rhs and solution x without allocate space for .Store */
6085d8b2efaSHong Zhang #if defined(PETSC_USE_COMPLEX)
6095d8b2efaSHong Zhang   zCreate_Dense_Matrix(&lu->B, m, 1, PETSC_NULL, m, SLU_DN, SLU_Z, SLU_GE);
6105d8b2efaSHong Zhang   zCreate_Dense_Matrix(&lu->X, m, 1, PETSC_NULL, m, SLU_DN, SLU_Z, SLU_GE);
6115d8b2efaSHong Zhang #else
6125d8b2efaSHong Zhang   dCreate_Dense_Matrix(&lu->B, m, 1, PETSC_NULL, m, SLU_DN, SLU_D, SLU_GE);
6135d8b2efaSHong Zhang   dCreate_Dense_Matrix(&lu->X, m, 1, PETSC_NULL, m, SLU_DN, SLU_D, SLU_GE);
6145d8b2efaSHong Zhang #endif
6155d8b2efaSHong Zhang 
61675af56d4SHong Zhang #ifdef SUPERLU2
6175c9eb25fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatCreateNull","MatCreateNull_SuperLU",(void(*)(void))MatCreateNull_SuperLU);CHKERRQ(ierr);
61875af56d4SHong Zhang #endif
61935bd34faSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_seqaij_superlu",MatFactorGetSolverPackage_seqaij_superlu);CHKERRQ(ierr);
6205ccb76cbSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSuperluSetILUDropTol_C","MatSuperluSetILUDropTol_SuperLU",MatSuperluSetILUDropTol_SuperLU);CHKERRQ(ierr);
6215c9eb25fSBarry Smith   B->spptr = lu;
622899d7b4fSKris Buschelman   *F = B;
62314b396bbSKris Buschelman   PetscFunctionReturn(0);
62414b396bbSKris Buschelman }
6255c9eb25fSBarry Smith EXTERN_C_END
626d954db57SHong Zhang 
627