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 /* 34245fece9SBarry Smith This is the data that defines the SuperLU factored 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; 52c147c726SHong Zhang GlobalLU_t Glu; 532366e350SStefano Zampini PetscBool needconversion; 5414b396bbSKris Buschelman 5514b396bbSKris Buschelman /* Flag to clean up (non-global) SuperLU objects during Destroy */ 56ace3abfcSBarry Smith PetscBool CleanUpSuperLU; 57f0c56d0fSKris Buschelman } Mat_SuperLU; 5814b396bbSKris Buschelman 595a46d3faSBarry Smith /* 605a46d3faSBarry Smith Utility function 615a46d3faSBarry Smith */ 62d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatView_Info_SuperLU(Mat A, PetscViewer viewer) 63d71ae5a4SJacob Faibussowitsch { 64245fece9SBarry Smith Mat_SuperLU *lu = (Mat_SuperLU *)A->data; 655a46d3faSBarry Smith superlu_options_t options; 665a46d3faSBarry Smith 675a46d3faSBarry Smith PetscFunctionBegin; 685a46d3faSBarry Smith options = lu->options; 692205254eSKarl Rupp 709566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "SuperLU run parameters:\n")); 719566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Equil: %s\n", (options.Equil != NO) ? "YES" : "NO")); 729566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " ColPerm: %" PetscInt_FMT "\n", options.ColPerm)); 739566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " IterRefine: %" PetscInt_FMT "\n", options.IterRefine)); 749566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " SymmetricMode: %s\n", (options.SymmetricMode != NO) ? "YES" : "NO")); 759566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " DiagPivotThresh: %g\n", options.DiagPivotThresh)); 769566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " PivotGrowth: %s\n", (options.PivotGrowth != NO) ? "YES" : "NO")); 779566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " ConditionNumber: %s\n", (options.ConditionNumber != NO) ? "YES" : "NO")); 789566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " RowPerm: %" PetscInt_FMT "\n", options.RowPerm)); 799566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " ReplaceTinyPivot: %s\n", (options.ReplaceTinyPivot != NO) ? "YES" : "NO")); 809566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " PrintStat: %s\n", (options.PrintStat != NO) ? "YES" : "NO")); 819566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " lwork: %" PetscInt_FMT "\n", lu->lwork)); 82d5f3da31SBarry Smith if (A->factortype == MAT_FACTOR_ILU) { 839566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " ILU_DropTol: %g\n", options.ILU_DropTol)); 849566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " ILU_FillTol: %g\n", options.ILU_FillTol)); 859566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " ILU_FillFactor: %g\n", options.ILU_FillFactor)); 869566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " ILU_DropRule: %" PetscInt_FMT "\n", options.ILU_DropRule)); 879566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " ILU_Norm: %" PetscInt_FMT "\n", options.ILU_Norm)); 889566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " ILU_MILU: %" PetscInt_FMT "\n", options.ILU_MILU)); 89cffbb591SHong Zhang } 905a46d3faSBarry Smith PetscFunctionReturn(0); 915a46d3faSBarry Smith } 925a46d3faSBarry Smith 93d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSolve_SuperLU_Private(Mat A, Vec b, Vec x) 94d71ae5a4SJacob Faibussowitsch { 95245fece9SBarry Smith Mat_SuperLU *lu = (Mat_SuperLU *)A->data; 96245fece9SBarry Smith const PetscScalar *barray; 97245fece9SBarry Smith PetscScalar *xarray; 98245fece9SBarry Smith PetscInt info, i, n; 99245fece9SBarry Smith PetscReal ferr, berr; 100245fece9SBarry Smith static PetscBool cite = PETSC_FALSE; 101245fece9SBarry Smith 102245fece9SBarry Smith PetscFunctionBegin; 103245fece9SBarry Smith if (lu->lwork == -1) PetscFunctionReturn(0); 1049371c9d4SSatish Balay PetscCall(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 " 1059371c9d4SSatish Balay "pivoting},\n journal = {SIAM J. Matrix Analysis and Applications},\n year = {1999},\n volume = {20},\n number = {3},\n pages = {720-755}\n}\n", 1069371c9d4SSatish Balay &cite)); 107245fece9SBarry Smith 1089566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(x, &n)); 109245fece9SBarry Smith lu->B.ncol = 1; /* Set the number of right-hand side */ 110245fece9SBarry Smith if (lu->options.Equil && !lu->rhs_dup) { 111245fece9SBarry Smith /* superlu overwrites b when Equil is used, thus create rhs_dup to keep user's b unchanged */ 1129566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &lu->rhs_dup)); 113245fece9SBarry Smith } 114245fece9SBarry Smith if (lu->options.Equil) { 115245fece9SBarry Smith /* Copy b into rsh_dup */ 1169566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(b, &barray)); 1179566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(lu->rhs_dup, barray, n)); 1189566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(b, &barray)); 119245fece9SBarry Smith barray = lu->rhs_dup; 120245fece9SBarry Smith } else { 1219566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(b, &barray)); 122245fece9SBarry Smith } 1239566063dSJacob Faibussowitsch PetscCall(VecGetArray(x, &xarray)); 124245fece9SBarry Smith 125245fece9SBarry Smith #if defined(PETSC_USE_COMPLEX) 126245fece9SBarry Smith #if defined(PETSC_USE_REAL_SINGLE) 127245fece9SBarry Smith ((DNformat *)lu->B.Store)->nzval = (singlecomplex *)barray; 128245fece9SBarry Smith ((DNformat *)lu->X.Store)->nzval = (singlecomplex *)xarray; 129245fece9SBarry Smith #else 130245fece9SBarry Smith ((DNformat *)lu->B.Store)->nzval = (doublecomplex *)barray; 131245fece9SBarry Smith ((DNformat *)lu->X.Store)->nzval = (doublecomplex *)xarray; 132245fece9SBarry Smith #endif 133245fece9SBarry Smith #else 134245fece9SBarry Smith ((DNformat *)lu->B.Store)->nzval = (void *)barray; 135245fece9SBarry Smith ((DNformat *)lu->X.Store)->nzval = xarray; 136245fece9SBarry Smith #endif 137245fece9SBarry Smith 138245fece9SBarry Smith lu->options.Fact = FACTORED; /* Indicate the factored form of A is supplied. */ 139245fece9SBarry Smith if (A->factortype == MAT_FACTOR_LU) { 140245fece9SBarry Smith #if defined(PETSC_USE_COMPLEX) 141245fece9SBarry Smith #if defined(PETSC_USE_REAL_SINGLE) 1429371c9d4SSatish Balay PetscStackCallExternalVoid("SuperLU:cgssvx", cgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, &lu->Glu, &lu->mem_usage, &lu->stat, &info)); 143245fece9SBarry Smith #else 1449371c9d4SSatish Balay PetscStackCallExternalVoid("SuperLU:zgssvx", zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, &lu->Glu, &lu->mem_usage, &lu->stat, &info)); 145245fece9SBarry Smith #endif 146245fece9SBarry Smith #else 147245fece9SBarry Smith #if defined(PETSC_USE_REAL_SINGLE) 1489371c9d4SSatish Balay PetscStackCallExternalVoid("SuperLU:sgssvx", sgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, &lu->Glu, &lu->mem_usage, &lu->stat, &info)); 149245fece9SBarry Smith #else 1509371c9d4SSatish Balay PetscStackCallExternalVoid("SuperLU:dgssvx", dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, &lu->Glu, &lu->mem_usage, &lu->stat, &info)); 151245fece9SBarry Smith #endif 152245fece9SBarry Smith #endif 153245fece9SBarry Smith } else if (A->factortype == MAT_FACTOR_ILU) { 154245fece9SBarry Smith #if defined(PETSC_USE_COMPLEX) 155245fece9SBarry Smith #if defined(PETSC_USE_REAL_SINGLE) 1569371c9d4SSatish Balay PetscStackCallExternalVoid("SuperLU:cgsisx", cgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &lu->Glu, &lu->mem_usage, &lu->stat, &info)); 157245fece9SBarry Smith #else 1589371c9d4SSatish Balay PetscStackCallExternalVoid("SuperLU:zgsisx", zgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &lu->Glu, &lu->mem_usage, &lu->stat, &info)); 159245fece9SBarry Smith #endif 160245fece9SBarry Smith #else 161245fece9SBarry Smith #if defined(PETSC_USE_REAL_SINGLE) 1629371c9d4SSatish Balay PetscStackCallExternalVoid("SuperLU:sgsisx", sgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &lu->Glu, &lu->mem_usage, &lu->stat, &info)); 163245fece9SBarry Smith #else 1649371c9d4SSatish Balay PetscStackCallExternalVoid("SuperLU:dgsisx", dgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &lu->Glu, &lu->mem_usage, &lu->stat, &info)); 165245fece9SBarry Smith #endif 166245fece9SBarry Smith #endif 167245fece9SBarry Smith } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Factor type not supported"); 16848a46eb9SPierre Jolivet if (!lu->options.Equil) PetscCall(VecRestoreArrayRead(b, &barray)); 1699566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(x, &xarray)); 170245fece9SBarry Smith 171245fece9SBarry Smith if (!info || info == lu->A.ncol + 1) { 172245fece9SBarry Smith if (lu->options.IterRefine) { 1739566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "Iterative Refinement:\n")); 1749566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " %8s%8s%16s%16s\n", "rhs", "Steps", "FERR", "BERR")); 1759566063dSJacob Faibussowitsch for (i = 0; i < 1; ++i) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %8d%8d%16e%16e\n", i + 1, lu->stat.RefineSteps, ferr, berr)); 176245fece9SBarry Smith } 177245fece9SBarry Smith } else if (info > 0) { 178245fece9SBarry Smith if (lu->lwork == -1) { 1799566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " ** Estimated memory: %" PetscInt_FMT " bytes\n", info - lu->A.ncol)); 180245fece9SBarry Smith } else { 1819566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " Warning: gssvx() returns info %" PetscInt_FMT "\n", info)); 182245fece9SBarry Smith } 1835f80ce2aSJacob Faibussowitsch } else PetscCheck(info >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "info = %" PetscInt_FMT ", the %" PetscInt_FMT "-th argument in gssvx() had an illegal value", info, -info); 184245fece9SBarry Smith 185245fece9SBarry Smith if (lu->options.PrintStat) { 1869566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "MatSolve__SuperLU():\n")); 187e77caa6dSBarry Smith PetscStackCallExternalVoid("SuperLU:StatPrint", StatPrint(&lu->stat)); 188245fece9SBarry Smith } 189245fece9SBarry Smith PetscFunctionReturn(0); 190245fece9SBarry Smith } 191245fece9SBarry Smith 192d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSolve_SuperLU(Mat A, Vec b, Vec x) 193d71ae5a4SJacob Faibussowitsch { 194245fece9SBarry Smith Mat_SuperLU *lu = (Mat_SuperLU *)A->data; 195245fece9SBarry Smith 196245fece9SBarry Smith PetscFunctionBegin; 197603e8f96SBarry Smith if (A->factorerrortype) { 1989566063dSJacob Faibussowitsch PetscCall(PetscInfo(A, "MatSolve is called with singular matrix factor, skip\n")); 1999566063dSJacob Faibussowitsch PetscCall(VecSetInf(x)); 200245fece9SBarry Smith PetscFunctionReturn(0); 201245fece9SBarry Smith } 202245fece9SBarry Smith 203245fece9SBarry Smith lu->options.Trans = TRANS; 2049566063dSJacob Faibussowitsch PetscCall(MatSolve_SuperLU_Private(A, b, x)); 205245fece9SBarry Smith PetscFunctionReturn(0); 206245fece9SBarry Smith } 207245fece9SBarry Smith 208d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSolveTranspose_SuperLU(Mat A, Vec b, Vec x) 209d71ae5a4SJacob Faibussowitsch { 210245fece9SBarry Smith Mat_SuperLU *lu = (Mat_SuperLU *)A->data; 211245fece9SBarry Smith 212245fece9SBarry Smith PetscFunctionBegin; 213603e8f96SBarry Smith if (A->factorerrortype) { 2149566063dSJacob Faibussowitsch PetscCall(PetscInfo(A, "MatSolve is called with singular matrix factor, skip\n")); 2159566063dSJacob Faibussowitsch PetscCall(VecSetInf(x)); 216245fece9SBarry Smith PetscFunctionReturn(0); 217245fece9SBarry Smith } 218245fece9SBarry Smith 219245fece9SBarry Smith lu->options.Trans = NOTRANS; 2209566063dSJacob Faibussowitsch PetscCall(MatSolve_SuperLU_Private(A, b, x)); 221245fece9SBarry Smith PetscFunctionReturn(0); 222245fece9SBarry Smith } 223245fece9SBarry Smith 224d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatLUFactorNumeric_SuperLU(Mat F, Mat A, const MatFactorInfo *info) 225d71ae5a4SJacob Faibussowitsch { 226245fece9SBarry Smith Mat_SuperLU *lu = (Mat_SuperLU *)F->data; 227cae5a23dSHong Zhang Mat_SeqAIJ *aa; 2285a46d3faSBarry Smith PetscInt sinfo; 2295a46d3faSBarry Smith PetscReal ferr, berr; 2305a46d3faSBarry Smith NCformat *Ustore; 2315a46d3faSBarry Smith SCformat *Lstore; 2325a46d3faSBarry Smith 2335a46d3faSBarry Smith PetscFunctionBegin; 2345a46d3faSBarry Smith if (lu->flg == SAME_NONZERO_PATTERN) { /* successing numerical factorization */ 2355a46d3faSBarry Smith lu->options.Fact = SamePattern; 2365a46d3faSBarry Smith /* Ref: ~SuperLU_3.0/EXAMPLE/dlinsolx2.c */ 2375a46d3faSBarry Smith Destroy_SuperMatrix_Store(&lu->A); 2381baa6e33SBarry Smith if (lu->A_dup) PetscCall(MatCopy_SeqAIJ(A, lu->A_dup, SAME_NONZERO_PATTERN)); 2392366e350SStefano Zampini 2405a46d3faSBarry Smith if (lu->lwork >= 0) { 241e77caa6dSBarry Smith PetscStackCallExternalVoid("SuperLU:Destroy_SuperNode_Matrix", Destroy_SuperNode_Matrix(&lu->L)); 242e77caa6dSBarry Smith PetscStackCallExternalVoid("SuperLU:Destroy_CompCol_Matrix", Destroy_CompCol_Matrix(&lu->U)); 2435a46d3faSBarry Smith lu->options.Fact = SamePattern; 2445a46d3faSBarry Smith } 2455a46d3faSBarry Smith } 2465a46d3faSBarry Smith 2475a46d3faSBarry Smith /* Create the SuperMatrix for lu->A=A^T: 2485a46d3faSBarry Smith Since SuperLU likes column-oriented matrices,we pass it the transpose, 2495a46d3faSBarry Smith and then solve A^T X = B in MatSolve(). */ 2502366e350SStefano Zampini if (lu->A_dup) { 251cae5a23dSHong Zhang aa = (Mat_SeqAIJ *)(lu->A_dup)->data; 252cae5a23dSHong Zhang } else { 253cae5a23dSHong Zhang aa = (Mat_SeqAIJ *)(A)->data; 254cae5a23dSHong Zhang } 2555a46d3faSBarry Smith #if defined(PETSC_USE_COMPLEX) 2563cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 257e77caa6dSBarry Smith PetscStackCallExternalVoid("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)); 2583cb270beSHong Zhang #else 259e77caa6dSBarry Smith PetscStackCallExternalVoid("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)); 2603cb270beSHong Zhang #endif 2613cb270beSHong Zhang #else 2623cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 263e77caa6dSBarry Smith PetscStackCallExternalVoid("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)); 2645a46d3faSBarry Smith #else 265e77caa6dSBarry Smith PetscStackCallExternalVoid("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)); 2665a46d3faSBarry Smith #endif 2673cb270beSHong Zhang #endif 2685a46d3faSBarry Smith 2695a46d3faSBarry Smith /* Numerical factorization */ 2705a46d3faSBarry Smith lu->B.ncol = 0; /* Indicate not to solve the system */ 271d5f3da31SBarry Smith if (F->factortype == MAT_FACTOR_LU) { 2725a46d3faSBarry Smith #if defined(PETSC_USE_COMPLEX) 2733cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 2749371c9d4SSatish Balay PetscStackCallExternalVoid("SuperLU:cgssvx", cgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, &lu->Glu, &lu->mem_usage, &lu->stat, &sinfo)); 2753cb270beSHong Zhang #else 2769371c9d4SSatish Balay PetscStackCallExternalVoid("SuperLU:zgssvx", zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, &lu->Glu, &lu->mem_usage, &lu->stat, &sinfo)); 2773cb270beSHong Zhang #endif 2783cb270beSHong Zhang #else 2793cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 2809371c9d4SSatish Balay PetscStackCallExternalVoid("SuperLU:sgssvx", sgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, &lu->Glu, &lu->mem_usage, &lu->stat, &sinfo)); 2815a46d3faSBarry Smith #else 2829371c9d4SSatish Balay PetscStackCallExternalVoid("SuperLU:dgssvx", dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, &lu->Glu, &lu->mem_usage, &lu->stat, &sinfo)); 2835a46d3faSBarry Smith #endif 2843cb270beSHong Zhang #endif 285d5f3da31SBarry Smith } else if (F->factortype == MAT_FACTOR_ILU) { 286cffbb591SHong Zhang /* Compute the incomplete factorization, condition number and pivot growth */ 287cffbb591SHong Zhang #if defined(PETSC_USE_COMPLEX) 2883cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 2899371c9d4SSatish Balay PetscStackCallExternalVoid("SuperLU:cgsisx", cgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &lu->Glu, &lu->mem_usage, &lu->stat, &sinfo)); 2903cb270beSHong Zhang #else 2919371c9d4SSatish Balay PetscStackCallExternalVoid("SuperLU:zgsisx", zgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &lu->Glu, &lu->mem_usage, &lu->stat, &sinfo)); 2923cb270beSHong Zhang #endif 2933cb270beSHong Zhang #else 2943cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 2959371c9d4SSatish Balay PetscStackCallExternalVoid("SuperLU:sgsisx", sgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &lu->Glu, &lu->mem_usage, &lu->stat, &sinfo)); 296cffbb591SHong Zhang #else 2979371c9d4SSatish Balay PetscStackCallExternalVoid("SuperLU:dgsisx", dgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &lu->Glu, &lu->mem_usage, &lu->stat, &sinfo)); 298cffbb591SHong Zhang #endif 2993cb270beSHong Zhang #endif 300f23aa3ddSBarry Smith } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Factor type not supported"); 3015a46d3faSBarry Smith if (!sinfo || sinfo == lu->A.ncol + 1) { 30248a46eb9SPierre Jolivet if (lu->options.PivotGrowth) PetscCall(PetscPrintf(PETSC_COMM_SELF, " Recip. pivot growth = %e\n", lu->rpg)); 30348a46eb9SPierre Jolivet if (lu->options.ConditionNumber) PetscCall(PetscPrintf(PETSC_COMM_SELF, " Recip. condition number = %e\n", lu->rcond)); 3045a46d3faSBarry Smith } else if (sinfo > 0) { 305675d1226SHong Zhang if (A->erroriffailure) { 30698921bdaSJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_MAT_LU_ZRPVT, "Zero pivot in row %" PetscInt_FMT, sinfo); 307675d1226SHong Zhang } else { 308675d1226SHong Zhang if (sinfo <= lu->A.ncol) { 309ad540459SPierre Jolivet if (lu->options.ILU_FillTol == 0.0) F->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 3109566063dSJacob Faibussowitsch PetscCall(PetscInfo(F, "Number of zero pivots %" PetscInt_FMT ", ILU_FillTol %g\n", sinfo, lu->options.ILU_FillTol)); 311675d1226SHong Zhang } else if (sinfo == lu->A.ncol + 1) { 312675d1226SHong Zhang /* 313675d1226SHong Zhang U is nonsingular, but RCOND is less than machine 314675d1226SHong Zhang precision, meaning that the matrix is singular to 315675d1226SHong Zhang working precision. Nevertheless, the solution and 316675d1226SHong Zhang error bounds are computed because there are a number 317675d1226SHong Zhang of situations where the computed solution can be more 318675d1226SHong Zhang accurate than the value of RCOND would suggest. 319675d1226SHong Zhang */ 3209566063dSJacob Faibussowitsch PetscCall(PetscInfo(F, "Matrix factor U is nonsingular, but is singular to working precision. The solution is computed. info %" PetscInt_FMT, sinfo)); 321675d1226SHong Zhang } else { /* sinfo > lu->A.ncol + 1 */ 322603e8f96SBarry Smith F->factorerrortype = MAT_FACTOR_OUTMEMORY; 3239566063dSJacob Faibussowitsch PetscCall(PetscInfo(F, "Number of bytes allocated when memory allocation fails %" PetscInt_FMT "\n", sinfo)); 324675d1226SHong Zhang } 325675d1226SHong Zhang } 32698921bdaSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "info = %" PetscInt_FMT ", the %" PetscInt_FMT "-th argument in gssvx() had an illegal value", sinfo, -sinfo); 3275a46d3faSBarry Smith 3285a46d3faSBarry Smith if (lu->options.PrintStat) { 3299566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "MatLUFactorNumeric_SuperLU():\n")); 330e77caa6dSBarry Smith PetscStackCallExternalVoid("SuperLU:StatPrint", StatPrint(&lu->stat)); 3315a46d3faSBarry Smith Lstore = (SCformat *)lu->L.Store; 3325a46d3faSBarry Smith Ustore = (NCformat *)lu->U.Store; 3339566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " No of nonzeros in factor L = %" PetscInt_FMT "\n", Lstore->nnz)); 3349566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " No of nonzeros in factor U = %" PetscInt_FMT "\n", Ustore->nnz)); 3359566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " No of nonzeros in L+U = %" PetscInt_FMT "\n", Lstore->nnz + Ustore->nnz - lu->A.ncol)); 3369566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " L\\U MB %.3f\ttotal MB needed %.3f\n", lu->mem_usage.for_lu / 1e6, lu->mem_usage.total_needed / 1e6)); 3375a46d3faSBarry Smith } 3385a46d3faSBarry Smith 3395a46d3faSBarry Smith lu->flg = SAME_NONZERO_PATTERN; 3401d5ca7f3SHong Zhang F->ops->solve = MatSolve_SuperLU; 3411d5ca7f3SHong Zhang F->ops->solvetranspose = MatSolveTranspose_SuperLU; 3421aef8b4cSStefano Zampini F->ops->matsolve = NULL; 3435a46d3faSBarry Smith PetscFunctionReturn(0); 3445a46d3faSBarry Smith } 3455a46d3faSBarry Smith 346d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatDestroy_SuperLU(Mat A) 347d71ae5a4SJacob Faibussowitsch { 348245fece9SBarry Smith Mat_SuperLU *lu = (Mat_SuperLU *)A->data; 34914b396bbSKris Buschelman 35014b396bbSKris Buschelman PetscFunctionBegin; 351245fece9SBarry Smith if (lu->CleanUpSuperLU) { /* Free the SuperLU datastructures */ 352e77caa6dSBarry Smith PetscStackCallExternalVoid("SuperLU:Destroy_SuperMatrix_Store", Destroy_SuperMatrix_Store(&lu->A)); 353e77caa6dSBarry Smith PetscStackCallExternalVoid("SuperLU:Destroy_SuperMatrix_Store", Destroy_SuperMatrix_Store(&lu->B)); 354e77caa6dSBarry Smith PetscStackCallExternalVoid("SuperLU:Destroy_SuperMatrix_Store", Destroy_SuperMatrix_Store(&lu->X)); 355e77caa6dSBarry Smith PetscStackCallExternalVoid("SuperLU:StatFree", StatFree(&lu->stat)); 3560e742b69SHong Zhang if (lu->lwork >= 0) { 357e77caa6dSBarry Smith PetscStackCallExternalVoid("SuperLU:Destroy_SuperNode_Matrix", Destroy_SuperNode_Matrix(&lu->L)); 358e77caa6dSBarry Smith PetscStackCallExternalVoid("SuperLU:Destroy_CompCol_Matrix", Destroy_CompCol_Matrix(&lu->U)); 3590e742b69SHong Zhang } 3600e742b69SHong Zhang } 3619566063dSJacob Faibussowitsch PetscCall(PetscFree(lu->etree)); 3629566063dSJacob Faibussowitsch PetscCall(PetscFree(lu->perm_r)); 3639566063dSJacob Faibussowitsch PetscCall(PetscFree(lu->perm_c)); 3649566063dSJacob Faibussowitsch PetscCall(PetscFree(lu->R)); 3659566063dSJacob Faibussowitsch PetscCall(PetscFree(lu->C)); 3669566063dSJacob Faibussowitsch PetscCall(PetscFree(lu->rhs_dup)); 3679566063dSJacob Faibussowitsch PetscCall(MatDestroy(&lu->A_dup)); 3689566063dSJacob Faibussowitsch PetscCall(PetscFree(A->data)); 3690e742b69SHong Zhang 370d954db57SHong Zhang /* clear composed functions */ 3719566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatFactorGetSolverType_C", NULL)); 3729566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSuperluSetILUDropTol_C", NULL)); 37314b396bbSKris Buschelman PetscFunctionReturn(0); 37414b396bbSKris Buschelman } 37514b396bbSKris Buschelman 376d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatView_SuperLU(Mat A, PetscViewer viewer) 377d71ae5a4SJacob Faibussowitsch { 378ace3abfcSBarry Smith PetscBool iascii; 37914b396bbSKris Buschelman PetscViewerFormat format; 38014b396bbSKris Buschelman 38114b396bbSKris Buschelman PetscFunctionBegin; 3829566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 38332077d6dSBarry Smith if (iascii) { 3849566063dSJacob Faibussowitsch PetscCall(PetscViewerGetFormat(viewer, &format)); 38548a46eb9SPierre Jolivet if (format == PETSC_VIEWER_ASCII_INFO) PetscCall(MatView_Info_SuperLU(A, viewer)); 38614b396bbSKris Buschelman } 38714b396bbSKris Buschelman PetscFunctionReturn(0); 38814b396bbSKris Buschelman } 38914b396bbSKris Buschelman 390d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMatSolve_SuperLU(Mat A, Mat B, Mat X) 391d71ae5a4SJacob Faibussowitsch { 392245fece9SBarry Smith Mat_SuperLU *lu = (Mat_SuperLU *)A->data; 393cd723cd1SBarry Smith PetscBool flg; 394e0b74bf9SHong Zhang 395e0b74bf9SHong Zhang PetscFunctionBegin; 3969566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompareAny((PetscObject)B, &flg, MATSEQDENSE, MATMPIDENSE, NULL)); 3975f80ce2aSJacob Faibussowitsch PetscCheck(flg, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Matrix B must be MATDENSE matrix"); 3989566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompareAny((PetscObject)X, &flg, MATSEQDENSE, MATMPIDENSE, NULL)); 3995f80ce2aSJacob Faibussowitsch PetscCheck(flg, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Matrix X must be MATDENSE matrix"); 4002205254eSKarl Rupp lu->options.Trans = TRANS; 401e0b74bf9SHong Zhang SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "MatMatSolve_SuperLU() is not implemented yet"); 402e0b74bf9SHong Zhang PetscFunctionReturn(0); 403e0b74bf9SHong Zhang } 404e0b74bf9SHong Zhang 405d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatLUFactorSymbolic_SuperLU(Mat F, Mat A, IS r, IS c, const MatFactorInfo *info) 406d71ae5a4SJacob Faibussowitsch { 407245fece9SBarry Smith Mat_SuperLU *lu = (Mat_SuperLU *)(F->data); 40826cc229bSBarry Smith PetscInt indx; 40926cc229bSBarry Smith PetscBool flg, set; 41026cc229bSBarry Smith PetscReal real_input; 41126cc229bSBarry Smith const char *colperm[] = {"NATURAL", "MMD_ATA", "MMD_AT_PLUS_A", "COLAMD"}; /* MY_PERMC - not supported by the petsc interface yet */ 41226cc229bSBarry Smith const char *iterrefine[] = {"NOREFINE", "SINGLE", "DOUBLE", "EXTRA"}; 41326cc229bSBarry Smith const char *rowperm[] = {"NOROWPERM", "LargeDiag"}; /* MY_PERMC - not supported by the petsc interface yet */ 414b24902e0SBarry Smith 415b24902e0SBarry Smith PetscFunctionBegin; 41626cc229bSBarry Smith /* Set options to F */ 41726cc229bSBarry Smith PetscOptionsBegin(PetscObjectComm((PetscObject)F), ((PetscObject)F)->prefix, "SuperLU Options", "Mat"); 41826cc229bSBarry Smith PetscCall(PetscOptionsBool("-mat_superlu_equil", "Equil", "None", (PetscBool)lu->options.Equil, (PetscBool *)&lu->options.Equil, NULL)); 41926cc229bSBarry Smith PetscCall(PetscOptionsEList("-mat_superlu_colperm", "ColPerm", "None", colperm, 4, colperm[3], &indx, &flg)); 42026cc229bSBarry Smith if (flg) lu->options.ColPerm = (colperm_t)indx; 42126cc229bSBarry Smith PetscCall(PetscOptionsEList("-mat_superlu_iterrefine", "IterRefine", "None", iterrefine, 4, iterrefine[0], &indx, &flg)); 42226cc229bSBarry Smith if (flg) lu->options.IterRefine = (IterRefine_t)indx; 42326cc229bSBarry Smith PetscCall(PetscOptionsBool("-mat_superlu_symmetricmode", "SymmetricMode", "None", (PetscBool)lu->options.SymmetricMode, &flg, &set)); 42426cc229bSBarry Smith if (set && flg) lu->options.SymmetricMode = YES; 42526cc229bSBarry Smith PetscCall(PetscOptionsReal("-mat_superlu_diagpivotthresh", "DiagPivotThresh", "None", lu->options.DiagPivotThresh, &real_input, &flg)); 42626cc229bSBarry Smith if (flg) lu->options.DiagPivotThresh = (double)real_input; 42726cc229bSBarry Smith PetscCall(PetscOptionsBool("-mat_superlu_pivotgrowth", "PivotGrowth", "None", (PetscBool)lu->options.PivotGrowth, &flg, &set)); 42826cc229bSBarry Smith if (set && flg) lu->options.PivotGrowth = YES; 42926cc229bSBarry Smith PetscCall(PetscOptionsBool("-mat_superlu_conditionnumber", "ConditionNumber", "None", (PetscBool)lu->options.ConditionNumber, &flg, &set)); 43026cc229bSBarry Smith if (set && flg) lu->options.ConditionNumber = YES; 43126cc229bSBarry Smith PetscCall(PetscOptionsEList("-mat_superlu_rowperm", "rowperm", "None", rowperm, 2, rowperm[lu->options.RowPerm], &indx, &flg)); 43226cc229bSBarry Smith if (flg) lu->options.RowPerm = (rowperm_t)indx; 43326cc229bSBarry Smith PetscCall(PetscOptionsBool("-mat_superlu_replacetinypivot", "ReplaceTinyPivot", "None", (PetscBool)lu->options.ReplaceTinyPivot, &flg, &set)); 43426cc229bSBarry Smith if (set && flg) lu->options.ReplaceTinyPivot = YES; 43526cc229bSBarry Smith PetscCall(PetscOptionsBool("-mat_superlu_printstat", "PrintStat", "None", (PetscBool)lu->options.PrintStat, &flg, &set)); 43626cc229bSBarry Smith if (set && flg) lu->options.PrintStat = YES; 43726cc229bSBarry Smith PetscCall(PetscOptionsInt("-mat_superlu_lwork", "size of work array in bytes used by factorization", "None", lu->lwork, &lu->lwork, NULL)); 43826cc229bSBarry Smith if (lu->lwork > 0) { 43926cc229bSBarry Smith /* lwork is in bytes, hence PetscMalloc() is used here, not PetscMalloc1()*/ 44026cc229bSBarry Smith PetscCall(PetscMalloc(lu->lwork, &lu->work)); 44126cc229bSBarry Smith } else if (lu->lwork != 0 && lu->lwork != -1) { 44226cc229bSBarry Smith PetscCall(PetscPrintf(PETSC_COMM_SELF, " Warning: lwork %" PetscInt_FMT " is not supported by SUPERLU. The default lwork=0 is used.\n", lu->lwork)); 44326cc229bSBarry Smith lu->lwork = 0; 44426cc229bSBarry Smith } 44526cc229bSBarry Smith /* ilu options */ 44626cc229bSBarry Smith PetscCall(PetscOptionsReal("-mat_superlu_ilu_droptol", "ILU_DropTol", "None", lu->options.ILU_DropTol, &real_input, &flg)); 44726cc229bSBarry Smith if (flg) lu->options.ILU_DropTol = (double)real_input; 44826cc229bSBarry Smith PetscCall(PetscOptionsReal("-mat_superlu_ilu_filltol", "ILU_FillTol", "None", lu->options.ILU_FillTol, &real_input, &flg)); 44926cc229bSBarry Smith if (flg) lu->options.ILU_FillTol = (double)real_input; 45026cc229bSBarry Smith PetscCall(PetscOptionsReal("-mat_superlu_ilu_fillfactor", "ILU_FillFactor", "None", lu->options.ILU_FillFactor, &real_input, &flg)); 45126cc229bSBarry Smith if (flg) lu->options.ILU_FillFactor = (double)real_input; 45226cc229bSBarry Smith PetscCall(PetscOptionsInt("-mat_superlu_ilu_droprull", "ILU_DropRule", "None", lu->options.ILU_DropRule, &lu->options.ILU_DropRule, NULL)); 45326cc229bSBarry Smith PetscCall(PetscOptionsInt("-mat_superlu_ilu_norm", "ILU_Norm", "None", lu->options.ILU_Norm, &indx, &flg)); 45426cc229bSBarry Smith if (flg) lu->options.ILU_Norm = (norm_t)indx; 45526cc229bSBarry Smith PetscCall(PetscOptionsInt("-mat_superlu_ilu_milu", "ILU_MILU", "None", lu->options.ILU_MILU, &indx, &flg)); 45626cc229bSBarry Smith if (flg) lu->options.ILU_MILU = (milu_t)indx; 45726cc229bSBarry Smith PetscOptionsEnd(); 45826cc229bSBarry Smith 459b24902e0SBarry Smith lu->flg = DIFFERENT_NONZERO_PATTERN; 460b24902e0SBarry Smith lu->CleanUpSuperLU = PETSC_TRUE; 4611d5ca7f3SHong Zhang F->ops->lufactornumeric = MatLUFactorNumeric_SuperLU; 4622366e350SStefano Zampini 4632366e350SStefano Zampini /* if we are here, the nonzero pattern has changed unless the user explicitly called MatLUFactorSymbolic */ 4649566063dSJacob Faibussowitsch PetscCall(MatDestroy(&lu->A_dup)); 46548a46eb9SPierre Jolivet if (lu->needconversion) PetscCall(MatConvert(A, MATSEQAIJ, MAT_INITIAL_MATRIX, &lu->A_dup)); 4662366e350SStefano Zampini if (lu->options.Equil == YES && !lu->A_dup) { /* superlu overwrites input matrix and rhs when Equil is used, thus create A_dup to keep user's A unchanged */ 4679566063dSJacob Faibussowitsch PetscCall(MatDuplicate_SeqAIJ(A, MAT_COPY_VALUES, &lu->A_dup)); 4682366e350SStefano Zampini } 469b24902e0SBarry Smith PetscFunctionReturn(0); 470b24902e0SBarry Smith } 471b24902e0SBarry Smith 472d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSuperluSetILUDropTol_SuperLU(Mat F, PetscReal dtol) 473d71ae5a4SJacob Faibussowitsch { 474245fece9SBarry Smith Mat_SuperLU *lu = (Mat_SuperLU *)F->data; 4755ccb76cbSHong Zhang 4765ccb76cbSHong Zhang PetscFunctionBegin; 4775ccb76cbSHong Zhang lu->options.ILU_DropTol = dtol; 4785ccb76cbSHong Zhang PetscFunctionReturn(0); 4795ccb76cbSHong Zhang } 4805ccb76cbSHong Zhang 4815ccb76cbSHong Zhang /*@ 4825ccb76cbSHong Zhang MatSuperluSetILUDropTol - Set SuperLU ILU drop tolerance 483147403d9SBarry Smith 484*c3339decSBarry Smith Logically Collective 4855ccb76cbSHong Zhang 4865ccb76cbSHong Zhang Input Parameters: 48711a5261eSBarry Smith + F - the factored matrix obtained by calling `MatGetFactor()` from PETSc-SuperLU interface 4885ccb76cbSHong Zhang - dtol - drop tolerance 4895ccb76cbSHong Zhang 4903c7db156SBarry Smith Options Database Key: 491147403d9SBarry Smith . -mat_superlu_ilu_droptol <dtol> - the drop tolerance 4925ccb76cbSHong Zhang 4935ccb76cbSHong Zhang Level: beginner 4945ccb76cbSHong Zhang 49596a0c994SBarry Smith References: 496606c0280SSatish Balay . * - SuperLU Users' Guide 4975ccb76cbSHong Zhang 498db781477SPatrick Sanan .seealso: `MatGetFactor()` 4995ccb76cbSHong Zhang @*/ 500d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSuperluSetILUDropTol(Mat F, PetscReal dtol) 501d71ae5a4SJacob Faibussowitsch { 5025ccb76cbSHong Zhang PetscFunctionBegin; 5035ccb76cbSHong Zhang PetscValidHeaderSpecific(F, MAT_CLASSID, 1); 50469fbec6eSBarry Smith PetscValidLogicalCollectiveReal(F, dtol, 2); 505cac4c232SBarry Smith PetscTryMethod(F, "MatSuperluSetILUDropTol_C", (Mat, PetscReal), (F, dtol)); 5065ccb76cbSHong Zhang PetscFunctionReturn(0); 5075ccb76cbSHong Zhang } 5085ccb76cbSHong Zhang 509d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFactorGetSolverType_seqaij_superlu(Mat A, MatSolverType *type) 510d71ae5a4SJacob Faibussowitsch { 51135bd34faSBarry Smith PetscFunctionBegin; 5122692d6eeSBarry Smith *type = MATSOLVERSUPERLU; 51335bd34faSBarry Smith PetscFunctionReturn(0); 51435bd34faSBarry Smith } 51535bd34faSBarry Smith 516b24902e0SBarry Smith /*MC 517ba992d64SSatish Balay MATSOLVERSUPERLU = "superlu" - A solver package providing solvers LU and ILU for sequential matrices 518b24902e0SBarry Smith via the external package SuperLU. 519b24902e0SBarry Smith 520e2e64c6bSBarry Smith Use ./configure --download-superlu to have PETSc installed with SuperLU 521b24902e0SBarry Smith 5223ca39a21SBarry Smith Use -pc_type lu -pc_factor_mat_solver_type superlu to use this direct solver 523c2b89b5dSBarry Smith 524b24902e0SBarry Smith Options Database Keys: 525e08999f5SMatthew G Knepley + -mat_superlu_equil <FALSE> - Equil (None) 526e08999f5SMatthew G Knepley . -mat_superlu_colperm <COLAMD> - (choose one of) NATURAL MMD_ATA MMD_AT_PLUS_A COLAMD 527e08999f5SMatthew G Knepley . -mat_superlu_iterrefine <NOREFINE> - (choose one of) NOREFINE SINGLE DOUBLE EXTRA 528e08999f5SMatthew G Knepley . -mat_superlu_symmetricmode: <FALSE> - SymmetricMode (None) 529e08999f5SMatthew G Knepley . -mat_superlu_diagpivotthresh <1> - DiagPivotThresh (None) 530e08999f5SMatthew G Knepley . -mat_superlu_pivotgrowth <FALSE> - PivotGrowth (None) 531e08999f5SMatthew G Knepley . -mat_superlu_conditionnumber <FALSE> - ConditionNumber (None) 532e08999f5SMatthew G Knepley . -mat_superlu_rowperm <NOROWPERM> - (choose one of) NOROWPERM LargeDiag 533e08999f5SMatthew G Knepley . -mat_superlu_replacetinypivot <FALSE> - ReplaceTinyPivot (None) 534e08999f5SMatthew G Knepley . -mat_superlu_printstat <FALSE> - PrintStat (None) 535e08999f5SMatthew G Knepley . -mat_superlu_lwork <0> - size of work array in bytes used by factorization (None) 536e08999f5SMatthew G Knepley . -mat_superlu_ilu_droptol <0> - ILU_DropTol (None) 537e08999f5SMatthew G Knepley . -mat_superlu_ilu_filltol <0> - ILU_FillTol (None) 538e08999f5SMatthew G Knepley . -mat_superlu_ilu_fillfactor <0> - ILU_FillFactor (None) 539e08999f5SMatthew G Knepley . -mat_superlu_ilu_droprull <0> - ILU_DropRule (None) 540e08999f5SMatthew G Knepley . -mat_superlu_ilu_norm <0> - ILU_Norm (None) 541e08999f5SMatthew G Knepley - -mat_superlu_ilu_milu <0> - ILU_MILU (None) 542b24902e0SBarry Smith 54395452b02SPatrick Sanan Notes: 54411a5261eSBarry Smith Do not confuse this with `MATSOLVERSUPERLU_DIST` which is for parallel sparse solves 5455c9eb25fSBarry Smith 54611a5261eSBarry Smith Cannot use ordering provided by PETSc, provides its own. 5472c7c0729SBarry Smith 548b24902e0SBarry Smith Level: beginner 549b24902e0SBarry Smith 550db781477SPatrick Sanan .seealso: `PCLU`, `PCILU`, `MATSOLVERSUPERLU_DIST`, `MATSOLVERMUMPS`, `PCFactorSetMatSolverType()`, `MatSolverType` 551b24902e0SBarry Smith M*/ 552b24902e0SBarry Smith 553d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGetFactor_seqaij_superlu(Mat A, MatFactorType ftype, Mat *F) 554d71ae5a4SJacob Faibussowitsch { 55514b396bbSKris Buschelman Mat B; 556f0c56d0fSKris Buschelman Mat_SuperLU *lu; 55726cc229bSBarry Smith PetscInt m = A->rmap->n, n = A->cmap->n; 55814b396bbSKris Buschelman 55914b396bbSKris Buschelman PetscFunctionBegin; 5609566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B)); 5619566063dSJacob Faibussowitsch PetscCall(MatSetSizes(B, A->rmap->n, A->cmap->n, PETSC_DETERMINE, PETSC_DETERMINE)); 5629566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy("superlu", &((PetscObject)B)->type_name)); 5639566063dSJacob Faibussowitsch PetscCall(MatSetUp(B)); 56466e17bc3SBarry Smith B->trivialsymbolic = PETSC_TRUE; 565cffbb591SHong Zhang if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU) { 566b24902e0SBarry Smith B->ops->lufactorsymbolic = MatLUFactorSymbolic_SuperLU; 567cffbb591SHong Zhang B->ops->ilufactorsymbolic = MatLUFactorSymbolic_SuperLU; 568b3a44c85SBarry Smith } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Factor type not supported"); 569cffbb591SHong Zhang 5709566063dSJacob Faibussowitsch PetscCall(PetscFree(B->solvertype)); 5719566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(MATSOLVERSUPERLU, &B->solvertype)); 57200c67f3bSHong Zhang 573b9c12af5SBarry Smith B->ops->getinfo = MatGetInfo_External; 574b24902e0SBarry Smith B->ops->destroy = MatDestroy_SuperLU; 5753519fcdcSHong Zhang B->ops->view = MatView_SuperLU; 576d5f3da31SBarry Smith B->factortype = ftype; 57794b7f48cSBarry Smith B->assembled = PETSC_TRUE; /* required by -ksp_view */ 5785c9eb25fSBarry Smith B->preallocated = PETSC_TRUE; 57914b396bbSKris Buschelman 5804dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&lu)); 581cae5a23dSHong Zhang 582cffbb591SHong Zhang if (ftype == MAT_FACTOR_LU) { 5839ce50919SHong Zhang set_default_options(&lu->options); 5843d6581fbSHong Zhang /* Comments from SuperLU_4.0/SRC/dgssvx.c: 5853d6581fbSHong Zhang "Whether or not the system will be equilibrated depends on the 5863d6581fbSHong Zhang scaling of the matrix A, but if equilibration is used, A is 5873d6581fbSHong Zhang overwritten by diag(R)*A*diag(C) and B by diag(R)*B 5883d6581fbSHong Zhang (if options->Trans=NOTRANS) or diag(C)*B (if options->Trans = TRANS or CONJ)." 5893d6581fbSHong Zhang We set 'options.Equil = NO' as default because additional space is needed for it. 5903d6581fbSHong Zhang */ 5913d6581fbSHong Zhang lu->options.Equil = NO; 592cffbb591SHong Zhang } else if (ftype == MAT_FACTOR_ILU) { 5930c74a584SJed Brown /* Set the default input options of ilu: */ 594e77caa6dSBarry Smith PetscStackCallExternalVoid("SuperLU:ilu_set_default_options", ilu_set_default_options(&lu->options)); 595cffbb591SHong Zhang } 5969ce50919SHong Zhang lu->options.PrintStat = NO; 5971a2e2f44SHong Zhang 5985d8b2efaSHong Zhang /* Initialize the statistics variables. */ 599e77caa6dSBarry Smith PetscStackCallExternalVoid("SuperLU:StatInit", StatInit(&lu->stat)); 6008914a3f7SHong Zhang lu->lwork = 0; /* allocate space internally by system malloc */ 6019ce50919SHong Zhang 6025d8b2efaSHong Zhang /* Allocate spaces (notice sizes are for the transpose) */ 6039566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m, &lu->etree)); 6049566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &lu->perm_r)); 6059566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m, &lu->perm_c)); 6069566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &lu->R)); 6079566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m, &lu->C)); 6085d8b2efaSHong Zhang 6095d8b2efaSHong Zhang /* create rhs and solution x without allocate space for .Store */ 6105d8b2efaSHong Zhang #if defined(PETSC_USE_COMPLEX) 6113cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 612e77caa6dSBarry Smith PetscStackCallExternalVoid("SuperLU:cCreate_Dense_Matrix(", cCreate_Dense_Matrix(&lu->B, m, 1, NULL, m, SLU_DN, SLU_C, SLU_GE)); 613e77caa6dSBarry Smith PetscStackCallExternalVoid("SuperLU:cCreate_Dense_Matrix(", cCreate_Dense_Matrix(&lu->X, m, 1, NULL, m, SLU_DN, SLU_C, SLU_GE)); 6143cb270beSHong Zhang #else 615e77caa6dSBarry Smith PetscStackCallExternalVoid("SuperLU:zCreate_Dense_Matrix", zCreate_Dense_Matrix(&lu->B, m, 1, NULL, m, SLU_DN, SLU_Z, SLU_GE)); 616e77caa6dSBarry Smith PetscStackCallExternalVoid("SuperLU:zCreate_Dense_Matrix", zCreate_Dense_Matrix(&lu->X, m, 1, NULL, m, SLU_DN, SLU_Z, SLU_GE)); 6173cb270beSHong Zhang #endif 6183cb270beSHong Zhang #else 6193cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 620e77caa6dSBarry Smith PetscStackCallExternalVoid("SuperLU:sCreate_Dense_Matrix", sCreate_Dense_Matrix(&lu->B, m, 1, NULL, m, SLU_DN, SLU_S, SLU_GE)); 621e77caa6dSBarry Smith PetscStackCallExternalVoid("SuperLU:sCreate_Dense_Matrix", sCreate_Dense_Matrix(&lu->X, m, 1, NULL, m, SLU_DN, SLU_S, SLU_GE)); 6225d8b2efaSHong Zhang #else 623e77caa6dSBarry Smith PetscStackCallExternalVoid("SuperLU:dCreate_Dense_Matrix", dCreate_Dense_Matrix(&lu->B, m, 1, NULL, m, SLU_DN, SLU_D, SLU_GE)); 624e77caa6dSBarry Smith PetscStackCallExternalVoid("SuperLU:dCreate_Dense_Matrix", dCreate_Dense_Matrix(&lu->X, m, 1, NULL, m, SLU_DN, SLU_D, SLU_GE)); 6255d8b2efaSHong Zhang #endif 6263cb270beSHong Zhang #endif 6275d8b2efaSHong Zhang 6289566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatFactorGetSolverType_C", MatFactorGetSolverType_seqaij_superlu)); 6299566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSuperluSetILUDropTol_C", MatSuperluSetILUDropTol_SuperLU)); 630245fece9SBarry Smith B->data = lu; 63120be8e61SHong Zhang 632899d7b4fSKris Buschelman *F = B; 63314b396bbSKris Buschelman PetscFunctionReturn(0); 63414b396bbSKris Buschelman } 635d954db57SHong Zhang 636d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGetFactor_seqsell_superlu(Mat A, MatFactorType ftype, Mat *F) 637d71ae5a4SJacob Faibussowitsch { 6382366e350SStefano Zampini Mat_SuperLU *lu; 6392366e350SStefano Zampini 6402366e350SStefano Zampini PetscFunctionBegin; 6419566063dSJacob Faibussowitsch PetscCall(MatGetFactor_seqaij_superlu(A, ftype, F)); 6422366e350SStefano Zampini lu = (Mat_SuperLU *)((*F)->data); 6432366e350SStefano Zampini lu->needconversion = PETSC_TRUE; 6442366e350SStefano Zampini PetscFunctionReturn(0); 6452366e350SStefano Zampini } 6462366e350SStefano Zampini 647d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatSolverTypeRegister_SuperLU(void) 648d71ae5a4SJacob Faibussowitsch { 64942c9c57cSBarry Smith PetscFunctionBegin; 6509566063dSJacob Faibussowitsch PetscCall(MatSolverTypeRegister(MATSOLVERSUPERLU, MATSEQAIJ, MAT_FACTOR_LU, MatGetFactor_seqaij_superlu)); 6519566063dSJacob Faibussowitsch PetscCall(MatSolverTypeRegister(MATSOLVERSUPERLU, MATSEQAIJ, MAT_FACTOR_ILU, MatGetFactor_seqaij_superlu)); 6529566063dSJacob Faibussowitsch PetscCall(MatSolverTypeRegister(MATSOLVERSUPERLU, MATSEQSELL, MAT_FACTOR_LU, MatGetFactor_seqsell_superlu)); 6539566063dSJacob Faibussowitsch PetscCall(MatSolverTypeRegister(MATSOLVERSUPERLU, MATSEQSELL, MAT_FACTOR_ILU, MatGetFactor_seqsell_superlu)); 65442c9c57cSBarry Smith PetscFunctionReturn(0); 65542c9c57cSBarry Smith } 656