xref: /petsc/src/mat/impls/aij/seq/superlu/superlu.c (revision 519f805a543c2a7f195bcba21173bb2cdfadb956)
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
213cffbb591SHong Zhang   } else {
214e32f2f54SBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported");
215cffbb591SHong Zhang   }
2165a46d3faSBarry Smith   if ( !sinfo || sinfo == lu->A.ncol+1 ) {
2175a46d3faSBarry Smith     if ( lu->options.PivotGrowth )
2185a46d3faSBarry Smith       ierr = PetscPrintf(PETSC_COMM_SELF,"  Recip. pivot growth = %e\n", lu->rpg);
2195a46d3faSBarry Smith     if ( lu->options.ConditionNumber )
2205a46d3faSBarry Smith       ierr = PetscPrintf(PETSC_COMM_SELF,"  Recip. condition number = %e\n", lu->rcond);
2215a46d3faSBarry Smith   } else if ( sinfo > 0 ) {
2225a46d3faSBarry Smith     if ( lu->lwork == -1 ) {
2235a46d3faSBarry Smith       ierr = PetscPrintf(PETSC_COMM_SELF,"  ** Estimated memory: %D bytes\n", sinfo - lu->A.ncol);
2249308c137SBarry Smith     } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot in row %D",sinfo);
2255a46d3faSBarry Smith   } else { /* sinfo < 0 */
226e32f2f54SBarry Smith     SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_LIB, "info = %D, the %D-th argument in gssvx() had an illegal value", sinfo,-sinfo);
2275a46d3faSBarry Smith   }
2285a46d3faSBarry Smith 
2295a46d3faSBarry Smith   if ( lu->options.PrintStat ) {
2305a46d3faSBarry Smith     ierr = PetscPrintf(PETSC_COMM_SELF,"MatLUFactorNumeric_SuperLU():\n");
2315d8b2efaSHong Zhang     StatPrint(&lu->stat);
2325a46d3faSBarry Smith     Lstore = (SCformat *) lu->L.Store;
2335a46d3faSBarry Smith     Ustore = (NCformat *) lu->U.Store;
2345a46d3faSBarry Smith     ierr = PetscPrintf(PETSC_COMM_SELF,"  No of nonzeros in factor L = %D\n", Lstore->nnz);
2355a46d3faSBarry Smith     ierr = PetscPrintf(PETSC_COMM_SELF,"  No of nonzeros in factor U = %D\n", Ustore->nnz);
2365a46d3faSBarry Smith     ierr = PetscPrintf(PETSC_COMM_SELF,"  No of nonzeros in L+U = %D\n", Lstore->nnz + Ustore->nnz - lu->A.ncol);
2376da386baSHong Zhang     ierr = PetscPrintf(PETSC_COMM_SELF,"  L\\U MB %.3f\ttotal MB needed %.3f\n",
2386da386baSHong Zhang     lu->mem_usage.for_lu/1e6, lu->mem_usage.total_needed/1e6);
2395a46d3faSBarry Smith   }
2405a46d3faSBarry Smith 
2415a46d3faSBarry Smith   lu->flg = SAME_NONZERO_PATTERN;
2421d5ca7f3SHong Zhang   F->ops->solve          = MatSolve_SuperLU;
2431d5ca7f3SHong Zhang   F->ops->solvetranspose = MatSolveTranspose_SuperLU;
244e0b74bf9SHong Zhang   F->ops->matsolve       = MatMatSolve_SuperLU;
2455a46d3faSBarry Smith   PetscFunctionReturn(0);
2465a46d3faSBarry Smith }
2475a46d3faSBarry Smith 
24814b396bbSKris Buschelman #undef __FUNCT__
249f0c56d0fSKris Buschelman #define __FUNCT__ "MatDestroy_SuperLU"
250dfbe8321SBarry Smith PetscErrorCode MatDestroy_SuperLU(Mat A)
25114b396bbSKris Buschelman {
252dfbe8321SBarry Smith   PetscErrorCode ierr;
253f0c56d0fSKris Buschelman   Mat_SuperLU    *lu=(Mat_SuperLU*)A->spptr;
25414b396bbSKris Buschelman 
25514b396bbSKris Buschelman   PetscFunctionBegin;
256bf0cc555SLisandro Dalcin   if (lu && lu->CleanUpSuperLU) { /* Free the SuperLU datastructures */
25775af56d4SHong Zhang     Destroy_SuperMatrix_Store(&lu->A);
25875af56d4SHong Zhang     Destroy_SuperMatrix_Store(&lu->B);
25975af56d4SHong Zhang     Destroy_SuperMatrix_Store(&lu->X);
2605d8b2efaSHong Zhang     StatFree(&lu->stat);
2610e742b69SHong Zhang     if (lu->lwork >= 0) {
2620e742b69SHong Zhang       Destroy_SuperNode_Matrix(&lu->L);
2630e742b69SHong Zhang       Destroy_CompCol_Matrix(&lu->U);
2640e742b69SHong Zhang     }
2650e742b69SHong Zhang   }
266bf0cc555SLisandro Dalcin   if (lu) {
2679ce50919SHong Zhang     ierr = PetscFree(lu->etree);CHKERRQ(ierr);
2689ce50919SHong Zhang     ierr = PetscFree(lu->perm_r);CHKERRQ(ierr);
2699ce50919SHong Zhang     ierr = PetscFree(lu->perm_c);CHKERRQ(ierr);
2709ce50919SHong Zhang     ierr = PetscFree(lu->R);CHKERRQ(ierr);
2719ce50919SHong Zhang     ierr = PetscFree(lu->C);CHKERRQ(ierr);
272bf0cc555SLisandro Dalcin     ierr = PetscFree(lu->rhs_dup);CHKERRQ(ierr);
273bf0cc555SLisandro Dalcin     ierr = MatDestroy(&lu->A_dup);CHKERRQ(ierr);
274bf0cc555SLisandro Dalcin   }
275bf0cc555SLisandro Dalcin   ierr = PetscFree(A->spptr);CHKERRQ(ierr);
2760e742b69SHong Zhang 
277d954db57SHong Zhang   /* clear composed functions */
278d954db57SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatFactorGetSolverPackage_C","",PETSC_NULL);CHKERRQ(ierr);
279790a96dfSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatSuperluSetILUDropTol_C","",PETSC_NULL);CHKERRQ(ierr);
280d954db57SHong Zhang 
281b24902e0SBarry Smith   ierr = MatDestroy_SeqAIJ(A);CHKERRQ(ierr);
28214b396bbSKris Buschelman   PetscFunctionReturn(0);
28314b396bbSKris Buschelman }
28414b396bbSKris Buschelman 
28514b396bbSKris Buschelman #undef __FUNCT__
286f0c56d0fSKris Buschelman #define __FUNCT__ "MatView_SuperLU"
287dfbe8321SBarry Smith PetscErrorCode MatView_SuperLU(Mat A,PetscViewer viewer)
28814b396bbSKris Buschelman {
289dfbe8321SBarry Smith   PetscErrorCode    ierr;
290ace3abfcSBarry Smith   PetscBool         iascii;
29114b396bbSKris Buschelman   PetscViewerFormat format;
29214b396bbSKris Buschelman 
29314b396bbSKris Buschelman   PetscFunctionBegin;
294251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
29532077d6dSBarry Smith   if (iascii) {
29614b396bbSKris Buschelman     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
2972f59403fSHong Zhang     if (format == PETSC_VIEWER_ASCII_INFO) {
298f0c56d0fSKris Buschelman       ierr = MatFactorInfo_SuperLU(A,viewer);CHKERRQ(ierr);
29914b396bbSKris Buschelman     }
30014b396bbSKris Buschelman   }
30114b396bbSKris Buschelman   PetscFunctionReturn(0);
30214b396bbSKris Buschelman }
30314b396bbSKris Buschelman 
30414b396bbSKris Buschelman 
30514b396bbSKris Buschelman #undef __FUNCT__
306c29e884eSHong Zhang #define __FUNCT__ "MatSolve_SuperLU_Private"
307c29e884eSHong Zhang PetscErrorCode MatSolve_SuperLU_Private(Mat A,Vec b,Vec x)
30814b396bbSKris Buschelman {
309f0c56d0fSKris Buschelman   Mat_SuperLU    *lu = (Mat_SuperLU*)A->spptr;
31075af56d4SHong Zhang   PetscScalar    *barray,*xarray;
311dfbe8321SBarry Smith   PetscErrorCode ierr;
312cae5a23dSHong Zhang   PetscInt       info,i,n=x->map->n;
313da7a1d00SHong Zhang   PetscReal      ferr,berr;
31414b396bbSKris Buschelman 
31514b396bbSKris Buschelman   PetscFunctionBegin;
316937a6b0eSHong Zhang   if ( lu->lwork == -1 ) {
317937a6b0eSHong Zhang     PetscFunctionReturn(0);
318937a6b0eSHong Zhang   }
319cae5a23dSHong Zhang 
32075af56d4SHong Zhang   lu->B.ncol = 1;   /* Set the number of right-hand side */
321cae5a23dSHong Zhang   if (lu->options.Equil && !lu->rhs_dup) {
322cae5a23dSHong Zhang     /* superlu overwrites b when Equil is used, thus create rhs_dup to keep user's b unchanged */
323cae5a23dSHong Zhang     ierr = PetscMalloc(n*sizeof(PetscScalar),&lu->rhs_dup);CHKERRQ(ierr);
324cae5a23dSHong Zhang   }
325cae5a23dSHong Zhang   if (lu->options.Equil) {
326cae5a23dSHong Zhang     /* Copy b into rsh_dup */
32775af56d4SHong Zhang     ierr = VecGetArray(b,&barray);CHKERRQ(ierr);
328cae5a23dSHong Zhang     ierr = PetscMemcpy(lu->rhs_dup,barray,n*sizeof(PetscScalar));CHKERRQ(ierr);
329cae5a23dSHong Zhang     ierr = VecRestoreArray(b,&barray);CHKERRQ(ierr);
330cae5a23dSHong Zhang     barray = lu->rhs_dup;
331cae5a23dSHong Zhang   } else {
332cae5a23dSHong Zhang     ierr = VecGetArray(b,&barray);CHKERRQ(ierr);
333cae5a23dSHong Zhang   }
33475af56d4SHong Zhang   ierr = VecGetArray(x,&xarray);CHKERRQ(ierr);
3355fe6150dSHong Zhang 
3365fe6150dSHong Zhang #if defined(PETSC_USE_COMPLEX)
3373cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
3383cb270beSHong Zhang   ((DNformat*)lu->B.Store)->nzval = (singlecomplex*)barray;
3393cb270beSHong Zhang   ((DNformat*)lu->X.Store)->nzval = (singlecomplex*)xarray;
3403cb270beSHong Zhang #else
3415fe6150dSHong Zhang   ((DNformat*)lu->B.Store)->nzval = (doublecomplex*)barray;
3425fe6150dSHong Zhang   ((DNformat*)lu->X.Store)->nzval = (doublecomplex*)xarray;
3433cb270beSHong Zhang #endif
3445fe6150dSHong Zhang #else
3455fe6150dSHong Zhang   ((DNformat*)lu->B.Store)->nzval = barray;
34675af56d4SHong Zhang   ((DNformat*)lu->X.Store)->nzval = xarray;
3475fe6150dSHong Zhang #endif
34875af56d4SHong Zhang 
34975af56d4SHong Zhang   lu->options.Fact = FACTORED; /* Indicate the factored form of A is supplied. */
350d5f3da31SBarry Smith   if (A->factortype == MAT_FACTOR_LU) {
3518914a3f7SHong Zhang #if defined(PETSC_USE_COMPLEX)
3523cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
3533cb270beSHong Zhang     cgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
3543cb270beSHong Zhang            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
3553cb270beSHong Zhang            &lu->mem_usage, &lu->stat, &info);
3563cb270beSHong Zhang #else
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,
3595d8b2efaSHong Zhang            &lu->mem_usage, &lu->stat, &info);
3603cb270beSHong Zhang #endif
3613cb270beSHong Zhang #else
3623cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
3633cb270beSHong Zhang     sgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
3643cb270beSHong Zhang            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
3653cb270beSHong Zhang            &lu->mem_usage, &lu->stat, &info);
3668914a3f7SHong Zhang #else
36775af56d4SHong Zhang     dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
36875af56d4SHong Zhang            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
3695d8b2efaSHong Zhang            &lu->mem_usage, &lu->stat, &info);
3708914a3f7SHong Zhang #endif
3713cb270beSHong Zhang #endif
372d5f3da31SBarry Smith   } else if (A->factortype == MAT_FACTOR_ILU) {
373cffbb591SHong Zhang #if defined(PETSC_USE_COMPLEX)
3743cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
3753cb270beSHong Zhang     cgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
3763cb270beSHong Zhang            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond,
3773cb270beSHong Zhang            &lu->mem_usage, &lu->stat, &info);
3783cb270beSHong Zhang #else
379cffbb591SHong Zhang     zgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
380cffbb591SHong Zhang            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond,
3815d8b2efaSHong Zhang            &lu->mem_usage, &lu->stat, &info);
3823cb270beSHong Zhang #endif
3833cb270beSHong Zhang #else
3843cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
3853cb270beSHong Zhang     sgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
3863cb270beSHong Zhang            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond,
3873cb270beSHong Zhang            &lu->mem_usage, &lu->stat, &info);
388cffbb591SHong Zhang #else
389cffbb591SHong Zhang     dgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
390cffbb591SHong Zhang            &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond,
3915d8b2efaSHong Zhang            &lu->mem_usage, &lu->stat, &info);
392cffbb591SHong Zhang #endif
3933cb270beSHong Zhang #endif
394cffbb591SHong Zhang   } else {
395e32f2f54SBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported");
396cffbb591SHong Zhang   }
397cae5a23dSHong Zhang   if (!lu->options.Equil) {
39875af56d4SHong Zhang     ierr = VecRestoreArray(b,&barray);CHKERRQ(ierr);
399cae5a23dSHong Zhang   }
40075af56d4SHong Zhang   ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
40175af56d4SHong Zhang 
402958c9bccSBarry Smith   if ( !info || info == lu->A.ncol+1 ) {
40375af56d4SHong Zhang     if ( lu->options.IterRefine ) {
4048914a3f7SHong Zhang       ierr = PetscPrintf(PETSC_COMM_SELF,"Iterative Refinement:\n");
4058914a3f7SHong Zhang       ierr = PetscPrintf(PETSC_COMM_SELF,"  %8s%8s%16s%16s\n", "rhs", "Steps", "FERR", "BERR");
40675af56d4SHong Zhang       for (i = 0; i < 1; ++i)
4075d8b2efaSHong Zhang         ierr = PetscPrintf(PETSC_COMM_SELF,"  %8d%8d%16e%16e\n", i+1, lu->stat.RefineSteps, ferr, berr);
40875af56d4SHong Zhang     }
4098914a3f7SHong Zhang   } else if ( info > 0 ) {
4108914a3f7SHong Zhang     if ( lu->lwork == -1 ) {
41177431f27SBarry Smith       ierr = PetscPrintf(PETSC_COMM_SELF,"  ** Estimated memory: %D bytes\n", info - lu->A.ncol);
4128914a3f7SHong Zhang     } else {
41377431f27SBarry Smith       ierr = PetscPrintf(PETSC_COMM_SELF,"  Warning: gssvx() returns info %D\n",info);
4148914a3f7SHong Zhang     }
4158914a3f7SHong Zhang   } else if (info < 0) {
416e32f2f54SBarry Smith     SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_LIB, "info = %D, the %D-th argument in gssvx() had an illegal value", info,-info);
41714b396bbSKris Buschelman   }
41814b396bbSKris Buschelman 
4198914a3f7SHong Zhang   if ( lu->options.PrintStat ) {
4208914a3f7SHong Zhang     ierr = PetscPrintf(PETSC_COMM_SELF,"MatSolve__SuperLU():\n");
4215d8b2efaSHong Zhang     StatPrint(&lu->stat);
4228914a3f7SHong Zhang   }
42375af56d4SHong Zhang   PetscFunctionReturn(0);
42475af56d4SHong Zhang }
42514b396bbSKris Buschelman 
42614b396bbSKris Buschelman #undef __FUNCT__
427c29e884eSHong Zhang #define __FUNCT__ "MatSolve_SuperLU"
428c29e884eSHong Zhang PetscErrorCode MatSolve_SuperLU(Mat A,Vec b,Vec x)
429c29e884eSHong Zhang {
430c29e884eSHong Zhang   Mat_SuperLU    *lu = (Mat_SuperLU*)A->spptr;
431c29e884eSHong Zhang   PetscErrorCode ierr;
432c29e884eSHong Zhang 
433c29e884eSHong Zhang   PetscFunctionBegin;
434c29e884eSHong Zhang   lu->options.Trans = TRANS;
435c29e884eSHong Zhang   ierr = MatSolve_SuperLU_Private(A,b,x);CHKERRQ(ierr);
436c29e884eSHong Zhang   PetscFunctionReturn(0);
437c29e884eSHong Zhang }
438c29e884eSHong Zhang 
439c29e884eSHong Zhang #undef __FUNCT__
440c7c1fe80SHong Zhang #define __FUNCT__ "MatSolveTranspose_SuperLU"
441c7c1fe80SHong Zhang PetscErrorCode MatSolveTranspose_SuperLU(Mat A,Vec b,Vec x)
442c7c1fe80SHong Zhang {
443c7c1fe80SHong Zhang   Mat_SuperLU    *lu = (Mat_SuperLU*)A->spptr;
444c7c1fe80SHong Zhang   PetscErrorCode ierr;
445c7c1fe80SHong Zhang 
446c7c1fe80SHong Zhang   PetscFunctionBegin;
447c7c1fe80SHong Zhang   lu->options.Trans = NOTRANS;
448c29e884eSHong Zhang   ierr = MatSolve_SuperLU_Private(A,b,x);CHKERRQ(ierr);
449c7c1fe80SHong Zhang   PetscFunctionReturn(0);
450c7c1fe80SHong Zhang }
451c7c1fe80SHong Zhang 
452e0b74bf9SHong Zhang #undef __FUNCT__
453e0b74bf9SHong Zhang #define __FUNCT__ "MatMatSolve_SuperLU"
454e0b74bf9SHong Zhang PetscErrorCode MatMatSolve_SuperLU(Mat A,Mat B,Mat X)
455e0b74bf9SHong Zhang {
456e0b74bf9SHong Zhang   Mat_SuperLU    *lu = (Mat_SuperLU*)A->spptr;
457cd723cd1SBarry Smith   PetscBool      flg;
458cd723cd1SBarry Smith   PetscErrorCode ierr;
459e0b74bf9SHong Zhang 
460e0b74bf9SHong Zhang   PetscFunctionBegin;
461251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompareAny((PetscObject)B,&flg,MATSEQDENSE,MATMPIDENSE,PETSC_NULL);CHKERRQ(ierr);
462cd723cd1SBarry Smith   if (!flg) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix");
463251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,PETSC_NULL);CHKERRQ(ierr);
464cd723cd1SBarry Smith   if (!flg) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix");  lu->options.Trans = TRANS;
465e0b74bf9SHong Zhang   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatMatSolve_SuperLU() is not implemented yet");
466e0b74bf9SHong Zhang   PetscFunctionReturn(0);
467e0b74bf9SHong Zhang }
468e0b74bf9SHong Zhang 
46914b396bbSKris Buschelman /*
47014b396bbSKris Buschelman    Note the r permutation is ignored
47114b396bbSKris Buschelman */
47214b396bbSKris Buschelman #undef __FUNCT__
473f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_SuperLU"
4740481f469SBarry Smith PetscErrorCode MatLUFactorSymbolic_SuperLU(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
47514b396bbSKris Buschelman {
4761d5ca7f3SHong Zhang   Mat_SuperLU    *lu = (Mat_SuperLU*)(F->spptr);
477b24902e0SBarry Smith 
478b24902e0SBarry Smith   PetscFunctionBegin;
479b24902e0SBarry Smith   lu->flg                 = DIFFERENT_NONZERO_PATTERN;
480b24902e0SBarry Smith   lu->CleanUpSuperLU      = PETSC_TRUE;
4811d5ca7f3SHong Zhang   F->ops->lufactornumeric = MatLUFactorNumeric_SuperLU;
482b24902e0SBarry Smith   PetscFunctionReturn(0);
483b24902e0SBarry Smith }
484b24902e0SBarry Smith 
48535bd34faSBarry Smith EXTERN_C_BEGIN
48635bd34faSBarry Smith #undef __FUNCT__
4875ccb76cbSHong Zhang #define __FUNCT__ "MatSuperluSetILUDropTol_SuperLU"
4885ccb76cbSHong Zhang PetscErrorCode MatSuperluSetILUDropTol_SuperLU(Mat F,PetscReal dtol)
4895ccb76cbSHong Zhang {
4905ccb76cbSHong Zhang   Mat_SuperLU *lu= (Mat_SuperLU*)F->spptr;
4915ccb76cbSHong Zhang 
4925ccb76cbSHong Zhang   PetscFunctionBegin;
4935ccb76cbSHong Zhang   lu->options.ILU_DropTol = dtol;
4945ccb76cbSHong Zhang   PetscFunctionReturn(0);
4955ccb76cbSHong Zhang }
4965ccb76cbSHong Zhang EXTERN_C_END
4975ccb76cbSHong Zhang 
4985ccb76cbSHong Zhang #undef __FUNCT__
4995ccb76cbSHong Zhang #define __FUNCT__ "MatSuperluSetILUDropTol"
5005ccb76cbSHong Zhang /*@
5015ccb76cbSHong Zhang   MatSuperluSetILUDropTol - Set SuperLU ILU drop tolerance
5025ccb76cbSHong Zhang    Logically Collective on Mat
5035ccb76cbSHong Zhang 
5045ccb76cbSHong Zhang    Input Parameters:
5055ccb76cbSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-SuperLU interface
5065ccb76cbSHong Zhang -  dtol - drop tolerance
5075ccb76cbSHong Zhang 
5085ccb76cbSHong Zhang   Options Database:
5095ccb76cbSHong Zhang .   -mat_superlu_ilu_droptol <dtol>
5105ccb76cbSHong Zhang 
5115ccb76cbSHong Zhang    Level: beginner
5125ccb76cbSHong Zhang 
5135ccb76cbSHong Zhang    References: SuperLU Users' Guide
5145ccb76cbSHong Zhang 
5155ccb76cbSHong Zhang .seealso: MatGetFactor()
5165ccb76cbSHong Zhang @*/
5175ccb76cbSHong Zhang PetscErrorCode MatSuperluSetILUDropTol(Mat F,PetscReal dtol)
5185ccb76cbSHong Zhang {
5195ccb76cbSHong Zhang   PetscErrorCode ierr;
5205ccb76cbSHong Zhang 
5215ccb76cbSHong Zhang   PetscFunctionBegin;
5225ccb76cbSHong Zhang   PetscValidHeaderSpecific(F,MAT_CLASSID,1);
5235ccb76cbSHong Zhang   PetscValidLogicalCollectiveInt(F,dtol,2);
5245ccb76cbSHong Zhang   ierr = PetscTryMethod(F,"MatSuperluSetILUDropTol_C",(Mat,PetscReal),(F,dtol));CHKERRQ(ierr);
5255ccb76cbSHong Zhang   PetscFunctionReturn(0);
5265ccb76cbSHong Zhang }
5275ccb76cbSHong Zhang 
5285ccb76cbSHong Zhang EXTERN_C_BEGIN
5295ccb76cbSHong Zhang #undef __FUNCT__
53035bd34faSBarry Smith #define __FUNCT__ "MatFactorGetSolverPackage_seqaij_superlu"
53135bd34faSBarry Smith PetscErrorCode MatFactorGetSolverPackage_seqaij_superlu(Mat A,const MatSolverPackage *type)
53235bd34faSBarry Smith {
53335bd34faSBarry Smith   PetscFunctionBegin;
5342692d6eeSBarry Smith   *type = MATSOLVERSUPERLU;
53535bd34faSBarry Smith   PetscFunctionReturn(0);
53635bd34faSBarry Smith }
53735bd34faSBarry Smith EXTERN_C_END
53835bd34faSBarry Smith 
539b24902e0SBarry Smith 
540b24902e0SBarry Smith /*MC
541ba992d64SSatish Balay   MATSOLVERSUPERLU = "superlu" - A solver package providing solvers LU and ILU for sequential matrices
542b24902e0SBarry Smith   via the external package SuperLU.
543b24902e0SBarry Smith 
544e2e64c6bSBarry Smith   Use ./configure --download-superlu to have PETSc installed with SuperLU
545b24902e0SBarry Smith 
546b24902e0SBarry Smith   Options Database Keys:
547e08999f5SMatthew G Knepley + -mat_superlu_equil <FALSE>            - Equil (None)
548e08999f5SMatthew G Knepley . -mat_superlu_colperm <COLAMD>         - (choose one of) NATURAL MMD_ATA MMD_AT_PLUS_A COLAMD
549e08999f5SMatthew G Knepley . -mat_superlu_iterrefine <NOREFINE>    - (choose one of) NOREFINE SINGLE DOUBLE EXTRA
550e08999f5SMatthew G Knepley . -mat_superlu_symmetricmode: <FALSE>   - SymmetricMode (None)
551e08999f5SMatthew G Knepley . -mat_superlu_diagpivotthresh <1>      - DiagPivotThresh (None)
552e08999f5SMatthew G Knepley . -mat_superlu_pivotgrowth <FALSE>      - PivotGrowth (None)
553e08999f5SMatthew G Knepley . -mat_superlu_conditionnumber <FALSE>  - ConditionNumber (None)
554e08999f5SMatthew G Knepley . -mat_superlu_rowperm <NOROWPERM>      - (choose one of) NOROWPERM LargeDiag
555e08999f5SMatthew G Knepley . -mat_superlu_replacetinypivot <FALSE> - ReplaceTinyPivot (None)
556e08999f5SMatthew G Knepley . -mat_superlu_printstat <FALSE>        - PrintStat (None)
557e08999f5SMatthew G Knepley . -mat_superlu_lwork <0>                - size of work array in bytes used by factorization (None)
558e08999f5SMatthew G Knepley . -mat_superlu_ilu_droptol <0>          - ILU_DropTol (None)
559e08999f5SMatthew G Knepley . -mat_superlu_ilu_filltol <0>          - ILU_FillTol (None)
560e08999f5SMatthew G Knepley . -mat_superlu_ilu_fillfactor <0>       - ILU_FillFactor (None)
561e08999f5SMatthew G Knepley . -mat_superlu_ilu_droprull <0>         - ILU_DropRule (None)
562e08999f5SMatthew G Knepley . -mat_superlu_ilu_norm <0>             - ILU_Norm (None)
563e08999f5SMatthew G Knepley - -mat_superlu_ilu_milu <0>             - ILU_MILU (None)
564b24902e0SBarry Smith 
5652692d6eeSBarry Smith    Notes: Do not confuse this with MATSOLVERSUPERLU_DIST which is for parallel sparse solves
5665c9eb25fSBarry Smith 
567b24902e0SBarry Smith    Level: beginner
568b24902e0SBarry Smith 
569d45987f3SHong Zhang .seealso: PCLU, PCILU, MATSOLVERSUPERLU_DIST, MATSOLVERMUMPS, PCFactorSetMatSolverPackage(), MatSolverPackage
570b24902e0SBarry Smith M*/
571b24902e0SBarry Smith 
5725c9eb25fSBarry Smith EXTERN_C_BEGIN
573b24902e0SBarry Smith #undef __FUNCT__
574b24902e0SBarry Smith #define __FUNCT__ "MatGetFactor_seqaij_superlu"
5755c9eb25fSBarry Smith PetscErrorCode MatGetFactor_seqaij_superlu(Mat A,MatFactorType ftype,Mat *F)
576b24902e0SBarry Smith {
57714b396bbSKris Buschelman   Mat            B;
578f0c56d0fSKris Buschelman   Mat_SuperLU    *lu;
5796849ba73SBarry Smith   PetscErrorCode ierr;
5805d8b2efaSHong Zhang   PetscInt       indx,m=A->rmap->n,n=A->cmap->n;
581ace3abfcSBarry Smith   PetscBool      flg;
5823cb270beSHong Zhang   PetscReal      real_input;
5835ca28756SHong Zhang   const char     *colperm[]={"NATURAL","MMD_ATA","MMD_AT_PLUS_A","COLAMD"}; /* MY_PERMC - not supported by the petsc interface yet */
5845ca28756SHong Zhang   const char     *iterrefine[]={"NOREFINE", "SINGLE", "DOUBLE", "EXTRA"};
5855ca28756SHong Zhang   const char     *rowperm[]={"NOROWPERM", "LargeDiag"}; /* MY_PERMC - not supported by the petsc interface yet */
58614b396bbSKris Buschelman 
58714b396bbSKris Buschelman   PetscFunctionBegin;
5887adad957SLisandro Dalcin   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
589d0f46423SBarry Smith   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
5907adad957SLisandro Dalcin   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
591899d7b4fSKris Buschelman   ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
592f0c56d0fSKris Buschelman 
593cffbb591SHong Zhang   if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU) {
594b24902e0SBarry Smith     B->ops->lufactorsymbolic  = MatLUFactorSymbolic_SuperLU;
595cffbb591SHong Zhang     B->ops->ilufactorsymbolic = MatLUFactorSymbolic_SuperLU;
596b3a44c85SBarry Smith   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported");
597cffbb591SHong Zhang 
598b24902e0SBarry Smith   B->ops->destroy          = MatDestroy_SuperLU;
5993519fcdcSHong Zhang   B->ops->view             = MatView_SuperLU;
600d5f3da31SBarry Smith   B->factortype            = ftype;
60194b7f48cSBarry Smith   B->assembled             = PETSC_TRUE;  /* required by -ksp_view */
6025c9eb25fSBarry Smith   B->preallocated          = PETSC_TRUE;
60314b396bbSKris Buschelman 
604b24902e0SBarry Smith   ierr = PetscNewLog(B,Mat_SuperLU,&lu);CHKERRQ(ierr);
605cae5a23dSHong Zhang 
606cffbb591SHong Zhang   if (ftype == MAT_FACTOR_LU) {
6079ce50919SHong Zhang     set_default_options(&lu->options);
6083d6581fbSHong Zhang     /* Comments from SuperLU_4.0/SRC/dgssvx.c:
6093d6581fbSHong Zhang       "Whether or not the system will be equilibrated depends on the
6103d6581fbSHong Zhang        scaling of the matrix A, but if equilibration is used, A is
6113d6581fbSHong Zhang        overwritten by diag(R)*A*diag(C) and B by diag(R)*B
6123d6581fbSHong Zhang        (if options->Trans=NOTRANS) or diag(C)*B (if options->Trans = TRANS or CONJ)."
6133d6581fbSHong Zhang      We set 'options.Equil = NO' as default because additional space is needed for it.
6143d6581fbSHong Zhang     */
6153d6581fbSHong Zhang     lu->options.Equil = NO;
616cffbb591SHong Zhang   } else if (ftype == MAT_FACTOR_ILU) {
6170c74a584SJed Brown     /* Set the default input options of ilu: */
618cffbb591SHong Zhang     ilu_set_default_options(&lu->options);
619cffbb591SHong Zhang   }
6209ce50919SHong Zhang   lu->options.PrintStat = NO;
6211a2e2f44SHong Zhang 
6225d8b2efaSHong Zhang   /* Initialize the statistics variables. */
6235d8b2efaSHong Zhang   StatInit(&lu->stat);
6248914a3f7SHong Zhang   lu->lwork = 0;   /* allocate space internally by system malloc */
6259ce50919SHong Zhang 
6267adad957SLisandro Dalcin   ierr = PetscOptionsBegin(((PetscObject)A)->comm,((PetscObject)A)->prefix,"SuperLU Options","Mat");CHKERRQ(ierr);
627acfcf0e5SJed Brown     ierr = PetscOptionsBool("-mat_superlu_equil","Equil","None",(PetscBool)lu->options.Equil,(PetscBool*)&lu->options.Equil,0);CHKERRQ(ierr);
6288914a3f7SHong Zhang     ierr = PetscOptionsEList("-mat_superlu_colperm","ColPerm","None",colperm,4,colperm[3],&indx,&flg);CHKERRQ(ierr);
6299ce50919SHong Zhang     if (flg) {lu->options.ColPerm = (colperm_t)indx;}
6308914a3f7SHong Zhang     ierr = PetscOptionsEList("-mat_superlu_iterrefine","IterRefine","None",iterrefine,4,iterrefine[0],&indx,&flg);CHKERRQ(ierr);
6319ce50919SHong Zhang     if (flg) { lu->options.IterRefine = (IterRefine_t)indx;}
632acfcf0e5SJed Brown     ierr = PetscOptionsBool("-mat_superlu_symmetricmode","SymmetricMode","None",(PetscBool)lu->options.SymmetricMode,&flg,0);CHKERRQ(ierr);
6339ce50919SHong Zhang     if (flg) lu->options.SymmetricMode = YES;
6343cb270beSHong Zhang     ierr = PetscOptionsReal("-mat_superlu_diagpivotthresh","DiagPivotThresh","None",lu->options.DiagPivotThresh,&real_input,&flg);CHKERRQ(ierr);
6353cb270beSHong Zhang     if (flg) lu->options.DiagPivotThresh = (double) real_input;
636acfcf0e5SJed Brown     ierr = PetscOptionsBool("-mat_superlu_pivotgrowth","PivotGrowth","None",(PetscBool)lu->options.PivotGrowth,&flg,0);CHKERRQ(ierr);
6379ce50919SHong Zhang     if (flg) lu->options.PivotGrowth = YES;
638acfcf0e5SJed Brown     ierr = PetscOptionsBool("-mat_superlu_conditionnumber","ConditionNumber","None",(PetscBool)lu->options.ConditionNumber,&flg,0);CHKERRQ(ierr);
6399ce50919SHong Zhang     if (flg) lu->options.ConditionNumber = YES;
640d7ebd59bSHong Zhang     ierr = PetscOptionsEList("-mat_superlu_rowperm","rowperm","None",rowperm,2,rowperm[lu->options.RowPerm],&indx,&flg);CHKERRQ(ierr);
6419ce50919SHong Zhang     if (flg) {lu->options.RowPerm = (rowperm_t)indx;}
642acfcf0e5SJed Brown     ierr = PetscOptionsBool("-mat_superlu_replacetinypivot","ReplaceTinyPivot","None",(PetscBool)lu->options.ReplaceTinyPivot,&flg,0);CHKERRQ(ierr);
6439ce50919SHong Zhang     if (flg) lu->options.ReplaceTinyPivot = YES;
644acfcf0e5SJed Brown     ierr = PetscOptionsBool("-mat_superlu_printstat","PrintStat","None",(PetscBool)lu->options.PrintStat,&flg,0);CHKERRQ(ierr);
6459ce50919SHong Zhang     if (flg) lu->options.PrintStat = YES;
6468914a3f7SHong Zhang     ierr = PetscOptionsInt("-mat_superlu_lwork","size of work array in bytes used by factorization","None",lu->lwork,&lu->lwork,PETSC_NULL);CHKERRQ(ierr);
6475fe6150dSHong Zhang     if (lu->lwork > 0 ) {
6485fe6150dSHong Zhang       ierr = PetscMalloc(lu->lwork,&lu->work);CHKERRQ(ierr);
6495fe6150dSHong Zhang     } else if (lu->lwork != 0 && lu->lwork != -1) {
65077431f27SBarry Smith       ierr = PetscPrintf(PETSC_COMM_SELF,"   Warning: lwork %D is not supported by SUPERLU. The default lwork=0 is used.\n",lu->lwork);
6518914a3f7SHong Zhang       lu->lwork = 0;
6528914a3f7SHong Zhang     }
653cffbb591SHong Zhang     /* ilu options */
6543cb270beSHong Zhang     ierr = PetscOptionsReal("-mat_superlu_ilu_droptol","ILU_DropTol","None",lu->options.ILU_DropTol,&real_input,&flg);CHKERRQ(ierr);
6553cb270beSHong Zhang     if (flg) lu->options.ILU_DropTol = (double) real_input;
6563cb270beSHong Zhang     ierr = PetscOptionsReal("-mat_superlu_ilu_filltol","ILU_FillTol","None",lu->options.ILU_FillTol,&real_input,&flg);CHKERRQ(ierr);
6573cb270beSHong Zhang     if (flg) lu->options.ILU_FillTol = (double) real_input;
6583cb270beSHong Zhang     ierr = PetscOptionsReal("-mat_superlu_ilu_fillfactor","ILU_FillFactor","None",lu->options.ILU_FillFactor,&real_input,&flg);CHKERRQ(ierr);
6593cb270beSHong Zhang     if (flg) lu->options.ILU_FillFactor = (double) real_input;
660cffbb591SHong Zhang     ierr = PetscOptionsInt("-mat_superlu_ilu_droprull","ILU_DropRule","None",lu->options.ILU_DropRule,&lu->options.ILU_DropRule,PETSC_NULL);CHKERRQ(ierr);
661cffbb591SHong Zhang     ierr = PetscOptionsInt("-mat_superlu_ilu_norm","ILU_Norm","None",lu->options.ILU_Norm,&indx,&flg);CHKERRQ(ierr);
662cffbb591SHong Zhang     if (flg) {
663cffbb591SHong Zhang       lu->options.ILU_Norm = (norm_t)indx;
664cffbb591SHong Zhang     }
665cffbb591SHong Zhang     ierr = PetscOptionsInt("-mat_superlu_ilu_milu","ILU_MILU","None",lu->options.ILU_MILU,&indx,&flg);CHKERRQ(ierr);
666cffbb591SHong Zhang     if (flg) {
667cffbb591SHong Zhang       lu->options.ILU_MILU = (milu_t)indx;
668cffbb591SHong Zhang     }
6699ce50919SHong Zhang   PetscOptionsEnd();
67038abfdaeSHong Zhang   if (lu->options.Equil == YES) {
67138abfdaeSHong Zhang     /* superlu overwrites input matrix and rhs when Equil is used, thus create A_dup to keep user's A unchanged */
67238abfdaeSHong Zhang     ierr = MatDuplicate_SeqAIJ(A,MAT_COPY_VALUES,&lu->A_dup);CHKERRQ(ierr);
67338abfdaeSHong Zhang   }
6749ce50919SHong Zhang 
6755d8b2efaSHong Zhang   /* Allocate spaces (notice sizes are for the transpose) */
6765d8b2efaSHong Zhang   ierr = PetscMalloc(m*sizeof(PetscInt),&lu->etree);CHKERRQ(ierr);
6775d8b2efaSHong Zhang   ierr = PetscMalloc(n*sizeof(PetscInt),&lu->perm_r);CHKERRQ(ierr);
6785d8b2efaSHong Zhang   ierr = PetscMalloc(m*sizeof(PetscInt),&lu->perm_c);CHKERRQ(ierr);
6795d8b2efaSHong Zhang   ierr = PetscMalloc(n*sizeof(PetscScalar),&lu->R);CHKERRQ(ierr);
6805d8b2efaSHong Zhang   ierr = PetscMalloc(m*sizeof(PetscScalar),&lu->C);CHKERRQ(ierr);
6815d8b2efaSHong Zhang 
6825d8b2efaSHong Zhang   /* create rhs and solution x without allocate space for .Store */
6835d8b2efaSHong Zhang #if defined(PETSC_USE_COMPLEX)
6843cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
6853cb270beSHong Zhang   cCreate_Dense_Matrix(&lu->B, m, 1, PETSC_NULL, m, SLU_DN, SLU_C, SLU_GE);
6863cb270beSHong Zhang   cCreate_Dense_Matrix(&lu->X, m, 1, PETSC_NULL, m, SLU_DN, SLU_C, SLU_GE);
6873cb270beSHong Zhang #else
6885d8b2efaSHong Zhang   zCreate_Dense_Matrix(&lu->B, m, 1, PETSC_NULL, m, SLU_DN, SLU_Z, SLU_GE);
6895d8b2efaSHong Zhang   zCreate_Dense_Matrix(&lu->X, m, 1, PETSC_NULL, m, SLU_DN, SLU_Z, SLU_GE);
6903cb270beSHong Zhang #endif
6913cb270beSHong Zhang #else
6923cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
6933cb270beSHong Zhang   sCreate_Dense_Matrix(&lu->B, m, 1, PETSC_NULL, m, SLU_DN, SLU_S, SLU_GE);
6943cb270beSHong Zhang   sCreate_Dense_Matrix(&lu->X, m, 1, PETSC_NULL, m, SLU_DN, SLU_S, SLU_GE);
6955d8b2efaSHong Zhang #else
6965d8b2efaSHong Zhang   dCreate_Dense_Matrix(&lu->B, m, 1, PETSC_NULL, m, SLU_DN, SLU_D, SLU_GE);
6975d8b2efaSHong Zhang   dCreate_Dense_Matrix(&lu->X, m, 1, PETSC_NULL, m, SLU_DN, SLU_D, SLU_GE);
6985d8b2efaSHong Zhang #endif
6993cb270beSHong Zhang #endif
7005d8b2efaSHong Zhang 
701*519f805aSKarl Rupp #if defined(SUPERLU2)
7025c9eb25fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatCreateNull","MatCreateNull_SuperLU",(void(*)(void))MatCreateNull_SuperLU);CHKERRQ(ierr);
70375af56d4SHong Zhang #endif
70435bd34faSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_seqaij_superlu",MatFactorGetSolverPackage_seqaij_superlu);CHKERRQ(ierr);
7055ccb76cbSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSuperluSetILUDropTol_C","MatSuperluSetILUDropTol_SuperLU",MatSuperluSetILUDropTol_SuperLU);CHKERRQ(ierr);
7065c9eb25fSBarry Smith   B->spptr = lu;
707899d7b4fSKris Buschelman   *F = B;
70814b396bbSKris Buschelman   PetscFunctionReturn(0);
70914b396bbSKris Buschelman }
7105c9eb25fSBarry Smith EXTERN_C_END
711d954db57SHong Zhang 
712