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; 90*2205254eSKarl 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) { 1395a46d3faSBarry Smith Destroy_SuperNode_Matrix(&lu->L); 1405a46d3faSBarry Smith 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) 155b0043f70SBarry Smith 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 157b0043f70SBarry Smith 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) 161b0043f70SBarry Smith 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 163b0043f70SBarry Smith 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) 1723cb270beSHong Zhang 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, 1743cb270beSHong Zhang &lu->mem_usage, &lu->stat, &sinfo); 1753cb270beSHong Zhang #else 1765a46d3faSBarry Smith 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, 1785d8b2efaSHong Zhang &lu->mem_usage, &lu->stat, &sinfo); 1793cb270beSHong Zhang #endif 1803cb270beSHong Zhang #else 1813cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 1823cb270beSHong Zhang 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, 1843cb270beSHong Zhang &lu->mem_usage, &lu->stat, &sinfo); 1855a46d3faSBarry Smith #else 1865a46d3faSBarry Smith 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, 1885d8b2efaSHong Zhang &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) 1953cb270beSHong Zhang 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, 1973cb270beSHong Zhang &lu->mem_usage, &lu->stat, &sinfo); 1983cb270beSHong Zhang #else 199cffbb591SHong Zhang 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, 2015d8b2efaSHong Zhang &lu->mem_usage, &lu->stat, &sinfo); 2023cb270beSHong Zhang #endif 2033cb270beSHong Zhang #else 2043cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 2053cb270beSHong Zhang 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, 2073cb270beSHong Zhang &lu->mem_usage, &lu->stat, &sinfo); 208cffbb591SHong Zhang #else 209cffbb591SHong Zhang 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, 2115d8b2efaSHong Zhang &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) { 216*2205254eSKarl Rupp if (lu->options.PivotGrowth) { 2175a46d3faSBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," Recip. pivot growth = %e\n", lu->rpg); 218*2205254eSKarl Rupp } 219*2205254eSKarl Rupp if (lu->options.ConditionNumber) { 2205a46d3faSBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," Recip. condition number = %e\n", lu->rcond); 221*2205254eSKarl 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"); 2305d8b2efaSHong Zhang 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 */ 25675af56d4SHong Zhang Destroy_SuperMatrix_Store(&lu->A); 25775af56d4SHong Zhang Destroy_SuperMatrix_Store(&lu->B); 25875af56d4SHong Zhang Destroy_SuperMatrix_Store(&lu->X); 2595d8b2efaSHong Zhang StatFree(&lu->stat); 2600e742b69SHong Zhang if (lu->lwork >= 0) { 2610e742b69SHong Zhang Destroy_SuperNode_Matrix(&lu->L); 2620e742b69SHong Zhang 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 */ 277d954db57SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatFactorGetSolverPackage_C","",PETSC_NULL);CHKERRQ(ierr); 278790a96dfSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatSuperluSetILUDropTol_C","",PETSC_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; 311cae5a23dSHong Zhang PetscInt info,i,n=x->map->n; 312da7a1d00SHong Zhang PetscReal ferr,berr; 31314b396bbSKris Buschelman 31414b396bbSKris Buschelman PetscFunctionBegin; 315*2205254eSKarl Rupp if (lu->lwork == -1) PetscFunctionReturn(0); 316cae5a23dSHong Zhang 31775af56d4SHong Zhang lu->B.ncol = 1; /* Set the number of right-hand side */ 318cae5a23dSHong Zhang if (lu->options.Equil && !lu->rhs_dup) { 319cae5a23dSHong Zhang /* superlu overwrites b when Equil is used, thus create rhs_dup to keep user's b unchanged */ 320cae5a23dSHong Zhang ierr = PetscMalloc(n*sizeof(PetscScalar),&lu->rhs_dup);CHKERRQ(ierr); 321cae5a23dSHong Zhang } 322cae5a23dSHong Zhang if (lu->options.Equil) { 323cae5a23dSHong Zhang /* Copy b into rsh_dup */ 32475af56d4SHong Zhang ierr = VecGetArray(b,&barray);CHKERRQ(ierr); 325cae5a23dSHong Zhang ierr = PetscMemcpy(lu->rhs_dup,barray,n*sizeof(PetscScalar));CHKERRQ(ierr); 326cae5a23dSHong Zhang ierr = VecRestoreArray(b,&barray);CHKERRQ(ierr); 327cae5a23dSHong Zhang barray = lu->rhs_dup; 328cae5a23dSHong Zhang } else { 329cae5a23dSHong Zhang ierr = VecGetArray(b,&barray);CHKERRQ(ierr); 330cae5a23dSHong Zhang } 33175af56d4SHong Zhang ierr = VecGetArray(x,&xarray);CHKERRQ(ierr); 3325fe6150dSHong Zhang 3335fe6150dSHong Zhang #if defined(PETSC_USE_COMPLEX) 3343cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 3353cb270beSHong Zhang ((DNformat*)lu->B.Store)->nzval = (singlecomplex*)barray; 3363cb270beSHong Zhang ((DNformat*)lu->X.Store)->nzval = (singlecomplex*)xarray; 3373cb270beSHong Zhang #else 3385fe6150dSHong Zhang ((DNformat*)lu->B.Store)->nzval = (doublecomplex*)barray; 3395fe6150dSHong Zhang ((DNformat*)lu->X.Store)->nzval = (doublecomplex*)xarray; 3403cb270beSHong Zhang #endif 3415fe6150dSHong Zhang #else 3425fe6150dSHong Zhang ((DNformat*)lu->B.Store)->nzval = barray; 34375af56d4SHong Zhang ((DNformat*)lu->X.Store)->nzval = xarray; 3445fe6150dSHong Zhang #endif 34575af56d4SHong Zhang 34675af56d4SHong Zhang lu->options.Fact = FACTORED; /* Indicate the factored form of A is supplied. */ 347d5f3da31SBarry Smith if (A->factortype == MAT_FACTOR_LU) { 3488914a3f7SHong Zhang #if defined(PETSC_USE_COMPLEX) 3493cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 3503cb270beSHong Zhang cgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 3513cb270beSHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 3523cb270beSHong Zhang &lu->mem_usage, &lu->stat, &info); 3533cb270beSHong Zhang #else 3548914a3f7SHong Zhang zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 3558914a3f7SHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 3565d8b2efaSHong Zhang &lu->mem_usage, &lu->stat, &info); 3573cb270beSHong Zhang #endif 3583cb270beSHong Zhang #else 3593cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 3603cb270beSHong Zhang sgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 3613cb270beSHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 3623cb270beSHong Zhang &lu->mem_usage, &lu->stat, &info); 3638914a3f7SHong Zhang #else 36475af56d4SHong Zhang dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 36575af56d4SHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 3665d8b2efaSHong Zhang &lu->mem_usage, &lu->stat, &info); 3678914a3f7SHong Zhang #endif 3683cb270beSHong Zhang #endif 369d5f3da31SBarry Smith } else if (A->factortype == MAT_FACTOR_ILU) { 370cffbb591SHong Zhang #if defined(PETSC_USE_COMPLEX) 3713cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 3723cb270beSHong Zhang cgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 3733cb270beSHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, 3743cb270beSHong Zhang &lu->mem_usage, &lu->stat, &info); 3753cb270beSHong Zhang #else 376cffbb591SHong Zhang zgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 377cffbb591SHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, 3785d8b2efaSHong Zhang &lu->mem_usage, &lu->stat, &info); 3793cb270beSHong Zhang #endif 3803cb270beSHong Zhang #else 3813cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 3823cb270beSHong Zhang sgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 3833cb270beSHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, 3843cb270beSHong Zhang &lu->mem_usage, &lu->stat, &info); 385cffbb591SHong Zhang #else 386cffbb591SHong Zhang dgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 387cffbb591SHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, 3885d8b2efaSHong Zhang &lu->mem_usage, &lu->stat, &info); 389cffbb591SHong Zhang #endif 3903cb270beSHong Zhang #endif 391f23aa3ddSBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported"); 392cae5a23dSHong Zhang if (!lu->options.Equil) { 39375af56d4SHong Zhang ierr = VecRestoreArray(b,&barray);CHKERRQ(ierr); 394cae5a23dSHong Zhang } 39575af56d4SHong Zhang ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr); 39675af56d4SHong Zhang 397958c9bccSBarry Smith if (!info || info == lu->A.ncol+1) { 39875af56d4SHong Zhang if (lu->options.IterRefine) { 3998914a3f7SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF,"Iterative Refinement:\n"); 4008914a3f7SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," %8s%8s%16s%16s\n", "rhs", "Steps", "FERR", "BERR"); 401*2205254eSKarl Rupp for (i = 0; i < 1; ++i) { 4025d8b2efaSHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," %8d%8d%16e%16e\n", i+1, lu->stat.RefineSteps, ferr, berr); 40375af56d4SHong Zhang } 404*2205254eSKarl Rupp } 4058914a3f7SHong Zhang } else if (info > 0) { 4068914a3f7SHong Zhang if (lu->lwork == -1) { 40777431f27SBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," ** Estimated memory: %D bytes\n", info - lu->A.ncol); 4088914a3f7SHong Zhang } else { 40977431f27SBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," Warning: gssvx() returns info %D\n",info); 4108914a3f7SHong Zhang } 411f23aa3ddSBarry 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); 41214b396bbSKris Buschelman 4138914a3f7SHong Zhang if (lu->options.PrintStat) { 4148914a3f7SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF,"MatSolve__SuperLU():\n"); 4155d8b2efaSHong Zhang StatPrint(&lu->stat); 4168914a3f7SHong Zhang } 41775af56d4SHong Zhang PetscFunctionReturn(0); 41875af56d4SHong Zhang } 41914b396bbSKris Buschelman 42014b396bbSKris Buschelman #undef __FUNCT__ 421c29e884eSHong Zhang #define __FUNCT__ "MatSolve_SuperLU" 422c29e884eSHong Zhang PetscErrorCode MatSolve_SuperLU(Mat A,Vec b,Vec x) 423c29e884eSHong Zhang { 424c29e884eSHong Zhang Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr; 425c29e884eSHong Zhang PetscErrorCode ierr; 426c29e884eSHong Zhang 427c29e884eSHong Zhang PetscFunctionBegin; 428c29e884eSHong Zhang lu->options.Trans = TRANS; 429*2205254eSKarl Rupp 430c29e884eSHong Zhang ierr = MatSolve_SuperLU_Private(A,b,x);CHKERRQ(ierr); 431c29e884eSHong Zhang PetscFunctionReturn(0); 432c29e884eSHong Zhang } 433c29e884eSHong Zhang 434c29e884eSHong Zhang #undef __FUNCT__ 435c7c1fe80SHong Zhang #define __FUNCT__ "MatSolveTranspose_SuperLU" 436c7c1fe80SHong Zhang PetscErrorCode MatSolveTranspose_SuperLU(Mat A,Vec b,Vec x) 437c7c1fe80SHong Zhang { 438c7c1fe80SHong Zhang Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr; 439c7c1fe80SHong Zhang PetscErrorCode ierr; 440c7c1fe80SHong Zhang 441c7c1fe80SHong Zhang PetscFunctionBegin; 442c7c1fe80SHong Zhang lu->options.Trans = NOTRANS; 443*2205254eSKarl Rupp 444c29e884eSHong Zhang ierr = MatSolve_SuperLU_Private(A,b,x);CHKERRQ(ierr); 445c7c1fe80SHong Zhang PetscFunctionReturn(0); 446c7c1fe80SHong Zhang } 447c7c1fe80SHong Zhang 448e0b74bf9SHong Zhang #undef __FUNCT__ 449e0b74bf9SHong Zhang #define __FUNCT__ "MatMatSolve_SuperLU" 450e0b74bf9SHong Zhang PetscErrorCode MatMatSolve_SuperLU(Mat A,Mat B,Mat X) 451e0b74bf9SHong Zhang { 452e0b74bf9SHong Zhang Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr; 453cd723cd1SBarry Smith PetscBool flg; 454cd723cd1SBarry Smith PetscErrorCode ierr; 455e0b74bf9SHong Zhang 456e0b74bf9SHong Zhang PetscFunctionBegin; 457251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompareAny((PetscObject)B,&flg,MATSEQDENSE,MATMPIDENSE,PETSC_NULL);CHKERRQ(ierr); 458cd723cd1SBarry Smith if (!flg) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix"); 459251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,PETSC_NULL);CHKERRQ(ierr); 460*2205254eSKarl Rupp if (!flg) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix"); 461*2205254eSKarl Rupp lu->options.Trans = TRANS; 462e0b74bf9SHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatMatSolve_SuperLU() is not implemented yet"); 463e0b74bf9SHong Zhang PetscFunctionReturn(0); 464e0b74bf9SHong Zhang } 465e0b74bf9SHong Zhang 46614b396bbSKris Buschelman /* 46714b396bbSKris Buschelman Note the r permutation is ignored 46814b396bbSKris Buschelman */ 46914b396bbSKris Buschelman #undef __FUNCT__ 470f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_SuperLU" 4710481f469SBarry Smith PetscErrorCode MatLUFactorSymbolic_SuperLU(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info) 47214b396bbSKris Buschelman { 4731d5ca7f3SHong Zhang Mat_SuperLU *lu = (Mat_SuperLU*)(F->spptr); 474b24902e0SBarry Smith 475b24902e0SBarry Smith PetscFunctionBegin; 476b24902e0SBarry Smith lu->flg = DIFFERENT_NONZERO_PATTERN; 477b24902e0SBarry Smith lu->CleanUpSuperLU = PETSC_TRUE; 4781d5ca7f3SHong Zhang F->ops->lufactornumeric = MatLUFactorNumeric_SuperLU; 479b24902e0SBarry Smith PetscFunctionReturn(0); 480b24902e0SBarry Smith } 481b24902e0SBarry Smith 48235bd34faSBarry Smith EXTERN_C_BEGIN 48335bd34faSBarry Smith #undef __FUNCT__ 4845ccb76cbSHong Zhang #define __FUNCT__ "MatSuperluSetILUDropTol_SuperLU" 4855ccb76cbSHong Zhang PetscErrorCode MatSuperluSetILUDropTol_SuperLU(Mat F,PetscReal dtol) 4865ccb76cbSHong Zhang { 4875ccb76cbSHong Zhang Mat_SuperLU *lu= (Mat_SuperLU*)F->spptr; 4885ccb76cbSHong Zhang 4895ccb76cbSHong Zhang PetscFunctionBegin; 4905ccb76cbSHong Zhang lu->options.ILU_DropTol = dtol; 4915ccb76cbSHong Zhang PetscFunctionReturn(0); 4925ccb76cbSHong Zhang } 4935ccb76cbSHong Zhang EXTERN_C_END 4945ccb76cbSHong Zhang 4955ccb76cbSHong Zhang #undef __FUNCT__ 4965ccb76cbSHong Zhang #define __FUNCT__ "MatSuperluSetILUDropTol" 4975ccb76cbSHong Zhang /*@ 4985ccb76cbSHong Zhang MatSuperluSetILUDropTol - Set SuperLU ILU drop tolerance 4995ccb76cbSHong Zhang Logically Collective on Mat 5005ccb76cbSHong Zhang 5015ccb76cbSHong Zhang Input Parameters: 5025ccb76cbSHong Zhang + F - the factored matrix obtained by calling MatGetFactor() from PETSc-SuperLU interface 5035ccb76cbSHong Zhang - dtol - drop tolerance 5045ccb76cbSHong Zhang 5055ccb76cbSHong Zhang Options Database: 5065ccb76cbSHong Zhang . -mat_superlu_ilu_droptol <dtol> 5075ccb76cbSHong Zhang 5085ccb76cbSHong Zhang Level: beginner 5095ccb76cbSHong Zhang 5105ccb76cbSHong Zhang References: SuperLU Users' Guide 5115ccb76cbSHong Zhang 5125ccb76cbSHong Zhang .seealso: MatGetFactor() 5135ccb76cbSHong Zhang @*/ 5145ccb76cbSHong Zhang PetscErrorCode MatSuperluSetILUDropTol(Mat F,PetscReal dtol) 5155ccb76cbSHong Zhang { 5165ccb76cbSHong Zhang PetscErrorCode ierr; 5175ccb76cbSHong Zhang 5185ccb76cbSHong Zhang PetscFunctionBegin; 5195ccb76cbSHong Zhang PetscValidHeaderSpecific(F,MAT_CLASSID,1); 5205ccb76cbSHong Zhang PetscValidLogicalCollectiveInt(F,dtol,2); 5215ccb76cbSHong Zhang ierr = PetscTryMethod(F,"MatSuperluSetILUDropTol_C",(Mat,PetscReal),(F,dtol));CHKERRQ(ierr); 5225ccb76cbSHong Zhang PetscFunctionReturn(0); 5235ccb76cbSHong Zhang } 5245ccb76cbSHong Zhang 5255ccb76cbSHong Zhang EXTERN_C_BEGIN 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 EXTERN_C_END 53535bd34faSBarry Smith 536b24902e0SBarry Smith 537b24902e0SBarry Smith /*MC 538ba992d64SSatish Balay MATSOLVERSUPERLU = "superlu" - A solver package providing solvers LU and ILU for sequential matrices 539b24902e0SBarry Smith via the external package SuperLU. 540b24902e0SBarry Smith 541e2e64c6bSBarry Smith Use ./configure --download-superlu to have PETSc installed with SuperLU 542b24902e0SBarry Smith 543b24902e0SBarry Smith Options Database Keys: 544e08999f5SMatthew G Knepley + -mat_superlu_equil <FALSE> - Equil (None) 545e08999f5SMatthew G Knepley . -mat_superlu_colperm <COLAMD> - (choose one of) NATURAL MMD_ATA MMD_AT_PLUS_A COLAMD 546e08999f5SMatthew G Knepley . -mat_superlu_iterrefine <NOREFINE> - (choose one of) NOREFINE SINGLE DOUBLE EXTRA 547e08999f5SMatthew G Knepley . -mat_superlu_symmetricmode: <FALSE> - SymmetricMode (None) 548e08999f5SMatthew G Knepley . -mat_superlu_diagpivotthresh <1> - DiagPivotThresh (None) 549e08999f5SMatthew G Knepley . -mat_superlu_pivotgrowth <FALSE> - PivotGrowth (None) 550e08999f5SMatthew G Knepley . -mat_superlu_conditionnumber <FALSE> - ConditionNumber (None) 551e08999f5SMatthew G Knepley . -mat_superlu_rowperm <NOROWPERM> - (choose one of) NOROWPERM LargeDiag 552e08999f5SMatthew G Knepley . -mat_superlu_replacetinypivot <FALSE> - ReplaceTinyPivot (None) 553e08999f5SMatthew G Knepley . -mat_superlu_printstat <FALSE> - PrintStat (None) 554e08999f5SMatthew G Knepley . -mat_superlu_lwork <0> - size of work array in bytes used by factorization (None) 555e08999f5SMatthew G Knepley . -mat_superlu_ilu_droptol <0> - ILU_DropTol (None) 556e08999f5SMatthew G Knepley . -mat_superlu_ilu_filltol <0> - ILU_FillTol (None) 557e08999f5SMatthew G Knepley . -mat_superlu_ilu_fillfactor <0> - ILU_FillFactor (None) 558e08999f5SMatthew G Knepley . -mat_superlu_ilu_droprull <0> - ILU_DropRule (None) 559e08999f5SMatthew G Knepley . -mat_superlu_ilu_norm <0> - ILU_Norm (None) 560e08999f5SMatthew G Knepley - -mat_superlu_ilu_milu <0> - ILU_MILU (None) 561b24902e0SBarry Smith 5622692d6eeSBarry Smith Notes: Do not confuse this with MATSOLVERSUPERLU_DIST which is for parallel sparse solves 5635c9eb25fSBarry Smith 564b24902e0SBarry Smith Level: beginner 565b24902e0SBarry Smith 566d45987f3SHong Zhang .seealso: PCLU, PCILU, MATSOLVERSUPERLU_DIST, MATSOLVERMUMPS, PCFactorSetMatSolverPackage(), MatSolverPackage 567b24902e0SBarry Smith M*/ 568b24902e0SBarry Smith 5695c9eb25fSBarry Smith EXTERN_C_BEGIN 570b24902e0SBarry Smith #undef __FUNCT__ 571b24902e0SBarry Smith #define __FUNCT__ "MatGetFactor_seqaij_superlu" 5725c9eb25fSBarry Smith PetscErrorCode MatGetFactor_seqaij_superlu(Mat A,MatFactorType ftype,Mat *F) 573b24902e0SBarry Smith { 57414b396bbSKris Buschelman Mat B; 575f0c56d0fSKris Buschelman Mat_SuperLU *lu; 5766849ba73SBarry Smith PetscErrorCode ierr; 5775d8b2efaSHong Zhang PetscInt indx,m=A->rmap->n,n=A->cmap->n; 578ace3abfcSBarry Smith PetscBool flg; 5793cb270beSHong Zhang PetscReal real_input; 5805ca28756SHong Zhang const char *colperm[] ={"NATURAL","MMD_ATA","MMD_AT_PLUS_A","COLAMD"}; /* MY_PERMC - not supported by the petsc interface yet */ 5815ca28756SHong Zhang const char *iterrefine[]={"NOREFINE", "SINGLE", "DOUBLE", "EXTRA"}; 5825ca28756SHong Zhang const char *rowperm[] ={"NOROWPERM", "LargeDiag"}; /* MY_PERMC - not supported by the petsc interface yet */ 58314b396bbSKris Buschelman 58414b396bbSKris Buschelman PetscFunctionBegin; 5857adad957SLisandro Dalcin ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr); 586d0f46423SBarry Smith ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 5877adad957SLisandro Dalcin ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 588899d7b4fSKris Buschelman ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr); 589f0c56d0fSKris Buschelman 590cffbb591SHong Zhang if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU) { 591b24902e0SBarry Smith B->ops->lufactorsymbolic = MatLUFactorSymbolic_SuperLU; 592cffbb591SHong Zhang B->ops->ilufactorsymbolic = MatLUFactorSymbolic_SuperLU; 593b3a44c85SBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported"); 594cffbb591SHong Zhang 595b24902e0SBarry Smith B->ops->destroy = MatDestroy_SuperLU; 5963519fcdcSHong Zhang B->ops->view = MatView_SuperLU; 597d5f3da31SBarry Smith B->factortype = ftype; 59894b7f48cSBarry Smith B->assembled = PETSC_TRUE; /* required by -ksp_view */ 5995c9eb25fSBarry Smith B->preallocated = PETSC_TRUE; 60014b396bbSKris Buschelman 601b24902e0SBarry Smith ierr = PetscNewLog(B,Mat_SuperLU,&lu);CHKERRQ(ierr); 602cae5a23dSHong Zhang 603cffbb591SHong Zhang if (ftype == MAT_FACTOR_LU) { 6049ce50919SHong Zhang set_default_options(&lu->options); 6053d6581fbSHong Zhang /* Comments from SuperLU_4.0/SRC/dgssvx.c: 6063d6581fbSHong Zhang "Whether or not the system will be equilibrated depends on the 6073d6581fbSHong Zhang scaling of the matrix A, but if equilibration is used, A is 6083d6581fbSHong Zhang overwritten by diag(R)*A*diag(C) and B by diag(R)*B 6093d6581fbSHong Zhang (if options->Trans=NOTRANS) or diag(C)*B (if options->Trans = TRANS or CONJ)." 6103d6581fbSHong Zhang We set 'options.Equil = NO' as default because additional space is needed for it. 6113d6581fbSHong Zhang */ 6123d6581fbSHong Zhang lu->options.Equil = NO; 613cffbb591SHong Zhang } else if (ftype == MAT_FACTOR_ILU) { 6140c74a584SJed Brown /* Set the default input options of ilu: */ 615cffbb591SHong Zhang ilu_set_default_options(&lu->options); 616cffbb591SHong Zhang } 6179ce50919SHong Zhang lu->options.PrintStat = NO; 6181a2e2f44SHong Zhang 6195d8b2efaSHong Zhang /* Initialize the statistics variables. */ 6205d8b2efaSHong Zhang StatInit(&lu->stat); 6218914a3f7SHong Zhang lu->lwork = 0; /* allocate space internally by system malloc */ 6229ce50919SHong Zhang 6237adad957SLisandro Dalcin ierr = PetscOptionsBegin(((PetscObject)A)->comm,((PetscObject)A)->prefix,"SuperLU Options","Mat");CHKERRQ(ierr); 624acfcf0e5SJed Brown ierr = PetscOptionsBool("-mat_superlu_equil","Equil","None",(PetscBool)lu->options.Equil,(PetscBool*)&lu->options.Equil,0);CHKERRQ(ierr); 6258914a3f7SHong Zhang ierr = PetscOptionsEList("-mat_superlu_colperm","ColPerm","None",colperm,4,colperm[3],&indx,&flg);CHKERRQ(ierr); 626*2205254eSKarl Rupp if (flg) lu->options.ColPerm = (colperm_t)indx; 6278914a3f7SHong Zhang ierr = PetscOptionsEList("-mat_superlu_iterrefine","IterRefine","None",iterrefine,4,iterrefine[0],&indx,&flg);CHKERRQ(ierr); 628*2205254eSKarl Rupp if (flg) lu->options.IterRefine = (IterRefine_t)indx; 629acfcf0e5SJed Brown ierr = PetscOptionsBool("-mat_superlu_symmetricmode","SymmetricMode","None",(PetscBool)lu->options.SymmetricMode,&flg,0);CHKERRQ(ierr); 6309ce50919SHong Zhang if (flg) lu->options.SymmetricMode = YES; 6313cb270beSHong Zhang ierr = PetscOptionsReal("-mat_superlu_diagpivotthresh","DiagPivotThresh","None",lu->options.DiagPivotThresh,&real_input,&flg);CHKERRQ(ierr); 6323cb270beSHong Zhang if (flg) lu->options.DiagPivotThresh = (double) real_input; 633acfcf0e5SJed Brown ierr = PetscOptionsBool("-mat_superlu_pivotgrowth","PivotGrowth","None",(PetscBool)lu->options.PivotGrowth,&flg,0);CHKERRQ(ierr); 6349ce50919SHong Zhang if (flg) lu->options.PivotGrowth = YES; 635acfcf0e5SJed Brown ierr = PetscOptionsBool("-mat_superlu_conditionnumber","ConditionNumber","None",(PetscBool)lu->options.ConditionNumber,&flg,0);CHKERRQ(ierr); 6369ce50919SHong Zhang if (flg) lu->options.ConditionNumber = YES; 637d7ebd59bSHong Zhang ierr = PetscOptionsEList("-mat_superlu_rowperm","rowperm","None",rowperm,2,rowperm[lu->options.RowPerm],&indx,&flg);CHKERRQ(ierr); 638*2205254eSKarl Rupp if (flg) lu->options.RowPerm = (rowperm_t)indx; 639acfcf0e5SJed Brown ierr = PetscOptionsBool("-mat_superlu_replacetinypivot","ReplaceTinyPivot","None",(PetscBool)lu->options.ReplaceTinyPivot,&flg,0);CHKERRQ(ierr); 6409ce50919SHong Zhang if (flg) lu->options.ReplaceTinyPivot = YES; 641acfcf0e5SJed Brown ierr = PetscOptionsBool("-mat_superlu_printstat","PrintStat","None",(PetscBool)lu->options.PrintStat,&flg,0);CHKERRQ(ierr); 6429ce50919SHong Zhang if (flg) lu->options.PrintStat = YES; 6438914a3f7SHong Zhang ierr = PetscOptionsInt("-mat_superlu_lwork","size of work array in bytes used by factorization","None",lu->lwork,&lu->lwork,PETSC_NULL);CHKERRQ(ierr); 6445fe6150dSHong Zhang if (lu->lwork > 0) { 6455fe6150dSHong Zhang ierr = PetscMalloc(lu->lwork,&lu->work);CHKERRQ(ierr); 6465fe6150dSHong Zhang } else if (lu->lwork != 0 && lu->lwork != -1) { 64777431f27SBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," Warning: lwork %D is not supported by SUPERLU. The default lwork=0 is used.\n",lu->lwork); 6488914a3f7SHong Zhang lu->lwork = 0; 6498914a3f7SHong Zhang } 650cffbb591SHong Zhang /* ilu options */ 6513cb270beSHong Zhang ierr = PetscOptionsReal("-mat_superlu_ilu_droptol","ILU_DropTol","None",lu->options.ILU_DropTol,&real_input,&flg);CHKERRQ(ierr); 6523cb270beSHong Zhang if (flg) lu->options.ILU_DropTol = (double) real_input; 6533cb270beSHong Zhang ierr = PetscOptionsReal("-mat_superlu_ilu_filltol","ILU_FillTol","None",lu->options.ILU_FillTol,&real_input,&flg);CHKERRQ(ierr); 6543cb270beSHong Zhang if (flg) lu->options.ILU_FillTol = (double) real_input; 6553cb270beSHong Zhang ierr = PetscOptionsReal("-mat_superlu_ilu_fillfactor","ILU_FillFactor","None",lu->options.ILU_FillFactor,&real_input,&flg);CHKERRQ(ierr); 6563cb270beSHong Zhang if (flg) lu->options.ILU_FillFactor = (double) real_input; 657cffbb591SHong Zhang ierr = PetscOptionsInt("-mat_superlu_ilu_droprull","ILU_DropRule","None",lu->options.ILU_DropRule,&lu->options.ILU_DropRule,PETSC_NULL);CHKERRQ(ierr); 658cffbb591SHong Zhang ierr = PetscOptionsInt("-mat_superlu_ilu_norm","ILU_Norm","None",lu->options.ILU_Norm,&indx,&flg);CHKERRQ(ierr); 659*2205254eSKarl Rupp if (flg) lu->options.ILU_Norm = (norm_t)indx; 660cffbb591SHong Zhang ierr = PetscOptionsInt("-mat_superlu_ilu_milu","ILU_MILU","None",lu->options.ILU_MILU,&indx,&flg);CHKERRQ(ierr); 661*2205254eSKarl Rupp if (flg) lu->options.ILU_MILU = (milu_t)indx; 6629ce50919SHong Zhang PetscOptionsEnd(); 66338abfdaeSHong Zhang if (lu->options.Equil == YES) { 66438abfdaeSHong Zhang /* superlu overwrites input matrix and rhs when Equil is used, thus create A_dup to keep user's A unchanged */ 66538abfdaeSHong Zhang ierr = MatDuplicate_SeqAIJ(A,MAT_COPY_VALUES,&lu->A_dup);CHKERRQ(ierr); 66638abfdaeSHong Zhang } 6679ce50919SHong Zhang 6685d8b2efaSHong Zhang /* Allocate spaces (notice sizes are for the transpose) */ 6695d8b2efaSHong Zhang ierr = PetscMalloc(m*sizeof(PetscInt),&lu->etree);CHKERRQ(ierr); 6705d8b2efaSHong Zhang ierr = PetscMalloc(n*sizeof(PetscInt),&lu->perm_r);CHKERRQ(ierr); 6715d8b2efaSHong Zhang ierr = PetscMalloc(m*sizeof(PetscInt),&lu->perm_c);CHKERRQ(ierr); 6725d8b2efaSHong Zhang ierr = PetscMalloc(n*sizeof(PetscScalar),&lu->R);CHKERRQ(ierr); 6735d8b2efaSHong Zhang ierr = PetscMalloc(m*sizeof(PetscScalar),&lu->C);CHKERRQ(ierr); 6745d8b2efaSHong Zhang 6755d8b2efaSHong Zhang /* create rhs and solution x without allocate space for .Store */ 6765d8b2efaSHong Zhang #if defined(PETSC_USE_COMPLEX) 6773cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 6783cb270beSHong Zhang cCreate_Dense_Matrix(&lu->B, m, 1, PETSC_NULL, m, SLU_DN, SLU_C, SLU_GE); 6793cb270beSHong Zhang cCreate_Dense_Matrix(&lu->X, m, 1, PETSC_NULL, m, SLU_DN, SLU_C, SLU_GE); 6803cb270beSHong Zhang #else 6815d8b2efaSHong Zhang zCreate_Dense_Matrix(&lu->B, m, 1, PETSC_NULL, m, SLU_DN, SLU_Z, SLU_GE); 6825d8b2efaSHong Zhang zCreate_Dense_Matrix(&lu->X, m, 1, PETSC_NULL, m, SLU_DN, SLU_Z, SLU_GE); 6833cb270beSHong Zhang #endif 6843cb270beSHong Zhang #else 6853cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 6863cb270beSHong Zhang sCreate_Dense_Matrix(&lu->B, m, 1, PETSC_NULL, m, SLU_DN, SLU_S, SLU_GE); 6873cb270beSHong Zhang sCreate_Dense_Matrix(&lu->X, m, 1, PETSC_NULL, m, SLU_DN, SLU_S, SLU_GE); 6885d8b2efaSHong Zhang #else 6895d8b2efaSHong Zhang dCreate_Dense_Matrix(&lu->B, m, 1, PETSC_NULL, m, SLU_DN, SLU_D, SLU_GE); 6905d8b2efaSHong Zhang dCreate_Dense_Matrix(&lu->X, m, 1, PETSC_NULL, m, SLU_DN, SLU_D, SLU_GE); 6915d8b2efaSHong Zhang #endif 6923cb270beSHong Zhang #endif 6935d8b2efaSHong Zhang 694519f805aSKarl Rupp #if defined(SUPERLU2) 6955c9eb25fSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatCreateNull","MatCreateNull_SuperLU",(void (*)(void))MatCreateNull_SuperLU);CHKERRQ(ierr); 69675af56d4SHong Zhang #endif 69735bd34faSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_seqaij_superlu",MatFactorGetSolverPackage_seqaij_superlu);CHKERRQ(ierr); 6985ccb76cbSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSuperluSetILUDropTol_C","MatSuperluSetILUDropTol_SuperLU",MatSuperluSetILUDropTol_SuperLU);CHKERRQ(ierr); 6995c9eb25fSBarry Smith B->spptr = lu; 700899d7b4fSKris Buschelman *F = B; 70114b396bbSKris Buschelman PetscFunctionReturn(0); 70214b396bbSKris Buschelman } 7035c9eb25fSBarry Smith EXTERN_C_END 704d954db57SHong Zhang 705