xref: /petsc/src/mat/impls/aij/seq/superlu/superlu.c (revision f23aa3dd738493dcb3a70a8c0c7f5454aa9150c2)
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 */
17c6db04a5SJed 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)
243cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
253cb270beSHong Zhang #include <slu_cdefs.h>
263cb270beSHong Zhang #else
27c6db04a5SJed Brown #include <slu_zdefs.h>
283cb270beSHong Zhang #endif
293cb270beSHong Zhang #else
303cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
313cb270beSHong Zhang #include <slu_sdefs.h>
3214b396bbSKris Buschelman #else
33c6db04a5SJed Brown #include <slu_ddefs.h>
3414b396bbSKris Buschelman #endif
353cb270beSHong Zhang #endif
36c6db04a5SJed Brown #include <slu_util.h>
3714b396bbSKris Buschelman EXTERN_C_END
3814b396bbSKris Buschelman 
395a46d3faSBarry Smith /*
405a46d3faSBarry Smith      This is the data we are "ADDING" to the SeqAIJ matrix type to get the SuperLU matrix type
415a46d3faSBarry Smith */
4214b396bbSKris Buschelman typedef struct {
4375af56d4SHong Zhang   SuperMatrix       A,L,U,B,X;
4475af56d4SHong Zhang   superlu_options_t options;
45da7a1d00SHong Zhang   PetscInt          *perm_c; /* column permutation vector */
46da7a1d00SHong Zhang   PetscInt          *perm_r; /* row permutations from partial pivoting */
47da7a1d00SHong Zhang   PetscInt          *etree;
48da7a1d00SHong Zhang   PetscReal         *R, *C;
4975af56d4SHong Zhang   char              equed[1];
50da7a1d00SHong Zhang   PetscInt          lwork;
5175af56d4SHong Zhang   void              *work;
52da7a1d00SHong Zhang   PetscReal         rpg, rcond;
5375af56d4SHong Zhang   mem_usage_t       mem_usage;
5475af56d4SHong Zhang   MatStructure      flg;
555d8b2efaSHong Zhang   SuperLUStat_t     stat;
56cae5a23dSHong Zhang   Mat               A_dup;
57cae5a23dSHong Zhang   PetscScalar       *rhs_dup;
5814b396bbSKris Buschelman 
5914b396bbSKris Buschelman   /* Flag to clean up (non-global) SuperLU objects during Destroy */
60ace3abfcSBarry Smith   PetscBool  CleanUpSuperLU;
61f0c56d0fSKris Buschelman } Mat_SuperLU;
6214b396bbSKris Buschelman 
63e0e586b9SHong Zhang extern PetscErrorCode MatFactorInfo_SuperLU(Mat,PetscViewer);
640481f469SBarry Smith extern PetscErrorCode MatLUFactorNumeric_SuperLU(Mat,Mat,const MatFactorInfo *);
65e0e586b9SHong Zhang extern PetscErrorCode MatDestroy_SuperLU(Mat);
66e0e586b9SHong Zhang extern PetscErrorCode MatView_SuperLU(Mat,PetscViewer);
67e0e586b9SHong Zhang extern PetscErrorCode MatAssemblyEnd_SuperLU(Mat,MatAssemblyType);
68e0e586b9SHong Zhang extern PetscErrorCode MatSolve_SuperLU(Mat,Vec,Vec);
69e0b74bf9SHong Zhang extern PetscErrorCode MatMatSolve_SuperLU(Mat,Mat,Mat);
70e0e586b9SHong Zhang extern PetscErrorCode MatSolveTranspose_SuperLU(Mat,Vec,Vec);
710481f469SBarry Smith extern PetscErrorCode MatLUFactorSymbolic_SuperLU(Mat,Mat,IS,IS,const MatFactorInfo*);
72e0e586b9SHong Zhang extern PetscErrorCode MatDuplicate_SuperLU(Mat, MatDuplicateOption, Mat *);
73e0e586b9SHong Zhang 
745a46d3faSBarry Smith /*
755a46d3faSBarry Smith     Utility function
765a46d3faSBarry Smith */
775a46d3faSBarry Smith #undef __FUNCT__
785a46d3faSBarry Smith #define __FUNCT__ "MatFactorInfo_SuperLU"
795a46d3faSBarry Smith PetscErrorCode MatFactorInfo_SuperLU(Mat A,PetscViewer viewer)
805a46d3faSBarry Smith {
815a46d3faSBarry Smith   Mat_SuperLU       *lu= (Mat_SuperLU*)A->spptr;
825a46d3faSBarry Smith   PetscErrorCode    ierr;
835a46d3faSBarry Smith   superlu_options_t options;
845a46d3faSBarry Smith 
855a46d3faSBarry Smith   PetscFunctionBegin;
865a46d3faSBarry Smith   /* check if matrix is superlu_dist type */
875a46d3faSBarry Smith   if (A->ops->solve != MatSolve_SuperLU) PetscFunctionReturn(0);
885a46d3faSBarry Smith 
895a46d3faSBarry Smith   options = lu->options;
905a46d3faSBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"SuperLU run parameters:\n");CHKERRQ(ierr);
915a46d3faSBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"  Equil: %s\n",(options.Equil != NO) ? "YES": "NO");CHKERRQ(ierr);
925a46d3faSBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"  ColPerm: %D\n",options.ColPerm);CHKERRQ(ierr);
935a46d3faSBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"  IterRefine: %D\n",options.IterRefine);CHKERRQ(ierr);
945a46d3faSBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"  SymmetricMode: %s\n",(options.SymmetricMode != NO) ? "YES": "NO");CHKERRQ(ierr);
955a46d3faSBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"  DiagPivotThresh: %g\n",options.DiagPivotThresh);CHKERRQ(ierr);
965a46d3faSBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"  PivotGrowth: %s\n",(options.PivotGrowth != NO) ? "YES": "NO");CHKERRQ(ierr);
975a46d3faSBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"  ConditionNumber: %s\n",(options.ConditionNumber != NO) ? "YES": "NO");CHKERRQ(ierr);
985a46d3faSBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"  RowPerm: %D\n",options.RowPerm);CHKERRQ(ierr);
995a46d3faSBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"  ReplaceTinyPivot: %s\n",(options.ReplaceTinyPivot != NO) ? "YES": "NO");CHKERRQ(ierr);
1005a46d3faSBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"  PrintStat: %s\n",(options.PrintStat != NO) ? "YES": "NO");CHKERRQ(ierr);
1015a46d3faSBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"  lwork: %D\n",lu->lwork);CHKERRQ(ierr);
102d5f3da31SBarry Smith   if (A->factortype == MAT_FACTOR_ILU) {
103cffbb591SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  ILU_DropTol: %g\n",options.ILU_DropTol);CHKERRQ(ierr);
104cffbb591SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  ILU_FillTol: %g\n",options.ILU_FillTol);CHKERRQ(ierr);
105cffbb591SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  ILU_FillFactor: %g\n",options.ILU_FillFactor);CHKERRQ(ierr);
106cffbb591SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  ILU_DropRule: %D\n",options.ILU_DropRule);CHKERRQ(ierr);
107cffbb591SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  ILU_Norm: %D\n",options.ILU_Norm);CHKERRQ(ierr);
108cffbb591SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  ILU_MILU: %D\n",options.ILU_MILU);CHKERRQ(ierr);
109cffbb591SHong Zhang   }
1105a46d3faSBarry Smith   PetscFunctionReturn(0);
1115a46d3faSBarry Smith }
1125a46d3faSBarry Smith 
1135a46d3faSBarry Smith /*
1145a46d3faSBarry Smith     These are the methods provided to REPLACE the corresponding methods of the
1155a46d3faSBarry Smith    SeqAIJ matrix class. Other methods of SeqAIJ are not replaced
1165a46d3faSBarry Smith */
1175a46d3faSBarry Smith #undef __FUNCT__
1185a46d3faSBarry Smith #define __FUNCT__ "MatLUFactorNumeric_SuperLU"
1190481f469SBarry Smith PetscErrorCode MatLUFactorNumeric_SuperLU(Mat F,Mat A,const MatFactorInfo *info)
1205a46d3faSBarry Smith {
1211d5ca7f3SHong Zhang   Mat_SuperLU    *lu = (Mat_SuperLU*)F->spptr;
122cae5a23dSHong Zhang   Mat_SeqAIJ     *aa;
1235a46d3faSBarry Smith   PetscErrorCode ierr;
1245a46d3faSBarry Smith   PetscInt       sinfo;
1255a46d3faSBarry Smith   PetscReal      ferr, berr;
1265a46d3faSBarry Smith   NCformat       *Ustore;
1275a46d3faSBarry Smith   SCformat       *Lstore;
1285a46d3faSBarry Smith 
1295a46d3faSBarry Smith   PetscFunctionBegin;
1305a46d3faSBarry Smith   if (lu->flg == SAME_NONZERO_PATTERN) { /* successing numerical factorization */
1315a46d3faSBarry Smith     lu->options.Fact = SamePattern;
1325a46d3faSBarry Smith     /* Ref: ~SuperLU_3.0/EXAMPLE/dlinsolx2.c */
1335a46d3faSBarry Smith     Destroy_SuperMatrix_Store(&lu->A);
134cae5a23dSHong Zhang     if (lu->options.Equil) {
135cae5a23dSHong Zhang       ierr = MatCopy_SeqAIJ(A,lu->A_dup,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
136cae5a23dSHong Zhang     }
1375a46d3faSBarry Smith     if (lu->lwork >= 0) {
1385a46d3faSBarry Smith       Destroy_SuperNode_Matrix(&lu->L);
1395a46d3faSBarry Smith       Destroy_CompCol_Matrix(&lu->U);
1405a46d3faSBarry Smith       lu->options.Fact = SamePattern;
1415a46d3faSBarry Smith     }
1425a46d3faSBarry Smith   }
1435a46d3faSBarry Smith 
1445a46d3faSBarry Smith   /* Create the SuperMatrix for lu->A=A^T:
1455a46d3faSBarry Smith        Since SuperLU likes column-oriented matrices,we pass it the transpose,
1465a46d3faSBarry Smith        and then solve A^T X = B in MatSolve(). */
147cae5a23dSHong Zhang   if (lu->options.Equil) {
148cae5a23dSHong Zhang     aa = (Mat_SeqAIJ*)(lu->A_dup)->data;
149cae5a23dSHong Zhang   } else {
150cae5a23dSHong Zhang     aa = (Mat_SeqAIJ*)(A)->data;
151cae5a23dSHong Zhang   }
1525a46d3faSBarry Smith #if defined(PETSC_USE_COMPLEX)
1533cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
154b0043f70SBarry Smith   cCreate_CompCol_Matrix(&lu->A,A->cmap->n,A->rmap->n,aa->nz,(singlecomplex*)aa->a,aa->j,aa->i,SLU_NC,SLU_C,SLU_GE);
1553cb270beSHong Zhang #else
156b0043f70SBarry Smith   zCreate_CompCol_Matrix(&lu->A,A->cmap->n,A->rmap->n,aa->nz,(doublecomplex*)aa->a,aa->j,aa->i,SLU_NC,SLU_Z,SLU_GE);
1573cb270beSHong Zhang #endif
1583cb270beSHong Zhang #else
1593cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
160b0043f70SBarry Smith   sCreate_CompCol_Matrix(&lu->A,A->cmap->n,A->rmap->n,aa->nz,aa->a,aa->j,aa->i,SLU_NC,SLU_S,SLU_GE);
1615a46d3faSBarry Smith #else
162b0043f70SBarry Smith   dCreate_CompCol_Matrix(&lu->A,A->cmap->n,A->rmap->n,aa->nz,aa->a,aa->j,aa->i,SLU_NC,SLU_D,SLU_GE);
1635a46d3faSBarry Smith #endif
1643cb270beSHong Zhang #endif
1655a46d3faSBarry Smith 
1665a46d3faSBarry Smith   /* Numerical factorization */
1675a46d3faSBarry Smith   lu->B.ncol = 0;  /* Indicate not to solve the system */
168d5f3da31SBarry Smith   if (F->factortype == MAT_FACTOR_LU) {
1695a46d3faSBarry Smith #if defined(PETSC_USE_COMPLEX)
1703cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
1713cb270beSHong Zhang     cgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
1723cb270beSHong Zhang            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
1733cb270beSHong Zhang            &lu->mem_usage, &lu->stat, &sinfo);
1743cb270beSHong Zhang #else
1755a46d3faSBarry Smith     zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
1765a46d3faSBarry Smith            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
1775d8b2efaSHong Zhang            &lu->mem_usage, &lu->stat, &sinfo);
1783cb270beSHong Zhang #endif
1793cb270beSHong Zhang #else
1803cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
1813cb270beSHong Zhang     sgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
1823cb270beSHong Zhang            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
1833cb270beSHong Zhang            &lu->mem_usage, &lu->stat, &sinfo);
1845a46d3faSBarry Smith #else
1855a46d3faSBarry Smith     dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
1865a46d3faSBarry Smith            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
1875d8b2efaSHong Zhang            &lu->mem_usage, &lu->stat, &sinfo);
1885a46d3faSBarry Smith #endif
1893cb270beSHong Zhang #endif
190d5f3da31SBarry Smith   } else if (F->factortype == MAT_FACTOR_ILU) {
191cffbb591SHong Zhang     /* Compute the incomplete factorization, condition number and pivot growth */
192cffbb591SHong Zhang #if defined(PETSC_USE_COMPLEX)
1933cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
1943cb270beSHong Zhang     cgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r,lu->etree, lu->equed, lu->R, lu->C,
1953cb270beSHong Zhang            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond,
1963cb270beSHong Zhang            &lu->mem_usage, &lu->stat, &sinfo);
1973cb270beSHong Zhang #else
198cffbb591SHong Zhang     zgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r,lu->etree, lu->equed, lu->R, lu->C,
199cffbb591SHong Zhang            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond,
2005d8b2efaSHong Zhang            &lu->mem_usage, &lu->stat, &sinfo);
2013cb270beSHong Zhang #endif
2023cb270beSHong Zhang #else
2033cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
2043cb270beSHong Zhang     sgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
2053cb270beSHong Zhang           &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond,
2063cb270beSHong Zhang           &lu->mem_usage, &lu->stat, &sinfo);
207cffbb591SHong Zhang #else
208cffbb591SHong Zhang     dgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
209cffbb591SHong Zhang           &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond,
2105d8b2efaSHong Zhang           &lu->mem_usage, &lu->stat, &sinfo);
211cffbb591SHong Zhang #endif
2123cb270beSHong Zhang #endif
213*f23aa3ddSBarry Smith   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported");
2145a46d3faSBarry Smith   if (!sinfo || sinfo == lu->A.ncol+1) {
2155a46d3faSBarry Smith     if (lu->options.PivotGrowth)
2165a46d3faSBarry Smith       ierr = PetscPrintf(PETSC_COMM_SELF,"  Recip. pivot growth = %e\n", lu->rpg);
2175a46d3faSBarry Smith     if (lu->options.ConditionNumber)
2185a46d3faSBarry Smith       ierr = PetscPrintf(PETSC_COMM_SELF,"  Recip. condition number = %e\n", lu->rcond);
2195a46d3faSBarry Smith   } else if (sinfo > 0) {
2205a46d3faSBarry Smith     if (lu->lwork == -1) {
2215a46d3faSBarry Smith       ierr = PetscPrintf(PETSC_COMM_SELF,"  ** Estimated memory: %D bytes\n", sinfo - lu->A.ncol);
2229308c137SBarry Smith     } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot in row %D",sinfo);
223*f23aa3ddSBarry Smith   } else SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_LIB, "info = %D, the %D-th argument in gssvx() had an illegal value", sinfo,-sinfo);
2245a46d3faSBarry Smith 
2255a46d3faSBarry Smith   if (lu->options.PrintStat) {
2265a46d3faSBarry Smith     ierr = PetscPrintf(PETSC_COMM_SELF,"MatLUFactorNumeric_SuperLU():\n");
2275d8b2efaSHong Zhang     StatPrint(&lu->stat);
2285a46d3faSBarry Smith     Lstore = (SCformat *) lu->L.Store;
2295a46d3faSBarry Smith     Ustore = (NCformat *) lu->U.Store;
2305a46d3faSBarry Smith     ierr = PetscPrintf(PETSC_COMM_SELF,"  No of nonzeros in factor L = %D\n", Lstore->nnz);
2315a46d3faSBarry Smith     ierr = PetscPrintf(PETSC_COMM_SELF,"  No of nonzeros in factor U = %D\n", Ustore->nnz);
2325a46d3faSBarry Smith     ierr = PetscPrintf(PETSC_COMM_SELF,"  No of nonzeros in L+U = %D\n", Lstore->nnz + Ustore->nnz - lu->A.ncol);
2336da386baSHong Zhang     ierr = PetscPrintf(PETSC_COMM_SELF,"  L\\U MB %.3f\ttotal MB needed %.3f\n",
2346da386baSHong Zhang     lu->mem_usage.for_lu/1e6, lu->mem_usage.total_needed/1e6);
2355a46d3faSBarry Smith   }
2365a46d3faSBarry Smith 
2375a46d3faSBarry Smith   lu->flg = SAME_NONZERO_PATTERN;
2381d5ca7f3SHong Zhang   F->ops->solve          = MatSolve_SuperLU;
2391d5ca7f3SHong Zhang   F->ops->solvetranspose = MatSolveTranspose_SuperLU;
240e0b74bf9SHong Zhang   F->ops->matsolve       = MatMatSolve_SuperLU;
2415a46d3faSBarry Smith   PetscFunctionReturn(0);
2425a46d3faSBarry Smith }
2435a46d3faSBarry Smith 
24414b396bbSKris Buschelman #undef __FUNCT__
245f0c56d0fSKris Buschelman #define __FUNCT__ "MatDestroy_SuperLU"
246dfbe8321SBarry Smith PetscErrorCode MatDestroy_SuperLU(Mat A)
24714b396bbSKris Buschelman {
248dfbe8321SBarry Smith   PetscErrorCode ierr;
249f0c56d0fSKris Buschelman   Mat_SuperLU    *lu=(Mat_SuperLU*)A->spptr;
25014b396bbSKris Buschelman 
25114b396bbSKris Buschelman   PetscFunctionBegin;
252bf0cc555SLisandro Dalcin   if (lu && lu->CleanUpSuperLU) { /* Free the SuperLU datastructures */
25375af56d4SHong Zhang     Destroy_SuperMatrix_Store(&lu->A);
25475af56d4SHong Zhang     Destroy_SuperMatrix_Store(&lu->B);
25575af56d4SHong Zhang     Destroy_SuperMatrix_Store(&lu->X);
2565d8b2efaSHong Zhang     StatFree(&lu->stat);
2570e742b69SHong Zhang     if (lu->lwork >= 0) {
2580e742b69SHong Zhang       Destroy_SuperNode_Matrix(&lu->L);
2590e742b69SHong Zhang       Destroy_CompCol_Matrix(&lu->U);
2600e742b69SHong Zhang     }
2610e742b69SHong Zhang   }
262bf0cc555SLisandro Dalcin   if (lu) {
2639ce50919SHong Zhang     ierr = PetscFree(lu->etree);CHKERRQ(ierr);
2649ce50919SHong Zhang     ierr = PetscFree(lu->perm_r);CHKERRQ(ierr);
2659ce50919SHong Zhang     ierr = PetscFree(lu->perm_c);CHKERRQ(ierr);
2669ce50919SHong Zhang     ierr = PetscFree(lu->R);CHKERRQ(ierr);
2679ce50919SHong Zhang     ierr = PetscFree(lu->C);CHKERRQ(ierr);
268bf0cc555SLisandro Dalcin     ierr = PetscFree(lu->rhs_dup);CHKERRQ(ierr);
269bf0cc555SLisandro Dalcin     ierr = MatDestroy(&lu->A_dup);CHKERRQ(ierr);
270bf0cc555SLisandro Dalcin   }
271bf0cc555SLisandro Dalcin   ierr = PetscFree(A->spptr);CHKERRQ(ierr);
2720e742b69SHong Zhang 
273d954db57SHong Zhang   /* clear composed functions */
274d954db57SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatFactorGetSolverPackage_C","",PETSC_NULL);CHKERRQ(ierr);
275790a96dfSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatSuperluSetILUDropTol_C","",PETSC_NULL);CHKERRQ(ierr);
276d954db57SHong Zhang 
277b24902e0SBarry Smith   ierr = MatDestroy_SeqAIJ(A);CHKERRQ(ierr);
27814b396bbSKris Buschelman   PetscFunctionReturn(0);
27914b396bbSKris Buschelman }
28014b396bbSKris Buschelman 
28114b396bbSKris Buschelman #undef __FUNCT__
282f0c56d0fSKris Buschelman #define __FUNCT__ "MatView_SuperLU"
283dfbe8321SBarry Smith PetscErrorCode MatView_SuperLU(Mat A,PetscViewer viewer)
28414b396bbSKris Buschelman {
285dfbe8321SBarry Smith   PetscErrorCode    ierr;
286ace3abfcSBarry Smith   PetscBool         iascii;
28714b396bbSKris Buschelman   PetscViewerFormat format;
28814b396bbSKris Buschelman 
28914b396bbSKris Buschelman   PetscFunctionBegin;
290251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
29132077d6dSBarry Smith   if (iascii) {
29214b396bbSKris Buschelman     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
2932f59403fSHong Zhang     if (format == PETSC_VIEWER_ASCII_INFO) {
294f0c56d0fSKris Buschelman       ierr = MatFactorInfo_SuperLU(A,viewer);CHKERRQ(ierr);
29514b396bbSKris Buschelman     }
29614b396bbSKris Buschelman   }
29714b396bbSKris Buschelman   PetscFunctionReturn(0);
29814b396bbSKris Buschelman }
29914b396bbSKris Buschelman 
30014b396bbSKris Buschelman 
30114b396bbSKris Buschelman #undef __FUNCT__
302c29e884eSHong Zhang #define __FUNCT__ "MatSolve_SuperLU_Private"
303c29e884eSHong Zhang PetscErrorCode MatSolve_SuperLU_Private(Mat A,Vec b,Vec x)
30414b396bbSKris Buschelman {
305f0c56d0fSKris Buschelman   Mat_SuperLU    *lu = (Mat_SuperLU*)A->spptr;
30675af56d4SHong Zhang   PetscScalar    *barray,*xarray;
307dfbe8321SBarry Smith   PetscErrorCode ierr;
308cae5a23dSHong Zhang   PetscInt       info,i,n=x->map->n;
309da7a1d00SHong Zhang   PetscReal      ferr,berr;
31014b396bbSKris Buschelman 
31114b396bbSKris Buschelman   PetscFunctionBegin;
312937a6b0eSHong Zhang   if (lu->lwork == -1) {
313937a6b0eSHong Zhang     PetscFunctionReturn(0);
314937a6b0eSHong Zhang   }
315cae5a23dSHong Zhang 
31675af56d4SHong Zhang   lu->B.ncol = 1;   /* Set the number of right-hand side */
317cae5a23dSHong Zhang   if (lu->options.Equil && !lu->rhs_dup) {
318cae5a23dSHong Zhang     /* superlu overwrites b when Equil is used, thus create rhs_dup to keep user's b unchanged */
319cae5a23dSHong Zhang     ierr = PetscMalloc(n*sizeof(PetscScalar),&lu->rhs_dup);CHKERRQ(ierr);
320cae5a23dSHong Zhang   }
321cae5a23dSHong Zhang   if (lu->options.Equil) {
322cae5a23dSHong Zhang     /* Copy b into rsh_dup */
32375af56d4SHong Zhang     ierr = VecGetArray(b,&barray);CHKERRQ(ierr);
324cae5a23dSHong Zhang     ierr = PetscMemcpy(lu->rhs_dup,barray,n*sizeof(PetscScalar));CHKERRQ(ierr);
325cae5a23dSHong Zhang     ierr = VecRestoreArray(b,&barray);CHKERRQ(ierr);
326cae5a23dSHong Zhang     barray = lu->rhs_dup;
327cae5a23dSHong Zhang   } else {
328cae5a23dSHong Zhang     ierr = VecGetArray(b,&barray);CHKERRQ(ierr);
329cae5a23dSHong Zhang   }
33075af56d4SHong Zhang   ierr = VecGetArray(x,&xarray);CHKERRQ(ierr);
3315fe6150dSHong Zhang 
3325fe6150dSHong Zhang #if defined(PETSC_USE_COMPLEX)
3333cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
3343cb270beSHong Zhang   ((DNformat*)lu->B.Store)->nzval = (singlecomplex*)barray;
3353cb270beSHong Zhang   ((DNformat*)lu->X.Store)->nzval = (singlecomplex*)xarray;
3363cb270beSHong Zhang #else
3375fe6150dSHong Zhang   ((DNformat*)lu->B.Store)->nzval = (doublecomplex*)barray;
3385fe6150dSHong Zhang   ((DNformat*)lu->X.Store)->nzval = (doublecomplex*)xarray;
3393cb270beSHong Zhang #endif
3405fe6150dSHong Zhang #else
3415fe6150dSHong Zhang   ((DNformat*)lu->B.Store)->nzval = barray;
34275af56d4SHong Zhang   ((DNformat*)lu->X.Store)->nzval = xarray;
3435fe6150dSHong Zhang #endif
34475af56d4SHong Zhang 
34575af56d4SHong Zhang   lu->options.Fact = FACTORED; /* Indicate the factored form of A is supplied. */
346d5f3da31SBarry Smith   if (A->factortype == MAT_FACTOR_LU) {
3478914a3f7SHong Zhang #if defined(PETSC_USE_COMPLEX)
3483cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
3493cb270beSHong Zhang     cgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
3503cb270beSHong Zhang            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
3513cb270beSHong Zhang            &lu->mem_usage, &lu->stat, &info);
3523cb270beSHong Zhang #else
3538914a3f7SHong Zhang     zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
3548914a3f7SHong Zhang            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
3555d8b2efaSHong Zhang            &lu->mem_usage, &lu->stat, &info);
3563cb270beSHong Zhang #endif
3573cb270beSHong Zhang #else
3583cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
3593cb270beSHong Zhang     sgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
3603cb270beSHong Zhang            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
3613cb270beSHong Zhang            &lu->mem_usage, &lu->stat, &info);
3628914a3f7SHong Zhang #else
36375af56d4SHong Zhang     dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
36475af56d4SHong Zhang            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
3655d8b2efaSHong Zhang            &lu->mem_usage, &lu->stat, &info);
3668914a3f7SHong Zhang #endif
3673cb270beSHong Zhang #endif
368d5f3da31SBarry Smith   } else if (A->factortype == MAT_FACTOR_ILU) {
369cffbb591SHong Zhang #if defined(PETSC_USE_COMPLEX)
3703cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
3713cb270beSHong Zhang     cgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
3723cb270beSHong Zhang            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond,
3733cb270beSHong Zhang            &lu->mem_usage, &lu->stat, &info);
3743cb270beSHong Zhang #else
375cffbb591SHong Zhang     zgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
376cffbb591SHong Zhang            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond,
3775d8b2efaSHong Zhang            &lu->mem_usage, &lu->stat, &info);
3783cb270beSHong Zhang #endif
3793cb270beSHong Zhang #else
3803cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
3813cb270beSHong Zhang     sgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
3823cb270beSHong Zhang            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond,
3833cb270beSHong Zhang            &lu->mem_usage, &lu->stat, &info);
384cffbb591SHong Zhang #else
385cffbb591SHong Zhang     dgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
386cffbb591SHong Zhang            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond,
3875d8b2efaSHong Zhang            &lu->mem_usage, &lu->stat, &info);
388cffbb591SHong Zhang #endif
3893cb270beSHong Zhang #endif
390*f23aa3ddSBarry Smith   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported");
391cae5a23dSHong Zhang   if (!lu->options.Equil) {
39275af56d4SHong Zhang     ierr = VecRestoreArray(b,&barray);CHKERRQ(ierr);
393cae5a23dSHong Zhang   }
39475af56d4SHong Zhang   ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
39575af56d4SHong Zhang 
396958c9bccSBarry Smith   if (!info || info == lu->A.ncol+1) {
39775af56d4SHong Zhang     if (lu->options.IterRefine) {
3988914a3f7SHong Zhang       ierr = PetscPrintf(PETSC_COMM_SELF,"Iterative Refinement:\n");
3998914a3f7SHong Zhang       ierr = PetscPrintf(PETSC_COMM_SELF,"  %8s%8s%16s%16s\n", "rhs", "Steps", "FERR", "BERR");
40075af56d4SHong Zhang       for (i = 0; i < 1; ++i)
4015d8b2efaSHong Zhang         ierr = PetscPrintf(PETSC_COMM_SELF,"  %8d%8d%16e%16e\n", i+1, lu->stat.RefineSteps, ferr, berr);
40275af56d4SHong Zhang     }
4038914a3f7SHong Zhang   } else if (info > 0) {
4048914a3f7SHong Zhang     if (lu->lwork == -1) {
40577431f27SBarry Smith       ierr = PetscPrintf(PETSC_COMM_SELF,"  ** Estimated memory: %D bytes\n", info - lu->A.ncol);
4068914a3f7SHong Zhang     } else {
40777431f27SBarry Smith       ierr = PetscPrintf(PETSC_COMM_SELF,"  Warning: gssvx() returns info %D\n",info);
4088914a3f7SHong Zhang     }
409*f23aa3ddSBarry Smith   } else if (info < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_LIB, "info = %D, the %D-th argument in gssvx() had an illegal value", info,-info);
41014b396bbSKris Buschelman 
4118914a3f7SHong Zhang   if (lu->options.PrintStat) {
4128914a3f7SHong Zhang     ierr = PetscPrintf(PETSC_COMM_SELF,"MatSolve__SuperLU():\n");
4135d8b2efaSHong Zhang     StatPrint(&lu->stat);
4148914a3f7SHong Zhang   }
41575af56d4SHong Zhang   PetscFunctionReturn(0);
41675af56d4SHong Zhang }
41714b396bbSKris Buschelman 
41814b396bbSKris Buschelman #undef __FUNCT__
419c29e884eSHong Zhang #define __FUNCT__ "MatSolve_SuperLU"
420c29e884eSHong Zhang PetscErrorCode MatSolve_SuperLU(Mat A,Vec b,Vec x)
421c29e884eSHong Zhang {
422c29e884eSHong Zhang   Mat_SuperLU    *lu = (Mat_SuperLU*)A->spptr;
423c29e884eSHong Zhang   PetscErrorCode ierr;
424c29e884eSHong Zhang 
425c29e884eSHong Zhang   PetscFunctionBegin;
426c29e884eSHong Zhang   lu->options.Trans = TRANS;
427c29e884eSHong Zhang   ierr = MatSolve_SuperLU_Private(A,b,x);CHKERRQ(ierr);
428c29e884eSHong Zhang   PetscFunctionReturn(0);
429c29e884eSHong Zhang }
430c29e884eSHong Zhang 
431c29e884eSHong Zhang #undef __FUNCT__
432c7c1fe80SHong Zhang #define __FUNCT__ "MatSolveTranspose_SuperLU"
433c7c1fe80SHong Zhang PetscErrorCode MatSolveTranspose_SuperLU(Mat A,Vec b,Vec x)
434c7c1fe80SHong Zhang {
435c7c1fe80SHong Zhang   Mat_SuperLU    *lu = (Mat_SuperLU*)A->spptr;
436c7c1fe80SHong Zhang   PetscErrorCode ierr;
437c7c1fe80SHong Zhang 
438c7c1fe80SHong Zhang   PetscFunctionBegin;
439c7c1fe80SHong Zhang   lu->options.Trans = NOTRANS;
440c29e884eSHong Zhang   ierr = MatSolve_SuperLU_Private(A,b,x);CHKERRQ(ierr);
441c7c1fe80SHong Zhang   PetscFunctionReturn(0);
442c7c1fe80SHong Zhang }
443c7c1fe80SHong Zhang 
444e0b74bf9SHong Zhang #undef __FUNCT__
445e0b74bf9SHong Zhang #define __FUNCT__ "MatMatSolve_SuperLU"
446e0b74bf9SHong Zhang PetscErrorCode MatMatSolve_SuperLU(Mat A,Mat B,Mat X)
447e0b74bf9SHong Zhang {
448e0b74bf9SHong Zhang   Mat_SuperLU    *lu = (Mat_SuperLU*)A->spptr;
449cd723cd1SBarry Smith   PetscBool      flg;
450cd723cd1SBarry Smith   PetscErrorCode ierr;
451e0b74bf9SHong Zhang 
452e0b74bf9SHong Zhang   PetscFunctionBegin;
453251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompareAny((PetscObject)B,&flg,MATSEQDENSE,MATMPIDENSE,PETSC_NULL);CHKERRQ(ierr);
454cd723cd1SBarry Smith   if (!flg) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix");
455251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,PETSC_NULL);CHKERRQ(ierr);
456cd723cd1SBarry Smith   if (!flg) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix");  lu->options.Trans = TRANS;
457e0b74bf9SHong Zhang   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatMatSolve_SuperLU() is not implemented yet");
458e0b74bf9SHong Zhang   PetscFunctionReturn(0);
459e0b74bf9SHong Zhang }
460e0b74bf9SHong Zhang 
46114b396bbSKris Buschelman /*
46214b396bbSKris Buschelman    Note the r permutation is ignored
46314b396bbSKris Buschelman */
46414b396bbSKris Buschelman #undef __FUNCT__
465f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_SuperLU"
4660481f469SBarry Smith PetscErrorCode MatLUFactorSymbolic_SuperLU(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
46714b396bbSKris Buschelman {
4681d5ca7f3SHong Zhang   Mat_SuperLU    *lu = (Mat_SuperLU*)(F->spptr);
469b24902e0SBarry Smith 
470b24902e0SBarry Smith   PetscFunctionBegin;
471b24902e0SBarry Smith   lu->flg                 = DIFFERENT_NONZERO_PATTERN;
472b24902e0SBarry Smith   lu->CleanUpSuperLU      = PETSC_TRUE;
4731d5ca7f3SHong Zhang   F->ops->lufactornumeric = MatLUFactorNumeric_SuperLU;
474b24902e0SBarry Smith   PetscFunctionReturn(0);
475b24902e0SBarry Smith }
476b24902e0SBarry Smith 
47735bd34faSBarry Smith EXTERN_C_BEGIN
47835bd34faSBarry Smith #undef __FUNCT__
4795ccb76cbSHong Zhang #define __FUNCT__ "MatSuperluSetILUDropTol_SuperLU"
4805ccb76cbSHong Zhang PetscErrorCode MatSuperluSetILUDropTol_SuperLU(Mat F,PetscReal dtol)
4815ccb76cbSHong Zhang {
4825ccb76cbSHong Zhang   Mat_SuperLU *lu= (Mat_SuperLU*)F->spptr;
4835ccb76cbSHong Zhang 
4845ccb76cbSHong Zhang   PetscFunctionBegin;
4855ccb76cbSHong Zhang   lu->options.ILU_DropTol = dtol;
4865ccb76cbSHong Zhang   PetscFunctionReturn(0);
4875ccb76cbSHong Zhang }
4885ccb76cbSHong Zhang EXTERN_C_END
4895ccb76cbSHong Zhang 
4905ccb76cbSHong Zhang #undef __FUNCT__
4915ccb76cbSHong Zhang #define __FUNCT__ "MatSuperluSetILUDropTol"
4925ccb76cbSHong Zhang /*@
4935ccb76cbSHong Zhang   MatSuperluSetILUDropTol - Set SuperLU ILU drop tolerance
4945ccb76cbSHong Zhang    Logically Collective on Mat
4955ccb76cbSHong Zhang 
4965ccb76cbSHong Zhang    Input Parameters:
4975ccb76cbSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-SuperLU interface
4985ccb76cbSHong Zhang -  dtol - drop tolerance
4995ccb76cbSHong Zhang 
5005ccb76cbSHong Zhang   Options Database:
5015ccb76cbSHong Zhang .   -mat_superlu_ilu_droptol <dtol>
5025ccb76cbSHong Zhang 
5035ccb76cbSHong Zhang    Level: beginner
5045ccb76cbSHong Zhang 
5055ccb76cbSHong Zhang    References: SuperLU Users' Guide
5065ccb76cbSHong Zhang 
5075ccb76cbSHong Zhang .seealso: MatGetFactor()
5085ccb76cbSHong Zhang @*/
5095ccb76cbSHong Zhang PetscErrorCode MatSuperluSetILUDropTol(Mat F,PetscReal dtol)
5105ccb76cbSHong Zhang {
5115ccb76cbSHong Zhang   PetscErrorCode ierr;
5125ccb76cbSHong Zhang 
5135ccb76cbSHong Zhang   PetscFunctionBegin;
5145ccb76cbSHong Zhang   PetscValidHeaderSpecific(F,MAT_CLASSID,1);
5155ccb76cbSHong Zhang   PetscValidLogicalCollectiveInt(F,dtol,2);
5165ccb76cbSHong Zhang   ierr = PetscTryMethod(F,"MatSuperluSetILUDropTol_C",(Mat,PetscReal),(F,dtol));CHKERRQ(ierr);
5175ccb76cbSHong Zhang   PetscFunctionReturn(0);
5185ccb76cbSHong Zhang }
5195ccb76cbSHong Zhang 
5205ccb76cbSHong Zhang EXTERN_C_BEGIN
5215ccb76cbSHong Zhang #undef __FUNCT__
52235bd34faSBarry Smith #define __FUNCT__ "MatFactorGetSolverPackage_seqaij_superlu"
52335bd34faSBarry Smith PetscErrorCode MatFactorGetSolverPackage_seqaij_superlu(Mat A,const MatSolverPackage *type)
52435bd34faSBarry Smith {
52535bd34faSBarry Smith   PetscFunctionBegin;
5262692d6eeSBarry Smith   *type = MATSOLVERSUPERLU;
52735bd34faSBarry Smith   PetscFunctionReturn(0);
52835bd34faSBarry Smith }
52935bd34faSBarry Smith EXTERN_C_END
53035bd34faSBarry Smith 
531b24902e0SBarry Smith 
532b24902e0SBarry Smith /*MC
533ba992d64SSatish Balay   MATSOLVERSUPERLU = "superlu" - A solver package providing solvers LU and ILU for sequential matrices
534b24902e0SBarry Smith   via the external package SuperLU.
535b24902e0SBarry Smith 
536e2e64c6bSBarry Smith   Use ./configure --download-superlu to have PETSc installed with SuperLU
537b24902e0SBarry Smith 
538b24902e0SBarry Smith   Options Database Keys:
539e08999f5SMatthew G Knepley + -mat_superlu_equil <FALSE>            - Equil (None)
540e08999f5SMatthew G Knepley . -mat_superlu_colperm <COLAMD>         - (choose one of) NATURAL MMD_ATA MMD_AT_PLUS_A COLAMD
541e08999f5SMatthew G Knepley . -mat_superlu_iterrefine <NOREFINE>    - (choose one of) NOREFINE SINGLE DOUBLE EXTRA
542e08999f5SMatthew G Knepley . -mat_superlu_symmetricmode: <FALSE>   - SymmetricMode (None)
543e08999f5SMatthew G Knepley . -mat_superlu_diagpivotthresh <1>      - DiagPivotThresh (None)
544e08999f5SMatthew G Knepley . -mat_superlu_pivotgrowth <FALSE>      - PivotGrowth (None)
545e08999f5SMatthew G Knepley . -mat_superlu_conditionnumber <FALSE>  - ConditionNumber (None)
546e08999f5SMatthew G Knepley . -mat_superlu_rowperm <NOROWPERM>      - (choose one of) NOROWPERM LargeDiag
547e08999f5SMatthew G Knepley . -mat_superlu_replacetinypivot <FALSE> - ReplaceTinyPivot (None)
548e08999f5SMatthew G Knepley . -mat_superlu_printstat <FALSE>        - PrintStat (None)
549e08999f5SMatthew G Knepley . -mat_superlu_lwork <0>                - size of work array in bytes used by factorization (None)
550e08999f5SMatthew G Knepley . -mat_superlu_ilu_droptol <0>          - ILU_DropTol (None)
551e08999f5SMatthew G Knepley . -mat_superlu_ilu_filltol <0>          - ILU_FillTol (None)
552e08999f5SMatthew G Knepley . -mat_superlu_ilu_fillfactor <0>       - ILU_FillFactor (None)
553e08999f5SMatthew G Knepley . -mat_superlu_ilu_droprull <0>         - ILU_DropRule (None)
554e08999f5SMatthew G Knepley . -mat_superlu_ilu_norm <0>             - ILU_Norm (None)
555e08999f5SMatthew G Knepley - -mat_superlu_ilu_milu <0>             - ILU_MILU (None)
556b24902e0SBarry Smith 
5572692d6eeSBarry Smith    Notes: Do not confuse this with MATSOLVERSUPERLU_DIST which is for parallel sparse solves
5585c9eb25fSBarry Smith 
559b24902e0SBarry Smith    Level: beginner
560b24902e0SBarry Smith 
561d45987f3SHong Zhang .seealso: PCLU, PCILU, MATSOLVERSUPERLU_DIST, MATSOLVERMUMPS, PCFactorSetMatSolverPackage(), MatSolverPackage
562b24902e0SBarry Smith M*/
563b24902e0SBarry Smith 
5645c9eb25fSBarry Smith EXTERN_C_BEGIN
565b24902e0SBarry Smith #undef __FUNCT__
566b24902e0SBarry Smith #define __FUNCT__ "MatGetFactor_seqaij_superlu"
5675c9eb25fSBarry Smith PetscErrorCode MatGetFactor_seqaij_superlu(Mat A,MatFactorType ftype,Mat *F)
568b24902e0SBarry Smith {
56914b396bbSKris Buschelman   Mat            B;
570f0c56d0fSKris Buschelman   Mat_SuperLU    *lu;
5716849ba73SBarry Smith   PetscErrorCode ierr;
5725d8b2efaSHong Zhang   PetscInt       indx,m=A->rmap->n,n=A->cmap->n;
573ace3abfcSBarry Smith   PetscBool      flg;
5743cb270beSHong Zhang   PetscReal      real_input;
5755ca28756SHong Zhang   const char     *colperm[]={"NATURAL","MMD_ATA","MMD_AT_PLUS_A","COLAMD"}; /* MY_PERMC - not supported by the petsc interface yet */
5765ca28756SHong Zhang   const char     *iterrefine[]={"NOREFINE", "SINGLE", "DOUBLE", "EXTRA"};
5775ca28756SHong Zhang   const char     *rowperm[]={"NOROWPERM", "LargeDiag"}; /* MY_PERMC - not supported by the petsc interface yet */
57814b396bbSKris Buschelman 
57914b396bbSKris Buschelman   PetscFunctionBegin;
5807adad957SLisandro Dalcin   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
581d0f46423SBarry Smith   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
5827adad957SLisandro Dalcin   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
583899d7b4fSKris Buschelman   ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
584f0c56d0fSKris Buschelman 
585cffbb591SHong Zhang   if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU) {
586b24902e0SBarry Smith     B->ops->lufactorsymbolic  = MatLUFactorSymbolic_SuperLU;
587cffbb591SHong Zhang     B->ops->ilufactorsymbolic = MatLUFactorSymbolic_SuperLU;
588b3a44c85SBarry Smith   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported");
589cffbb591SHong Zhang 
590b24902e0SBarry Smith   B->ops->destroy          = MatDestroy_SuperLU;
5913519fcdcSHong Zhang   B->ops->view             = MatView_SuperLU;
592d5f3da31SBarry Smith   B->factortype            = ftype;
59394b7f48cSBarry Smith   B->assembled             = PETSC_TRUE;  /* required by -ksp_view */
5945c9eb25fSBarry Smith   B->preallocated          = PETSC_TRUE;
59514b396bbSKris Buschelman 
596b24902e0SBarry Smith   ierr = PetscNewLog(B,Mat_SuperLU,&lu);CHKERRQ(ierr);
597cae5a23dSHong Zhang 
598cffbb591SHong Zhang   if (ftype == MAT_FACTOR_LU) {
5999ce50919SHong Zhang     set_default_options(&lu->options);
6003d6581fbSHong Zhang     /* Comments from SuperLU_4.0/SRC/dgssvx.c:
6013d6581fbSHong Zhang       "Whether or not the system will be equilibrated depends on the
6023d6581fbSHong Zhang        scaling of the matrix A, but if equilibration is used, A is
6033d6581fbSHong Zhang        overwritten by diag(R)*A*diag(C) and B by diag(R)*B
6043d6581fbSHong Zhang        (if options->Trans=NOTRANS) or diag(C)*B (if options->Trans = TRANS or CONJ)."
6053d6581fbSHong Zhang      We set 'options.Equil = NO' as default because additional space is needed for it.
6063d6581fbSHong Zhang     */
6073d6581fbSHong Zhang     lu->options.Equil = NO;
608cffbb591SHong Zhang   } else if (ftype == MAT_FACTOR_ILU) {
6090c74a584SJed Brown     /* Set the default input options of ilu: */
610cffbb591SHong Zhang     ilu_set_default_options(&lu->options);
611cffbb591SHong Zhang   }
6129ce50919SHong Zhang   lu->options.PrintStat = NO;
6131a2e2f44SHong Zhang 
6145d8b2efaSHong Zhang   /* Initialize the statistics variables. */
6155d8b2efaSHong Zhang   StatInit(&lu->stat);
6168914a3f7SHong Zhang   lu->lwork = 0;   /* allocate space internally by system malloc */
6179ce50919SHong Zhang 
6187adad957SLisandro Dalcin   ierr = PetscOptionsBegin(((PetscObject)A)->comm,((PetscObject)A)->prefix,"SuperLU Options","Mat");CHKERRQ(ierr);
619acfcf0e5SJed Brown     ierr = PetscOptionsBool("-mat_superlu_equil","Equil","None",(PetscBool)lu->options.Equil,(PetscBool*)&lu->options.Equil,0);CHKERRQ(ierr);
6208914a3f7SHong Zhang     ierr = PetscOptionsEList("-mat_superlu_colperm","ColPerm","None",colperm,4,colperm[3],&indx,&flg);CHKERRQ(ierr);
6219ce50919SHong Zhang     if (flg) {lu->options.ColPerm = (colperm_t)indx;}
6228914a3f7SHong Zhang     ierr = PetscOptionsEList("-mat_superlu_iterrefine","IterRefine","None",iterrefine,4,iterrefine[0],&indx,&flg);CHKERRQ(ierr);
6239ce50919SHong Zhang     if (flg) { lu->options.IterRefine = (IterRefine_t)indx;}
624acfcf0e5SJed Brown     ierr = PetscOptionsBool("-mat_superlu_symmetricmode","SymmetricMode","None",(PetscBool)lu->options.SymmetricMode,&flg,0);CHKERRQ(ierr);
6259ce50919SHong Zhang     if (flg) lu->options.SymmetricMode = YES;
6263cb270beSHong Zhang     ierr = PetscOptionsReal("-mat_superlu_diagpivotthresh","DiagPivotThresh","None",lu->options.DiagPivotThresh,&real_input,&flg);CHKERRQ(ierr);
6273cb270beSHong Zhang     if (flg) lu->options.DiagPivotThresh = (double) real_input;
628acfcf0e5SJed Brown     ierr = PetscOptionsBool("-mat_superlu_pivotgrowth","PivotGrowth","None",(PetscBool)lu->options.PivotGrowth,&flg,0);CHKERRQ(ierr);
6299ce50919SHong Zhang     if (flg) lu->options.PivotGrowth = YES;
630acfcf0e5SJed Brown     ierr = PetscOptionsBool("-mat_superlu_conditionnumber","ConditionNumber","None",(PetscBool)lu->options.ConditionNumber,&flg,0);CHKERRQ(ierr);
6319ce50919SHong Zhang     if (flg) lu->options.ConditionNumber = YES;
632d7ebd59bSHong Zhang     ierr = PetscOptionsEList("-mat_superlu_rowperm","rowperm","None",rowperm,2,rowperm[lu->options.RowPerm],&indx,&flg);CHKERRQ(ierr);
6339ce50919SHong Zhang     if (flg) {lu->options.RowPerm = (rowperm_t)indx;}
634acfcf0e5SJed Brown     ierr = PetscOptionsBool("-mat_superlu_replacetinypivot","ReplaceTinyPivot","None",(PetscBool)lu->options.ReplaceTinyPivot,&flg,0);CHKERRQ(ierr);
6359ce50919SHong Zhang     if (flg) lu->options.ReplaceTinyPivot = YES;
636acfcf0e5SJed Brown     ierr = PetscOptionsBool("-mat_superlu_printstat","PrintStat","None",(PetscBool)lu->options.PrintStat,&flg,0);CHKERRQ(ierr);
6379ce50919SHong Zhang     if (flg) lu->options.PrintStat = YES;
6388914a3f7SHong Zhang     ierr = PetscOptionsInt("-mat_superlu_lwork","size of work array in bytes used by factorization","None",lu->lwork,&lu->lwork,PETSC_NULL);CHKERRQ(ierr);
6395fe6150dSHong Zhang     if (lu->lwork > 0) {
6405fe6150dSHong Zhang       ierr = PetscMalloc(lu->lwork,&lu->work);CHKERRQ(ierr);
6415fe6150dSHong Zhang     } else if (lu->lwork != 0 && lu->lwork != -1) {
64277431f27SBarry Smith       ierr = PetscPrintf(PETSC_COMM_SELF,"   Warning: lwork %D is not supported by SUPERLU. The default lwork=0 is used.\n",lu->lwork);
6438914a3f7SHong Zhang       lu->lwork = 0;
6448914a3f7SHong Zhang     }
645cffbb591SHong Zhang     /* ilu options */
6463cb270beSHong Zhang     ierr = PetscOptionsReal("-mat_superlu_ilu_droptol","ILU_DropTol","None",lu->options.ILU_DropTol,&real_input,&flg);CHKERRQ(ierr);
6473cb270beSHong Zhang     if (flg) lu->options.ILU_DropTol = (double) real_input;
6483cb270beSHong Zhang     ierr = PetscOptionsReal("-mat_superlu_ilu_filltol","ILU_FillTol","None",lu->options.ILU_FillTol,&real_input,&flg);CHKERRQ(ierr);
6493cb270beSHong Zhang     if (flg) lu->options.ILU_FillTol = (double) real_input;
6503cb270beSHong Zhang     ierr = PetscOptionsReal("-mat_superlu_ilu_fillfactor","ILU_FillFactor","None",lu->options.ILU_FillFactor,&real_input,&flg);CHKERRQ(ierr);
6513cb270beSHong Zhang     if (flg) lu->options.ILU_FillFactor = (double) real_input;
652cffbb591SHong Zhang     ierr = PetscOptionsInt("-mat_superlu_ilu_droprull","ILU_DropRule","None",lu->options.ILU_DropRule,&lu->options.ILU_DropRule,PETSC_NULL);CHKERRQ(ierr);
653cffbb591SHong Zhang     ierr = PetscOptionsInt("-mat_superlu_ilu_norm","ILU_Norm","None",lu->options.ILU_Norm,&indx,&flg);CHKERRQ(ierr);
654cffbb591SHong Zhang     if (flg) {
655cffbb591SHong Zhang       lu->options.ILU_Norm = (norm_t)indx;
656cffbb591SHong Zhang     }
657cffbb591SHong Zhang     ierr = PetscOptionsInt("-mat_superlu_ilu_milu","ILU_MILU","None",lu->options.ILU_MILU,&indx,&flg);CHKERRQ(ierr);
658cffbb591SHong Zhang     if (flg) {
659cffbb591SHong Zhang       lu->options.ILU_MILU = (milu_t)indx;
660cffbb591SHong Zhang     }
6619ce50919SHong Zhang   PetscOptionsEnd();
66238abfdaeSHong Zhang   if (lu->options.Equil == YES) {
66338abfdaeSHong Zhang     /* superlu overwrites input matrix and rhs when Equil is used, thus create A_dup to keep user's A unchanged */
66438abfdaeSHong Zhang     ierr = MatDuplicate_SeqAIJ(A,MAT_COPY_VALUES,&lu->A_dup);CHKERRQ(ierr);
66538abfdaeSHong Zhang   }
6669ce50919SHong Zhang 
6675d8b2efaSHong Zhang   /* Allocate spaces (notice sizes are for the transpose) */
6685d8b2efaSHong Zhang   ierr = PetscMalloc(m*sizeof(PetscInt),&lu->etree);CHKERRQ(ierr);
6695d8b2efaSHong Zhang   ierr = PetscMalloc(n*sizeof(PetscInt),&lu->perm_r);CHKERRQ(ierr);
6705d8b2efaSHong Zhang   ierr = PetscMalloc(m*sizeof(PetscInt),&lu->perm_c);CHKERRQ(ierr);
6715d8b2efaSHong Zhang   ierr = PetscMalloc(n*sizeof(PetscScalar),&lu->R);CHKERRQ(ierr);
6725d8b2efaSHong Zhang   ierr = PetscMalloc(m*sizeof(PetscScalar),&lu->C);CHKERRQ(ierr);
6735d8b2efaSHong Zhang 
6745d8b2efaSHong Zhang   /* create rhs and solution x without allocate space for .Store */
6755d8b2efaSHong Zhang #if defined(PETSC_USE_COMPLEX)
6763cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
6773cb270beSHong Zhang   cCreate_Dense_Matrix(&lu->B, m, 1, PETSC_NULL, m, SLU_DN, SLU_C, SLU_GE);
6783cb270beSHong Zhang   cCreate_Dense_Matrix(&lu->X, m, 1, PETSC_NULL, m, SLU_DN, SLU_C, SLU_GE);
6793cb270beSHong Zhang #else
6805d8b2efaSHong Zhang   zCreate_Dense_Matrix(&lu->B, m, 1, PETSC_NULL, m, SLU_DN, SLU_Z, SLU_GE);
6815d8b2efaSHong Zhang   zCreate_Dense_Matrix(&lu->X, m, 1, PETSC_NULL, m, SLU_DN, SLU_Z, SLU_GE);
6823cb270beSHong Zhang #endif
6833cb270beSHong Zhang #else
6843cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
6853cb270beSHong Zhang   sCreate_Dense_Matrix(&lu->B, m, 1, PETSC_NULL, m, SLU_DN, SLU_S, SLU_GE);
6863cb270beSHong Zhang   sCreate_Dense_Matrix(&lu->X, m, 1, PETSC_NULL, m, SLU_DN, SLU_S, SLU_GE);
6875d8b2efaSHong Zhang #else
6885d8b2efaSHong Zhang   dCreate_Dense_Matrix(&lu->B, m, 1, PETSC_NULL, m, SLU_DN, SLU_D, SLU_GE);
6895d8b2efaSHong Zhang   dCreate_Dense_Matrix(&lu->X, m, 1, PETSC_NULL, m, SLU_DN, SLU_D, SLU_GE);
6905d8b2efaSHong Zhang #endif
6913cb270beSHong Zhang #endif
6925d8b2efaSHong Zhang 
693519f805aSKarl Rupp #if defined(SUPERLU2)
6945c9eb25fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatCreateNull","MatCreateNull_SuperLU",(void(*)(void))MatCreateNull_SuperLU);CHKERRQ(ierr);
69575af56d4SHong Zhang #endif
69635bd34faSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_seqaij_superlu",MatFactorGetSolverPackage_seqaij_superlu);CHKERRQ(ierr);
6975ccb76cbSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSuperluSetILUDropTol_C","MatSuperluSetILUDropTol_SuperLU",MatSuperluSetILUDropTol_SuperLU);CHKERRQ(ierr);
6985c9eb25fSBarry Smith   B->spptr = lu;
699899d7b4fSKris Buschelman   *F = B;
70014b396bbSKris Buschelman   PetscFunctionReturn(0);
70114b396bbSKris Buschelman }
7025c9eb25fSBarry Smith EXTERN_C_END
703d954db57SHong Zhang 
704