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 */ 62860c79edSBarry Smith static PetscErrorCode MatView_Info_SuperLU(Mat A,PetscViewer viewer) 635a46d3faSBarry Smith { 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 93245fece9SBarry Smith PetscErrorCode MatSolve_SuperLU_Private(Mat A,Vec b,Vec x) 94245fece9SBarry Smith { 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); 1049566063dSJacob Faibussowitsch 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 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)); 105245fece9SBarry Smith 1069566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(x,&n)); 107245fece9SBarry Smith lu->B.ncol = 1; /* Set the number of right-hand side */ 108245fece9SBarry Smith if (lu->options.Equil && !lu->rhs_dup) { 109245fece9SBarry Smith /* superlu overwrites b when Equil is used, thus create rhs_dup to keep user's b unchanged */ 1109566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n,&lu->rhs_dup)); 111245fece9SBarry Smith } 112245fece9SBarry Smith if (lu->options.Equil) { 113245fece9SBarry Smith /* Copy b into rsh_dup */ 1149566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(b,&barray)); 1159566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(lu->rhs_dup,barray,n)); 1169566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(b,&barray)); 117245fece9SBarry Smith barray = lu->rhs_dup; 118245fece9SBarry Smith } else { 1199566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(b,&barray)); 120245fece9SBarry Smith } 1219566063dSJacob Faibussowitsch PetscCall(VecGetArray(x,&xarray)); 122245fece9SBarry Smith 123245fece9SBarry Smith #if defined(PETSC_USE_COMPLEX) 124245fece9SBarry Smith #if defined(PETSC_USE_REAL_SINGLE) 125245fece9SBarry Smith ((DNformat*)lu->B.Store)->nzval = (singlecomplex*)barray; 126245fece9SBarry Smith ((DNformat*)lu->X.Store)->nzval = (singlecomplex*)xarray; 127245fece9SBarry Smith #else 128245fece9SBarry Smith ((DNformat*)lu->B.Store)->nzval = (doublecomplex*)barray; 129245fece9SBarry Smith ((DNformat*)lu->X.Store)->nzval = (doublecomplex*)xarray; 130245fece9SBarry Smith #endif 131245fece9SBarry Smith #else 132245fece9SBarry Smith ((DNformat*)lu->B.Store)->nzval = (void*)barray; 133245fece9SBarry Smith ((DNformat*)lu->X.Store)->nzval = xarray; 134245fece9SBarry Smith #endif 135245fece9SBarry Smith 136245fece9SBarry Smith lu->options.Fact = FACTORED; /* Indicate the factored form of A is supplied. */ 137245fece9SBarry Smith if (A->factortype == MAT_FACTOR_LU) { 138245fece9SBarry Smith #if defined(PETSC_USE_COMPLEX) 139245fece9SBarry Smith #if defined(PETSC_USE_REAL_SINGLE) 140*e77caa6dSBarry Smith PetscStackCallExternalVoid("SuperLU:cgssvx",cgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 141245fece9SBarry Smith &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 142245fece9SBarry Smith &lu->Glu, &lu->mem_usage, &lu->stat, &info)); 143245fece9SBarry Smith #else 144*e77caa6dSBarry Smith PetscStackCallExternalVoid("SuperLU:zgssvx",zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 145245fece9SBarry Smith &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 146245fece9SBarry Smith &lu->Glu, &lu->mem_usage, &lu->stat, &info)); 147245fece9SBarry Smith #endif 148245fece9SBarry Smith #else 149245fece9SBarry Smith #if defined(PETSC_USE_REAL_SINGLE) 150*e77caa6dSBarry Smith PetscStackCallExternalVoid("SuperLU:sgssvx",sgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 151245fece9SBarry Smith &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 152245fece9SBarry Smith &lu->Glu, &lu->mem_usage, &lu->stat, &info)); 153245fece9SBarry Smith #else 154*e77caa6dSBarry Smith PetscStackCallExternalVoid("SuperLU:dgssvx",dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 155245fece9SBarry Smith &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 156245fece9SBarry Smith &lu->Glu,&lu->mem_usage, &lu->stat, &info)); 157245fece9SBarry Smith #endif 158245fece9SBarry Smith #endif 159245fece9SBarry Smith } else if (A->factortype == MAT_FACTOR_ILU) { 160245fece9SBarry Smith #if defined(PETSC_USE_COMPLEX) 161245fece9SBarry Smith #if defined(PETSC_USE_REAL_SINGLE) 162*e77caa6dSBarry Smith PetscStackCallExternalVoid("SuperLU:cgsisx",cgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 163245fece9SBarry Smith &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, 164245fece9SBarry Smith &lu->Glu, &lu->mem_usage, &lu->stat, &info)); 165245fece9SBarry Smith #else 166*e77caa6dSBarry Smith PetscStackCallExternalVoid("SuperLU:zgsisx",zgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 167245fece9SBarry Smith &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, 168245fece9SBarry Smith &lu->Glu, &lu->mem_usage, &lu->stat, &info)); 169245fece9SBarry Smith #endif 170245fece9SBarry Smith #else 171245fece9SBarry Smith #if defined(PETSC_USE_REAL_SINGLE) 172*e77caa6dSBarry Smith PetscStackCallExternalVoid("SuperLU:sgsisx",sgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 173245fece9SBarry Smith &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, 174245fece9SBarry Smith &lu->Glu, &lu->mem_usage, &lu->stat, &info)); 175245fece9SBarry Smith #else 176*e77caa6dSBarry Smith PetscStackCallExternalVoid("SuperLU:dgsisx",dgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 177245fece9SBarry Smith &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, 178245fece9SBarry Smith &lu->Glu, &lu->mem_usage, &lu->stat, &info)); 179245fece9SBarry Smith #endif 180245fece9SBarry Smith #endif 181245fece9SBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported"); 182245fece9SBarry Smith if (!lu->options.Equil) { 1839566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(b,&barray)); 184245fece9SBarry Smith } 1859566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(x,&xarray)); 186245fece9SBarry Smith 187245fece9SBarry Smith if (!info || info == lu->A.ncol+1) { 188245fece9SBarry Smith if (lu->options.IterRefine) { 1899566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Iterative Refinement:\n")); 1909566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF," %8s%8s%16s%16s\n", "rhs", "Steps", "FERR", "BERR")); 1919566063dSJacob Faibussowitsch for (i = 0; i < 1; ++i) PetscCall(PetscPrintf(PETSC_COMM_SELF," %8d%8d%16e%16e\n", i+1, lu->stat.RefineSteps, ferr, berr)); 192245fece9SBarry Smith } 193245fece9SBarry Smith } else if (info > 0) { 194245fece9SBarry Smith if (lu->lwork == -1) { 1959566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF," ** Estimated memory: %" PetscInt_FMT " bytes\n", info - lu->A.ncol)); 196245fece9SBarry Smith } else { 1979566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF," Warning: gssvx() returns info %" PetscInt_FMT "\n",info)); 198245fece9SBarry Smith } 1995f80ce2aSJacob 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); 200245fece9SBarry Smith 201245fece9SBarry Smith if (lu->options.PrintStat) { 2029566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"MatSolve__SuperLU():\n")); 203*e77caa6dSBarry Smith PetscStackCallExternalVoid("SuperLU:StatPrint",StatPrint(&lu->stat)); 204245fece9SBarry Smith } 205245fece9SBarry Smith PetscFunctionReturn(0); 206245fece9SBarry Smith } 207245fece9SBarry Smith 208245fece9SBarry Smith PetscErrorCode MatSolve_SuperLU(Mat A,Vec b,Vec x) 209245fece9SBarry Smith { 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 = TRANS; 2209566063dSJacob Faibussowitsch PetscCall(MatSolve_SuperLU_Private(A,b,x)); 221245fece9SBarry Smith PetscFunctionReturn(0); 222245fece9SBarry Smith } 223245fece9SBarry Smith 224245fece9SBarry Smith PetscErrorCode MatSolveTranspose_SuperLU(Mat A,Vec b,Vec x) 225245fece9SBarry Smith { 226245fece9SBarry Smith Mat_SuperLU *lu = (Mat_SuperLU*)A->data; 227245fece9SBarry Smith 228245fece9SBarry Smith PetscFunctionBegin; 229603e8f96SBarry Smith if (A->factorerrortype) { 2309566063dSJacob Faibussowitsch PetscCall(PetscInfo(A,"MatSolve is called with singular matrix factor, skip\n")); 2319566063dSJacob Faibussowitsch PetscCall(VecSetInf(x)); 232245fece9SBarry Smith PetscFunctionReturn(0); 233245fece9SBarry Smith } 234245fece9SBarry Smith 235245fece9SBarry Smith lu->options.Trans = NOTRANS; 2369566063dSJacob Faibussowitsch PetscCall(MatSolve_SuperLU_Private(A,b,x)); 237245fece9SBarry Smith PetscFunctionReturn(0); 238245fece9SBarry Smith } 239245fece9SBarry Smith 240245fece9SBarry Smith static PetscErrorCode MatLUFactorNumeric_SuperLU(Mat F,Mat A,const MatFactorInfo *info) 2415a46d3faSBarry Smith { 242245fece9SBarry Smith Mat_SuperLU *lu = (Mat_SuperLU*)F->data; 243cae5a23dSHong Zhang Mat_SeqAIJ *aa; 2445a46d3faSBarry Smith PetscInt sinfo; 2455a46d3faSBarry Smith PetscReal ferr, berr; 2465a46d3faSBarry Smith NCformat *Ustore; 2475a46d3faSBarry Smith SCformat *Lstore; 2485a46d3faSBarry Smith 2495a46d3faSBarry Smith PetscFunctionBegin; 2505a46d3faSBarry Smith if (lu->flg == SAME_NONZERO_PATTERN) { /* successing numerical factorization */ 2515a46d3faSBarry Smith lu->options.Fact = SamePattern; 2525a46d3faSBarry Smith /* Ref: ~SuperLU_3.0/EXAMPLE/dlinsolx2.c */ 2535a46d3faSBarry Smith Destroy_SuperMatrix_Store(&lu->A); 2541baa6e33SBarry Smith if (lu->A_dup) PetscCall(MatCopy_SeqAIJ(A,lu->A_dup,SAME_NONZERO_PATTERN)); 2552366e350SStefano Zampini 2565a46d3faSBarry Smith if (lu->lwork >= 0) { 257*e77caa6dSBarry Smith PetscStackCallExternalVoid("SuperLU:Destroy_SuperNode_Matrix",Destroy_SuperNode_Matrix(&lu->L)); 258*e77caa6dSBarry Smith PetscStackCallExternalVoid("SuperLU:Destroy_CompCol_Matrix",Destroy_CompCol_Matrix(&lu->U)); 2595a46d3faSBarry Smith lu->options.Fact = SamePattern; 2605a46d3faSBarry Smith } 2615a46d3faSBarry Smith } 2625a46d3faSBarry Smith 2635a46d3faSBarry Smith /* Create the SuperMatrix for lu->A=A^T: 2645a46d3faSBarry Smith Since SuperLU likes column-oriented matrices,we pass it the transpose, 2655a46d3faSBarry Smith and then solve A^T X = B in MatSolve(). */ 2662366e350SStefano Zampini if (lu->A_dup) { 267cae5a23dSHong Zhang aa = (Mat_SeqAIJ*)(lu->A_dup)->data; 268cae5a23dSHong Zhang } else { 269cae5a23dSHong Zhang aa = (Mat_SeqAIJ*)(A)->data; 270cae5a23dSHong Zhang } 2715a46d3faSBarry Smith #if defined(PETSC_USE_COMPLEX) 2723cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 273*e77caa6dSBarry 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)); 2743cb270beSHong Zhang #else 275*e77caa6dSBarry 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)); 2763cb270beSHong Zhang #endif 2773cb270beSHong Zhang #else 2783cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 279*e77caa6dSBarry 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)); 2805a46d3faSBarry Smith #else 281*e77caa6dSBarry 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)); 2825a46d3faSBarry Smith #endif 2833cb270beSHong Zhang #endif 2845a46d3faSBarry Smith 2855a46d3faSBarry Smith /* Numerical factorization */ 2865a46d3faSBarry Smith lu->B.ncol = 0; /* Indicate not to solve the system */ 287d5f3da31SBarry Smith if (F->factortype == MAT_FACTOR_LU) { 2885a46d3faSBarry Smith #if defined(PETSC_USE_COMPLEX) 2893cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 290*e77caa6dSBarry Smith PetscStackCallExternalVoid("SuperLU:cgssvx",cgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 2913cb270beSHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 2925db81dd2SSatish Balay &lu->Glu, &lu->mem_usage, &lu->stat, &sinfo)); 2933cb270beSHong Zhang #else 294*e77caa6dSBarry Smith PetscStackCallExternalVoid("SuperLU:zgssvx",zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 2955a46d3faSBarry Smith &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 2965db81dd2SSatish Balay &lu->Glu, &lu->mem_usage, &lu->stat, &sinfo)); 2973cb270beSHong Zhang #endif 2983cb270beSHong Zhang #else 2993cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 300*e77caa6dSBarry Smith PetscStackCallExternalVoid("SuperLU:sgssvx",sgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 3013cb270beSHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 3025db81dd2SSatish Balay &lu->Glu, &lu->mem_usage, &lu->stat, &sinfo)); 3035a46d3faSBarry Smith #else 304*e77caa6dSBarry Smith PetscStackCallExternalVoid("SuperLU:dgssvx",dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 3055a46d3faSBarry Smith &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 306c147c726SHong Zhang &lu->Glu,&lu->mem_usage, &lu->stat, &sinfo)); 3075a46d3faSBarry Smith #endif 3083cb270beSHong Zhang #endif 309d5f3da31SBarry Smith } else if (F->factortype == MAT_FACTOR_ILU) { 310cffbb591SHong Zhang /* Compute the incomplete factorization, condition number and pivot growth */ 311cffbb591SHong Zhang #if defined(PETSC_USE_COMPLEX) 3123cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 313*e77caa6dSBarry Smith PetscStackCallExternalVoid("SuperLU:cgsisx",cgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r,lu->etree, lu->equed, lu->R, lu->C, 3143cb270beSHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, 3155db81dd2SSatish Balay &lu->Glu, &lu->mem_usage, &lu->stat, &sinfo)); 3163cb270beSHong Zhang #else 317*e77caa6dSBarry Smith PetscStackCallExternalVoid("SuperLU:zgsisx",zgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r,lu->etree, lu->equed, lu->R, lu->C, 318cffbb591SHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, 3195db81dd2SSatish Balay &lu->Glu, &lu->mem_usage, &lu->stat, &sinfo)); 3203cb270beSHong Zhang #endif 3213cb270beSHong Zhang #else 3223cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 323*e77caa6dSBarry Smith PetscStackCallExternalVoid("SuperLU:sgsisx",sgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 3243cb270beSHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, 3255db81dd2SSatish Balay &lu->Glu, &lu->mem_usage, &lu->stat, &sinfo)); 326cffbb591SHong Zhang #else 327*e77caa6dSBarry Smith PetscStackCallExternalVoid("SuperLU:dgsisx",dgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 328cffbb591SHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, 329c147c726SHong Zhang &lu->Glu, &lu->mem_usage, &lu->stat, &sinfo)); 330cffbb591SHong Zhang #endif 3313cb270beSHong Zhang #endif 332f23aa3ddSBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported"); 3335a46d3faSBarry Smith if (!sinfo || sinfo == lu->A.ncol+1) { 3342205254eSKarl Rupp if (lu->options.PivotGrowth) { 3359566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF," Recip. pivot growth = %e\n", lu->rpg)); 3362205254eSKarl Rupp } 3372205254eSKarl Rupp if (lu->options.ConditionNumber) { 3389566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF," Recip. condition number = %e\n", lu->rcond)); 3392205254eSKarl Rupp } 3405a46d3faSBarry Smith } else if (sinfo > 0) { 341675d1226SHong Zhang if (A->erroriffailure) { 34298921bdaSJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot in row %" PetscInt_FMT,sinfo); 343675d1226SHong Zhang } else { 344675d1226SHong Zhang if (sinfo <= lu->A.ncol) { 345675d1226SHong Zhang if (lu->options.ILU_FillTol == 0.0) { 346603e8f96SBarry Smith F->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 347675d1226SHong Zhang } 3489566063dSJacob Faibussowitsch PetscCall(PetscInfo(F,"Number of zero pivots %" PetscInt_FMT ", ILU_FillTol %g\n",sinfo,lu->options.ILU_FillTol)); 349675d1226SHong Zhang } else if (sinfo == lu->A.ncol + 1) { 350675d1226SHong Zhang /* 351675d1226SHong Zhang U is nonsingular, but RCOND is less than machine 352675d1226SHong Zhang precision, meaning that the matrix is singular to 353675d1226SHong Zhang working precision. Nevertheless, the solution and 354675d1226SHong Zhang error bounds are computed because there are a number 355675d1226SHong Zhang of situations where the computed solution can be more 356675d1226SHong Zhang accurate than the value of RCOND would suggest. 357675d1226SHong Zhang */ 3589566063dSJacob Faibussowitsch PetscCall(PetscInfo(F,"Matrix factor U is nonsingular, but is singular to working precision. The solution is computed. info %" PetscInt_FMT,sinfo)); 359675d1226SHong Zhang } else { /* sinfo > lu->A.ncol + 1 */ 360603e8f96SBarry Smith F->factorerrortype = MAT_FACTOR_OUTMEMORY; 3619566063dSJacob Faibussowitsch PetscCall(PetscInfo(F,"Number of bytes allocated when memory allocation fails %" PetscInt_FMT "\n",sinfo)); 362675d1226SHong Zhang } 363675d1226SHong Zhang } 36498921bdaSJacob 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); 3655a46d3faSBarry Smith 3665a46d3faSBarry Smith if (lu->options.PrintStat) { 3679566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"MatLUFactorNumeric_SuperLU():\n")); 368*e77caa6dSBarry Smith PetscStackCallExternalVoid("SuperLU:StatPrint",StatPrint(&lu->stat)); 3695a46d3faSBarry Smith Lstore = (SCformat*) lu->L.Store; 3705a46d3faSBarry Smith Ustore = (NCformat*) lu->U.Store; 3719566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF," No of nonzeros in factor L = %" PetscInt_FMT "\n", Lstore->nnz)); 3729566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF," No of nonzeros in factor U = %" PetscInt_FMT "\n", Ustore->nnz)); 3739566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF," No of nonzeros in L+U = %" PetscInt_FMT "\n", Lstore->nnz + Ustore->nnz - lu->A.ncol)); 3749566063dSJacob 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)); 3755a46d3faSBarry Smith } 3765a46d3faSBarry Smith 3775a46d3faSBarry Smith lu->flg = SAME_NONZERO_PATTERN; 3781d5ca7f3SHong Zhang F->ops->solve = MatSolve_SuperLU; 3791d5ca7f3SHong Zhang F->ops->solvetranspose = MatSolveTranspose_SuperLU; 3801aef8b4cSStefano Zampini F->ops->matsolve = NULL; 3815a46d3faSBarry Smith PetscFunctionReturn(0); 3825a46d3faSBarry Smith } 3835a46d3faSBarry Smith 384245fece9SBarry Smith static PetscErrorCode MatDestroy_SuperLU(Mat A) 38514b396bbSKris Buschelman { 386245fece9SBarry Smith Mat_SuperLU *lu=(Mat_SuperLU*)A->data; 38714b396bbSKris Buschelman 38814b396bbSKris Buschelman PetscFunctionBegin; 389245fece9SBarry Smith if (lu->CleanUpSuperLU) { /* Free the SuperLU datastructures */ 390*e77caa6dSBarry Smith PetscStackCallExternalVoid("SuperLU:Destroy_SuperMatrix_Store",Destroy_SuperMatrix_Store(&lu->A)); 391*e77caa6dSBarry Smith PetscStackCallExternalVoid("SuperLU:Destroy_SuperMatrix_Store",Destroy_SuperMatrix_Store(&lu->B)); 392*e77caa6dSBarry Smith PetscStackCallExternalVoid("SuperLU:Destroy_SuperMatrix_Store",Destroy_SuperMatrix_Store(&lu->X)); 393*e77caa6dSBarry Smith PetscStackCallExternalVoid("SuperLU:StatFree",StatFree(&lu->stat)); 3940e742b69SHong Zhang if (lu->lwork >= 0) { 395*e77caa6dSBarry Smith PetscStackCallExternalVoid("SuperLU:Destroy_SuperNode_Matrix",Destroy_SuperNode_Matrix(&lu->L)); 396*e77caa6dSBarry Smith PetscStackCallExternalVoid("SuperLU:Destroy_CompCol_Matrix",Destroy_CompCol_Matrix(&lu->U)); 3970e742b69SHong Zhang } 3980e742b69SHong Zhang } 3999566063dSJacob Faibussowitsch PetscCall(PetscFree(lu->etree)); 4009566063dSJacob Faibussowitsch PetscCall(PetscFree(lu->perm_r)); 4019566063dSJacob Faibussowitsch PetscCall(PetscFree(lu->perm_c)); 4029566063dSJacob Faibussowitsch PetscCall(PetscFree(lu->R)); 4039566063dSJacob Faibussowitsch PetscCall(PetscFree(lu->C)); 4049566063dSJacob Faibussowitsch PetscCall(PetscFree(lu->rhs_dup)); 4059566063dSJacob Faibussowitsch PetscCall(MatDestroy(&lu->A_dup)); 4069566063dSJacob Faibussowitsch PetscCall(PetscFree(A->data)); 4070e742b69SHong Zhang 408d954db57SHong Zhang /* clear composed functions */ 4099566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatFactorGetSolverType_C",NULL)); 4109566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatSuperluSetILUDropTol_C",NULL)); 41114b396bbSKris Buschelman PetscFunctionReturn(0); 41214b396bbSKris Buschelman } 41314b396bbSKris Buschelman 414245fece9SBarry Smith static PetscErrorCode MatView_SuperLU(Mat A,PetscViewer viewer) 41514b396bbSKris Buschelman { 416ace3abfcSBarry Smith PetscBool iascii; 41714b396bbSKris Buschelman PetscViewerFormat format; 41814b396bbSKris Buschelman 41914b396bbSKris Buschelman PetscFunctionBegin; 4209566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii)); 42132077d6dSBarry Smith if (iascii) { 4229566063dSJacob Faibussowitsch PetscCall(PetscViewerGetFormat(viewer,&format)); 4232f59403fSHong Zhang if (format == PETSC_VIEWER_ASCII_INFO) { 4249566063dSJacob Faibussowitsch PetscCall(MatView_Info_SuperLU(A,viewer)); 42514b396bbSKris Buschelman } 42614b396bbSKris Buschelman } 42714b396bbSKris Buschelman PetscFunctionReturn(0); 42814b396bbSKris Buschelman } 42914b396bbSKris Buschelman 430e0b74bf9SHong Zhang PetscErrorCode MatMatSolve_SuperLU(Mat A,Mat B,Mat X) 431e0b74bf9SHong Zhang { 432245fece9SBarry Smith Mat_SuperLU *lu = (Mat_SuperLU*)A->data; 433cd723cd1SBarry Smith PetscBool flg; 434e0b74bf9SHong Zhang 435e0b74bf9SHong Zhang PetscFunctionBegin; 4369566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompareAny((PetscObject)B,&flg,MATSEQDENSE,MATMPIDENSE,NULL)); 4375f80ce2aSJacob Faibussowitsch PetscCheck(flg,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix"); 4389566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,NULL)); 4395f80ce2aSJacob Faibussowitsch PetscCheck(flg,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix"); 4402205254eSKarl Rupp lu->options.Trans = TRANS; 441e0b74bf9SHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatMatSolve_SuperLU() is not implemented yet"); 442e0b74bf9SHong Zhang PetscFunctionReturn(0); 443e0b74bf9SHong Zhang } 444e0b74bf9SHong Zhang 445245fece9SBarry Smith static PetscErrorCode MatLUFactorSymbolic_SuperLU(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info) 44614b396bbSKris Buschelman { 447245fece9SBarry Smith Mat_SuperLU *lu = (Mat_SuperLU*)(F->data); 44826cc229bSBarry Smith PetscInt indx; 44926cc229bSBarry Smith PetscBool flg,set; 45026cc229bSBarry Smith PetscReal real_input; 45126cc229bSBarry Smith const char *colperm[] ={"NATURAL","MMD_ATA","MMD_AT_PLUS_A","COLAMD"}; /* MY_PERMC - not supported by the petsc interface yet */ 45226cc229bSBarry Smith const char *iterrefine[]={"NOREFINE", "SINGLE", "DOUBLE", "EXTRA"}; 45326cc229bSBarry Smith const char *rowperm[] ={"NOROWPERM", "LargeDiag"}; /* MY_PERMC - not supported by the petsc interface yet */ 454b24902e0SBarry Smith 455b24902e0SBarry Smith PetscFunctionBegin; 45626cc229bSBarry Smith /* Set options to F */ 45726cc229bSBarry Smith PetscOptionsBegin(PetscObjectComm((PetscObject)F),((PetscObject)F)->prefix,"SuperLU Options","Mat"); 45826cc229bSBarry Smith PetscCall(PetscOptionsBool("-mat_superlu_equil","Equil","None",(PetscBool)lu->options.Equil,(PetscBool*)&lu->options.Equil,NULL)); 45926cc229bSBarry Smith PetscCall(PetscOptionsEList("-mat_superlu_colperm","ColPerm","None",colperm,4,colperm[3],&indx,&flg)); 46026cc229bSBarry Smith if (flg) lu->options.ColPerm = (colperm_t)indx; 46126cc229bSBarry Smith PetscCall(PetscOptionsEList("-mat_superlu_iterrefine","IterRefine","None",iterrefine,4,iterrefine[0],&indx,&flg)); 46226cc229bSBarry Smith if (flg) lu->options.IterRefine = (IterRefine_t)indx; 46326cc229bSBarry Smith PetscCall(PetscOptionsBool("-mat_superlu_symmetricmode","SymmetricMode","None",(PetscBool)lu->options.SymmetricMode,&flg,&set)); 46426cc229bSBarry Smith if (set && flg) lu->options.SymmetricMode = YES; 46526cc229bSBarry Smith PetscCall(PetscOptionsReal("-mat_superlu_diagpivotthresh","DiagPivotThresh","None",lu->options.DiagPivotThresh,&real_input,&flg)); 46626cc229bSBarry Smith if (flg) lu->options.DiagPivotThresh = (double) real_input; 46726cc229bSBarry Smith PetscCall(PetscOptionsBool("-mat_superlu_pivotgrowth","PivotGrowth","None",(PetscBool)lu->options.PivotGrowth,&flg,&set)); 46826cc229bSBarry Smith if (set && flg) lu->options.PivotGrowth = YES; 46926cc229bSBarry Smith PetscCall(PetscOptionsBool("-mat_superlu_conditionnumber","ConditionNumber","None",(PetscBool)lu->options.ConditionNumber,&flg,&set)); 47026cc229bSBarry Smith if (set && flg) lu->options.ConditionNumber = YES; 47126cc229bSBarry Smith PetscCall(PetscOptionsEList("-mat_superlu_rowperm","rowperm","None",rowperm,2,rowperm[lu->options.RowPerm],&indx,&flg)); 47226cc229bSBarry Smith if (flg) lu->options.RowPerm = (rowperm_t)indx; 47326cc229bSBarry Smith PetscCall(PetscOptionsBool("-mat_superlu_replacetinypivot","ReplaceTinyPivot","None",(PetscBool)lu->options.ReplaceTinyPivot,&flg,&set)); 47426cc229bSBarry Smith if (set && flg) lu->options.ReplaceTinyPivot = YES; 47526cc229bSBarry Smith PetscCall(PetscOptionsBool("-mat_superlu_printstat","PrintStat","None",(PetscBool)lu->options.PrintStat,&flg,&set)); 47626cc229bSBarry Smith if (set && flg) lu->options.PrintStat = YES; 47726cc229bSBarry Smith PetscCall(PetscOptionsInt("-mat_superlu_lwork","size of work array in bytes used by factorization","None",lu->lwork,&lu->lwork,NULL)); 47826cc229bSBarry Smith if (lu->lwork > 0) { 47926cc229bSBarry Smith /* lwork is in bytes, hence PetscMalloc() is used here, not PetscMalloc1()*/ 48026cc229bSBarry Smith PetscCall(PetscMalloc(lu->lwork,&lu->work)); 48126cc229bSBarry Smith } else if (lu->lwork != 0 && lu->lwork != -1) { 48226cc229bSBarry Smith PetscCall(PetscPrintf(PETSC_COMM_SELF," Warning: lwork %" PetscInt_FMT " is not supported by SUPERLU. The default lwork=0 is used.\n",lu->lwork)); 48326cc229bSBarry Smith lu->lwork = 0; 48426cc229bSBarry Smith } 48526cc229bSBarry Smith /* ilu options */ 48626cc229bSBarry Smith PetscCall(PetscOptionsReal("-mat_superlu_ilu_droptol","ILU_DropTol","None",lu->options.ILU_DropTol,&real_input,&flg)); 48726cc229bSBarry Smith if (flg) lu->options.ILU_DropTol = (double) real_input; 48826cc229bSBarry Smith PetscCall(PetscOptionsReal("-mat_superlu_ilu_filltol","ILU_FillTol","None",lu->options.ILU_FillTol,&real_input,&flg)); 48926cc229bSBarry Smith if (flg) lu->options.ILU_FillTol = (double) real_input; 49026cc229bSBarry Smith PetscCall(PetscOptionsReal("-mat_superlu_ilu_fillfactor","ILU_FillFactor","None",lu->options.ILU_FillFactor,&real_input,&flg)); 49126cc229bSBarry Smith if (flg) lu->options.ILU_FillFactor = (double) real_input; 49226cc229bSBarry Smith PetscCall(PetscOptionsInt("-mat_superlu_ilu_droprull","ILU_DropRule","None",lu->options.ILU_DropRule,&lu->options.ILU_DropRule,NULL)); 49326cc229bSBarry Smith PetscCall(PetscOptionsInt("-mat_superlu_ilu_norm","ILU_Norm","None",lu->options.ILU_Norm,&indx,&flg)); 49426cc229bSBarry Smith if (flg) lu->options.ILU_Norm = (norm_t)indx; 49526cc229bSBarry Smith PetscCall(PetscOptionsInt("-mat_superlu_ilu_milu","ILU_MILU","None",lu->options.ILU_MILU,&indx,&flg)); 49626cc229bSBarry Smith if (flg) lu->options.ILU_MILU = (milu_t)indx; 49726cc229bSBarry Smith PetscOptionsEnd(); 49826cc229bSBarry Smith 499b24902e0SBarry Smith lu->flg = DIFFERENT_NONZERO_PATTERN; 500b24902e0SBarry Smith lu->CleanUpSuperLU = PETSC_TRUE; 5011d5ca7f3SHong Zhang F->ops->lufactornumeric = MatLUFactorNumeric_SuperLU; 5022366e350SStefano Zampini 5032366e350SStefano Zampini /* if we are here, the nonzero pattern has changed unless the user explicitly called MatLUFactorSymbolic */ 5049566063dSJacob Faibussowitsch PetscCall(MatDestroy(&lu->A_dup)); 5052366e350SStefano Zampini if (lu->needconversion) { 5069566063dSJacob Faibussowitsch PetscCall(MatConvert(A,MATSEQAIJ,MAT_INITIAL_MATRIX,&lu->A_dup)); 5072366e350SStefano Zampini } 5082366e350SStefano 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 */ 5099566063dSJacob Faibussowitsch PetscCall(MatDuplicate_SeqAIJ(A,MAT_COPY_VALUES,&lu->A_dup)); 5102366e350SStefano Zampini } 511b24902e0SBarry Smith PetscFunctionReturn(0); 512b24902e0SBarry Smith } 513b24902e0SBarry Smith 514b2573a8aSBarry Smith static PetscErrorCode MatSuperluSetILUDropTol_SuperLU(Mat F,PetscReal dtol) 5155ccb76cbSHong Zhang { 516245fece9SBarry Smith Mat_SuperLU *lu= (Mat_SuperLU*)F->data; 5175ccb76cbSHong Zhang 5185ccb76cbSHong Zhang PetscFunctionBegin; 5195ccb76cbSHong Zhang lu->options.ILU_DropTol = dtol; 5205ccb76cbSHong Zhang PetscFunctionReturn(0); 5215ccb76cbSHong Zhang } 5225ccb76cbSHong Zhang 5235ccb76cbSHong Zhang /*@ 5245ccb76cbSHong Zhang MatSuperluSetILUDropTol - Set SuperLU ILU drop tolerance 525147403d9SBarry Smith 5265ccb76cbSHong Zhang Logically Collective on Mat 5275ccb76cbSHong Zhang 5285ccb76cbSHong Zhang Input Parameters: 5295ccb76cbSHong Zhang + F - the factored matrix obtained by calling MatGetFactor() from PETSc-SuperLU interface 5305ccb76cbSHong Zhang - dtol - drop tolerance 5315ccb76cbSHong Zhang 5325ccb76cbSHong Zhang Options Database: 533147403d9SBarry Smith . -mat_superlu_ilu_droptol <dtol> - the drop tolerance 5345ccb76cbSHong Zhang 5355ccb76cbSHong Zhang Level: beginner 5365ccb76cbSHong Zhang 53796a0c994SBarry Smith References: 538606c0280SSatish Balay . * - SuperLU Users' Guide 5395ccb76cbSHong Zhang 540db781477SPatrick Sanan .seealso: `MatGetFactor()` 5415ccb76cbSHong Zhang @*/ 5425ccb76cbSHong Zhang PetscErrorCode MatSuperluSetILUDropTol(Mat F,PetscReal dtol) 5435ccb76cbSHong Zhang { 5445ccb76cbSHong Zhang PetscFunctionBegin; 5455ccb76cbSHong Zhang PetscValidHeaderSpecific(F,MAT_CLASSID,1); 54669fbec6eSBarry Smith PetscValidLogicalCollectiveReal(F,dtol,2); 547cac4c232SBarry Smith PetscTryMethod(F,"MatSuperluSetILUDropTol_C",(Mat,PetscReal),(F,dtol)); 5485ccb76cbSHong Zhang PetscFunctionReturn(0); 5495ccb76cbSHong Zhang } 5505ccb76cbSHong Zhang 551ea799195SBarry Smith PetscErrorCode MatFactorGetSolverType_seqaij_superlu(Mat A,MatSolverType *type) 55235bd34faSBarry Smith { 55335bd34faSBarry Smith PetscFunctionBegin; 5542692d6eeSBarry Smith *type = MATSOLVERSUPERLU; 55535bd34faSBarry Smith PetscFunctionReturn(0); 55635bd34faSBarry Smith } 55735bd34faSBarry Smith 558b24902e0SBarry Smith /*MC 559ba992d64SSatish Balay MATSOLVERSUPERLU = "superlu" - A solver package providing solvers LU and ILU for sequential matrices 560b24902e0SBarry Smith via the external package SuperLU. 561b24902e0SBarry Smith 562e2e64c6bSBarry Smith Use ./configure --download-superlu to have PETSc installed with SuperLU 563b24902e0SBarry Smith 5643ca39a21SBarry Smith Use -pc_type lu -pc_factor_mat_solver_type superlu to use this direct solver 565c2b89b5dSBarry Smith 566b24902e0SBarry Smith Options Database Keys: 567e08999f5SMatthew G Knepley + -mat_superlu_equil <FALSE> - Equil (None) 568e08999f5SMatthew G Knepley . -mat_superlu_colperm <COLAMD> - (choose one of) NATURAL MMD_ATA MMD_AT_PLUS_A COLAMD 569e08999f5SMatthew G Knepley . -mat_superlu_iterrefine <NOREFINE> - (choose one of) NOREFINE SINGLE DOUBLE EXTRA 570e08999f5SMatthew G Knepley . -mat_superlu_symmetricmode: <FALSE> - SymmetricMode (None) 571e08999f5SMatthew G Knepley . -mat_superlu_diagpivotthresh <1> - DiagPivotThresh (None) 572e08999f5SMatthew G Knepley . -mat_superlu_pivotgrowth <FALSE> - PivotGrowth (None) 573e08999f5SMatthew G Knepley . -mat_superlu_conditionnumber <FALSE> - ConditionNumber (None) 574e08999f5SMatthew G Knepley . -mat_superlu_rowperm <NOROWPERM> - (choose one of) NOROWPERM LargeDiag 575e08999f5SMatthew G Knepley . -mat_superlu_replacetinypivot <FALSE> - ReplaceTinyPivot (None) 576e08999f5SMatthew G Knepley . -mat_superlu_printstat <FALSE> - PrintStat (None) 577e08999f5SMatthew G Knepley . -mat_superlu_lwork <0> - size of work array in bytes used by factorization (None) 578e08999f5SMatthew G Knepley . -mat_superlu_ilu_droptol <0> - ILU_DropTol (None) 579e08999f5SMatthew G Knepley . -mat_superlu_ilu_filltol <0> - ILU_FillTol (None) 580e08999f5SMatthew G Knepley . -mat_superlu_ilu_fillfactor <0> - ILU_FillFactor (None) 581e08999f5SMatthew G Knepley . -mat_superlu_ilu_droprull <0> - ILU_DropRule (None) 582e08999f5SMatthew G Knepley . -mat_superlu_ilu_norm <0> - ILU_Norm (None) 583e08999f5SMatthew G Knepley - -mat_superlu_ilu_milu <0> - ILU_MILU (None) 584b24902e0SBarry Smith 58595452b02SPatrick Sanan Notes: 58695452b02SPatrick Sanan Do not confuse this with MATSOLVERSUPERLU_DIST which is for parallel sparse solves 5875c9eb25fSBarry Smith 5882c7c0729SBarry Smith Cannot currently use ordering provided by PETSc. 5892c7c0729SBarry Smith 590b24902e0SBarry Smith Level: beginner 591b24902e0SBarry Smith 592db781477SPatrick Sanan .seealso: `PCLU`, `PCILU`, `MATSOLVERSUPERLU_DIST`, `MATSOLVERMUMPS`, `PCFactorSetMatSolverType()`, `MatSolverType` 593b24902e0SBarry Smith M*/ 594b24902e0SBarry Smith 595db87b0f2SBarry Smith static PetscErrorCode MatGetFactor_seqaij_superlu(Mat A,MatFactorType ftype,Mat *F) 596b24902e0SBarry Smith { 59714b396bbSKris Buschelman Mat B; 598f0c56d0fSKris Buschelman Mat_SuperLU *lu; 59926cc229bSBarry Smith PetscInt m=A->rmap->n,n=A->cmap->n; 60014b396bbSKris Buschelman 60114b396bbSKris Buschelman PetscFunctionBegin; 6029566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A),&B)); 6039566063dSJacob Faibussowitsch PetscCall(MatSetSizes(B,A->rmap->n,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE)); 6049566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy("superlu",&((PetscObject)B)->type_name)); 6059566063dSJacob Faibussowitsch PetscCall(MatSetUp(B)); 60666e17bc3SBarry Smith B->trivialsymbolic = PETSC_TRUE; 607cffbb591SHong Zhang if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU) { 608b24902e0SBarry Smith B->ops->lufactorsymbolic = MatLUFactorSymbolic_SuperLU; 609cffbb591SHong Zhang B->ops->ilufactorsymbolic = MatLUFactorSymbolic_SuperLU; 610b3a44c85SBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported"); 611cffbb591SHong Zhang 6129566063dSJacob Faibussowitsch PetscCall(PetscFree(B->solvertype)); 6139566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(MATSOLVERSUPERLU,&B->solvertype)); 61400c67f3bSHong Zhang 615b9c12af5SBarry Smith B->ops->getinfo = MatGetInfo_External; 616b24902e0SBarry Smith B->ops->destroy = MatDestroy_SuperLU; 6173519fcdcSHong Zhang B->ops->view = MatView_SuperLU; 618d5f3da31SBarry Smith B->factortype = ftype; 61994b7f48cSBarry Smith B->assembled = PETSC_TRUE; /* required by -ksp_view */ 6205c9eb25fSBarry Smith B->preallocated = PETSC_TRUE; 62114b396bbSKris Buschelman 6229566063dSJacob Faibussowitsch PetscCall(PetscNewLog(B,&lu)); 623cae5a23dSHong Zhang 624cffbb591SHong Zhang if (ftype == MAT_FACTOR_LU) { 6259ce50919SHong Zhang set_default_options(&lu->options); 6263d6581fbSHong Zhang /* Comments from SuperLU_4.0/SRC/dgssvx.c: 6273d6581fbSHong Zhang "Whether or not the system will be equilibrated depends on the 6283d6581fbSHong Zhang scaling of the matrix A, but if equilibration is used, A is 6293d6581fbSHong Zhang overwritten by diag(R)*A*diag(C) and B by diag(R)*B 6303d6581fbSHong Zhang (if options->Trans=NOTRANS) or diag(C)*B (if options->Trans = TRANS or CONJ)." 6313d6581fbSHong Zhang We set 'options.Equil = NO' as default because additional space is needed for it. 6323d6581fbSHong Zhang */ 6333d6581fbSHong Zhang lu->options.Equil = NO; 634cffbb591SHong Zhang } else if (ftype == MAT_FACTOR_ILU) { 6350c74a584SJed Brown /* Set the default input options of ilu: */ 636*e77caa6dSBarry Smith PetscStackCallExternalVoid("SuperLU:ilu_set_default_options",ilu_set_default_options(&lu->options)); 637cffbb591SHong Zhang } 6389ce50919SHong Zhang lu->options.PrintStat = NO; 6391a2e2f44SHong Zhang 6405d8b2efaSHong Zhang /* Initialize the statistics variables. */ 641*e77caa6dSBarry Smith PetscStackCallExternalVoid("SuperLU:StatInit",StatInit(&lu->stat)); 6428914a3f7SHong Zhang lu->lwork = 0; /* allocate space internally by system malloc */ 6439ce50919SHong Zhang 6445d8b2efaSHong Zhang /* Allocate spaces (notice sizes are for the transpose) */ 6459566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m,&lu->etree)); 6469566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n,&lu->perm_r)); 6479566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m,&lu->perm_c)); 6489566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n,&lu->R)); 6499566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m,&lu->C)); 6505d8b2efaSHong Zhang 6515d8b2efaSHong Zhang /* create rhs and solution x without allocate space for .Store */ 6525d8b2efaSHong Zhang #if defined(PETSC_USE_COMPLEX) 6533cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 654*e77caa6dSBarry Smith PetscStackCallExternalVoid("SuperLU:cCreate_Dense_Matrix(",cCreate_Dense_Matrix(&lu->B, m, 1, NULL, m, SLU_DN, SLU_C, SLU_GE)); 655*e77caa6dSBarry Smith PetscStackCallExternalVoid("SuperLU:cCreate_Dense_Matrix(",cCreate_Dense_Matrix(&lu->X, m, 1, NULL, m, SLU_DN, SLU_C, SLU_GE)); 6563cb270beSHong Zhang #else 657*e77caa6dSBarry Smith PetscStackCallExternalVoid("SuperLU:zCreate_Dense_Matrix",zCreate_Dense_Matrix(&lu->B, m, 1, NULL, m, SLU_DN, SLU_Z, SLU_GE)); 658*e77caa6dSBarry Smith PetscStackCallExternalVoid("SuperLU:zCreate_Dense_Matrix",zCreate_Dense_Matrix(&lu->X, m, 1, NULL, m, SLU_DN, SLU_Z, SLU_GE)); 6593cb270beSHong Zhang #endif 6603cb270beSHong Zhang #else 6613cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 662*e77caa6dSBarry Smith PetscStackCallExternalVoid("SuperLU:sCreate_Dense_Matrix",sCreate_Dense_Matrix(&lu->B, m, 1, NULL, m, SLU_DN, SLU_S, SLU_GE)); 663*e77caa6dSBarry Smith PetscStackCallExternalVoid("SuperLU:sCreate_Dense_Matrix",sCreate_Dense_Matrix(&lu->X, m, 1, NULL, m, SLU_DN, SLU_S, SLU_GE)); 6645d8b2efaSHong Zhang #else 665*e77caa6dSBarry Smith PetscStackCallExternalVoid("SuperLU:dCreate_Dense_Matrix",dCreate_Dense_Matrix(&lu->B, m, 1, NULL, m, SLU_DN, SLU_D, SLU_GE)); 666*e77caa6dSBarry Smith PetscStackCallExternalVoid("SuperLU:dCreate_Dense_Matrix",dCreate_Dense_Matrix(&lu->X, m, 1, NULL, m, SLU_DN, SLU_D, SLU_GE)); 6675d8b2efaSHong Zhang #endif 6683cb270beSHong Zhang #endif 6695d8b2efaSHong Zhang 6709566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverType_C",MatFactorGetSolverType_seqaij_superlu)); 6719566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatSuperluSetILUDropTol_C",MatSuperluSetILUDropTol_SuperLU)); 672245fece9SBarry Smith B->data = lu; 67320be8e61SHong Zhang 674899d7b4fSKris Buschelman *F = B; 67514b396bbSKris Buschelman PetscFunctionReturn(0); 67614b396bbSKris Buschelman } 677d954db57SHong Zhang 6782366e350SStefano Zampini static PetscErrorCode MatGetFactor_seqsell_superlu(Mat A,MatFactorType ftype,Mat *F) 6792366e350SStefano Zampini { 6802366e350SStefano Zampini Mat_SuperLU *lu; 6812366e350SStefano Zampini 6822366e350SStefano Zampini PetscFunctionBegin; 6839566063dSJacob Faibussowitsch PetscCall(MatGetFactor_seqaij_superlu(A,ftype,F)); 6842366e350SStefano Zampini lu = (Mat_SuperLU*)((*F)->data); 6852366e350SStefano Zampini lu->needconversion = PETSC_TRUE; 6862366e350SStefano Zampini PetscFunctionReturn(0); 6872366e350SStefano Zampini } 6882366e350SStefano Zampini 6893ca39a21SBarry Smith PETSC_EXTERN PetscErrorCode MatSolverTypeRegister_SuperLU(void) 69042c9c57cSBarry Smith { 69142c9c57cSBarry Smith PetscFunctionBegin; 6929566063dSJacob Faibussowitsch PetscCall(MatSolverTypeRegister(MATSOLVERSUPERLU,MATSEQAIJ,MAT_FACTOR_LU,MatGetFactor_seqaij_superlu)); 6939566063dSJacob Faibussowitsch PetscCall(MatSolverTypeRegister(MATSOLVERSUPERLU,MATSEQAIJ,MAT_FACTOR_ILU,MatGetFactor_seqaij_superlu)); 6949566063dSJacob Faibussowitsch PetscCall(MatSolverTypeRegister(MATSOLVERSUPERLU,MATSEQSELL,MAT_FACTOR_LU,MatGetFactor_seqsell_superlu)); 6959566063dSJacob Faibussowitsch PetscCall(MatSolverTypeRegister(MATSOLVERSUPERLU,MATSEQSELL,MAT_FACTOR_ILU,MatGetFactor_seqsell_superlu)); 69642c9c57cSBarry Smith PetscFunctionReturn(0); 69742c9c57cSBarry Smith } 698