114b396bbSKris Buschelman 25a46d3faSBarry Smith /* -------------------------------------------------------------------- 35a46d3faSBarry Smith 45a46d3faSBarry Smith This file implements a subclass of the SeqAIJ matrix class that uses 55d8b2efaSHong Zhang the SuperLU sparse solver. You can use this as a starting point for 65a46d3faSBarry Smith implementing your own subclass of a PETSc matrix class. 75a46d3faSBarry Smith 85a46d3faSBarry Smith This demonstrates a way to make an implementation inheritence of a PETSc 95a46d3faSBarry Smith matrix type. This means constructing a new matrix type (SuperLU) by changing some 105a46d3faSBarry Smith of the methods of the previous type (SeqAIJ), adding additional data, and possibly 115a46d3faSBarry Smith additional method. (See any book on object oriented programming). 1214b396bbSKris Buschelman */ 1314b396bbSKris Buschelman 145a46d3faSBarry Smith /* 155a46d3faSBarry Smith Defines the data structure for the base matrix type (SeqAIJ) 165a46d3faSBarry Smith */ 17c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h> /*I "petscmat.h" I*/ 1814b396bbSKris Buschelman 195a46d3faSBarry Smith /* 205a46d3faSBarry Smith SuperLU include files 215a46d3faSBarry Smith */ 2214b396bbSKris Buschelman EXTERN_C_BEGIN 2314b396bbSKris Buschelman #if defined(PETSC_USE_COMPLEX) 24*3cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 25*3cb270beSHong Zhang #include <slu_cdefs.h> 26*3cb270beSHong Zhang #else 27c6db04a5SJed Brown #include <slu_zdefs.h> 28*3cb270beSHong Zhang #endif 29*3cb270beSHong Zhang #else 30*3cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 31*3cb270beSHong Zhang #include <slu_sdefs.h> 3214b396bbSKris Buschelman #else 33c6db04a5SJed Brown #include <slu_ddefs.h> 3414b396bbSKris Buschelman #endif 35*3cb270beSHong Zhang #endif 36c6db04a5SJed Brown #include <slu_util.h> 3714b396bbSKris Buschelman EXTERN_C_END 3814b396bbSKris Buschelman 395a46d3faSBarry Smith /* 405a46d3faSBarry Smith This is the data we are "ADDING" to the SeqAIJ matrix type to get the SuperLU matrix type 415a46d3faSBarry Smith */ 4214b396bbSKris Buschelman typedef struct { 4375af56d4SHong Zhang SuperMatrix A,L,U,B,X; 4475af56d4SHong Zhang superlu_options_t options; 45da7a1d00SHong Zhang PetscInt *perm_c; /* column permutation vector */ 46da7a1d00SHong Zhang PetscInt *perm_r; /* row permutations from partial pivoting */ 47da7a1d00SHong Zhang PetscInt *etree; 48da7a1d00SHong Zhang PetscReal *R, *C; 4975af56d4SHong Zhang char equed[1]; 50da7a1d00SHong Zhang PetscInt lwork; 5175af56d4SHong Zhang void *work; 52da7a1d00SHong Zhang PetscReal rpg, rcond; 5375af56d4SHong Zhang mem_usage_t mem_usage; 5475af56d4SHong Zhang MatStructure flg; 555d8b2efaSHong Zhang SuperLUStat_t stat; 56cae5a23dSHong Zhang Mat A_dup; 57cae5a23dSHong Zhang PetscScalar *rhs_dup; 5814b396bbSKris Buschelman 5914b396bbSKris Buschelman /* Flag to clean up (non-global) SuperLU objects during Destroy */ 60ace3abfcSBarry Smith PetscBool CleanUpSuperLU; 61f0c56d0fSKris Buschelman } Mat_SuperLU; 6214b396bbSKris Buschelman 63e0e586b9SHong Zhang extern PetscErrorCode MatFactorInfo_SuperLU(Mat,PetscViewer); 640481f469SBarry Smith extern PetscErrorCode MatLUFactorNumeric_SuperLU(Mat,Mat,const MatFactorInfo *); 65e0e586b9SHong Zhang extern PetscErrorCode MatDestroy_SuperLU(Mat); 66e0e586b9SHong Zhang extern PetscErrorCode MatView_SuperLU(Mat,PetscViewer); 67e0e586b9SHong Zhang extern PetscErrorCode MatAssemblyEnd_SuperLU(Mat,MatAssemblyType); 68e0e586b9SHong Zhang extern PetscErrorCode MatSolve_SuperLU(Mat,Vec,Vec); 69e0b74bf9SHong Zhang extern PetscErrorCode MatMatSolve_SuperLU(Mat,Mat,Mat); 70e0e586b9SHong Zhang extern PetscErrorCode MatSolveTranspose_SuperLU(Mat,Vec,Vec); 710481f469SBarry Smith extern PetscErrorCode MatLUFactorSymbolic_SuperLU(Mat,Mat,IS,IS,const MatFactorInfo*); 72e0e586b9SHong Zhang extern PetscErrorCode MatDuplicate_SuperLU(Mat, MatDuplicateOption, Mat *); 73e0e586b9SHong Zhang 745a46d3faSBarry Smith /* 755a46d3faSBarry Smith Utility function 765a46d3faSBarry Smith */ 775a46d3faSBarry Smith #undef __FUNCT__ 785a46d3faSBarry Smith #define __FUNCT__ "MatFactorInfo_SuperLU" 795a46d3faSBarry Smith PetscErrorCode MatFactorInfo_SuperLU(Mat A,PetscViewer viewer) 805a46d3faSBarry Smith { 815a46d3faSBarry Smith Mat_SuperLU *lu= (Mat_SuperLU*)A->spptr; 825a46d3faSBarry Smith PetscErrorCode ierr; 835a46d3faSBarry Smith superlu_options_t options; 845a46d3faSBarry Smith 855a46d3faSBarry Smith PetscFunctionBegin; 865a46d3faSBarry Smith /* check if matrix is superlu_dist type */ 875a46d3faSBarry Smith if (A->ops->solve != MatSolve_SuperLU) PetscFunctionReturn(0); 885a46d3faSBarry Smith 895a46d3faSBarry Smith options = lu->options; 905a46d3faSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"SuperLU run parameters:\n");CHKERRQ(ierr); 915a46d3faSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Equil: %s\n",(options.Equil != NO) ? "YES": "NO");CHKERRQ(ierr); 925a46d3faSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," ColPerm: %D\n",options.ColPerm);CHKERRQ(ierr); 935a46d3faSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," IterRefine: %D\n",options.IterRefine);CHKERRQ(ierr); 945a46d3faSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," SymmetricMode: %s\n",(options.SymmetricMode != NO) ? "YES": "NO");CHKERRQ(ierr); 955a46d3faSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," DiagPivotThresh: %g\n",options.DiagPivotThresh);CHKERRQ(ierr); 965a46d3faSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," PivotGrowth: %s\n",(options.PivotGrowth != NO) ? "YES": "NO");CHKERRQ(ierr); 975a46d3faSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," ConditionNumber: %s\n",(options.ConditionNumber != NO) ? "YES": "NO");CHKERRQ(ierr); 985a46d3faSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," RowPerm: %D\n",options.RowPerm);CHKERRQ(ierr); 995a46d3faSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," ReplaceTinyPivot: %s\n",(options.ReplaceTinyPivot != NO) ? "YES": "NO");CHKERRQ(ierr); 1005a46d3faSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," PrintStat: %s\n",(options.PrintStat != NO) ? "YES": "NO");CHKERRQ(ierr); 1015a46d3faSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," lwork: %D\n",lu->lwork);CHKERRQ(ierr); 102d5f3da31SBarry Smith if (A->factortype == MAT_FACTOR_ILU){ 103cffbb591SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ILU_DropTol: %g\n",options.ILU_DropTol);CHKERRQ(ierr); 104cffbb591SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ILU_FillTol: %g\n",options.ILU_FillTol);CHKERRQ(ierr); 105cffbb591SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ILU_FillFactor: %g\n",options.ILU_FillFactor);CHKERRQ(ierr); 106cffbb591SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ILU_DropRule: %D\n",options.ILU_DropRule);CHKERRQ(ierr); 107cffbb591SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ILU_Norm: %D\n",options.ILU_Norm);CHKERRQ(ierr); 108cffbb591SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ILU_MILU: %D\n",options.ILU_MILU);CHKERRQ(ierr); 109cffbb591SHong Zhang } 1105a46d3faSBarry Smith PetscFunctionReturn(0); 1115a46d3faSBarry Smith } 1125a46d3faSBarry Smith 1135a46d3faSBarry Smith /* 1145a46d3faSBarry Smith These are the methods provided to REPLACE the corresponding methods of the 1155a46d3faSBarry Smith SeqAIJ matrix class. Other methods of SeqAIJ are not replaced 1165a46d3faSBarry Smith */ 1175a46d3faSBarry Smith #undef __FUNCT__ 1185a46d3faSBarry Smith #define __FUNCT__ "MatLUFactorNumeric_SuperLU" 1190481f469SBarry Smith PetscErrorCode MatLUFactorNumeric_SuperLU(Mat F,Mat A,const MatFactorInfo *info) 1205a46d3faSBarry Smith { 1211d5ca7f3SHong Zhang Mat_SuperLU *lu = (Mat_SuperLU*)F->spptr; 122cae5a23dSHong Zhang Mat_SeqAIJ *aa; 1235a46d3faSBarry Smith PetscErrorCode ierr; 1245a46d3faSBarry Smith PetscInt sinfo; 1255a46d3faSBarry Smith PetscReal ferr, berr; 1265a46d3faSBarry Smith NCformat *Ustore; 1275a46d3faSBarry Smith SCformat *Lstore; 1285a46d3faSBarry Smith 1295a46d3faSBarry Smith PetscFunctionBegin; 1305a46d3faSBarry Smith if (lu->flg == SAME_NONZERO_PATTERN){ /* successing numerical factorization */ 1315a46d3faSBarry Smith lu->options.Fact = SamePattern; 1325a46d3faSBarry Smith /* Ref: ~SuperLU_3.0/EXAMPLE/dlinsolx2.c */ 1335a46d3faSBarry Smith Destroy_SuperMatrix_Store(&lu->A); 134cae5a23dSHong Zhang if (lu->options.Equil){ 135cae5a23dSHong Zhang ierr = MatCopy_SeqAIJ(A,lu->A_dup,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 136cae5a23dSHong Zhang } 1375a46d3faSBarry Smith if ( lu->lwork >= 0 ) { 1385a46d3faSBarry Smith Destroy_SuperNode_Matrix(&lu->L); 1395a46d3faSBarry Smith Destroy_CompCol_Matrix(&lu->U); 1405a46d3faSBarry Smith lu->options.Fact = SamePattern; 1415a46d3faSBarry Smith } 1425a46d3faSBarry Smith } 1435a46d3faSBarry Smith 1445a46d3faSBarry Smith /* Create the SuperMatrix for lu->A=A^T: 1455a46d3faSBarry Smith Since SuperLU likes column-oriented matrices,we pass it the transpose, 1465a46d3faSBarry Smith and then solve A^T X = B in MatSolve(). */ 147cae5a23dSHong Zhang if (lu->options.Equil){ 148cae5a23dSHong Zhang aa = (Mat_SeqAIJ*)(lu->A_dup)->data; 149cae5a23dSHong Zhang } else { 150cae5a23dSHong Zhang aa = (Mat_SeqAIJ*)(A)->data; 151cae5a23dSHong Zhang } 1525a46d3faSBarry Smith #if defined(PETSC_USE_COMPLEX) 153*3cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 154*3cb270beSHong Zhang cCreate_CompCol_Matrix(&lu->A,A->cmap->n,A->rmap->n,aa->nz,(singlecomplex*)aa->a,aa->j,aa->i, 155*3cb270beSHong Zhang SLU_NC,SLU_C,SLU_GE); 156*3cb270beSHong Zhang #else 157d0f46423SBarry Smith zCreate_CompCol_Matrix(&lu->A,A->cmap->n,A->rmap->n,aa->nz,(doublecomplex*)aa->a,aa->j,aa->i, 1585a46d3faSBarry Smith SLU_NC,SLU_Z,SLU_GE); 159*3cb270beSHong Zhang #endif 160*3cb270beSHong Zhang #else 161*3cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 162*3cb270beSHong Zhang sCreate_CompCol_Matrix(&lu->A,A->cmap->n,A->rmap->n,aa->nz,aa->a,aa->j,aa->i, 163*3cb270beSHong Zhang SLU_NC,SLU_S,SLU_GE); 1645a46d3faSBarry Smith #else 165d0f46423SBarry Smith dCreate_CompCol_Matrix(&lu->A,A->cmap->n,A->rmap->n,aa->nz,aa->a,aa->j,aa->i, 1665a46d3faSBarry Smith SLU_NC,SLU_D,SLU_GE); 1675a46d3faSBarry Smith #endif 168*3cb270beSHong Zhang #endif 1695a46d3faSBarry Smith 1705a46d3faSBarry Smith /* Numerical factorization */ 1715a46d3faSBarry Smith lu->B.ncol = 0; /* Indicate not to solve the system */ 172d5f3da31SBarry Smith if (F->factortype == MAT_FACTOR_LU){ 1735a46d3faSBarry Smith #if defined(PETSC_USE_COMPLEX) 174*3cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 175*3cb270beSHong Zhang cgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 176*3cb270beSHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 177*3cb270beSHong Zhang &lu->mem_usage, &lu->stat, &sinfo); 178*3cb270beSHong Zhang #else 1795a46d3faSBarry Smith zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 1805a46d3faSBarry Smith &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 1815d8b2efaSHong Zhang &lu->mem_usage, &lu->stat, &sinfo); 182*3cb270beSHong Zhang #endif 183*3cb270beSHong Zhang #else 184*3cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 185*3cb270beSHong Zhang sgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 186*3cb270beSHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 187*3cb270beSHong Zhang &lu->mem_usage, &lu->stat, &sinfo); 1885a46d3faSBarry Smith #else 1895a46d3faSBarry Smith dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 1905a46d3faSBarry Smith &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 1915d8b2efaSHong Zhang &lu->mem_usage, &lu->stat, &sinfo); 1925a46d3faSBarry Smith #endif 193*3cb270beSHong Zhang #endif 194d5f3da31SBarry Smith } else if (F->factortype == MAT_FACTOR_ILU){ 195cffbb591SHong Zhang /* Compute the incomplete factorization, condition number and pivot growth */ 196cffbb591SHong Zhang #if defined(PETSC_USE_COMPLEX) 197*3cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 198*3cb270beSHong Zhang cgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r,lu->etree, lu->equed, lu->R, lu->C, 199*3cb270beSHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, 200*3cb270beSHong Zhang &lu->mem_usage, &lu->stat, &sinfo); 201*3cb270beSHong Zhang #else 202cffbb591SHong Zhang zgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r,lu->etree, lu->equed, lu->R, lu->C, 203cffbb591SHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, 2045d8b2efaSHong Zhang &lu->mem_usage, &lu->stat, &sinfo); 205*3cb270beSHong Zhang #endif 206*3cb270beSHong Zhang #else 207*3cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 208*3cb270beSHong Zhang sgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 209*3cb270beSHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, 210*3cb270beSHong Zhang &lu->mem_usage, &lu->stat, &sinfo); 211cffbb591SHong Zhang #else 212cffbb591SHong Zhang dgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 213cffbb591SHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, 2145d8b2efaSHong Zhang &lu->mem_usage, &lu->stat, &sinfo); 215cffbb591SHong Zhang #endif 216*3cb270beSHong Zhang #endif 217cffbb591SHong Zhang } else { 218e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported"); 219cffbb591SHong Zhang } 2205a46d3faSBarry Smith if ( !sinfo || sinfo == lu->A.ncol+1 ) { 2215a46d3faSBarry Smith if ( lu->options.PivotGrowth ) 2225a46d3faSBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," Recip. pivot growth = %e\n", lu->rpg); 2235a46d3faSBarry Smith if ( lu->options.ConditionNumber ) 2245a46d3faSBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," Recip. condition number = %e\n", lu->rcond); 2255a46d3faSBarry Smith } else if ( sinfo > 0 ){ 2265a46d3faSBarry Smith if ( lu->lwork == -1 ) { 2275a46d3faSBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," ** Estimated memory: %D bytes\n", sinfo - lu->A.ncol); 2289308c137SBarry Smith } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot in row %D",sinfo); 2295a46d3faSBarry Smith } else { /* sinfo < 0 */ 230e32f2f54SBarry Smith SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_LIB, "info = %D, the %D-th argument in gssvx() had an illegal value", sinfo,-sinfo); 2315a46d3faSBarry Smith } 2325a46d3faSBarry Smith 2335a46d3faSBarry Smith if ( lu->options.PrintStat ) { 2345a46d3faSBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF,"MatLUFactorNumeric_SuperLU():\n"); 2355d8b2efaSHong Zhang StatPrint(&lu->stat); 2365a46d3faSBarry Smith Lstore = (SCformat *) lu->L.Store; 2375a46d3faSBarry Smith Ustore = (NCformat *) lu->U.Store; 2385a46d3faSBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," No of nonzeros in factor L = %D\n", Lstore->nnz); 2395a46d3faSBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," No of nonzeros in factor U = %D\n", Ustore->nnz); 2405a46d3faSBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," No of nonzeros in L+U = %D\n", Lstore->nnz + Ustore->nnz - lu->A.ncol); 2416da386baSHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," L\\U MB %.3f\ttotal MB needed %.3f\n", 2426da386baSHong Zhang lu->mem_usage.for_lu/1e6, lu->mem_usage.total_needed/1e6); 2435a46d3faSBarry Smith } 2445a46d3faSBarry Smith 2455a46d3faSBarry Smith lu->flg = SAME_NONZERO_PATTERN; 2461d5ca7f3SHong Zhang F->ops->solve = MatSolve_SuperLU; 2471d5ca7f3SHong Zhang F->ops->solvetranspose = MatSolveTranspose_SuperLU; 248e0b74bf9SHong Zhang F->ops->matsolve = MatMatSolve_SuperLU; 2495a46d3faSBarry Smith PetscFunctionReturn(0); 2505a46d3faSBarry Smith } 2515a46d3faSBarry Smith 25214b396bbSKris Buschelman #undef __FUNCT__ 253f0c56d0fSKris Buschelman #define __FUNCT__ "MatDestroy_SuperLU" 254dfbe8321SBarry Smith PetscErrorCode MatDestroy_SuperLU(Mat A) 25514b396bbSKris Buschelman { 256dfbe8321SBarry Smith PetscErrorCode ierr; 257f0c56d0fSKris Buschelman Mat_SuperLU *lu=(Mat_SuperLU*)A->spptr; 25814b396bbSKris Buschelman 25914b396bbSKris Buschelman PetscFunctionBegin; 260bf0cc555SLisandro Dalcin if (lu && lu->CleanUpSuperLU) { /* Free the SuperLU datastructures */ 26175af56d4SHong Zhang Destroy_SuperMatrix_Store(&lu->A); 26275af56d4SHong Zhang Destroy_SuperMatrix_Store(&lu->B); 26375af56d4SHong Zhang Destroy_SuperMatrix_Store(&lu->X); 2645d8b2efaSHong Zhang StatFree(&lu->stat); 2650e742b69SHong Zhang if (lu->lwork >= 0) { 2660e742b69SHong Zhang Destroy_SuperNode_Matrix(&lu->L); 2670e742b69SHong Zhang Destroy_CompCol_Matrix(&lu->U); 2680e742b69SHong Zhang } 2690e742b69SHong Zhang } 270bf0cc555SLisandro Dalcin if (lu) { 2719ce50919SHong Zhang ierr = PetscFree(lu->etree);CHKERRQ(ierr); 2729ce50919SHong Zhang ierr = PetscFree(lu->perm_r);CHKERRQ(ierr); 2739ce50919SHong Zhang ierr = PetscFree(lu->perm_c);CHKERRQ(ierr); 2749ce50919SHong Zhang ierr = PetscFree(lu->R);CHKERRQ(ierr); 2759ce50919SHong Zhang ierr = PetscFree(lu->C);CHKERRQ(ierr); 276bf0cc555SLisandro Dalcin ierr = PetscFree(lu->rhs_dup);CHKERRQ(ierr); 277bf0cc555SLisandro Dalcin ierr = MatDestroy(&lu->A_dup);CHKERRQ(ierr); 278bf0cc555SLisandro Dalcin } 279bf0cc555SLisandro Dalcin ierr = PetscFree(A->spptr);CHKERRQ(ierr); 2800e742b69SHong Zhang 281d954db57SHong Zhang /* clear composed functions */ 282d954db57SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatFactorGetSolverPackage_C","",PETSC_NULL);CHKERRQ(ierr); 283790a96dfSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatSuperluSetILUDropTol_C","",PETSC_NULL);CHKERRQ(ierr); 284d954db57SHong Zhang 285b24902e0SBarry Smith ierr = MatDestroy_SeqAIJ(A);CHKERRQ(ierr); 28614b396bbSKris Buschelman PetscFunctionReturn(0); 28714b396bbSKris Buschelman } 28814b396bbSKris Buschelman 28914b396bbSKris Buschelman #undef __FUNCT__ 290f0c56d0fSKris Buschelman #define __FUNCT__ "MatView_SuperLU" 291dfbe8321SBarry Smith PetscErrorCode MatView_SuperLU(Mat A,PetscViewer viewer) 29214b396bbSKris Buschelman { 293dfbe8321SBarry Smith PetscErrorCode ierr; 294ace3abfcSBarry Smith PetscBool iascii; 29514b396bbSKris Buschelman PetscViewerFormat format; 29614b396bbSKris Buschelman 29714b396bbSKris Buschelman PetscFunctionBegin; 298251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 29932077d6dSBarry Smith if (iascii) { 30014b396bbSKris Buschelman ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 3012f59403fSHong Zhang if (format == PETSC_VIEWER_ASCII_INFO) { 302f0c56d0fSKris Buschelman ierr = MatFactorInfo_SuperLU(A,viewer);CHKERRQ(ierr); 30314b396bbSKris Buschelman } 30414b396bbSKris Buschelman } 30514b396bbSKris Buschelman PetscFunctionReturn(0); 30614b396bbSKris Buschelman } 30714b396bbSKris Buschelman 30814b396bbSKris Buschelman 30914b396bbSKris Buschelman #undef __FUNCT__ 310c29e884eSHong Zhang #define __FUNCT__ "MatSolve_SuperLU_Private" 311c29e884eSHong Zhang PetscErrorCode MatSolve_SuperLU_Private(Mat A,Vec b,Vec x) 31214b396bbSKris Buschelman { 313f0c56d0fSKris Buschelman Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr; 31475af56d4SHong Zhang PetscScalar *barray,*xarray; 315dfbe8321SBarry Smith PetscErrorCode ierr; 316cae5a23dSHong Zhang PetscInt info,i,n=x->map->n; 317da7a1d00SHong Zhang PetscReal ferr,berr; 31814b396bbSKris Buschelman 31914b396bbSKris Buschelman PetscFunctionBegin; 320937a6b0eSHong Zhang if ( lu->lwork == -1 ) { 321937a6b0eSHong Zhang PetscFunctionReturn(0); 322937a6b0eSHong Zhang } 323cae5a23dSHong Zhang 32475af56d4SHong Zhang lu->B.ncol = 1; /* Set the number of right-hand side */ 325cae5a23dSHong Zhang if (lu->options.Equil && !lu->rhs_dup){ 326cae5a23dSHong Zhang /* superlu overwrites b when Equil is used, thus create rhs_dup to keep user's b unchanged */ 327cae5a23dSHong Zhang ierr = PetscMalloc(n*sizeof(PetscScalar),&lu->rhs_dup);CHKERRQ(ierr); 328cae5a23dSHong Zhang } 329cae5a23dSHong Zhang if (lu->options.Equil){ 330cae5a23dSHong Zhang /* Copy b into rsh_dup */ 33175af56d4SHong Zhang ierr = VecGetArray(b,&barray);CHKERRQ(ierr); 332cae5a23dSHong Zhang ierr = PetscMemcpy(lu->rhs_dup,barray,n*sizeof(PetscScalar));CHKERRQ(ierr); 333cae5a23dSHong Zhang ierr = VecRestoreArray(b,&barray);CHKERRQ(ierr); 334cae5a23dSHong Zhang barray = lu->rhs_dup; 335cae5a23dSHong Zhang } else { 336cae5a23dSHong Zhang ierr = VecGetArray(b,&barray);CHKERRQ(ierr); 337cae5a23dSHong Zhang } 33875af56d4SHong Zhang ierr = VecGetArray(x,&xarray);CHKERRQ(ierr); 3395fe6150dSHong Zhang 3405fe6150dSHong Zhang #if defined(PETSC_USE_COMPLEX) 341*3cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 342*3cb270beSHong Zhang ((DNformat*)lu->B.Store)->nzval = (singlecomplex*)barray; 343*3cb270beSHong Zhang ((DNformat*)lu->X.Store)->nzval = (singlecomplex*)xarray; 344*3cb270beSHong Zhang #else 3455fe6150dSHong Zhang ((DNformat*)lu->B.Store)->nzval = (doublecomplex*)barray; 3465fe6150dSHong Zhang ((DNformat*)lu->X.Store)->nzval = (doublecomplex*)xarray; 347*3cb270beSHong Zhang #endif 3485fe6150dSHong Zhang #else 3495fe6150dSHong Zhang ((DNformat*)lu->B.Store)->nzval = barray; 35075af56d4SHong Zhang ((DNformat*)lu->X.Store)->nzval = xarray; 3515fe6150dSHong Zhang #endif 35275af56d4SHong Zhang 35375af56d4SHong Zhang lu->options.Fact = FACTORED; /* Indicate the factored form of A is supplied. */ 354d5f3da31SBarry Smith if (A->factortype == MAT_FACTOR_LU){ 3558914a3f7SHong Zhang #if defined(PETSC_USE_COMPLEX) 356*3cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 357*3cb270beSHong Zhang cgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 358*3cb270beSHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 359*3cb270beSHong Zhang &lu->mem_usage, &lu->stat, &info); 360*3cb270beSHong Zhang #else 3618914a3f7SHong Zhang zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 3628914a3f7SHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 3635d8b2efaSHong Zhang &lu->mem_usage, &lu->stat, &info); 364*3cb270beSHong Zhang #endif 365*3cb270beSHong Zhang #else 366*3cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 367*3cb270beSHong Zhang sgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 368*3cb270beSHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 369*3cb270beSHong Zhang &lu->mem_usage, &lu->stat, &info); 3708914a3f7SHong Zhang #else 37175af56d4SHong Zhang dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 37275af56d4SHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 3735d8b2efaSHong Zhang &lu->mem_usage, &lu->stat, &info); 3748914a3f7SHong Zhang #endif 375*3cb270beSHong Zhang #endif 376d5f3da31SBarry Smith } else if (A->factortype == MAT_FACTOR_ILU){ 377cffbb591SHong Zhang #if defined(PETSC_USE_COMPLEX) 378*3cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 379*3cb270beSHong Zhang cgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 380*3cb270beSHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, 381*3cb270beSHong Zhang &lu->mem_usage, &lu->stat, &info); 382*3cb270beSHong Zhang #else 383cffbb591SHong Zhang zgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 384cffbb591SHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, 3855d8b2efaSHong Zhang &lu->mem_usage, &lu->stat, &info); 386*3cb270beSHong Zhang #endif 387*3cb270beSHong Zhang #else 388*3cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 389*3cb270beSHong Zhang sgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 390*3cb270beSHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, 391*3cb270beSHong Zhang &lu->mem_usage, &lu->stat, &info); 392cffbb591SHong Zhang #else 393cffbb591SHong Zhang dgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 394cffbb591SHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, 3955d8b2efaSHong Zhang &lu->mem_usage, &lu->stat, &info); 396cffbb591SHong Zhang #endif 397*3cb270beSHong Zhang #endif 398cffbb591SHong Zhang } else { 399e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported"); 400cffbb591SHong Zhang } 401cae5a23dSHong Zhang if (!lu->options.Equil){ 40275af56d4SHong Zhang ierr = VecRestoreArray(b,&barray);CHKERRQ(ierr); 403cae5a23dSHong Zhang } 40475af56d4SHong Zhang ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr); 40575af56d4SHong Zhang 406958c9bccSBarry Smith if ( !info || info == lu->A.ncol+1 ) { 40775af56d4SHong Zhang if ( lu->options.IterRefine ) { 4088914a3f7SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF,"Iterative Refinement:\n"); 4098914a3f7SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," %8s%8s%16s%16s\n", "rhs", "Steps", "FERR", "BERR"); 41075af56d4SHong Zhang for (i = 0; i < 1; ++i) 4115d8b2efaSHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," %8d%8d%16e%16e\n", i+1, lu->stat.RefineSteps, ferr, berr); 41275af56d4SHong Zhang } 4138914a3f7SHong Zhang } else if ( info > 0 ){ 4148914a3f7SHong Zhang if ( lu->lwork == -1 ) { 41577431f27SBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," ** Estimated memory: %D bytes\n", info - lu->A.ncol); 4168914a3f7SHong Zhang } else { 41777431f27SBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," Warning: gssvx() returns info %D\n",info); 4188914a3f7SHong Zhang } 4198914a3f7SHong Zhang } else if (info < 0){ 420e32f2f54SBarry Smith SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_LIB, "info = %D, the %D-th argument in gssvx() had an illegal value", info,-info); 42114b396bbSKris Buschelman } 42214b396bbSKris Buschelman 4238914a3f7SHong Zhang if ( lu->options.PrintStat ) { 4248914a3f7SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF,"MatSolve__SuperLU():\n"); 4255d8b2efaSHong Zhang StatPrint(&lu->stat); 4268914a3f7SHong Zhang } 42775af56d4SHong Zhang PetscFunctionReturn(0); 42875af56d4SHong Zhang } 42914b396bbSKris Buschelman 43014b396bbSKris Buschelman #undef __FUNCT__ 431c29e884eSHong Zhang #define __FUNCT__ "MatSolve_SuperLU" 432c29e884eSHong Zhang PetscErrorCode MatSolve_SuperLU(Mat A,Vec b,Vec x) 433c29e884eSHong Zhang { 434c29e884eSHong Zhang Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr; 435c29e884eSHong Zhang PetscErrorCode ierr; 436c29e884eSHong Zhang 437c29e884eSHong Zhang PetscFunctionBegin; 438c29e884eSHong Zhang lu->options.Trans = TRANS; 439c29e884eSHong Zhang ierr = MatSolve_SuperLU_Private(A,b,x);CHKERRQ(ierr); 440c29e884eSHong Zhang PetscFunctionReturn(0); 441c29e884eSHong Zhang } 442c29e884eSHong Zhang 443c29e884eSHong Zhang #undef __FUNCT__ 444c7c1fe80SHong Zhang #define __FUNCT__ "MatSolveTranspose_SuperLU" 445c7c1fe80SHong Zhang PetscErrorCode MatSolveTranspose_SuperLU(Mat A,Vec b,Vec x) 446c7c1fe80SHong Zhang { 447c7c1fe80SHong Zhang Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr; 448c7c1fe80SHong Zhang PetscErrorCode ierr; 449c7c1fe80SHong Zhang 450c7c1fe80SHong Zhang PetscFunctionBegin; 451c7c1fe80SHong Zhang lu->options.Trans = NOTRANS; 452c29e884eSHong Zhang ierr = MatSolve_SuperLU_Private(A,b,x);CHKERRQ(ierr); 453c7c1fe80SHong Zhang PetscFunctionReturn(0); 454c7c1fe80SHong Zhang } 455c7c1fe80SHong Zhang 456e0b74bf9SHong Zhang #undef __FUNCT__ 457e0b74bf9SHong Zhang #define __FUNCT__ "MatMatSolve_SuperLU" 458e0b74bf9SHong Zhang PetscErrorCode MatMatSolve_SuperLU(Mat A,Mat B,Mat X) 459e0b74bf9SHong Zhang { 460e0b74bf9SHong Zhang Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr; 461cd723cd1SBarry Smith PetscBool flg; 462cd723cd1SBarry Smith PetscErrorCode ierr; 463e0b74bf9SHong Zhang 464e0b74bf9SHong Zhang PetscFunctionBegin; 465251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompareAny((PetscObject)B,&flg,MATSEQDENSE,MATMPIDENSE,PETSC_NULL);CHKERRQ(ierr); 466cd723cd1SBarry Smith if (!flg) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix"); 467251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,PETSC_NULL);CHKERRQ(ierr); 468cd723cd1SBarry Smith if (!flg) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix"); lu->options.Trans = TRANS; 469e0b74bf9SHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatMatSolve_SuperLU() is not implemented yet"); 470e0b74bf9SHong Zhang PetscFunctionReturn(0); 471e0b74bf9SHong Zhang } 472e0b74bf9SHong Zhang 47314b396bbSKris Buschelman /* 47414b396bbSKris Buschelman Note the r permutation is ignored 47514b396bbSKris Buschelman */ 47614b396bbSKris Buschelman #undef __FUNCT__ 477f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_SuperLU" 4780481f469SBarry Smith PetscErrorCode MatLUFactorSymbolic_SuperLU(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info) 47914b396bbSKris Buschelman { 4801d5ca7f3SHong Zhang Mat_SuperLU *lu = (Mat_SuperLU*)(F->spptr); 481b24902e0SBarry Smith 482b24902e0SBarry Smith PetscFunctionBegin; 483b24902e0SBarry Smith lu->flg = DIFFERENT_NONZERO_PATTERN; 484b24902e0SBarry Smith lu->CleanUpSuperLU = PETSC_TRUE; 4851d5ca7f3SHong Zhang F->ops->lufactornumeric = MatLUFactorNumeric_SuperLU; 486b24902e0SBarry Smith PetscFunctionReturn(0); 487b24902e0SBarry Smith } 488b24902e0SBarry Smith 48935bd34faSBarry Smith EXTERN_C_BEGIN 49035bd34faSBarry Smith #undef __FUNCT__ 4915ccb76cbSHong Zhang #define __FUNCT__ "MatSuperluSetILUDropTol_SuperLU" 4925ccb76cbSHong Zhang PetscErrorCode MatSuperluSetILUDropTol_SuperLU(Mat F,PetscReal dtol) 4935ccb76cbSHong Zhang { 4945ccb76cbSHong Zhang Mat_SuperLU *lu= (Mat_SuperLU*)F->spptr; 4955ccb76cbSHong Zhang 4965ccb76cbSHong Zhang PetscFunctionBegin; 4975ccb76cbSHong Zhang lu->options.ILU_DropTol = dtol; 4985ccb76cbSHong Zhang PetscFunctionReturn(0); 4995ccb76cbSHong Zhang } 5005ccb76cbSHong Zhang EXTERN_C_END 5015ccb76cbSHong Zhang 5025ccb76cbSHong Zhang #undef __FUNCT__ 5035ccb76cbSHong Zhang #define __FUNCT__ "MatSuperluSetILUDropTol" 5045ccb76cbSHong Zhang /*@ 5055ccb76cbSHong Zhang MatSuperluSetILUDropTol - Set SuperLU ILU drop tolerance 5065ccb76cbSHong Zhang Logically Collective on Mat 5075ccb76cbSHong Zhang 5085ccb76cbSHong Zhang Input Parameters: 5095ccb76cbSHong Zhang + F - the factored matrix obtained by calling MatGetFactor() from PETSc-SuperLU interface 5105ccb76cbSHong Zhang - dtol - drop tolerance 5115ccb76cbSHong Zhang 5125ccb76cbSHong Zhang Options Database: 5135ccb76cbSHong Zhang . -mat_superlu_ilu_droptol <dtol> 5145ccb76cbSHong Zhang 5155ccb76cbSHong Zhang Level: beginner 5165ccb76cbSHong Zhang 5175ccb76cbSHong Zhang References: SuperLU Users' Guide 5185ccb76cbSHong Zhang 5195ccb76cbSHong Zhang .seealso: MatGetFactor() 5205ccb76cbSHong Zhang @*/ 5215ccb76cbSHong Zhang PetscErrorCode MatSuperluSetILUDropTol(Mat F,PetscReal dtol) 5225ccb76cbSHong Zhang { 5235ccb76cbSHong Zhang PetscErrorCode ierr; 5245ccb76cbSHong Zhang 5255ccb76cbSHong Zhang PetscFunctionBegin; 5265ccb76cbSHong Zhang PetscValidHeaderSpecific(F,MAT_CLASSID,1); 5275ccb76cbSHong Zhang PetscValidLogicalCollectiveInt(F,dtol,2); 5285ccb76cbSHong Zhang ierr = PetscTryMethod(F,"MatSuperluSetILUDropTol_C",(Mat,PetscReal),(F,dtol));CHKERRQ(ierr); 5295ccb76cbSHong Zhang PetscFunctionReturn(0); 5305ccb76cbSHong Zhang } 5315ccb76cbSHong Zhang 5325ccb76cbSHong Zhang EXTERN_C_BEGIN 5335ccb76cbSHong Zhang #undef __FUNCT__ 53435bd34faSBarry Smith #define __FUNCT__ "MatFactorGetSolverPackage_seqaij_superlu" 53535bd34faSBarry Smith PetscErrorCode MatFactorGetSolverPackage_seqaij_superlu(Mat A,const MatSolverPackage *type) 53635bd34faSBarry Smith { 53735bd34faSBarry Smith PetscFunctionBegin; 5382692d6eeSBarry Smith *type = MATSOLVERSUPERLU; 53935bd34faSBarry Smith PetscFunctionReturn(0); 54035bd34faSBarry Smith } 54135bd34faSBarry Smith EXTERN_C_END 54235bd34faSBarry Smith 543b24902e0SBarry Smith 544b24902e0SBarry Smith /*MC 545ba992d64SSatish Balay MATSOLVERSUPERLU = "superlu" - A solver package providing solvers LU and ILU for sequential matrices 546b24902e0SBarry Smith via the external package SuperLU. 547b24902e0SBarry Smith 548e2e64c6bSBarry Smith Use ./configure --download-superlu to have PETSc installed with SuperLU 549b24902e0SBarry Smith 550b24902e0SBarry Smith Options Database Keys: 551e08999f5SMatthew G Knepley + -mat_superlu_equil <FALSE> - Equil (None) 552e08999f5SMatthew G Knepley . -mat_superlu_colperm <COLAMD> - (choose one of) NATURAL MMD_ATA MMD_AT_PLUS_A COLAMD 553e08999f5SMatthew G Knepley . -mat_superlu_iterrefine <NOREFINE> - (choose one of) NOREFINE SINGLE DOUBLE EXTRA 554e08999f5SMatthew G Knepley . -mat_superlu_symmetricmode: <FALSE> - SymmetricMode (None) 555e08999f5SMatthew G Knepley . -mat_superlu_diagpivotthresh <1> - DiagPivotThresh (None) 556e08999f5SMatthew G Knepley . -mat_superlu_pivotgrowth <FALSE> - PivotGrowth (None) 557e08999f5SMatthew G Knepley . -mat_superlu_conditionnumber <FALSE> - ConditionNumber (None) 558e08999f5SMatthew G Knepley . -mat_superlu_rowperm <NOROWPERM> - (choose one of) NOROWPERM LargeDiag 559e08999f5SMatthew G Knepley . -mat_superlu_replacetinypivot <FALSE> - ReplaceTinyPivot (None) 560e08999f5SMatthew G Knepley . -mat_superlu_printstat <FALSE> - PrintStat (None) 561e08999f5SMatthew G Knepley . -mat_superlu_lwork <0> - size of work array in bytes used by factorization (None) 562e08999f5SMatthew G Knepley . -mat_superlu_ilu_droptol <0> - ILU_DropTol (None) 563e08999f5SMatthew G Knepley . -mat_superlu_ilu_filltol <0> - ILU_FillTol (None) 564e08999f5SMatthew G Knepley . -mat_superlu_ilu_fillfactor <0> - ILU_FillFactor (None) 565e08999f5SMatthew G Knepley . -mat_superlu_ilu_droprull <0> - ILU_DropRule (None) 566e08999f5SMatthew G Knepley . -mat_superlu_ilu_norm <0> - ILU_Norm (None) 567e08999f5SMatthew G Knepley - -mat_superlu_ilu_milu <0> - ILU_MILU (None) 568b24902e0SBarry Smith 5692692d6eeSBarry Smith Notes: Do not confuse this with MATSOLVERSUPERLU_DIST which is for parallel sparse solves 5705c9eb25fSBarry Smith 571b24902e0SBarry Smith Level: beginner 572b24902e0SBarry Smith 573d45987f3SHong Zhang .seealso: PCLU, PCILU, MATSOLVERSUPERLU_DIST, MATSOLVERMUMPS, PCFactorSetMatSolverPackage(), MatSolverPackage 574b24902e0SBarry Smith M*/ 575b24902e0SBarry Smith 5765c9eb25fSBarry Smith EXTERN_C_BEGIN 577b24902e0SBarry Smith #undef __FUNCT__ 578b24902e0SBarry Smith #define __FUNCT__ "MatGetFactor_seqaij_superlu" 5795c9eb25fSBarry Smith PetscErrorCode MatGetFactor_seqaij_superlu(Mat A,MatFactorType ftype,Mat *F) 580b24902e0SBarry Smith { 58114b396bbSKris Buschelman Mat B; 582f0c56d0fSKris Buschelman Mat_SuperLU *lu; 5836849ba73SBarry Smith PetscErrorCode ierr; 5845d8b2efaSHong Zhang PetscInt indx,m=A->rmap->n,n=A->cmap->n; 585ace3abfcSBarry Smith PetscBool flg; 586*3cb270beSHong Zhang PetscReal real_input; 5875ca28756SHong Zhang const char *colperm[]={"NATURAL","MMD_ATA","MMD_AT_PLUS_A","COLAMD"}; /* MY_PERMC - not supported by the petsc interface yet */ 5885ca28756SHong Zhang const char *iterrefine[]={"NOREFINE", "SINGLE", "DOUBLE", "EXTRA"}; 5895ca28756SHong Zhang const char *rowperm[]={"NOROWPERM", "LargeDiag"}; /* MY_PERMC - not supported by the petsc interface yet */ 59014b396bbSKris Buschelman 59114b396bbSKris Buschelman PetscFunctionBegin; 5927adad957SLisandro Dalcin ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr); 593d0f46423SBarry Smith ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 5947adad957SLisandro Dalcin ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 595899d7b4fSKris Buschelman ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr); 596f0c56d0fSKris Buschelman 597cffbb591SHong Zhang if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU){ 598b24902e0SBarry Smith B->ops->lufactorsymbolic = MatLUFactorSymbolic_SuperLU; 599cffbb591SHong Zhang B->ops->ilufactorsymbolic = MatLUFactorSymbolic_SuperLU; 600b3a44c85SBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported"); 601cffbb591SHong Zhang 602b24902e0SBarry Smith B->ops->destroy = MatDestroy_SuperLU; 6033519fcdcSHong Zhang B->ops->view = MatView_SuperLU; 604d5f3da31SBarry Smith B->factortype = ftype; 60594b7f48cSBarry Smith B->assembled = PETSC_TRUE; /* required by -ksp_view */ 6065c9eb25fSBarry Smith B->preallocated = PETSC_TRUE; 60714b396bbSKris Buschelman 608b24902e0SBarry Smith ierr = PetscNewLog(B,Mat_SuperLU,&lu);CHKERRQ(ierr); 609cae5a23dSHong Zhang 610cffbb591SHong Zhang if (ftype == MAT_FACTOR_LU){ 6119ce50919SHong Zhang set_default_options(&lu->options); 6123d6581fbSHong Zhang /* Comments from SuperLU_4.0/SRC/dgssvx.c: 6133d6581fbSHong Zhang "Whether or not the system will be equilibrated depends on the 6143d6581fbSHong Zhang scaling of the matrix A, but if equilibration is used, A is 6153d6581fbSHong Zhang overwritten by diag(R)*A*diag(C) and B by diag(R)*B 6163d6581fbSHong Zhang (if options->Trans=NOTRANS) or diag(C)*B (if options->Trans = TRANS or CONJ)." 6173d6581fbSHong Zhang We set 'options.Equil = NO' as default because additional space is needed for it. 6183d6581fbSHong Zhang */ 6193d6581fbSHong Zhang lu->options.Equil = NO; 620cffbb591SHong Zhang } else if (ftype == MAT_FACTOR_ILU){ 6210c74a584SJed Brown /* Set the default input options of ilu: */ 622cffbb591SHong Zhang ilu_set_default_options(&lu->options); 623cffbb591SHong Zhang } 6249ce50919SHong Zhang lu->options.PrintStat = NO; 6251a2e2f44SHong Zhang 6265d8b2efaSHong Zhang /* Initialize the statistics variables. */ 6275d8b2efaSHong Zhang StatInit(&lu->stat); 6288914a3f7SHong Zhang lu->lwork = 0; /* allocate space internally by system malloc */ 6299ce50919SHong Zhang 6307adad957SLisandro Dalcin ierr = PetscOptionsBegin(((PetscObject)A)->comm,((PetscObject)A)->prefix,"SuperLU Options","Mat");CHKERRQ(ierr); 631acfcf0e5SJed Brown ierr = PetscOptionsBool("-mat_superlu_equil","Equil","None",(PetscBool)lu->options.Equil,(PetscBool*)&lu->options.Equil,0);CHKERRQ(ierr); 6328914a3f7SHong Zhang ierr = PetscOptionsEList("-mat_superlu_colperm","ColPerm","None",colperm,4,colperm[3],&indx,&flg);CHKERRQ(ierr); 6339ce50919SHong Zhang if (flg) {lu->options.ColPerm = (colperm_t)indx;} 6348914a3f7SHong Zhang ierr = PetscOptionsEList("-mat_superlu_iterrefine","IterRefine","None",iterrefine,4,iterrefine[0],&indx,&flg);CHKERRQ(ierr); 6359ce50919SHong Zhang if (flg) { lu->options.IterRefine = (IterRefine_t)indx;} 636acfcf0e5SJed Brown ierr = PetscOptionsBool("-mat_superlu_symmetricmode","SymmetricMode","None",(PetscBool)lu->options.SymmetricMode,&flg,0);CHKERRQ(ierr); 6379ce50919SHong Zhang if (flg) lu->options.SymmetricMode = YES; 638*3cb270beSHong Zhang ierr = PetscOptionsReal("-mat_superlu_diagpivotthresh","DiagPivotThresh","None",lu->options.DiagPivotThresh,&real_input,&flg);CHKERRQ(ierr); 639*3cb270beSHong Zhang if (flg) lu->options.DiagPivotThresh = (double) real_input; 640acfcf0e5SJed Brown ierr = PetscOptionsBool("-mat_superlu_pivotgrowth","PivotGrowth","None",(PetscBool)lu->options.PivotGrowth,&flg,0);CHKERRQ(ierr); 6419ce50919SHong Zhang if (flg) lu->options.PivotGrowth = YES; 642acfcf0e5SJed Brown ierr = PetscOptionsBool("-mat_superlu_conditionnumber","ConditionNumber","None",(PetscBool)lu->options.ConditionNumber,&flg,0);CHKERRQ(ierr); 6439ce50919SHong Zhang if (flg) lu->options.ConditionNumber = YES; 644d7ebd59bSHong Zhang ierr = PetscOptionsEList("-mat_superlu_rowperm","rowperm","None",rowperm,2,rowperm[lu->options.RowPerm],&indx,&flg);CHKERRQ(ierr); 6459ce50919SHong Zhang if (flg) {lu->options.RowPerm = (rowperm_t)indx;} 646acfcf0e5SJed Brown ierr = PetscOptionsBool("-mat_superlu_replacetinypivot","ReplaceTinyPivot","None",(PetscBool)lu->options.ReplaceTinyPivot,&flg,0);CHKERRQ(ierr); 6479ce50919SHong Zhang if (flg) lu->options.ReplaceTinyPivot = YES; 648acfcf0e5SJed Brown ierr = PetscOptionsBool("-mat_superlu_printstat","PrintStat","None",(PetscBool)lu->options.PrintStat,&flg,0);CHKERRQ(ierr); 6499ce50919SHong Zhang if (flg) lu->options.PrintStat = YES; 6508914a3f7SHong Zhang ierr = PetscOptionsInt("-mat_superlu_lwork","size of work array in bytes used by factorization","None",lu->lwork,&lu->lwork,PETSC_NULL);CHKERRQ(ierr); 6515fe6150dSHong Zhang if (lu->lwork > 0 ){ 6525fe6150dSHong Zhang ierr = PetscMalloc(lu->lwork,&lu->work);CHKERRQ(ierr); 6535fe6150dSHong Zhang } else if (lu->lwork != 0 && lu->lwork != -1){ 65477431f27SBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," Warning: lwork %D is not supported by SUPERLU. The default lwork=0 is used.\n",lu->lwork); 6558914a3f7SHong Zhang lu->lwork = 0; 6568914a3f7SHong Zhang } 657cffbb591SHong Zhang /* ilu options */ 658*3cb270beSHong Zhang ierr = PetscOptionsReal("-mat_superlu_ilu_droptol","ILU_DropTol","None",lu->options.ILU_DropTol,&real_input,&flg);CHKERRQ(ierr); 659*3cb270beSHong Zhang if (flg) lu->options.ILU_DropTol = (double) real_input; 660*3cb270beSHong Zhang ierr = PetscOptionsReal("-mat_superlu_ilu_filltol","ILU_FillTol","None",lu->options.ILU_FillTol,&real_input,&flg);CHKERRQ(ierr); 661*3cb270beSHong Zhang if (flg) lu->options.ILU_FillTol = (double) real_input; 662*3cb270beSHong Zhang ierr = PetscOptionsReal("-mat_superlu_ilu_fillfactor","ILU_FillFactor","None",lu->options.ILU_FillFactor,&real_input,&flg);CHKERRQ(ierr); 663*3cb270beSHong Zhang if (flg) lu->options.ILU_FillFactor = (double) real_input; 664cffbb591SHong Zhang ierr = PetscOptionsInt("-mat_superlu_ilu_droprull","ILU_DropRule","None",lu->options.ILU_DropRule,&lu->options.ILU_DropRule,PETSC_NULL);CHKERRQ(ierr); 665cffbb591SHong Zhang ierr = PetscOptionsInt("-mat_superlu_ilu_norm","ILU_Norm","None",lu->options.ILU_Norm,&indx,&flg);CHKERRQ(ierr); 666cffbb591SHong Zhang if (flg){ 667cffbb591SHong Zhang lu->options.ILU_Norm = (norm_t)indx; 668cffbb591SHong Zhang } 669cffbb591SHong Zhang ierr = PetscOptionsInt("-mat_superlu_ilu_milu","ILU_MILU","None",lu->options.ILU_MILU,&indx,&flg);CHKERRQ(ierr); 670cffbb591SHong Zhang if (flg){ 671cffbb591SHong Zhang lu->options.ILU_MILU = (milu_t)indx; 672cffbb591SHong Zhang } 6739ce50919SHong Zhang PetscOptionsEnd(); 67438abfdaeSHong Zhang if (lu->options.Equil == YES) { 67538abfdaeSHong Zhang /* superlu overwrites input matrix and rhs when Equil is used, thus create A_dup to keep user's A unchanged */ 67638abfdaeSHong Zhang ierr = MatDuplicate_SeqAIJ(A,MAT_COPY_VALUES,&lu->A_dup);CHKERRQ(ierr); 67738abfdaeSHong Zhang } 6789ce50919SHong Zhang 6795d8b2efaSHong Zhang /* Allocate spaces (notice sizes are for the transpose) */ 6805d8b2efaSHong Zhang ierr = PetscMalloc(m*sizeof(PetscInt),&lu->etree);CHKERRQ(ierr); 6815d8b2efaSHong Zhang ierr = PetscMalloc(n*sizeof(PetscInt),&lu->perm_r);CHKERRQ(ierr); 6825d8b2efaSHong Zhang ierr = PetscMalloc(m*sizeof(PetscInt),&lu->perm_c);CHKERRQ(ierr); 6835d8b2efaSHong Zhang ierr = PetscMalloc(n*sizeof(PetscScalar),&lu->R);CHKERRQ(ierr); 6845d8b2efaSHong Zhang ierr = PetscMalloc(m*sizeof(PetscScalar),&lu->C);CHKERRQ(ierr); 6855d8b2efaSHong Zhang 6865d8b2efaSHong Zhang /* create rhs and solution x without allocate space for .Store */ 6875d8b2efaSHong Zhang #if defined(PETSC_USE_COMPLEX) 688*3cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 689*3cb270beSHong Zhang cCreate_Dense_Matrix(&lu->B, m, 1, PETSC_NULL, m, SLU_DN, SLU_C, SLU_GE); 690*3cb270beSHong Zhang cCreate_Dense_Matrix(&lu->X, m, 1, PETSC_NULL, m, SLU_DN, SLU_C, SLU_GE); 691*3cb270beSHong Zhang #else 6925d8b2efaSHong Zhang zCreate_Dense_Matrix(&lu->B, m, 1, PETSC_NULL, m, SLU_DN, SLU_Z, SLU_GE); 6935d8b2efaSHong Zhang zCreate_Dense_Matrix(&lu->X, m, 1, PETSC_NULL, m, SLU_DN, SLU_Z, SLU_GE); 694*3cb270beSHong Zhang #endif 695*3cb270beSHong Zhang #else 696*3cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 697*3cb270beSHong Zhang sCreate_Dense_Matrix(&lu->B, m, 1, PETSC_NULL, m, SLU_DN, SLU_S, SLU_GE); 698*3cb270beSHong Zhang sCreate_Dense_Matrix(&lu->X, m, 1, PETSC_NULL, m, SLU_DN, SLU_S, SLU_GE); 6995d8b2efaSHong Zhang #else 7005d8b2efaSHong Zhang dCreate_Dense_Matrix(&lu->B, m, 1, PETSC_NULL, m, SLU_DN, SLU_D, SLU_GE); 7015d8b2efaSHong Zhang dCreate_Dense_Matrix(&lu->X, m, 1, PETSC_NULL, m, SLU_DN, SLU_D, SLU_GE); 7025d8b2efaSHong Zhang #endif 703*3cb270beSHong Zhang #endif 7045d8b2efaSHong Zhang 70575af56d4SHong Zhang #ifdef SUPERLU2 7065c9eb25fSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatCreateNull","MatCreateNull_SuperLU",(void(*)(void))MatCreateNull_SuperLU);CHKERRQ(ierr); 70775af56d4SHong Zhang #endif 70835bd34faSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_seqaij_superlu",MatFactorGetSolverPackage_seqaij_superlu);CHKERRQ(ierr); 7095ccb76cbSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSuperluSetILUDropTol_C","MatSuperluSetILUDropTol_SuperLU",MatSuperluSetILUDropTol_SuperLU);CHKERRQ(ierr); 7105c9eb25fSBarry Smith B->spptr = lu; 711899d7b4fSKris Buschelman *F = B; 71214b396bbSKris Buschelman PetscFunctionReturn(0); 71314b396bbSKris Buschelman } 7145c9eb25fSBarry Smith EXTERN_C_END 715d954db57SHong Zhang 716