114b396bbSKris Buschelman 25a46d3faSBarry Smith /* -------------------------------------------------------------------- 35a46d3faSBarry Smith 45a46d3faSBarry Smith This file implements a subclass of the SeqAIJ matrix class that uses 5c2b89b5dSBarry Smith the SuperLU sparse solver. 614b396bbSKris Buschelman */ 714b396bbSKris Buschelman 85a46d3faSBarry Smith /* 95a46d3faSBarry Smith Defines the data structure for the base matrix type (SeqAIJ) 105a46d3faSBarry Smith */ 11c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h> /*I "petscmat.h" I*/ 1214b396bbSKris Buschelman 135a46d3faSBarry Smith /* 145a46d3faSBarry Smith SuperLU include files 155a46d3faSBarry Smith */ 1614b396bbSKris Buschelman EXTERN_C_BEGIN 1714b396bbSKris Buschelman #if defined(PETSC_USE_COMPLEX) 183cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 193cb270beSHong Zhang #include <slu_cdefs.h> 203cb270beSHong Zhang #else 21c6db04a5SJed Brown #include <slu_zdefs.h> 223cb270beSHong Zhang #endif 233cb270beSHong Zhang #else 243cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 253cb270beSHong Zhang #include <slu_sdefs.h> 2614b396bbSKris Buschelman #else 27c6db04a5SJed Brown #include <slu_ddefs.h> 2814b396bbSKris Buschelman #endif 293cb270beSHong Zhang #endif 30c6db04a5SJed Brown #include <slu_util.h> 3114b396bbSKris Buschelman EXTERN_C_END 3214b396bbSKris Buschelman 335a46d3faSBarry Smith /* 345a46d3faSBarry Smith This is the data we are "ADDING" to the SeqAIJ matrix type to get the SuperLU matrix type 355a46d3faSBarry Smith */ 3614b396bbSKris Buschelman typedef struct { 3775af56d4SHong Zhang SuperMatrix A,L,U,B,X; 3875af56d4SHong Zhang superlu_options_t options; 39da7a1d00SHong Zhang PetscInt *perm_c; /* column permutation vector */ 40da7a1d00SHong Zhang PetscInt *perm_r; /* row permutations from partial pivoting */ 41da7a1d00SHong Zhang PetscInt *etree; 42da7a1d00SHong Zhang PetscReal *R, *C; 4375af56d4SHong Zhang char equed[1]; 44da7a1d00SHong Zhang PetscInt lwork; 4575af56d4SHong Zhang void *work; 46da7a1d00SHong Zhang PetscReal rpg, rcond; 4775af56d4SHong Zhang mem_usage_t mem_usage; 4875af56d4SHong Zhang MatStructure flg; 495d8b2efaSHong Zhang SuperLUStat_t stat; 50cae5a23dSHong Zhang Mat A_dup; 51cae5a23dSHong Zhang PetscScalar *rhs_dup; 5214b396bbSKris Buschelman 5314b396bbSKris Buschelman /* Flag to clean up (non-global) SuperLU objects during Destroy */ 54ace3abfcSBarry Smith PetscBool CleanUpSuperLU; 55f0c56d0fSKris Buschelman } Mat_SuperLU; 5614b396bbSKris Buschelman 57e0e586b9SHong Zhang extern PetscErrorCode MatFactorInfo_SuperLU(Mat,PetscViewer); 580481f469SBarry Smith extern PetscErrorCode MatLUFactorNumeric_SuperLU(Mat,Mat,const MatFactorInfo*); 59e0e586b9SHong Zhang extern PetscErrorCode MatDestroy_SuperLU(Mat); 60e0e586b9SHong Zhang extern PetscErrorCode MatView_SuperLU(Mat,PetscViewer); 61e0e586b9SHong Zhang extern PetscErrorCode MatAssemblyEnd_SuperLU(Mat,MatAssemblyType); 62e0e586b9SHong Zhang extern PetscErrorCode MatSolve_SuperLU(Mat,Vec,Vec); 63e0b74bf9SHong Zhang extern PetscErrorCode MatMatSolve_SuperLU(Mat,Mat,Mat); 64e0e586b9SHong Zhang extern PetscErrorCode MatSolveTranspose_SuperLU(Mat,Vec,Vec); 650481f469SBarry Smith extern PetscErrorCode MatLUFactorSymbolic_SuperLU(Mat,Mat,IS,IS,const MatFactorInfo*); 66e0e586b9SHong Zhang extern PetscErrorCode MatDuplicate_SuperLU(Mat, MatDuplicateOption, Mat*); 67e0e586b9SHong Zhang 685a46d3faSBarry Smith /* 695a46d3faSBarry Smith Utility function 705a46d3faSBarry Smith */ 715a46d3faSBarry Smith #undef __FUNCT__ 725a46d3faSBarry Smith #define __FUNCT__ "MatFactorInfo_SuperLU" 735a46d3faSBarry Smith PetscErrorCode MatFactorInfo_SuperLU(Mat A,PetscViewer viewer) 745a46d3faSBarry Smith { 755a46d3faSBarry Smith Mat_SuperLU *lu= (Mat_SuperLU*)A->spptr; 765a46d3faSBarry Smith PetscErrorCode ierr; 775a46d3faSBarry Smith superlu_options_t options; 785a46d3faSBarry Smith 795a46d3faSBarry Smith PetscFunctionBegin; 805a46d3faSBarry Smith /* check if matrix is superlu_dist type */ 815a46d3faSBarry Smith if (A->ops->solve != MatSolve_SuperLU) PetscFunctionReturn(0); 825a46d3faSBarry Smith 835a46d3faSBarry Smith options = lu->options; 842205254eSKarl Rupp 855a46d3faSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"SuperLU run parameters:\n");CHKERRQ(ierr); 865a46d3faSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Equil: %s\n",(options.Equil != NO) ? "YES" : "NO");CHKERRQ(ierr); 875a46d3faSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," ColPerm: %D\n",options.ColPerm);CHKERRQ(ierr); 885a46d3faSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," IterRefine: %D\n",options.IterRefine);CHKERRQ(ierr); 895a46d3faSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," SymmetricMode: %s\n",(options.SymmetricMode != NO) ? "YES" : "NO");CHKERRQ(ierr); 905a46d3faSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," DiagPivotThresh: %g\n",options.DiagPivotThresh);CHKERRQ(ierr); 915a46d3faSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," PivotGrowth: %s\n",(options.PivotGrowth != NO) ? "YES" : "NO");CHKERRQ(ierr); 925a46d3faSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," ConditionNumber: %s\n",(options.ConditionNumber != NO) ? "YES" : "NO");CHKERRQ(ierr); 935a46d3faSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," RowPerm: %D\n",options.RowPerm);CHKERRQ(ierr); 945a46d3faSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," ReplaceTinyPivot: %s\n",(options.ReplaceTinyPivot != NO) ? "YES" : "NO");CHKERRQ(ierr); 955a46d3faSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," PrintStat: %s\n",(options.PrintStat != NO) ? "YES" : "NO");CHKERRQ(ierr); 965a46d3faSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," lwork: %D\n",lu->lwork);CHKERRQ(ierr); 97d5f3da31SBarry Smith if (A->factortype == MAT_FACTOR_ILU) { 98cffbb591SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ILU_DropTol: %g\n",options.ILU_DropTol);CHKERRQ(ierr); 99cffbb591SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ILU_FillTol: %g\n",options.ILU_FillTol);CHKERRQ(ierr); 100cffbb591SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ILU_FillFactor: %g\n",options.ILU_FillFactor);CHKERRQ(ierr); 101cffbb591SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ILU_DropRule: %D\n",options.ILU_DropRule);CHKERRQ(ierr); 102cffbb591SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ILU_Norm: %D\n",options.ILU_Norm);CHKERRQ(ierr); 103cffbb591SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ILU_MILU: %D\n",options.ILU_MILU);CHKERRQ(ierr); 104cffbb591SHong Zhang } 1055a46d3faSBarry Smith PetscFunctionReturn(0); 1065a46d3faSBarry Smith } 1075a46d3faSBarry Smith 1085a46d3faSBarry Smith /* 1095a46d3faSBarry Smith These are the methods provided to REPLACE the corresponding methods of the 1105a46d3faSBarry Smith SeqAIJ matrix class. Other methods of SeqAIJ are not replaced 1115a46d3faSBarry Smith */ 1125a46d3faSBarry Smith #undef __FUNCT__ 1135a46d3faSBarry Smith #define __FUNCT__ "MatLUFactorNumeric_SuperLU" 1140481f469SBarry Smith PetscErrorCode MatLUFactorNumeric_SuperLU(Mat F,Mat A,const MatFactorInfo *info) 1155a46d3faSBarry Smith { 1161d5ca7f3SHong Zhang Mat_SuperLU *lu = (Mat_SuperLU*)F->spptr; 117cae5a23dSHong Zhang Mat_SeqAIJ *aa; 1185a46d3faSBarry Smith PetscErrorCode ierr; 1195a46d3faSBarry Smith PetscInt sinfo; 1205a46d3faSBarry Smith PetscReal ferr, berr; 1215a46d3faSBarry Smith NCformat *Ustore; 1225a46d3faSBarry Smith SCformat *Lstore; 1235a46d3faSBarry Smith 1245a46d3faSBarry Smith PetscFunctionBegin; 1255a46d3faSBarry Smith if (lu->flg == SAME_NONZERO_PATTERN) { /* successing numerical factorization */ 1265a46d3faSBarry Smith lu->options.Fact = SamePattern; 1275a46d3faSBarry Smith /* Ref: ~SuperLU_3.0/EXAMPLE/dlinsolx2.c */ 1285a46d3faSBarry Smith Destroy_SuperMatrix_Store(&lu->A); 129cae5a23dSHong Zhang if (lu->options.Equil) { 130cae5a23dSHong Zhang ierr = MatCopy_SeqAIJ(A,lu->A_dup,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 131cae5a23dSHong Zhang } 1325a46d3faSBarry Smith if (lu->lwork >= 0) { 133d387c056SBarry Smith PetscStackCall("SuperLU:Destroy_SuperNode_Matrix",Destroy_SuperNode_Matrix(&lu->L)); 134d387c056SBarry Smith PetscStackCall("SuperLU:Destroy_CompCol_Matrix",Destroy_CompCol_Matrix(&lu->U)); 1355a46d3faSBarry Smith lu->options.Fact = SamePattern; 1365a46d3faSBarry Smith } 1375a46d3faSBarry Smith } 1385a46d3faSBarry Smith 1395a46d3faSBarry Smith /* Create the SuperMatrix for lu->A=A^T: 1405a46d3faSBarry Smith Since SuperLU likes column-oriented matrices,we pass it the transpose, 1415a46d3faSBarry Smith and then solve A^T X = B in MatSolve(). */ 142cae5a23dSHong Zhang if (lu->options.Equil) { 143cae5a23dSHong Zhang aa = (Mat_SeqAIJ*)(lu->A_dup)->data; 144cae5a23dSHong Zhang } else { 145cae5a23dSHong Zhang aa = (Mat_SeqAIJ*)(A)->data; 146cae5a23dSHong Zhang } 1475a46d3faSBarry Smith #if defined(PETSC_USE_COMPLEX) 1483cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 149d387c056SBarry Smith PetscStackCall("SuperLU:cCreate_CompCol_Matrix",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)); 1503cb270beSHong Zhang #else 151d387c056SBarry Smith PetscStackCall("SuperLU:zCreate_CompCol_Matrix",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)); 1523cb270beSHong Zhang #endif 1533cb270beSHong Zhang #else 1543cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 155d387c056SBarry Smith PetscStackCall("SuperLU:sCreate_CompCol_Matrix",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)); 1565a46d3faSBarry Smith #else 157d387c056SBarry Smith PetscStackCall("SuperLU:dCreate_CompCol_Matrix",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)); 1585a46d3faSBarry Smith #endif 1593cb270beSHong Zhang #endif 1605a46d3faSBarry Smith 1615a46d3faSBarry Smith /* Numerical factorization */ 1625a46d3faSBarry Smith lu->B.ncol = 0; /* Indicate not to solve the system */ 163d5f3da31SBarry Smith if (F->factortype == MAT_FACTOR_LU) { 1645a46d3faSBarry Smith #if defined(PETSC_USE_COMPLEX) 1653cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 166d387c056SBarry Smith PetscStackCall("SuperLU:cgssvx",cgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 1673cb270beSHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 168d387c056SBarry Smith &lu->mem_usage, &lu->stat, &sinfo)); 1693cb270beSHong Zhang #else 170d387c056SBarry Smith PetscStackCall("SuperLU:zgssvx",zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 1715a46d3faSBarry Smith &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 172d387c056SBarry Smith &lu->mem_usage, &lu->stat, &sinfo)); 1733cb270beSHong Zhang #endif 1743cb270beSHong Zhang #else 1753cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 176d387c056SBarry Smith PetscStackCall("SuperLU:sgssvx",sgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 1773cb270beSHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 178d387c056SBarry Smith &lu->mem_usage, &lu->stat, &sinfo)); 1795a46d3faSBarry Smith #else 180d387c056SBarry Smith PetscStackCall("SuperLU:dgssvx",dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 1815a46d3faSBarry Smith &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 182d387c056SBarry Smith &lu->mem_usage, &lu->stat, &sinfo)); 1835a46d3faSBarry Smith #endif 1843cb270beSHong Zhang #endif 185d5f3da31SBarry Smith } else if (F->factortype == MAT_FACTOR_ILU) { 186cffbb591SHong Zhang /* Compute the incomplete factorization, condition number and pivot growth */ 187cffbb591SHong Zhang #if defined(PETSC_USE_COMPLEX) 1883cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 189d387c056SBarry Smith PetscStackCall("SuperLU:cgsisx",cgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r,lu->etree, lu->equed, lu->R, lu->C, 1903cb270beSHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, 191d387c056SBarry Smith &lu->mem_usage, &lu->stat, &sinfo)); 1923cb270beSHong Zhang #else 193d387c056SBarry Smith PetscStackCall("SuperLU:zgsisx",zgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r,lu->etree, lu->equed, lu->R, lu->C, 194cffbb591SHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, 195d387c056SBarry Smith &lu->mem_usage, &lu->stat, &sinfo)); 1963cb270beSHong Zhang #endif 1973cb270beSHong Zhang #else 1983cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 199d387c056SBarry Smith PetscStackCall("SuperLU:sgsisx",sgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 2003cb270beSHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, 201d387c056SBarry Smith &lu->mem_usage, &lu->stat, &sinfo)); 202cffbb591SHong Zhang #else 203d387c056SBarry Smith PetscStackCall("SuperLU:dgsisx",dgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 204cffbb591SHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, 205d387c056SBarry Smith &lu->mem_usage, &lu->stat, &sinfo)); 206cffbb591SHong Zhang #endif 2073cb270beSHong Zhang #endif 208f23aa3ddSBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported"); 2095a46d3faSBarry Smith if (!sinfo || sinfo == lu->A.ncol+1) { 2102205254eSKarl Rupp if (lu->options.PivotGrowth) { 2115a46d3faSBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," Recip. pivot growth = %e\n", lu->rpg); 2122205254eSKarl Rupp } 2132205254eSKarl Rupp if (lu->options.ConditionNumber) { 2145a46d3faSBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," Recip. condition number = %e\n", lu->rcond); 2152205254eSKarl Rupp } 2165a46d3faSBarry Smith } else if (sinfo > 0) { 2175a46d3faSBarry Smith if (lu->lwork == -1) { 2185a46d3faSBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," ** Estimated memory: %D bytes\n", sinfo - lu->A.ncol); 2199308c137SBarry Smith } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot in row %D",sinfo); 220f23aa3ddSBarry Smith } else SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_LIB, "info = %D, the %D-th argument in gssvx() had an illegal value", sinfo,-sinfo); 2215a46d3faSBarry Smith 2225a46d3faSBarry Smith if (lu->options.PrintStat) { 2235a46d3faSBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF,"MatLUFactorNumeric_SuperLU():\n"); 224d387c056SBarry Smith PetscStackCall("SuperLU:StatPrint",StatPrint(&lu->stat)); 2255a46d3faSBarry Smith Lstore = (SCformat*) lu->L.Store; 2265a46d3faSBarry Smith Ustore = (NCformat*) lu->U.Store; 2275a46d3faSBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," No of nonzeros in factor L = %D\n", Lstore->nnz); 2285a46d3faSBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," No of nonzeros in factor U = %D\n", Ustore->nnz); 2295a46d3faSBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," No of nonzeros in L+U = %D\n", Lstore->nnz + Ustore->nnz - lu->A.ncol); 2306da386baSHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," L\\U MB %.3f\ttotal MB needed %.3f\n", 2316da386baSHong Zhang lu->mem_usage.for_lu/1e6, lu->mem_usage.total_needed/1e6); 2325a46d3faSBarry Smith } 2335a46d3faSBarry Smith 2345a46d3faSBarry Smith lu->flg = SAME_NONZERO_PATTERN; 2351d5ca7f3SHong Zhang F->ops->solve = MatSolve_SuperLU; 2361d5ca7f3SHong Zhang F->ops->solvetranspose = MatSolveTranspose_SuperLU; 237*1aef8b4cSStefano Zampini /* Better to use MatMatSolve_Basic instead of raising an error */ 238*1aef8b4cSStefano Zampini /* F->ops->matsolve = MatMatSolve_SuperLU; */ 239*1aef8b4cSStefano Zampini F->ops->matsolve = NULL; 2405a46d3faSBarry Smith PetscFunctionReturn(0); 2415a46d3faSBarry Smith } 2425a46d3faSBarry Smith 24314b396bbSKris Buschelman #undef __FUNCT__ 24420be8e61SHong Zhang #define __FUNCT__ "MatGetDiagonal_SuperLU" 24520be8e61SHong Zhang PetscErrorCode MatGetDiagonal_SuperLU(Mat A,Vec v) 24620be8e61SHong Zhang { 24720be8e61SHong Zhang PetscFunctionBegin; 24820be8e61SHong Zhang SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Mat type: SuperLU factor"); 24920be8e61SHong Zhang PetscFunctionReturn(0); 25020be8e61SHong Zhang } 25120be8e61SHong Zhang 25220be8e61SHong Zhang #undef __FUNCT__ 253f0c56d0fSKris Buschelman #define __FUNCT__ "MatDestroy_SuperLU" 254dfbe8321SBarry Smith PetscErrorCode MatDestroy_SuperLU(Mat A) 25514b396bbSKris Buschelman { 256dfbe8321SBarry Smith PetscErrorCode ierr; 257f0c56d0fSKris Buschelman Mat_SuperLU *lu=(Mat_SuperLU*)A->spptr; 25814b396bbSKris Buschelman 25914b396bbSKris Buschelman PetscFunctionBegin; 260bf0cc555SLisandro Dalcin if (lu && lu->CleanUpSuperLU) { /* Free the SuperLU datastructures */ 261d387c056SBarry Smith PetscStackCall("SuperLU:Destroy_SuperMatrix_Store",Destroy_SuperMatrix_Store(&lu->A)); 262d387c056SBarry Smith PetscStackCall("SuperLU:Destroy_SuperMatrix_Store",Destroy_SuperMatrix_Store(&lu->B)); 263d387c056SBarry Smith PetscStackCall("SuperLU:Destroy_SuperMatrix_Store",Destroy_SuperMatrix_Store(&lu->X)); 264d387c056SBarry Smith PetscStackCall("SuperLU:StatFree",StatFree(&lu->stat)); 2650e742b69SHong Zhang if (lu->lwork >= 0) { 266d387c056SBarry Smith PetscStackCall("SuperLU:Destroy_SuperNode_Matrix",Destroy_SuperNode_Matrix(&lu->L)); 267d387c056SBarry Smith PetscStackCall("SuperLU:Destroy_CompCol_Matrix",Destroy_CompCol_Matrix(&lu->U)); 2680e742b69SHong Zhang } 2690e742b69SHong Zhang } 270bf0cc555SLisandro Dalcin if (lu) { 2719ce50919SHong Zhang ierr = PetscFree(lu->etree);CHKERRQ(ierr); 2729ce50919SHong Zhang ierr = PetscFree(lu->perm_r);CHKERRQ(ierr); 2739ce50919SHong Zhang ierr = PetscFree(lu->perm_c);CHKERRQ(ierr); 2749ce50919SHong Zhang ierr = PetscFree(lu->R);CHKERRQ(ierr); 2759ce50919SHong Zhang ierr = PetscFree(lu->C);CHKERRQ(ierr); 276bf0cc555SLisandro Dalcin ierr = PetscFree(lu->rhs_dup);CHKERRQ(ierr); 277bf0cc555SLisandro Dalcin ierr = MatDestroy(&lu->A_dup);CHKERRQ(ierr); 278bf0cc555SLisandro Dalcin } 279bf0cc555SLisandro Dalcin ierr = PetscFree(A->spptr);CHKERRQ(ierr); 2800e742b69SHong Zhang 281d954db57SHong Zhang /* clear composed functions */ 282bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorGetSolverPackage_C",NULL);CHKERRQ(ierr); 283bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatSuperluSetILUDropTol_C",NULL);CHKERRQ(ierr); 284d954db57SHong Zhang 285b24902e0SBarry Smith ierr = MatDestroy_SeqAIJ(A);CHKERRQ(ierr); 28614b396bbSKris Buschelman PetscFunctionReturn(0); 28714b396bbSKris Buschelman } 28814b396bbSKris Buschelman 28914b396bbSKris Buschelman #undef __FUNCT__ 290f0c56d0fSKris Buschelman #define __FUNCT__ "MatView_SuperLU" 291dfbe8321SBarry Smith PetscErrorCode MatView_SuperLU(Mat A,PetscViewer viewer) 29214b396bbSKris Buschelman { 293dfbe8321SBarry Smith PetscErrorCode ierr; 294ace3abfcSBarry Smith PetscBool iascii; 29514b396bbSKris Buschelman PetscViewerFormat format; 29614b396bbSKris Buschelman 29714b396bbSKris Buschelman PetscFunctionBegin; 298251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 29932077d6dSBarry Smith if (iascii) { 30014b396bbSKris Buschelman ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 3012f59403fSHong Zhang if (format == PETSC_VIEWER_ASCII_INFO) { 302f0c56d0fSKris Buschelman ierr = MatFactorInfo_SuperLU(A,viewer);CHKERRQ(ierr); 30314b396bbSKris Buschelman } 30414b396bbSKris Buschelman } 30514b396bbSKris Buschelman PetscFunctionReturn(0); 30614b396bbSKris Buschelman } 30714b396bbSKris Buschelman 30814b396bbSKris Buschelman 30914b396bbSKris Buschelman #undef __FUNCT__ 310c29e884eSHong Zhang #define __FUNCT__ "MatSolve_SuperLU_Private" 311c29e884eSHong Zhang PetscErrorCode MatSolve_SuperLU_Private(Mat A,Vec b,Vec x) 31214b396bbSKris Buschelman { 313f0c56d0fSKris Buschelman Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr; 314d9ca1df4SBarry Smith const PetscScalar *barray; 315d9ca1df4SBarry Smith PetscScalar *xarray; 316dfbe8321SBarry Smith PetscErrorCode ierr; 317077aedafSJed Brown PetscInt info,i,n; 318da7a1d00SHong Zhang PetscReal ferr,berr; 319dff31646SBarry Smith static PetscBool cite = PETSC_FALSE; 32014b396bbSKris Buschelman 32114b396bbSKris Buschelman PetscFunctionBegin; 3222205254eSKarl Rupp if (lu->lwork == -1) PetscFunctionReturn(0); 323dff31646SBarry Smith ierr = PetscCitationsRegister("@article{superlu99,\n author = {James W. Demmel and Stanley C. Eisenstat and\n John R. Gilbert and Xiaoye S. Li and Joseph W. H. Liu},\n title = {A supernodal approach to sparse partial pivoting},\n journal = {SIAM J. Matrix Analysis and Applications},\n year = {1999},\n volume = {20},\n number = {3},\n pages = {720-755}\n}\n",&cite);CHKERRQ(ierr); 324cae5a23dSHong Zhang 325077aedafSJed Brown ierr = VecGetLocalSize(x,&n);CHKERRQ(ierr); 32675af56d4SHong Zhang lu->B.ncol = 1; /* Set the number of right-hand side */ 327cae5a23dSHong Zhang if (lu->options.Equil && !lu->rhs_dup) { 328cae5a23dSHong Zhang /* superlu overwrites b when Equil is used, thus create rhs_dup to keep user's b unchanged */ 329785e854fSJed Brown ierr = PetscMalloc1(n,&lu->rhs_dup);CHKERRQ(ierr); 330cae5a23dSHong Zhang } 331cae5a23dSHong Zhang if (lu->options.Equil) { 332cae5a23dSHong Zhang /* Copy b into rsh_dup */ 333d9ca1df4SBarry Smith ierr = VecGetArrayRead(b,&barray);CHKERRQ(ierr); 334cae5a23dSHong Zhang ierr = PetscMemcpy(lu->rhs_dup,barray,n*sizeof(PetscScalar));CHKERRQ(ierr); 335d9ca1df4SBarry Smith ierr = VecRestoreArrayRead(b,&barray);CHKERRQ(ierr); 336cae5a23dSHong Zhang barray = lu->rhs_dup; 337cae5a23dSHong Zhang } else { 338d9ca1df4SBarry Smith ierr = VecGetArrayRead(b,&barray);CHKERRQ(ierr); 339cae5a23dSHong Zhang } 34075af56d4SHong Zhang ierr = VecGetArray(x,&xarray);CHKERRQ(ierr); 3415fe6150dSHong Zhang 3425fe6150dSHong Zhang #if defined(PETSC_USE_COMPLEX) 3433cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 3443cb270beSHong Zhang ((DNformat*)lu->B.Store)->nzval = (singlecomplex*)barray; 3453cb270beSHong Zhang ((DNformat*)lu->X.Store)->nzval = (singlecomplex*)xarray; 3463cb270beSHong Zhang #else 3475fe6150dSHong Zhang ((DNformat*)lu->B.Store)->nzval = (doublecomplex*)barray; 3485fe6150dSHong Zhang ((DNformat*)lu->X.Store)->nzval = (doublecomplex*)xarray; 3493cb270beSHong Zhang #endif 3505fe6150dSHong Zhang #else 351d9ca1df4SBarry Smith ((DNformat*)lu->B.Store)->nzval = (void*)barray; 35275af56d4SHong Zhang ((DNformat*)lu->X.Store)->nzval = xarray; 3535fe6150dSHong Zhang #endif 35475af56d4SHong Zhang 35575af56d4SHong Zhang lu->options.Fact = FACTORED; /* Indicate the factored form of A is supplied. */ 356d5f3da31SBarry Smith if (A->factortype == MAT_FACTOR_LU) { 3578914a3f7SHong Zhang #if defined(PETSC_USE_COMPLEX) 3583cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 359d387c056SBarry Smith PetscStackCall("SuperLU:cgssvx",cgssvx(&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, 361d387c056SBarry Smith &lu->mem_usage, &lu->stat, &info)); 3623cb270beSHong Zhang #else 363d387c056SBarry Smith PetscStackCall("SuperLU:zgssvx",zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 3648914a3f7SHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 365d387c056SBarry Smith &lu->mem_usage, &lu->stat, &info)); 3663cb270beSHong Zhang #endif 3673cb270beSHong Zhang #else 3683cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 369d387c056SBarry Smith PetscStackCall("SuperLU:sgssvx",sgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 3703cb270beSHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 371d387c056SBarry Smith &lu->mem_usage, &lu->stat, &info)); 3728914a3f7SHong Zhang #else 373d387c056SBarry Smith PetscStackCall("SuperLU:dgssvx",dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 37475af56d4SHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 375d387c056SBarry Smith &lu->mem_usage, &lu->stat, &info)); 3768914a3f7SHong Zhang #endif 3773cb270beSHong Zhang #endif 378d5f3da31SBarry Smith } else if (A->factortype == MAT_FACTOR_ILU) { 379cffbb591SHong Zhang #if defined(PETSC_USE_COMPLEX) 3803cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 381d387c056SBarry Smith PetscStackCall("SuperLU:cgsisx",cgsisx(&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, 383d387c056SBarry Smith &lu->mem_usage, &lu->stat, &info)); 3843cb270beSHong Zhang #else 385d387c056SBarry Smith PetscStackCall("SuperLU:zgsisx",zgsisx(&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, 387d387c056SBarry Smith &lu->mem_usage, &lu->stat, &info)); 3883cb270beSHong Zhang #endif 3893cb270beSHong Zhang #else 3903cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 391d387c056SBarry Smith PetscStackCall("SuperLU:sgsisx",sgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 3923cb270beSHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, 393d387c056SBarry Smith &lu->mem_usage, &lu->stat, &info)); 394cffbb591SHong Zhang #else 395d387c056SBarry Smith PetscStackCall("SuperLU:dgsisx",dgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 396cffbb591SHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, 397d387c056SBarry Smith &lu->mem_usage, &lu->stat, &info)); 398cffbb591SHong Zhang #endif 3993cb270beSHong Zhang #endif 400f23aa3ddSBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported"); 401cae5a23dSHong Zhang if (!lu->options.Equil) { 402d9ca1df4SBarry Smith ierr = VecRestoreArrayRead(b,&barray);CHKERRQ(ierr); 403cae5a23dSHong Zhang } 40475af56d4SHong Zhang ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr); 40575af56d4SHong Zhang 406958c9bccSBarry Smith if (!info || info == lu->A.ncol+1) { 40775af56d4SHong Zhang if (lu->options.IterRefine) { 4088914a3f7SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF,"Iterative Refinement:\n"); 4098914a3f7SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," %8s%8s%16s%16s\n", "rhs", "Steps", "FERR", "BERR"); 4102205254eSKarl Rupp for (i = 0; i < 1; ++i) { 4115d8b2efaSHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," %8d%8d%16e%16e\n", i+1, lu->stat.RefineSteps, ferr, berr); 41275af56d4SHong Zhang } 4132205254eSKarl Rupp } 4148914a3f7SHong Zhang } else if (info > 0) { 4158914a3f7SHong Zhang if (lu->lwork == -1) { 41677431f27SBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," ** Estimated memory: %D bytes\n", info - lu->A.ncol); 4178914a3f7SHong Zhang } else { 41877431f27SBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," Warning: gssvx() returns info %D\n",info); 4198914a3f7SHong Zhang } 420f23aa3ddSBarry 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); 42114b396bbSKris Buschelman 4228914a3f7SHong Zhang if (lu->options.PrintStat) { 4238914a3f7SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF,"MatSolve__SuperLU():\n"); 424d387c056SBarry Smith PetscStackCall("SuperLU:StatPrint",StatPrint(&lu->stat)); 4258914a3f7SHong Zhang } 42675af56d4SHong Zhang PetscFunctionReturn(0); 42775af56d4SHong Zhang } 42814b396bbSKris Buschelman 42914b396bbSKris Buschelman #undef __FUNCT__ 430c29e884eSHong Zhang #define __FUNCT__ "MatSolve_SuperLU" 431c29e884eSHong Zhang PetscErrorCode MatSolve_SuperLU(Mat A,Vec b,Vec x) 432c29e884eSHong Zhang { 433c29e884eSHong Zhang Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr; 434c29e884eSHong Zhang PetscErrorCode ierr; 435c29e884eSHong Zhang 436c29e884eSHong Zhang PetscFunctionBegin; 437c29e884eSHong Zhang lu->options.Trans = TRANS; 4382205254eSKarl Rupp 439c29e884eSHong Zhang ierr = MatSolve_SuperLU_Private(A,b,x);CHKERRQ(ierr); 440c29e884eSHong Zhang PetscFunctionReturn(0); 441c29e884eSHong Zhang } 442c29e884eSHong Zhang 443c29e884eSHong Zhang #undef __FUNCT__ 444c7c1fe80SHong Zhang #define __FUNCT__ "MatSolveTranspose_SuperLU" 445c7c1fe80SHong Zhang PetscErrorCode MatSolveTranspose_SuperLU(Mat A,Vec b,Vec x) 446c7c1fe80SHong Zhang { 447c7c1fe80SHong Zhang Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr; 448c7c1fe80SHong Zhang PetscErrorCode ierr; 449c7c1fe80SHong Zhang 450c7c1fe80SHong Zhang PetscFunctionBegin; 451c7c1fe80SHong Zhang lu->options.Trans = NOTRANS; 4522205254eSKarl Rupp 453c29e884eSHong Zhang ierr = MatSolve_SuperLU_Private(A,b,x);CHKERRQ(ierr); 454c7c1fe80SHong Zhang PetscFunctionReturn(0); 455c7c1fe80SHong Zhang } 456c7c1fe80SHong Zhang 457e0b74bf9SHong Zhang #undef __FUNCT__ 458e0b74bf9SHong Zhang #define __FUNCT__ "MatMatSolve_SuperLU" 459e0b74bf9SHong Zhang PetscErrorCode MatMatSolve_SuperLU(Mat A,Mat B,Mat X) 460e0b74bf9SHong Zhang { 461e0b74bf9SHong Zhang Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr; 462cd723cd1SBarry Smith PetscBool flg; 463cd723cd1SBarry Smith PetscErrorCode ierr; 464e0b74bf9SHong Zhang 465e0b74bf9SHong Zhang PetscFunctionBegin; 4660298fd71SBarry Smith ierr = PetscObjectTypeCompareAny((PetscObject)B,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr); 467ce94432eSBarry Smith if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix"); 4680298fd71SBarry Smith ierr = PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr); 469ce94432eSBarry Smith if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix"); 4702205254eSKarl Rupp lu->options.Trans = TRANS; 471e0b74bf9SHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatMatSolve_SuperLU() is not implemented yet"); 472e0b74bf9SHong Zhang PetscFunctionReturn(0); 473e0b74bf9SHong Zhang } 474e0b74bf9SHong Zhang 47514b396bbSKris Buschelman /* 47614b396bbSKris Buschelman Note the r permutation is ignored 47714b396bbSKris Buschelman */ 47814b396bbSKris Buschelman #undef __FUNCT__ 479f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_SuperLU" 4800481f469SBarry Smith PetscErrorCode MatLUFactorSymbolic_SuperLU(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info) 48114b396bbSKris Buschelman { 4821d5ca7f3SHong Zhang Mat_SuperLU *lu = (Mat_SuperLU*)(F->spptr); 483b24902e0SBarry Smith 484b24902e0SBarry Smith PetscFunctionBegin; 485b24902e0SBarry Smith lu->flg = DIFFERENT_NONZERO_PATTERN; 486b24902e0SBarry Smith lu->CleanUpSuperLU = PETSC_TRUE; 4871d5ca7f3SHong Zhang F->ops->lufactornumeric = MatLUFactorNumeric_SuperLU; 488b24902e0SBarry Smith PetscFunctionReturn(0); 489b24902e0SBarry Smith } 490b24902e0SBarry Smith 49135bd34faSBarry Smith #undef __FUNCT__ 4925ccb76cbSHong Zhang #define __FUNCT__ "MatSuperluSetILUDropTol_SuperLU" 493b2573a8aSBarry Smith static PetscErrorCode MatSuperluSetILUDropTol_SuperLU(Mat F,PetscReal dtol) 4945ccb76cbSHong Zhang { 4955ccb76cbSHong Zhang Mat_SuperLU *lu= (Mat_SuperLU*)F->spptr; 4965ccb76cbSHong Zhang 4975ccb76cbSHong Zhang PetscFunctionBegin; 4985ccb76cbSHong Zhang lu->options.ILU_DropTol = dtol; 4995ccb76cbSHong Zhang PetscFunctionReturn(0); 5005ccb76cbSHong Zhang } 5015ccb76cbSHong Zhang 5025ccb76cbSHong Zhang #undef __FUNCT__ 5035ccb76cbSHong Zhang #define __FUNCT__ "MatSuperluSetILUDropTol" 5045ccb76cbSHong Zhang /*@ 5055ccb76cbSHong Zhang MatSuperluSetILUDropTol - Set SuperLU ILU drop tolerance 5065ccb76cbSHong Zhang Logically Collective on Mat 5075ccb76cbSHong Zhang 5085ccb76cbSHong Zhang Input Parameters: 5095ccb76cbSHong Zhang + F - the factored matrix obtained by calling MatGetFactor() from PETSc-SuperLU interface 5105ccb76cbSHong Zhang - dtol - drop tolerance 5115ccb76cbSHong Zhang 5125ccb76cbSHong Zhang Options Database: 5135ccb76cbSHong Zhang . -mat_superlu_ilu_droptol <dtol> 5145ccb76cbSHong Zhang 5155ccb76cbSHong Zhang Level: beginner 5165ccb76cbSHong Zhang 5175ccb76cbSHong Zhang References: SuperLU Users' Guide 5185ccb76cbSHong Zhang 5195ccb76cbSHong Zhang .seealso: MatGetFactor() 5205ccb76cbSHong Zhang @*/ 5215ccb76cbSHong Zhang PetscErrorCode MatSuperluSetILUDropTol(Mat F,PetscReal dtol) 5225ccb76cbSHong Zhang { 5235ccb76cbSHong Zhang PetscErrorCode ierr; 5245ccb76cbSHong Zhang 5255ccb76cbSHong Zhang PetscFunctionBegin; 5265ccb76cbSHong Zhang PetscValidHeaderSpecific(F,MAT_CLASSID,1); 5275ccb76cbSHong Zhang PetscValidLogicalCollectiveInt(F,dtol,2); 5285ccb76cbSHong Zhang ierr = PetscTryMethod(F,"MatSuperluSetILUDropTol_C",(Mat,PetscReal),(F,dtol));CHKERRQ(ierr); 5295ccb76cbSHong Zhang PetscFunctionReturn(0); 5305ccb76cbSHong Zhang } 5315ccb76cbSHong Zhang 5325ccb76cbSHong Zhang #undef __FUNCT__ 53335bd34faSBarry Smith #define __FUNCT__ "MatFactorGetSolverPackage_seqaij_superlu" 53435bd34faSBarry Smith PetscErrorCode MatFactorGetSolverPackage_seqaij_superlu(Mat A,const MatSolverPackage *type) 53535bd34faSBarry Smith { 53635bd34faSBarry Smith PetscFunctionBegin; 5372692d6eeSBarry Smith *type = MATSOLVERSUPERLU; 53835bd34faSBarry Smith PetscFunctionReturn(0); 53935bd34faSBarry Smith } 54035bd34faSBarry Smith 541b24902e0SBarry Smith 542b24902e0SBarry Smith /*MC 543ba992d64SSatish Balay MATSOLVERSUPERLU = "superlu" - A solver package providing solvers LU and ILU for sequential matrices 544b24902e0SBarry Smith via the external package SuperLU. 545b24902e0SBarry Smith 546e2e64c6bSBarry Smith Use ./configure --download-superlu to have PETSc installed with SuperLU 547b24902e0SBarry Smith 548c2b89b5dSBarry Smith Use -pc_type lu -pc_factor_mat_solver_package superlu to us this direct solver 549c2b89b5dSBarry Smith 550b24902e0SBarry Smith Options Database Keys: 551e08999f5SMatthew G Knepley + -mat_superlu_equil <FALSE> - Equil (None) 552e08999f5SMatthew G Knepley . -mat_superlu_colperm <COLAMD> - (choose one of) NATURAL MMD_ATA MMD_AT_PLUS_A COLAMD 553e08999f5SMatthew G Knepley . -mat_superlu_iterrefine <NOREFINE> - (choose one of) NOREFINE SINGLE DOUBLE EXTRA 554e08999f5SMatthew G Knepley . -mat_superlu_symmetricmode: <FALSE> - SymmetricMode (None) 555e08999f5SMatthew G Knepley . -mat_superlu_diagpivotthresh <1> - DiagPivotThresh (None) 556e08999f5SMatthew G Knepley . -mat_superlu_pivotgrowth <FALSE> - PivotGrowth (None) 557e08999f5SMatthew G Knepley . -mat_superlu_conditionnumber <FALSE> - ConditionNumber (None) 558e08999f5SMatthew G Knepley . -mat_superlu_rowperm <NOROWPERM> - (choose one of) NOROWPERM LargeDiag 559e08999f5SMatthew G Knepley . -mat_superlu_replacetinypivot <FALSE> - ReplaceTinyPivot (None) 560e08999f5SMatthew G Knepley . -mat_superlu_printstat <FALSE> - PrintStat (None) 561e08999f5SMatthew G Knepley . -mat_superlu_lwork <0> - size of work array in bytes used by factorization (None) 562e08999f5SMatthew G Knepley . -mat_superlu_ilu_droptol <0> - ILU_DropTol (None) 563e08999f5SMatthew G Knepley . -mat_superlu_ilu_filltol <0> - ILU_FillTol (None) 564e08999f5SMatthew G Knepley . -mat_superlu_ilu_fillfactor <0> - ILU_FillFactor (None) 565e08999f5SMatthew G Knepley . -mat_superlu_ilu_droprull <0> - ILU_DropRule (None) 566e08999f5SMatthew G Knepley . -mat_superlu_ilu_norm <0> - ILU_Norm (None) 567e08999f5SMatthew G Knepley - -mat_superlu_ilu_milu <0> - ILU_MILU (None) 568b24902e0SBarry Smith 5692692d6eeSBarry Smith Notes: Do not confuse this with MATSOLVERSUPERLU_DIST which is for parallel sparse solves 5705c9eb25fSBarry Smith 571b24902e0SBarry Smith Level: beginner 572b24902e0SBarry Smith 573d45987f3SHong Zhang .seealso: PCLU, PCILU, MATSOLVERSUPERLU_DIST, MATSOLVERMUMPS, PCFactorSetMatSolverPackage(), MatSolverPackage 574b24902e0SBarry Smith M*/ 575b24902e0SBarry Smith 576b24902e0SBarry Smith #undef __FUNCT__ 577b24902e0SBarry Smith #define __FUNCT__ "MatGetFactor_seqaij_superlu" 5788cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactor_seqaij_superlu(Mat A,MatFactorType ftype,Mat *F) 579b24902e0SBarry Smith { 58014b396bbSKris Buschelman Mat B; 581f0c56d0fSKris Buschelman Mat_SuperLU *lu; 5826849ba73SBarry Smith PetscErrorCode ierr; 5835d8b2efaSHong Zhang PetscInt indx,m=A->rmap->n,n=A->cmap->n; 5848afaa268SBarry Smith PetscBool flg,set; 5853cb270beSHong Zhang PetscReal real_input; 5865ca28756SHong Zhang const char *colperm[] ={"NATURAL","MMD_ATA","MMD_AT_PLUS_A","COLAMD"}; /* MY_PERMC - not supported by the petsc interface yet */ 5875ca28756SHong Zhang const char *iterrefine[]={"NOREFINE", "SINGLE", "DOUBLE", "EXTRA"}; 5885ca28756SHong Zhang const char *rowperm[] ={"NOROWPERM", "LargeDiag"}; /* MY_PERMC - not supported by the petsc interface yet */ 58914b396bbSKris Buschelman 59014b396bbSKris Buschelman PetscFunctionBegin; 591ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr); 592d0f46423SBarry Smith ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 5937adad957SLisandro Dalcin ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 5940298fd71SBarry Smith ierr = MatSeqAIJSetPreallocation(B,0,NULL);CHKERRQ(ierr); 595f0c56d0fSKris Buschelman 596cffbb591SHong Zhang if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU) { 597b24902e0SBarry Smith B->ops->lufactorsymbolic = MatLUFactorSymbolic_SuperLU; 598cffbb591SHong Zhang B->ops->ilufactorsymbolic = MatLUFactorSymbolic_SuperLU; 599b3a44c85SBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported"); 600cffbb591SHong Zhang 601b24902e0SBarry Smith B->ops->destroy = MatDestroy_SuperLU; 6023519fcdcSHong Zhang B->ops->view = MatView_SuperLU; 60320be8e61SHong Zhang B->ops->getdiagonal = MatGetDiagonal_SuperLU; 604d5f3da31SBarry Smith B->factortype = ftype; 60594b7f48cSBarry Smith B->assembled = PETSC_TRUE; /* required by -ksp_view */ 6065c9eb25fSBarry Smith B->preallocated = PETSC_TRUE; 60714b396bbSKris Buschelman 608b00a9115SJed Brown ierr = PetscNewLog(B,&lu);CHKERRQ(ierr); 609cae5a23dSHong Zhang 610cffbb591SHong Zhang if (ftype == MAT_FACTOR_LU) { 6119ce50919SHong Zhang set_default_options(&lu->options); 6123d6581fbSHong Zhang /* Comments from SuperLU_4.0/SRC/dgssvx.c: 6133d6581fbSHong Zhang "Whether or not the system will be equilibrated depends on the 6143d6581fbSHong Zhang scaling of the matrix A, but if equilibration is used, A is 6153d6581fbSHong Zhang overwritten by diag(R)*A*diag(C) and B by diag(R)*B 6163d6581fbSHong Zhang (if options->Trans=NOTRANS) or diag(C)*B (if options->Trans = TRANS or CONJ)." 6173d6581fbSHong Zhang We set 'options.Equil = NO' as default because additional space is needed for it. 6183d6581fbSHong Zhang */ 6193d6581fbSHong Zhang lu->options.Equil = NO; 620cffbb591SHong Zhang } else if (ftype == MAT_FACTOR_ILU) { 6210c74a584SJed Brown /* Set the default input options of ilu: */ 622d387c056SBarry Smith PetscStackCall("SuperLU:ilu_set_default_options",ilu_set_default_options(&lu->options)); 623cffbb591SHong Zhang } 6249ce50919SHong Zhang lu->options.PrintStat = NO; 6251a2e2f44SHong Zhang 6265d8b2efaSHong Zhang /* Initialize the statistics variables. */ 627d387c056SBarry Smith PetscStackCall("SuperLU:StatInit",StatInit(&lu->stat)); 6288914a3f7SHong Zhang lu->lwork = 0; /* allocate space internally by system malloc */ 6299ce50919SHong Zhang 630ce94432eSBarry Smith ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"SuperLU Options","Mat");CHKERRQ(ierr); 6318afaa268SBarry Smith ierr = PetscOptionsBool("-mat_superlu_equil","Equil","None",(PetscBool)lu->options.Equil,(PetscBool*)&lu->options.Equil,NULL);CHKERRQ(ierr); 6328914a3f7SHong Zhang ierr = PetscOptionsEList("-mat_superlu_colperm","ColPerm","None",colperm,4,colperm[3],&indx,&flg);CHKERRQ(ierr); 6332205254eSKarl Rupp if (flg) lu->options.ColPerm = (colperm_t)indx; 6348914a3f7SHong Zhang ierr = PetscOptionsEList("-mat_superlu_iterrefine","IterRefine","None",iterrefine,4,iterrefine[0],&indx,&flg);CHKERRQ(ierr); 6352205254eSKarl Rupp if (flg) lu->options.IterRefine = (IterRefine_t)indx; 6368afaa268SBarry Smith ierr = PetscOptionsBool("-mat_superlu_symmetricmode","SymmetricMode","None",(PetscBool)lu->options.SymmetricMode,&flg,&set);CHKERRQ(ierr); 6378afaa268SBarry Smith if (set && flg) lu->options.SymmetricMode = YES; 6383cb270beSHong Zhang ierr = PetscOptionsReal("-mat_superlu_diagpivotthresh","DiagPivotThresh","None",lu->options.DiagPivotThresh,&real_input,&flg);CHKERRQ(ierr); 6393cb270beSHong Zhang if (flg) lu->options.DiagPivotThresh = (double) real_input; 6408afaa268SBarry Smith ierr = PetscOptionsBool("-mat_superlu_pivotgrowth","PivotGrowth","None",(PetscBool)lu->options.PivotGrowth,&flg,&set);CHKERRQ(ierr); 6418afaa268SBarry Smith if (set && flg) lu->options.PivotGrowth = YES; 6428afaa268SBarry Smith ierr = PetscOptionsBool("-mat_superlu_conditionnumber","ConditionNumber","None",(PetscBool)lu->options.ConditionNumber,&flg,&set);CHKERRQ(ierr); 6438afaa268SBarry Smith if (set && flg) lu->options.ConditionNumber = YES; 644d7ebd59bSHong Zhang ierr = PetscOptionsEList("-mat_superlu_rowperm","rowperm","None",rowperm,2,rowperm[lu->options.RowPerm],&indx,&flg);CHKERRQ(ierr); 6452205254eSKarl Rupp if (flg) lu->options.RowPerm = (rowperm_t)indx; 6468afaa268SBarry Smith ierr = PetscOptionsBool("-mat_superlu_replacetinypivot","ReplaceTinyPivot","None",(PetscBool)lu->options.ReplaceTinyPivot,&flg,&set);CHKERRQ(ierr); 6478afaa268SBarry Smith if (set && flg) lu->options.ReplaceTinyPivot = YES; 6488afaa268SBarry Smith ierr = PetscOptionsBool("-mat_superlu_printstat","PrintStat","None",(PetscBool)lu->options.PrintStat,&flg,&set);CHKERRQ(ierr); 6498afaa268SBarry Smith if (set && flg) lu->options.PrintStat = YES; 6500298fd71SBarry Smith ierr = PetscOptionsInt("-mat_superlu_lwork","size of work array in bytes used by factorization","None",lu->lwork,&lu->lwork,NULL);CHKERRQ(ierr); 6515fe6150dSHong Zhang if (lu->lwork > 0) { 652d87de817SBarry Smith /* lwork is in bytes, hence PetscMalloc() is used here, not PetscMalloc1()*/ 6535fe6150dSHong Zhang ierr = PetscMalloc(lu->lwork,&lu->work);CHKERRQ(ierr); 6545fe6150dSHong Zhang } else if (lu->lwork != 0 && lu->lwork != -1) { 65577431f27SBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," Warning: lwork %D is not supported by SUPERLU. The default lwork=0 is used.\n",lu->lwork); 6568914a3f7SHong Zhang lu->lwork = 0; 6578914a3f7SHong Zhang } 658cffbb591SHong Zhang /* ilu options */ 6593cb270beSHong Zhang ierr = PetscOptionsReal("-mat_superlu_ilu_droptol","ILU_DropTol","None",lu->options.ILU_DropTol,&real_input,&flg);CHKERRQ(ierr); 6603cb270beSHong Zhang if (flg) lu->options.ILU_DropTol = (double) real_input; 6613cb270beSHong Zhang ierr = PetscOptionsReal("-mat_superlu_ilu_filltol","ILU_FillTol","None",lu->options.ILU_FillTol,&real_input,&flg);CHKERRQ(ierr); 6623cb270beSHong Zhang if (flg) lu->options.ILU_FillTol = (double) real_input; 6633cb270beSHong Zhang ierr = PetscOptionsReal("-mat_superlu_ilu_fillfactor","ILU_FillFactor","None",lu->options.ILU_FillFactor,&real_input,&flg);CHKERRQ(ierr); 6643cb270beSHong Zhang if (flg) lu->options.ILU_FillFactor = (double) real_input; 6650298fd71SBarry Smith ierr = PetscOptionsInt("-mat_superlu_ilu_droprull","ILU_DropRule","None",lu->options.ILU_DropRule,&lu->options.ILU_DropRule,NULL);CHKERRQ(ierr); 666cffbb591SHong Zhang ierr = PetscOptionsInt("-mat_superlu_ilu_norm","ILU_Norm","None",lu->options.ILU_Norm,&indx,&flg);CHKERRQ(ierr); 6672205254eSKarl Rupp if (flg) lu->options.ILU_Norm = (norm_t)indx; 668cffbb591SHong Zhang ierr = PetscOptionsInt("-mat_superlu_ilu_milu","ILU_MILU","None",lu->options.ILU_MILU,&indx,&flg);CHKERRQ(ierr); 6692205254eSKarl Rupp if (flg) lu->options.ILU_MILU = (milu_t)indx; 6709ce50919SHong Zhang PetscOptionsEnd(); 67138abfdaeSHong Zhang if (lu->options.Equil == YES) { 67238abfdaeSHong Zhang /* superlu overwrites input matrix and rhs when Equil is used, thus create A_dup to keep user's A unchanged */ 67338abfdaeSHong Zhang ierr = MatDuplicate_SeqAIJ(A,MAT_COPY_VALUES,&lu->A_dup);CHKERRQ(ierr); 67438abfdaeSHong Zhang } 6759ce50919SHong Zhang 6765d8b2efaSHong Zhang /* Allocate spaces (notice sizes are for the transpose) */ 677785e854fSJed Brown ierr = PetscMalloc1(m,&lu->etree);CHKERRQ(ierr); 678785e854fSJed Brown ierr = PetscMalloc1(n,&lu->perm_r);CHKERRQ(ierr); 679785e854fSJed Brown ierr = PetscMalloc1(m,&lu->perm_c);CHKERRQ(ierr); 680785e854fSJed Brown ierr = PetscMalloc1(n,&lu->R);CHKERRQ(ierr); 681785e854fSJed Brown ierr = PetscMalloc1(m,&lu->C);CHKERRQ(ierr); 6825d8b2efaSHong Zhang 6835d8b2efaSHong Zhang /* create rhs and solution x without allocate space for .Store */ 6845d8b2efaSHong Zhang #if defined(PETSC_USE_COMPLEX) 6853cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 686d387c056SBarry Smith PetscStackCall("SuperLU:cCreate_Dense_Matrix(",cCreate_Dense_Matrix(&lu->B, m, 1, NULL, m, SLU_DN, SLU_C, SLU_GE)); 687d387c056SBarry Smith PetscStackCall("SuperLU:cCreate_Dense_Matrix(",cCreate_Dense_Matrix(&lu->X, m, 1, NULL, m, SLU_DN, SLU_C, SLU_GE)); 6883cb270beSHong Zhang #else 689d387c056SBarry Smith PetscStackCall("SuperLU:zCreate_Dense_Matrix",zCreate_Dense_Matrix(&lu->B, m, 1, NULL, m, SLU_DN, SLU_Z, SLU_GE)); 690d387c056SBarry Smith PetscStackCall("SuperLU:zCreate_Dense_Matrix",zCreate_Dense_Matrix(&lu->X, m, 1, NULL, m, SLU_DN, SLU_Z, SLU_GE)); 6913cb270beSHong Zhang #endif 6923cb270beSHong Zhang #else 6933cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 694d387c056SBarry Smith PetscStackCall("SuperLU:sCreate_Dense_Matrix",sCreate_Dense_Matrix(&lu->B, m, 1, NULL, m, SLU_DN, SLU_S, SLU_GE)); 695d387c056SBarry Smith PetscStackCall("SuperLU:sCreate_Dense_Matrix",sCreate_Dense_Matrix(&lu->X, m, 1, NULL, m, SLU_DN, SLU_S, SLU_GE)); 6965d8b2efaSHong Zhang #else 697d387c056SBarry Smith PetscStackCall("SuperLU:dCreate_Dense_Matrix",dCreate_Dense_Matrix(&lu->B, m, 1, NULL, m, SLU_DN, SLU_D, SLU_GE)); 698d387c056SBarry Smith PetscStackCall("SuperLU:dCreate_Dense_Matrix",dCreate_Dense_Matrix(&lu->X, m, 1, NULL, m, SLU_DN, SLU_D, SLU_GE)); 6995d8b2efaSHong Zhang #endif 7003cb270beSHong Zhang #endif 7015d8b2efaSHong Zhang 702bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_seqaij_superlu);CHKERRQ(ierr); 703bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatSuperluSetILUDropTol_C",MatSuperluSetILUDropTol_SuperLU);CHKERRQ(ierr); 7045c9eb25fSBarry Smith B->spptr = lu; 70520be8e61SHong Zhang 706899d7b4fSKris Buschelman *F = B; 70714b396bbSKris Buschelman PetscFunctionReturn(0); 70814b396bbSKris Buschelman } 709d954db57SHong Zhang 71042c9c57cSBarry Smith #undef __FUNCT__ 71142c9c57cSBarry Smith #define __FUNCT__ "MatSolverPackageRegister_SuperLU" 71229b38603SBarry Smith PETSC_EXTERN PetscErrorCode MatSolverPackageRegister_SuperLU(void) 71342c9c57cSBarry Smith { 71442c9c57cSBarry Smith PetscErrorCode ierr; 71542c9c57cSBarry Smith 71642c9c57cSBarry Smith PetscFunctionBegin; 71742c9c57cSBarry Smith ierr = MatSolverPackageRegister(MATSOLVERSUPERLU,MATSEQAIJ, MAT_FACTOR_LU,MatGetFactor_seqaij_superlu);CHKERRQ(ierr); 71842c9c57cSBarry Smith ierr = MatSolverPackageRegister(MATSOLVERSUPERLU,MATSEQAIJ, MAT_FACTOR_ILU,MatGetFactor_seqaij_superlu);CHKERRQ(ierr); 71942c9c57cSBarry Smith PetscFunctionReturn(0); 72042c9c57cSBarry Smith } 721