114b396bbSKris Buschelman 25a46d3faSBarry Smith /* -------------------------------------------------------------------- 35a46d3faSBarry Smith 45a46d3faSBarry Smith This file implements a subclass of the SeqAIJ matrix class that uses 55d8b2efaSHong Zhang the SuperLU sparse solver. You can use this as a starting point for 65a46d3faSBarry Smith implementing your own subclass of a PETSc matrix class. 75a46d3faSBarry Smith 85a46d3faSBarry Smith This demonstrates a way to make an implementation inheritence of a PETSc 95a46d3faSBarry Smith matrix type. This means constructing a new matrix type (SuperLU) by changing some 105a46d3faSBarry Smith of the methods of the previous type (SeqAIJ), adding additional data, and possibly 115a46d3faSBarry Smith additional method. (See any book on object oriented programming). 1214b396bbSKris Buschelman */ 1314b396bbSKris Buschelman 145a46d3faSBarry Smith /* 155a46d3faSBarry Smith Defines the data structure for the base matrix type (SeqAIJ) 165a46d3faSBarry Smith */ 17c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h> /*I "petscmat.h" I*/ 1814b396bbSKris Buschelman 195a46d3faSBarry Smith /* 205a46d3faSBarry Smith SuperLU include files 215a46d3faSBarry Smith */ 2214b396bbSKris Buschelman EXTERN_C_BEGIN 2314b396bbSKris Buschelman #if defined(PETSC_USE_COMPLEX) 243cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 253cb270beSHong Zhang #include <slu_cdefs.h> 263cb270beSHong Zhang #else 27c6db04a5SJed Brown #include <slu_zdefs.h> 283cb270beSHong Zhang #endif 293cb270beSHong Zhang #else 303cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 313cb270beSHong Zhang #include <slu_sdefs.h> 3214b396bbSKris Buschelman #else 33c6db04a5SJed Brown #include <slu_ddefs.h> 3414b396bbSKris Buschelman #endif 353cb270beSHong Zhang #endif 36c6db04a5SJed Brown #include <slu_util.h> 3714b396bbSKris Buschelman EXTERN_C_END 3814b396bbSKris Buschelman 395a46d3faSBarry Smith /* 405a46d3faSBarry Smith This is the data we are "ADDING" to the SeqAIJ matrix type to get the SuperLU matrix type 415a46d3faSBarry Smith */ 4214b396bbSKris Buschelman typedef struct { 4375af56d4SHong Zhang SuperMatrix A,L,U,B,X; 4475af56d4SHong Zhang superlu_options_t options; 45da7a1d00SHong Zhang PetscInt *perm_c; /* column permutation vector */ 46da7a1d00SHong Zhang PetscInt *perm_r; /* row permutations from partial pivoting */ 47da7a1d00SHong Zhang PetscInt *etree; 48da7a1d00SHong Zhang PetscReal *R, *C; 4975af56d4SHong Zhang char equed[1]; 50da7a1d00SHong Zhang PetscInt lwork; 5175af56d4SHong Zhang void *work; 52da7a1d00SHong Zhang PetscReal rpg, rcond; 5375af56d4SHong Zhang mem_usage_t mem_usage; 5475af56d4SHong Zhang MatStructure flg; 555d8b2efaSHong Zhang SuperLUStat_t stat; 56cae5a23dSHong Zhang Mat A_dup; 57cae5a23dSHong Zhang PetscScalar *rhs_dup; 5814b396bbSKris Buschelman 5914b396bbSKris Buschelman /* Flag to clean up (non-global) SuperLU objects during Destroy */ 60ace3abfcSBarry Smith PetscBool CleanUpSuperLU; 61f0c56d0fSKris Buschelman } Mat_SuperLU; 6214b396bbSKris Buschelman 63e0e586b9SHong Zhang extern PetscErrorCode MatFactorInfo_SuperLU(Mat,PetscViewer); 640481f469SBarry Smith extern PetscErrorCode MatLUFactorNumeric_SuperLU(Mat,Mat,const MatFactorInfo*); 65e0e586b9SHong Zhang extern PetscErrorCode MatDestroy_SuperLU(Mat); 66e0e586b9SHong Zhang extern PetscErrorCode MatView_SuperLU(Mat,PetscViewer); 67e0e586b9SHong Zhang extern PetscErrorCode MatAssemblyEnd_SuperLU(Mat,MatAssemblyType); 68e0e586b9SHong Zhang extern PetscErrorCode MatSolve_SuperLU(Mat,Vec,Vec); 69e0b74bf9SHong Zhang extern PetscErrorCode MatMatSolve_SuperLU(Mat,Mat,Mat); 70e0e586b9SHong Zhang extern PetscErrorCode MatSolveTranspose_SuperLU(Mat,Vec,Vec); 710481f469SBarry Smith extern PetscErrorCode MatLUFactorSymbolic_SuperLU(Mat,Mat,IS,IS,const MatFactorInfo*); 72e0e586b9SHong Zhang extern PetscErrorCode MatDuplicate_SuperLU(Mat, MatDuplicateOption, Mat*); 73e0e586b9SHong Zhang 745a46d3faSBarry Smith /* 755a46d3faSBarry Smith Utility function 765a46d3faSBarry Smith */ 775a46d3faSBarry Smith #undef __FUNCT__ 785a46d3faSBarry Smith #define __FUNCT__ "MatFactorInfo_SuperLU" 795a46d3faSBarry Smith PetscErrorCode MatFactorInfo_SuperLU(Mat A,PetscViewer viewer) 805a46d3faSBarry Smith { 815a46d3faSBarry Smith Mat_SuperLU *lu= (Mat_SuperLU*)A->spptr; 825a46d3faSBarry Smith PetscErrorCode ierr; 835a46d3faSBarry Smith superlu_options_t options; 845a46d3faSBarry Smith 855a46d3faSBarry Smith PetscFunctionBegin; 865a46d3faSBarry Smith /* check if matrix is superlu_dist type */ 875a46d3faSBarry Smith if (A->ops->solve != MatSolve_SuperLU) PetscFunctionReturn(0); 885a46d3faSBarry Smith 895a46d3faSBarry Smith options = lu->options; 902205254eSKarl Rupp 915a46d3faSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"SuperLU run parameters:\n");CHKERRQ(ierr); 925a46d3faSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Equil: %s\n",(options.Equil != NO) ? "YES" : "NO");CHKERRQ(ierr); 935a46d3faSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," ColPerm: %D\n",options.ColPerm);CHKERRQ(ierr); 945a46d3faSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," IterRefine: %D\n",options.IterRefine);CHKERRQ(ierr); 955a46d3faSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," SymmetricMode: %s\n",(options.SymmetricMode != NO) ? "YES" : "NO");CHKERRQ(ierr); 965a46d3faSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," DiagPivotThresh: %g\n",options.DiagPivotThresh);CHKERRQ(ierr); 975a46d3faSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," PivotGrowth: %s\n",(options.PivotGrowth != NO) ? "YES" : "NO");CHKERRQ(ierr); 985a46d3faSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," ConditionNumber: %s\n",(options.ConditionNumber != NO) ? "YES" : "NO");CHKERRQ(ierr); 995a46d3faSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," RowPerm: %D\n",options.RowPerm);CHKERRQ(ierr); 1005a46d3faSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," ReplaceTinyPivot: %s\n",(options.ReplaceTinyPivot != NO) ? "YES" : "NO");CHKERRQ(ierr); 1015a46d3faSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," PrintStat: %s\n",(options.PrintStat != NO) ? "YES" : "NO");CHKERRQ(ierr); 1025a46d3faSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," lwork: %D\n",lu->lwork);CHKERRQ(ierr); 103d5f3da31SBarry Smith if (A->factortype == MAT_FACTOR_ILU) { 104cffbb591SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ILU_DropTol: %g\n",options.ILU_DropTol);CHKERRQ(ierr); 105cffbb591SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ILU_FillTol: %g\n",options.ILU_FillTol);CHKERRQ(ierr); 106cffbb591SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ILU_FillFactor: %g\n",options.ILU_FillFactor);CHKERRQ(ierr); 107cffbb591SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ILU_DropRule: %D\n",options.ILU_DropRule);CHKERRQ(ierr); 108cffbb591SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ILU_Norm: %D\n",options.ILU_Norm);CHKERRQ(ierr); 109cffbb591SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ILU_MILU: %D\n",options.ILU_MILU);CHKERRQ(ierr); 110cffbb591SHong Zhang } 1115a46d3faSBarry Smith PetscFunctionReturn(0); 1125a46d3faSBarry Smith } 1135a46d3faSBarry Smith 1145a46d3faSBarry Smith /* 1155a46d3faSBarry Smith These are the methods provided to REPLACE the corresponding methods of the 1165a46d3faSBarry Smith SeqAIJ matrix class. Other methods of SeqAIJ are not replaced 1175a46d3faSBarry Smith */ 1185a46d3faSBarry Smith #undef __FUNCT__ 1195a46d3faSBarry Smith #define __FUNCT__ "MatLUFactorNumeric_SuperLU" 1200481f469SBarry Smith PetscErrorCode MatLUFactorNumeric_SuperLU(Mat F,Mat A,const MatFactorInfo *info) 1215a46d3faSBarry Smith { 1221d5ca7f3SHong Zhang Mat_SuperLU *lu = (Mat_SuperLU*)F->spptr; 123cae5a23dSHong Zhang Mat_SeqAIJ *aa; 1245a46d3faSBarry Smith PetscErrorCode ierr; 1255a46d3faSBarry Smith PetscInt sinfo; 1265a46d3faSBarry Smith PetscReal ferr, berr; 1275a46d3faSBarry Smith NCformat *Ustore; 1285a46d3faSBarry Smith SCformat *Lstore; 1295a46d3faSBarry Smith 1305a46d3faSBarry Smith PetscFunctionBegin; 1315a46d3faSBarry Smith if (lu->flg == SAME_NONZERO_PATTERN) { /* successing numerical factorization */ 1325a46d3faSBarry Smith lu->options.Fact = SamePattern; 1335a46d3faSBarry Smith /* Ref: ~SuperLU_3.0/EXAMPLE/dlinsolx2.c */ 1345a46d3faSBarry Smith Destroy_SuperMatrix_Store(&lu->A); 135cae5a23dSHong Zhang if (lu->options.Equil) { 136cae5a23dSHong Zhang ierr = MatCopy_SeqAIJ(A,lu->A_dup,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 137cae5a23dSHong Zhang } 1385a46d3faSBarry Smith if (lu->lwork >= 0) { 139d387c056SBarry Smith PetscStackCall("SuperLU:Destroy_SuperNode_Matrix",Destroy_SuperNode_Matrix(&lu->L)); 140d387c056SBarry Smith PetscStackCall("SuperLU:Destroy_CompCol_Matrix",Destroy_CompCol_Matrix(&lu->U)); 1415a46d3faSBarry Smith lu->options.Fact = SamePattern; 1425a46d3faSBarry Smith } 1435a46d3faSBarry Smith } 1445a46d3faSBarry Smith 1455a46d3faSBarry Smith /* Create the SuperMatrix for lu->A=A^T: 1465a46d3faSBarry Smith Since SuperLU likes column-oriented matrices,we pass it the transpose, 1475a46d3faSBarry Smith and then solve A^T X = B in MatSolve(). */ 148cae5a23dSHong Zhang if (lu->options.Equil) { 149cae5a23dSHong Zhang aa = (Mat_SeqAIJ*)(lu->A_dup)->data; 150cae5a23dSHong Zhang } else { 151cae5a23dSHong Zhang aa = (Mat_SeqAIJ*)(A)->data; 152cae5a23dSHong Zhang } 1535a46d3faSBarry Smith #if defined(PETSC_USE_COMPLEX) 1543cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 155d387c056SBarry Smith PetscStackCall("SuperLU:cCreate_CompCol_Matrix",cCreate_CompCol_Matrix(&lu->A,A->cmap->n,A->rmap->n,aa->nz,(singlecomplex*)aa->a,aa->j,aa->i,SLU_NC,SLU_C,SLU_GE)); 1563cb270beSHong Zhang #else 157d387c056SBarry Smith PetscStackCall("SuperLU:zCreate_CompCol_Matrix",zCreate_CompCol_Matrix(&lu->A,A->cmap->n,A->rmap->n,aa->nz,(doublecomplex*)aa->a,aa->j,aa->i,SLU_NC,SLU_Z,SLU_GE)); 1583cb270beSHong Zhang #endif 1593cb270beSHong Zhang #else 1603cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 161d387c056SBarry Smith PetscStackCall("SuperLU:sCreate_CompCol_Matrix",sCreate_CompCol_Matrix(&lu->A,A->cmap->n,A->rmap->n,aa->nz,aa->a,aa->j,aa->i,SLU_NC,SLU_S,SLU_GE)); 1625a46d3faSBarry Smith #else 163d387c056SBarry Smith PetscStackCall("SuperLU:dCreate_CompCol_Matrix",dCreate_CompCol_Matrix(&lu->A,A->cmap->n,A->rmap->n,aa->nz,aa->a,aa->j,aa->i,SLU_NC,SLU_D,SLU_GE)); 1645a46d3faSBarry Smith #endif 1653cb270beSHong Zhang #endif 1665a46d3faSBarry Smith 1675a46d3faSBarry Smith /* Numerical factorization */ 1685a46d3faSBarry Smith lu->B.ncol = 0; /* Indicate not to solve the system */ 169d5f3da31SBarry Smith if (F->factortype == MAT_FACTOR_LU) { 1705a46d3faSBarry Smith #if defined(PETSC_USE_COMPLEX) 1713cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 172d387c056SBarry Smith PetscStackCall("SuperLU:cgssvx",cgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 1733cb270beSHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 174d387c056SBarry Smith &lu->mem_usage, &lu->stat, &sinfo)); 1753cb270beSHong Zhang #else 176d387c056SBarry Smith PetscStackCall("SuperLU:zgssvx",zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 1775a46d3faSBarry Smith &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 178d387c056SBarry Smith &lu->mem_usage, &lu->stat, &sinfo)); 1793cb270beSHong Zhang #endif 1803cb270beSHong Zhang #else 1813cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 182d387c056SBarry Smith PetscStackCall("SuperLU:sgssvx",sgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 1833cb270beSHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 184d387c056SBarry Smith &lu->mem_usage, &lu->stat, &sinfo)); 1855a46d3faSBarry Smith #else 186d387c056SBarry Smith PetscStackCall("SuperLU:dgssvx",dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 1875a46d3faSBarry Smith &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 188d387c056SBarry Smith &lu->mem_usage, &lu->stat, &sinfo)); 1895a46d3faSBarry Smith #endif 1903cb270beSHong Zhang #endif 191d5f3da31SBarry Smith } else if (F->factortype == MAT_FACTOR_ILU) { 192cffbb591SHong Zhang /* Compute the incomplete factorization, condition number and pivot growth */ 193cffbb591SHong Zhang #if defined(PETSC_USE_COMPLEX) 1943cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 195d387c056SBarry Smith PetscStackCall("SuperLU:cgsisx",cgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r,lu->etree, lu->equed, lu->R, lu->C, 1963cb270beSHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, 197d387c056SBarry Smith &lu->mem_usage, &lu->stat, &sinfo)); 1983cb270beSHong Zhang #else 199d387c056SBarry Smith PetscStackCall("SuperLU:zgsisx",zgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r,lu->etree, lu->equed, lu->R, lu->C, 200cffbb591SHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, 201d387c056SBarry Smith &lu->mem_usage, &lu->stat, &sinfo)); 2023cb270beSHong Zhang #endif 2033cb270beSHong Zhang #else 2043cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 205d387c056SBarry Smith PetscStackCall("SuperLU:sgsisx",sgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 2063cb270beSHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, 207d387c056SBarry Smith &lu->mem_usage, &lu->stat, &sinfo)); 208cffbb591SHong Zhang #else 209d387c056SBarry Smith PetscStackCall("SuperLU:dgsisx",dgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 210cffbb591SHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, 211d387c056SBarry Smith &lu->mem_usage, &lu->stat, &sinfo)); 212cffbb591SHong Zhang #endif 2133cb270beSHong Zhang #endif 214f23aa3ddSBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported"); 2155a46d3faSBarry Smith if (!sinfo || sinfo == lu->A.ncol+1) { 2162205254eSKarl Rupp if (lu->options.PivotGrowth) { 2175a46d3faSBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," Recip. pivot growth = %e\n", lu->rpg); 2182205254eSKarl Rupp } 2192205254eSKarl Rupp if (lu->options.ConditionNumber) { 2205a46d3faSBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," Recip. condition number = %e\n", lu->rcond); 2212205254eSKarl Rupp } 2225a46d3faSBarry Smith } else if (sinfo > 0) { 2235a46d3faSBarry Smith if (lu->lwork == -1) { 2245a46d3faSBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," ** Estimated memory: %D bytes\n", sinfo - lu->A.ncol); 2259308c137SBarry Smith } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot in row %D",sinfo); 226f23aa3ddSBarry Smith } else SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_LIB, "info = %D, the %D-th argument in gssvx() had an illegal value", sinfo,-sinfo); 2275a46d3faSBarry Smith 2285a46d3faSBarry Smith if (lu->options.PrintStat) { 2295a46d3faSBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF,"MatLUFactorNumeric_SuperLU():\n"); 230d387c056SBarry Smith PetscStackCall("SuperLU:StatPrint",StatPrint(&lu->stat)); 2315a46d3faSBarry Smith Lstore = (SCformat*) lu->L.Store; 2325a46d3faSBarry Smith Ustore = (NCformat*) lu->U.Store; 2335a46d3faSBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," No of nonzeros in factor L = %D\n", Lstore->nnz); 2345a46d3faSBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," No of nonzeros in factor U = %D\n", Ustore->nnz); 2355a46d3faSBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," No of nonzeros in L+U = %D\n", Lstore->nnz + Ustore->nnz - lu->A.ncol); 2366da386baSHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," L\\U MB %.3f\ttotal MB needed %.3f\n", 2376da386baSHong Zhang lu->mem_usage.for_lu/1e6, lu->mem_usage.total_needed/1e6); 2385a46d3faSBarry Smith } 2395a46d3faSBarry Smith 2405a46d3faSBarry Smith lu->flg = SAME_NONZERO_PATTERN; 2411d5ca7f3SHong Zhang F->ops->solve = MatSolve_SuperLU; 2421d5ca7f3SHong Zhang F->ops->solvetranspose = MatSolveTranspose_SuperLU; 243e0b74bf9SHong Zhang F->ops->matsolve = MatMatSolve_SuperLU; 2445a46d3faSBarry Smith PetscFunctionReturn(0); 2455a46d3faSBarry Smith } 2465a46d3faSBarry Smith 24714b396bbSKris Buschelman #undef __FUNCT__ 248f0c56d0fSKris Buschelman #define __FUNCT__ "MatDestroy_SuperLU" 249dfbe8321SBarry Smith PetscErrorCode MatDestroy_SuperLU(Mat A) 25014b396bbSKris Buschelman { 251dfbe8321SBarry Smith PetscErrorCode ierr; 252f0c56d0fSKris Buschelman Mat_SuperLU *lu=(Mat_SuperLU*)A->spptr; 25314b396bbSKris Buschelman 25414b396bbSKris Buschelman PetscFunctionBegin; 255bf0cc555SLisandro Dalcin if (lu && lu->CleanUpSuperLU) { /* Free the SuperLU datastructures */ 256d387c056SBarry Smith PetscStackCall("SuperLU:Destroy_SuperMatrix_Store",Destroy_SuperMatrix_Store(&lu->A)); 257d387c056SBarry Smith PetscStackCall("SuperLU:Destroy_SuperMatrix_Store",Destroy_SuperMatrix_Store(&lu->B)); 258d387c056SBarry Smith PetscStackCall("SuperLU:Destroy_SuperMatrix_Store",Destroy_SuperMatrix_Store(&lu->X)); 259d387c056SBarry Smith PetscStackCall("SuperLU:StatFree",StatFree(&lu->stat)); 2600e742b69SHong Zhang if (lu->lwork >= 0) { 261d387c056SBarry Smith PetscStackCall("SuperLU:Destroy_SuperNode_Matrix",Destroy_SuperNode_Matrix(&lu->L)); 262d387c056SBarry Smith PetscStackCall("SuperLU:Destroy_CompCol_Matrix",Destroy_CompCol_Matrix(&lu->U)); 2630e742b69SHong Zhang } 2640e742b69SHong Zhang } 265bf0cc555SLisandro Dalcin if (lu) { 2669ce50919SHong Zhang ierr = PetscFree(lu->etree);CHKERRQ(ierr); 2679ce50919SHong Zhang ierr = PetscFree(lu->perm_r);CHKERRQ(ierr); 2689ce50919SHong Zhang ierr = PetscFree(lu->perm_c);CHKERRQ(ierr); 2699ce50919SHong Zhang ierr = PetscFree(lu->R);CHKERRQ(ierr); 2709ce50919SHong Zhang ierr = PetscFree(lu->C);CHKERRQ(ierr); 271bf0cc555SLisandro Dalcin ierr = PetscFree(lu->rhs_dup);CHKERRQ(ierr); 272bf0cc555SLisandro Dalcin ierr = MatDestroy(&lu->A_dup);CHKERRQ(ierr); 273bf0cc555SLisandro Dalcin } 274bf0cc555SLisandro Dalcin ierr = PetscFree(A->spptr);CHKERRQ(ierr); 2750e742b69SHong Zhang 276d954db57SHong Zhang /* clear composed functions */ 277bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorGetSolverPackage_C",NULL);CHKERRQ(ierr); 278bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatSuperluSetILUDropTol_C",NULL);CHKERRQ(ierr); 279d954db57SHong Zhang 280b24902e0SBarry Smith ierr = MatDestroy_SeqAIJ(A);CHKERRQ(ierr); 28114b396bbSKris Buschelman PetscFunctionReturn(0); 28214b396bbSKris Buschelman } 28314b396bbSKris Buschelman 28414b396bbSKris Buschelman #undef __FUNCT__ 285f0c56d0fSKris Buschelman #define __FUNCT__ "MatView_SuperLU" 286dfbe8321SBarry Smith PetscErrorCode MatView_SuperLU(Mat A,PetscViewer viewer) 28714b396bbSKris Buschelman { 288dfbe8321SBarry Smith PetscErrorCode ierr; 289ace3abfcSBarry Smith PetscBool iascii; 29014b396bbSKris Buschelman PetscViewerFormat format; 29114b396bbSKris Buschelman 29214b396bbSKris Buschelman PetscFunctionBegin; 293251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 29432077d6dSBarry Smith if (iascii) { 29514b396bbSKris Buschelman ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 2962f59403fSHong Zhang if (format == PETSC_VIEWER_ASCII_INFO) { 297f0c56d0fSKris Buschelman ierr = MatFactorInfo_SuperLU(A,viewer);CHKERRQ(ierr); 29814b396bbSKris Buschelman } 29914b396bbSKris Buschelman } 30014b396bbSKris Buschelman PetscFunctionReturn(0); 30114b396bbSKris Buschelman } 30214b396bbSKris Buschelman 30314b396bbSKris Buschelman 30414b396bbSKris Buschelman #undef __FUNCT__ 305c29e884eSHong Zhang #define __FUNCT__ "MatSolve_SuperLU_Private" 306c29e884eSHong Zhang PetscErrorCode MatSolve_SuperLU_Private(Mat A,Vec b,Vec x) 30714b396bbSKris Buschelman { 308f0c56d0fSKris Buschelman Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr; 30975af56d4SHong Zhang PetscScalar *barray,*xarray; 310dfbe8321SBarry Smith PetscErrorCode ierr; 311077aedafSJed Brown PetscInt info,i,n; 312da7a1d00SHong Zhang PetscReal ferr,berr; 313dff31646SBarry Smith static PetscBool cite = PETSC_FALSE; 31414b396bbSKris Buschelman 31514b396bbSKris Buschelman PetscFunctionBegin; 3162205254eSKarl Rupp if (lu->lwork == -1) PetscFunctionReturn(0); 317dff31646SBarry Smith ierr = PetscCitationsRegister("@article{superlu99,\n author = {James W. Demmel and Stanley C. Eisenstat and\n John R. Gilbert and Xiaoye S. Li and Joseph W. H. Liu},\n title = {A supernodal approach to sparse partial pivoting},\n journal = {SIAM J. Matrix Analysis and Applications},\n year = {1999},\n volume = {20},\n number = {3},\n pages = {720-755}\n}\n",&cite);CHKERRQ(ierr); 318cae5a23dSHong Zhang 319077aedafSJed Brown ierr = VecGetLocalSize(x,&n);CHKERRQ(ierr); 32075af56d4SHong Zhang lu->B.ncol = 1; /* Set the number of right-hand side */ 321cae5a23dSHong Zhang if (lu->options.Equil && !lu->rhs_dup) { 322cae5a23dSHong Zhang /* superlu overwrites b when Equil is used, thus create rhs_dup to keep user's b unchanged */ 323785e854fSJed Brown ierr = PetscMalloc1(n,&lu->rhs_dup);CHKERRQ(ierr); 324cae5a23dSHong Zhang } 325cae5a23dSHong Zhang if (lu->options.Equil) { 326cae5a23dSHong Zhang /* Copy b into rsh_dup */ 32775af56d4SHong Zhang ierr = VecGetArray(b,&barray);CHKERRQ(ierr); 328cae5a23dSHong Zhang ierr = PetscMemcpy(lu->rhs_dup,barray,n*sizeof(PetscScalar));CHKERRQ(ierr); 329cae5a23dSHong Zhang ierr = VecRestoreArray(b,&barray);CHKERRQ(ierr); 330cae5a23dSHong Zhang barray = lu->rhs_dup; 331cae5a23dSHong Zhang } else { 332cae5a23dSHong Zhang ierr = VecGetArray(b,&barray);CHKERRQ(ierr); 333cae5a23dSHong Zhang } 33475af56d4SHong Zhang ierr = VecGetArray(x,&xarray);CHKERRQ(ierr); 3355fe6150dSHong Zhang 3365fe6150dSHong Zhang #if defined(PETSC_USE_COMPLEX) 3373cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 3383cb270beSHong Zhang ((DNformat*)lu->B.Store)->nzval = (singlecomplex*)barray; 3393cb270beSHong Zhang ((DNformat*)lu->X.Store)->nzval = (singlecomplex*)xarray; 3403cb270beSHong Zhang #else 3415fe6150dSHong Zhang ((DNformat*)lu->B.Store)->nzval = (doublecomplex*)barray; 3425fe6150dSHong Zhang ((DNformat*)lu->X.Store)->nzval = (doublecomplex*)xarray; 3433cb270beSHong Zhang #endif 3445fe6150dSHong Zhang #else 3455fe6150dSHong Zhang ((DNformat*)lu->B.Store)->nzval = barray; 34675af56d4SHong Zhang ((DNformat*)lu->X.Store)->nzval = xarray; 3475fe6150dSHong Zhang #endif 34875af56d4SHong Zhang 34975af56d4SHong Zhang lu->options.Fact = FACTORED; /* Indicate the factored form of A is supplied. */ 350d5f3da31SBarry Smith if (A->factortype == MAT_FACTOR_LU) { 3518914a3f7SHong Zhang #if defined(PETSC_USE_COMPLEX) 3523cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 353d387c056SBarry Smith PetscStackCall("SuperLU:cgssvx",cgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 3543cb270beSHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 355d387c056SBarry Smith &lu->mem_usage, &lu->stat, &info)); 3563cb270beSHong Zhang #else 357d387c056SBarry Smith PetscStackCall("SuperLU:zgssvx",zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 3588914a3f7SHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 359d387c056SBarry Smith &lu->mem_usage, &lu->stat, &info)); 3603cb270beSHong Zhang #endif 3613cb270beSHong Zhang #else 3623cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 363d387c056SBarry Smith PetscStackCall("SuperLU:sgssvx",sgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 3643cb270beSHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 365d387c056SBarry Smith &lu->mem_usage, &lu->stat, &info)); 3668914a3f7SHong Zhang #else 367d387c056SBarry Smith PetscStackCall("SuperLU:dgssvx",dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 36875af56d4SHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 369d387c056SBarry Smith &lu->mem_usage, &lu->stat, &info)); 3708914a3f7SHong Zhang #endif 3713cb270beSHong Zhang #endif 372d5f3da31SBarry Smith } else if (A->factortype == MAT_FACTOR_ILU) { 373cffbb591SHong Zhang #if defined(PETSC_USE_COMPLEX) 3743cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 375d387c056SBarry Smith PetscStackCall("SuperLU:cgsisx",cgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 3763cb270beSHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, 377d387c056SBarry Smith &lu->mem_usage, &lu->stat, &info)); 3783cb270beSHong Zhang #else 379d387c056SBarry Smith PetscStackCall("SuperLU:zgsisx",zgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 380cffbb591SHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, 381d387c056SBarry Smith &lu->mem_usage, &lu->stat, &info)); 3823cb270beSHong Zhang #endif 3833cb270beSHong Zhang #else 3843cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 385d387c056SBarry Smith PetscStackCall("SuperLU:sgsisx",sgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 3863cb270beSHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, 387d387c056SBarry Smith &lu->mem_usage, &lu->stat, &info)); 388cffbb591SHong Zhang #else 389d387c056SBarry Smith PetscStackCall("SuperLU:dgsisx",dgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 390cffbb591SHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, 391d387c056SBarry Smith &lu->mem_usage, &lu->stat, &info)); 392cffbb591SHong Zhang #endif 3933cb270beSHong Zhang #endif 394f23aa3ddSBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported"); 395cae5a23dSHong Zhang if (!lu->options.Equil) { 39675af56d4SHong Zhang ierr = VecRestoreArray(b,&barray);CHKERRQ(ierr); 397cae5a23dSHong Zhang } 39875af56d4SHong Zhang ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr); 39975af56d4SHong Zhang 400958c9bccSBarry Smith if (!info || info == lu->A.ncol+1) { 40175af56d4SHong Zhang if (lu->options.IterRefine) { 4028914a3f7SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF,"Iterative Refinement:\n"); 4038914a3f7SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," %8s%8s%16s%16s\n", "rhs", "Steps", "FERR", "BERR"); 4042205254eSKarl Rupp for (i = 0; i < 1; ++i) { 4055d8b2efaSHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," %8d%8d%16e%16e\n", i+1, lu->stat.RefineSteps, ferr, berr); 40675af56d4SHong Zhang } 4072205254eSKarl Rupp } 4088914a3f7SHong Zhang } else if (info > 0) { 4098914a3f7SHong Zhang if (lu->lwork == -1) { 41077431f27SBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," ** Estimated memory: %D bytes\n", info - lu->A.ncol); 4118914a3f7SHong Zhang } else { 41277431f27SBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," Warning: gssvx() returns info %D\n",info); 4138914a3f7SHong Zhang } 414f23aa3ddSBarry Smith } else if (info < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_LIB, "info = %D, the %D-th argument in gssvx() had an illegal value", info,-info); 41514b396bbSKris Buschelman 4168914a3f7SHong Zhang if (lu->options.PrintStat) { 4178914a3f7SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF,"MatSolve__SuperLU():\n"); 418d387c056SBarry Smith PetscStackCall("SuperLU:StatPrint",StatPrint(&lu->stat)); 4198914a3f7SHong Zhang } 42075af56d4SHong Zhang PetscFunctionReturn(0); 42175af56d4SHong Zhang } 42214b396bbSKris Buschelman 42314b396bbSKris Buschelman #undef __FUNCT__ 424c29e884eSHong Zhang #define __FUNCT__ "MatSolve_SuperLU" 425c29e884eSHong Zhang PetscErrorCode MatSolve_SuperLU(Mat A,Vec b,Vec x) 426c29e884eSHong Zhang { 427c29e884eSHong Zhang Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr; 428c29e884eSHong Zhang PetscErrorCode ierr; 429c29e884eSHong Zhang 430c29e884eSHong Zhang PetscFunctionBegin; 431c29e884eSHong Zhang lu->options.Trans = TRANS; 4322205254eSKarl Rupp 433c29e884eSHong Zhang ierr = MatSolve_SuperLU_Private(A,b,x);CHKERRQ(ierr); 434c29e884eSHong Zhang PetscFunctionReturn(0); 435c29e884eSHong Zhang } 436c29e884eSHong Zhang 437c29e884eSHong Zhang #undef __FUNCT__ 438c7c1fe80SHong Zhang #define __FUNCT__ "MatSolveTranspose_SuperLU" 439c7c1fe80SHong Zhang PetscErrorCode MatSolveTranspose_SuperLU(Mat A,Vec b,Vec x) 440c7c1fe80SHong Zhang { 441c7c1fe80SHong Zhang Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr; 442c7c1fe80SHong Zhang PetscErrorCode ierr; 443c7c1fe80SHong Zhang 444c7c1fe80SHong Zhang PetscFunctionBegin; 445c7c1fe80SHong Zhang lu->options.Trans = NOTRANS; 4462205254eSKarl Rupp 447c29e884eSHong Zhang ierr = MatSolve_SuperLU_Private(A,b,x);CHKERRQ(ierr); 448c7c1fe80SHong Zhang PetscFunctionReturn(0); 449c7c1fe80SHong Zhang } 450c7c1fe80SHong Zhang 451e0b74bf9SHong Zhang #undef __FUNCT__ 452e0b74bf9SHong Zhang #define __FUNCT__ "MatMatSolve_SuperLU" 453e0b74bf9SHong Zhang PetscErrorCode MatMatSolve_SuperLU(Mat A,Mat B,Mat X) 454e0b74bf9SHong Zhang { 455e0b74bf9SHong Zhang Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr; 456cd723cd1SBarry Smith PetscBool flg; 457cd723cd1SBarry Smith PetscErrorCode ierr; 458e0b74bf9SHong Zhang 459e0b74bf9SHong Zhang PetscFunctionBegin; 4600298fd71SBarry Smith ierr = PetscObjectTypeCompareAny((PetscObject)B,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr); 461ce94432eSBarry Smith if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix"); 4620298fd71SBarry Smith ierr = PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr); 463ce94432eSBarry Smith if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix"); 4642205254eSKarl Rupp lu->options.Trans = TRANS; 465e0b74bf9SHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatMatSolve_SuperLU() is not implemented yet"); 466e0b74bf9SHong Zhang PetscFunctionReturn(0); 467e0b74bf9SHong Zhang } 468e0b74bf9SHong Zhang 46914b396bbSKris Buschelman /* 47014b396bbSKris Buschelman Note the r permutation is ignored 47114b396bbSKris Buschelman */ 47214b396bbSKris Buschelman #undef __FUNCT__ 473f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_SuperLU" 4740481f469SBarry Smith PetscErrorCode MatLUFactorSymbolic_SuperLU(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info) 47514b396bbSKris Buschelman { 4761d5ca7f3SHong Zhang Mat_SuperLU *lu = (Mat_SuperLU*)(F->spptr); 477b24902e0SBarry Smith 478b24902e0SBarry Smith PetscFunctionBegin; 479b24902e0SBarry Smith lu->flg = DIFFERENT_NONZERO_PATTERN; 480b24902e0SBarry Smith lu->CleanUpSuperLU = PETSC_TRUE; 4811d5ca7f3SHong Zhang F->ops->lufactornumeric = MatLUFactorNumeric_SuperLU; 482b24902e0SBarry Smith PetscFunctionReturn(0); 483b24902e0SBarry Smith } 484b24902e0SBarry Smith 48535bd34faSBarry Smith #undef __FUNCT__ 4865ccb76cbSHong Zhang #define __FUNCT__ "MatSuperluSetILUDropTol_SuperLU" 487b2573a8aSBarry Smith static PetscErrorCode MatSuperluSetILUDropTol_SuperLU(Mat F,PetscReal dtol) 4885ccb76cbSHong Zhang { 4895ccb76cbSHong Zhang Mat_SuperLU *lu= (Mat_SuperLU*)F->spptr; 4905ccb76cbSHong Zhang 4915ccb76cbSHong Zhang PetscFunctionBegin; 4925ccb76cbSHong Zhang lu->options.ILU_DropTol = dtol; 4935ccb76cbSHong Zhang PetscFunctionReturn(0); 4945ccb76cbSHong Zhang } 4955ccb76cbSHong Zhang 4965ccb76cbSHong Zhang #undef __FUNCT__ 4975ccb76cbSHong Zhang #define __FUNCT__ "MatSuperluSetILUDropTol" 4985ccb76cbSHong Zhang /*@ 4995ccb76cbSHong Zhang MatSuperluSetILUDropTol - Set SuperLU ILU drop tolerance 5005ccb76cbSHong Zhang Logically Collective on Mat 5015ccb76cbSHong Zhang 5025ccb76cbSHong Zhang Input Parameters: 5035ccb76cbSHong Zhang + F - the factored matrix obtained by calling MatGetFactor() from PETSc-SuperLU interface 5045ccb76cbSHong Zhang - dtol - drop tolerance 5055ccb76cbSHong Zhang 5065ccb76cbSHong Zhang Options Database: 5075ccb76cbSHong Zhang . -mat_superlu_ilu_droptol <dtol> 5085ccb76cbSHong Zhang 5095ccb76cbSHong Zhang Level: beginner 5105ccb76cbSHong Zhang 5115ccb76cbSHong Zhang References: SuperLU Users' Guide 5125ccb76cbSHong Zhang 5135ccb76cbSHong Zhang .seealso: MatGetFactor() 5145ccb76cbSHong Zhang @*/ 5155ccb76cbSHong Zhang PetscErrorCode MatSuperluSetILUDropTol(Mat F,PetscReal dtol) 5165ccb76cbSHong Zhang { 5175ccb76cbSHong Zhang PetscErrorCode ierr; 5185ccb76cbSHong Zhang 5195ccb76cbSHong Zhang PetscFunctionBegin; 5205ccb76cbSHong Zhang PetscValidHeaderSpecific(F,MAT_CLASSID,1); 5215ccb76cbSHong Zhang PetscValidLogicalCollectiveInt(F,dtol,2); 5225ccb76cbSHong Zhang ierr = PetscTryMethod(F,"MatSuperluSetILUDropTol_C",(Mat,PetscReal),(F,dtol));CHKERRQ(ierr); 5235ccb76cbSHong Zhang PetscFunctionReturn(0); 5245ccb76cbSHong Zhang } 5255ccb76cbSHong Zhang 5265ccb76cbSHong Zhang #undef __FUNCT__ 52735bd34faSBarry Smith #define __FUNCT__ "MatFactorGetSolverPackage_seqaij_superlu" 52835bd34faSBarry Smith PetscErrorCode MatFactorGetSolverPackage_seqaij_superlu(Mat A,const MatSolverPackage *type) 52935bd34faSBarry Smith { 53035bd34faSBarry Smith PetscFunctionBegin; 5312692d6eeSBarry Smith *type = MATSOLVERSUPERLU; 53235bd34faSBarry Smith PetscFunctionReturn(0); 53335bd34faSBarry Smith } 53435bd34faSBarry Smith 535b24902e0SBarry Smith 536b24902e0SBarry Smith /*MC 537ba992d64SSatish Balay MATSOLVERSUPERLU = "superlu" - A solver package providing solvers LU and ILU for sequential matrices 538b24902e0SBarry Smith via the external package SuperLU. 539b24902e0SBarry Smith 540e2e64c6bSBarry Smith Use ./configure --download-superlu to have PETSc installed with SuperLU 541b24902e0SBarry Smith 542b24902e0SBarry Smith Options Database Keys: 543e08999f5SMatthew G Knepley + -mat_superlu_equil <FALSE> - Equil (None) 544e08999f5SMatthew G Knepley . -mat_superlu_colperm <COLAMD> - (choose one of) NATURAL MMD_ATA MMD_AT_PLUS_A COLAMD 545e08999f5SMatthew G Knepley . -mat_superlu_iterrefine <NOREFINE> - (choose one of) NOREFINE SINGLE DOUBLE EXTRA 546e08999f5SMatthew G Knepley . -mat_superlu_symmetricmode: <FALSE> - SymmetricMode (None) 547e08999f5SMatthew G Knepley . -mat_superlu_diagpivotthresh <1> - DiagPivotThresh (None) 548e08999f5SMatthew G Knepley . -mat_superlu_pivotgrowth <FALSE> - PivotGrowth (None) 549e08999f5SMatthew G Knepley . -mat_superlu_conditionnumber <FALSE> - ConditionNumber (None) 550e08999f5SMatthew G Knepley . -mat_superlu_rowperm <NOROWPERM> - (choose one of) NOROWPERM LargeDiag 551e08999f5SMatthew G Knepley . -mat_superlu_replacetinypivot <FALSE> - ReplaceTinyPivot (None) 552e08999f5SMatthew G Knepley . -mat_superlu_printstat <FALSE> - PrintStat (None) 553e08999f5SMatthew G Knepley . -mat_superlu_lwork <0> - size of work array in bytes used by factorization (None) 554e08999f5SMatthew G Knepley . -mat_superlu_ilu_droptol <0> - ILU_DropTol (None) 555e08999f5SMatthew G Knepley . -mat_superlu_ilu_filltol <0> - ILU_FillTol (None) 556e08999f5SMatthew G Knepley . -mat_superlu_ilu_fillfactor <0> - ILU_FillFactor (None) 557e08999f5SMatthew G Knepley . -mat_superlu_ilu_droprull <0> - ILU_DropRule (None) 558e08999f5SMatthew G Knepley . -mat_superlu_ilu_norm <0> - ILU_Norm (None) 559e08999f5SMatthew G Knepley - -mat_superlu_ilu_milu <0> - ILU_MILU (None) 560b24902e0SBarry Smith 5612692d6eeSBarry Smith Notes: Do not confuse this with MATSOLVERSUPERLU_DIST which is for parallel sparse solves 5625c9eb25fSBarry Smith 563b24902e0SBarry Smith Level: beginner 564b24902e0SBarry Smith 565d45987f3SHong Zhang .seealso: PCLU, PCILU, MATSOLVERSUPERLU_DIST, MATSOLVERMUMPS, PCFactorSetMatSolverPackage(), MatSolverPackage 566b24902e0SBarry Smith M*/ 567b24902e0SBarry Smith 568b24902e0SBarry Smith #undef __FUNCT__ 569b24902e0SBarry Smith #define __FUNCT__ "MatGetFactor_seqaij_superlu" 5708cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactor_seqaij_superlu(Mat A,MatFactorType ftype,Mat *F) 571b24902e0SBarry Smith { 57214b396bbSKris Buschelman Mat B; 573f0c56d0fSKris Buschelman Mat_SuperLU *lu; 5746849ba73SBarry Smith PetscErrorCode ierr; 5755d8b2efaSHong Zhang PetscInt indx,m=A->rmap->n,n=A->cmap->n; 576ace3abfcSBarry Smith PetscBool flg; 5773cb270beSHong Zhang PetscReal real_input; 5785ca28756SHong Zhang const char *colperm[] ={"NATURAL","MMD_ATA","MMD_AT_PLUS_A","COLAMD"}; /* MY_PERMC - not supported by the petsc interface yet */ 5795ca28756SHong Zhang const char *iterrefine[]={"NOREFINE", "SINGLE", "DOUBLE", "EXTRA"}; 5805ca28756SHong Zhang const char *rowperm[] ={"NOROWPERM", "LargeDiag"}; /* MY_PERMC - not supported by the petsc interface yet */ 58114b396bbSKris Buschelman 58214b396bbSKris Buschelman PetscFunctionBegin; 583ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr); 584d0f46423SBarry Smith ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 5857adad957SLisandro Dalcin ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 5860298fd71SBarry Smith ierr = MatSeqAIJSetPreallocation(B,0,NULL);CHKERRQ(ierr); 587f0c56d0fSKris Buschelman 588cffbb591SHong Zhang if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU) { 589b24902e0SBarry Smith B->ops->lufactorsymbolic = MatLUFactorSymbolic_SuperLU; 590cffbb591SHong Zhang B->ops->ilufactorsymbolic = MatLUFactorSymbolic_SuperLU; 591b3a44c85SBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported"); 592cffbb591SHong Zhang 593b24902e0SBarry Smith B->ops->destroy = MatDestroy_SuperLU; 5943519fcdcSHong Zhang B->ops->view = MatView_SuperLU; 595d5f3da31SBarry Smith B->factortype = ftype; 59694b7f48cSBarry Smith B->assembled = PETSC_TRUE; /* required by -ksp_view */ 5975c9eb25fSBarry Smith B->preallocated = PETSC_TRUE; 59814b396bbSKris Buschelman 599*b00a9115SJed Brown ierr = PetscNewLog(B,&lu);CHKERRQ(ierr); 600cae5a23dSHong Zhang 601cffbb591SHong Zhang if (ftype == MAT_FACTOR_LU) { 6029ce50919SHong Zhang set_default_options(&lu->options); 6033d6581fbSHong Zhang /* Comments from SuperLU_4.0/SRC/dgssvx.c: 6043d6581fbSHong Zhang "Whether or not the system will be equilibrated depends on the 6053d6581fbSHong Zhang scaling of the matrix A, but if equilibration is used, A is 6063d6581fbSHong Zhang overwritten by diag(R)*A*diag(C) and B by diag(R)*B 6073d6581fbSHong Zhang (if options->Trans=NOTRANS) or diag(C)*B (if options->Trans = TRANS or CONJ)." 6083d6581fbSHong Zhang We set 'options.Equil = NO' as default because additional space is needed for it. 6093d6581fbSHong Zhang */ 6103d6581fbSHong Zhang lu->options.Equil = NO; 611cffbb591SHong Zhang } else if (ftype == MAT_FACTOR_ILU) { 6120c74a584SJed Brown /* Set the default input options of ilu: */ 613d387c056SBarry Smith PetscStackCall("SuperLU:ilu_set_default_options",ilu_set_default_options(&lu->options)); 614cffbb591SHong Zhang } 6159ce50919SHong Zhang lu->options.PrintStat = NO; 6161a2e2f44SHong Zhang 6175d8b2efaSHong Zhang /* Initialize the statistics variables. */ 618d387c056SBarry Smith PetscStackCall("SuperLU:StatInit",StatInit(&lu->stat)); 6198914a3f7SHong Zhang lu->lwork = 0; /* allocate space internally by system malloc */ 6209ce50919SHong Zhang 621ce94432eSBarry Smith ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"SuperLU Options","Mat");CHKERRQ(ierr); 622acfcf0e5SJed Brown ierr = PetscOptionsBool("-mat_superlu_equil","Equil","None",(PetscBool)lu->options.Equil,(PetscBool*)&lu->options.Equil,0);CHKERRQ(ierr); 6238914a3f7SHong Zhang ierr = PetscOptionsEList("-mat_superlu_colperm","ColPerm","None",colperm,4,colperm[3],&indx,&flg);CHKERRQ(ierr); 6242205254eSKarl Rupp if (flg) lu->options.ColPerm = (colperm_t)indx; 6258914a3f7SHong Zhang ierr = PetscOptionsEList("-mat_superlu_iterrefine","IterRefine","None",iterrefine,4,iterrefine[0],&indx,&flg);CHKERRQ(ierr); 6262205254eSKarl Rupp if (flg) lu->options.IterRefine = (IterRefine_t)indx; 627acfcf0e5SJed Brown ierr = PetscOptionsBool("-mat_superlu_symmetricmode","SymmetricMode","None",(PetscBool)lu->options.SymmetricMode,&flg,0);CHKERRQ(ierr); 6289ce50919SHong Zhang if (flg) lu->options.SymmetricMode = YES; 6293cb270beSHong Zhang ierr = PetscOptionsReal("-mat_superlu_diagpivotthresh","DiagPivotThresh","None",lu->options.DiagPivotThresh,&real_input,&flg);CHKERRQ(ierr); 6303cb270beSHong Zhang if (flg) lu->options.DiagPivotThresh = (double) real_input; 631acfcf0e5SJed Brown ierr = PetscOptionsBool("-mat_superlu_pivotgrowth","PivotGrowth","None",(PetscBool)lu->options.PivotGrowth,&flg,0);CHKERRQ(ierr); 6329ce50919SHong Zhang if (flg) lu->options.PivotGrowth = YES; 633acfcf0e5SJed Brown ierr = PetscOptionsBool("-mat_superlu_conditionnumber","ConditionNumber","None",(PetscBool)lu->options.ConditionNumber,&flg,0);CHKERRQ(ierr); 6349ce50919SHong Zhang if (flg) lu->options.ConditionNumber = YES; 635d7ebd59bSHong Zhang ierr = PetscOptionsEList("-mat_superlu_rowperm","rowperm","None",rowperm,2,rowperm[lu->options.RowPerm],&indx,&flg);CHKERRQ(ierr); 6362205254eSKarl Rupp if (flg) lu->options.RowPerm = (rowperm_t)indx; 637acfcf0e5SJed Brown ierr = PetscOptionsBool("-mat_superlu_replacetinypivot","ReplaceTinyPivot","None",(PetscBool)lu->options.ReplaceTinyPivot,&flg,0);CHKERRQ(ierr); 6389ce50919SHong Zhang if (flg) lu->options.ReplaceTinyPivot = YES; 639acfcf0e5SJed Brown ierr = PetscOptionsBool("-mat_superlu_printstat","PrintStat","None",(PetscBool)lu->options.PrintStat,&flg,0);CHKERRQ(ierr); 6409ce50919SHong Zhang if (flg) lu->options.PrintStat = YES; 6410298fd71SBarry Smith ierr = PetscOptionsInt("-mat_superlu_lwork","size of work array in bytes used by factorization","None",lu->lwork,&lu->lwork,NULL);CHKERRQ(ierr); 6425fe6150dSHong Zhang if (lu->lwork > 0) { 6435fe6150dSHong Zhang ierr = PetscMalloc(lu->lwork,&lu->work);CHKERRQ(ierr); 6445fe6150dSHong Zhang } else if (lu->lwork != 0 && lu->lwork != -1) { 64577431f27SBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," Warning: lwork %D is not supported by SUPERLU. The default lwork=0 is used.\n",lu->lwork); 6468914a3f7SHong Zhang lu->lwork = 0; 6478914a3f7SHong Zhang } 648cffbb591SHong Zhang /* ilu options */ 6493cb270beSHong Zhang ierr = PetscOptionsReal("-mat_superlu_ilu_droptol","ILU_DropTol","None",lu->options.ILU_DropTol,&real_input,&flg);CHKERRQ(ierr); 6503cb270beSHong Zhang if (flg) lu->options.ILU_DropTol = (double) real_input; 6513cb270beSHong Zhang ierr = PetscOptionsReal("-mat_superlu_ilu_filltol","ILU_FillTol","None",lu->options.ILU_FillTol,&real_input,&flg);CHKERRQ(ierr); 6523cb270beSHong Zhang if (flg) lu->options.ILU_FillTol = (double) real_input; 6533cb270beSHong Zhang ierr = PetscOptionsReal("-mat_superlu_ilu_fillfactor","ILU_FillFactor","None",lu->options.ILU_FillFactor,&real_input,&flg);CHKERRQ(ierr); 6543cb270beSHong Zhang if (flg) lu->options.ILU_FillFactor = (double) real_input; 6550298fd71SBarry Smith ierr = PetscOptionsInt("-mat_superlu_ilu_droprull","ILU_DropRule","None",lu->options.ILU_DropRule,&lu->options.ILU_DropRule,NULL);CHKERRQ(ierr); 656cffbb591SHong Zhang ierr = PetscOptionsInt("-mat_superlu_ilu_norm","ILU_Norm","None",lu->options.ILU_Norm,&indx,&flg);CHKERRQ(ierr); 6572205254eSKarl Rupp if (flg) lu->options.ILU_Norm = (norm_t)indx; 658cffbb591SHong Zhang ierr = PetscOptionsInt("-mat_superlu_ilu_milu","ILU_MILU","None",lu->options.ILU_MILU,&indx,&flg);CHKERRQ(ierr); 6592205254eSKarl Rupp if (flg) lu->options.ILU_MILU = (milu_t)indx; 6609ce50919SHong Zhang PetscOptionsEnd(); 66138abfdaeSHong Zhang if (lu->options.Equil == YES) { 66238abfdaeSHong Zhang /* superlu overwrites input matrix and rhs when Equil is used, thus create A_dup to keep user's A unchanged */ 66338abfdaeSHong Zhang ierr = MatDuplicate_SeqAIJ(A,MAT_COPY_VALUES,&lu->A_dup);CHKERRQ(ierr); 66438abfdaeSHong Zhang } 6659ce50919SHong Zhang 6665d8b2efaSHong Zhang /* Allocate spaces (notice sizes are for the transpose) */ 667785e854fSJed Brown ierr = PetscMalloc1(m,&lu->etree);CHKERRQ(ierr); 668785e854fSJed Brown ierr = PetscMalloc1(n,&lu->perm_r);CHKERRQ(ierr); 669785e854fSJed Brown ierr = PetscMalloc1(m,&lu->perm_c);CHKERRQ(ierr); 670785e854fSJed Brown ierr = PetscMalloc1(n,&lu->R);CHKERRQ(ierr); 671785e854fSJed Brown ierr = PetscMalloc1(m,&lu->C);CHKERRQ(ierr); 6725d8b2efaSHong Zhang 6735d8b2efaSHong Zhang /* create rhs and solution x without allocate space for .Store */ 6745d8b2efaSHong Zhang #if defined(PETSC_USE_COMPLEX) 6753cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 676d387c056SBarry Smith PetscStackCall("SuperLU:cCreate_Dense_Matrix(",cCreate_Dense_Matrix(&lu->B, m, 1, NULL, m, SLU_DN, SLU_C, SLU_GE)); 677d387c056SBarry Smith PetscStackCall("SuperLU:cCreate_Dense_Matrix(",cCreate_Dense_Matrix(&lu->X, m, 1, NULL, m, SLU_DN, SLU_C, SLU_GE)); 6783cb270beSHong Zhang #else 679d387c056SBarry Smith PetscStackCall("SuperLU:zCreate_Dense_Matrix",zCreate_Dense_Matrix(&lu->B, m, 1, NULL, m, SLU_DN, SLU_Z, SLU_GE)); 680d387c056SBarry Smith PetscStackCall("SuperLU:zCreate_Dense_Matrix",zCreate_Dense_Matrix(&lu->X, m, 1, NULL, m, SLU_DN, SLU_Z, SLU_GE)); 6813cb270beSHong Zhang #endif 6823cb270beSHong Zhang #else 6833cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 684d387c056SBarry Smith PetscStackCall("SuperLU:sCreate_Dense_Matrix",sCreate_Dense_Matrix(&lu->B, m, 1, NULL, m, SLU_DN, SLU_S, SLU_GE)); 685d387c056SBarry Smith PetscStackCall("SuperLU:sCreate_Dense_Matrix",sCreate_Dense_Matrix(&lu->X, m, 1, NULL, m, SLU_DN, SLU_S, SLU_GE)); 6865d8b2efaSHong Zhang #else 687d387c056SBarry Smith PetscStackCall("SuperLU:dCreate_Dense_Matrix",dCreate_Dense_Matrix(&lu->B, m, 1, NULL, m, SLU_DN, SLU_D, SLU_GE)); 688d387c056SBarry Smith PetscStackCall("SuperLU:dCreate_Dense_Matrix",dCreate_Dense_Matrix(&lu->X, m, 1, NULL, m, SLU_DN, SLU_D, SLU_GE)); 6895d8b2efaSHong Zhang #endif 6903cb270beSHong Zhang #endif 6915d8b2efaSHong Zhang 692bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_seqaij_superlu);CHKERRQ(ierr); 693bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatSuperluSetILUDropTol_C",MatSuperluSetILUDropTol_SuperLU);CHKERRQ(ierr); 6945c9eb25fSBarry Smith B->spptr = lu; 695899d7b4fSKris Buschelman *F = B; 69614b396bbSKris Buschelman PetscFunctionReturn(0); 69714b396bbSKris Buschelman } 698d954db57SHong Zhang 699