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) 24c6db04a5SJed Brown #include <slu_zdefs.h> 2514b396bbSKris Buschelman #else 26c6db04a5SJed Brown #include <slu_ddefs.h> 2714b396bbSKris Buschelman #endif 28c6db04a5SJed Brown #include <slu_util.h> 2914b396bbSKris Buschelman EXTERN_C_END 3014b396bbSKris Buschelman 315a46d3faSBarry Smith /* 325a46d3faSBarry Smith This is the data we are "ADDING" to the SeqAIJ matrix type to get the SuperLU matrix type 335a46d3faSBarry Smith */ 3414b396bbSKris Buschelman typedef struct { 3575af56d4SHong Zhang SuperMatrix A,L,U,B,X; 3675af56d4SHong Zhang superlu_options_t options; 37da7a1d00SHong Zhang PetscInt *perm_c; /* column permutation vector */ 38da7a1d00SHong Zhang PetscInt *perm_r; /* row permutations from partial pivoting */ 39da7a1d00SHong Zhang PetscInt *etree; 40da7a1d00SHong Zhang PetscReal *R, *C; 4175af56d4SHong Zhang char equed[1]; 42da7a1d00SHong Zhang PetscInt lwork; 4375af56d4SHong Zhang void *work; 44da7a1d00SHong Zhang PetscReal rpg, rcond; 4575af56d4SHong Zhang mem_usage_t mem_usage; 4675af56d4SHong Zhang MatStructure flg; 475d8b2efaSHong Zhang SuperLUStat_t stat; 48cae5a23dSHong Zhang Mat A_dup; 49cae5a23dSHong Zhang PetscScalar *rhs_dup; 5014b396bbSKris Buschelman 5114b396bbSKris Buschelman /* Flag to clean up (non-global) SuperLU objects during Destroy */ 52ace3abfcSBarry Smith PetscBool CleanUpSuperLU; 53f0c56d0fSKris Buschelman } Mat_SuperLU; 5414b396bbSKris Buschelman 55e0e586b9SHong Zhang extern PetscErrorCode MatFactorInfo_SuperLU(Mat,PetscViewer); 560481f469SBarry Smith extern PetscErrorCode MatLUFactorNumeric_SuperLU(Mat,Mat,const MatFactorInfo *); 57e0e586b9SHong Zhang extern PetscErrorCode MatDestroy_SuperLU(Mat); 58e0e586b9SHong Zhang extern PetscErrorCode MatView_SuperLU(Mat,PetscViewer); 59e0e586b9SHong Zhang extern PetscErrorCode MatAssemblyEnd_SuperLU(Mat,MatAssemblyType); 60e0e586b9SHong Zhang extern PetscErrorCode MatSolve_SuperLU(Mat,Vec,Vec); 61*e0b74bf9SHong Zhang extern PetscErrorCode MatMatSolve_SuperLU(Mat,Mat,Mat); 62e0e586b9SHong Zhang extern PetscErrorCode MatSolveTranspose_SuperLU(Mat,Vec,Vec); 630481f469SBarry Smith extern PetscErrorCode MatLUFactorSymbolic_SuperLU(Mat,Mat,IS,IS,const MatFactorInfo*); 64e0e586b9SHong Zhang extern PetscErrorCode MatDuplicate_SuperLU(Mat, MatDuplicateOption, Mat *); 65e0e586b9SHong Zhang 665a46d3faSBarry Smith /* 675a46d3faSBarry Smith Utility function 685a46d3faSBarry Smith */ 695a46d3faSBarry Smith #undef __FUNCT__ 705a46d3faSBarry Smith #define __FUNCT__ "MatFactorInfo_SuperLU" 715a46d3faSBarry Smith PetscErrorCode MatFactorInfo_SuperLU(Mat A,PetscViewer viewer) 725a46d3faSBarry Smith { 735a46d3faSBarry Smith Mat_SuperLU *lu= (Mat_SuperLU*)A->spptr; 745a46d3faSBarry Smith PetscErrorCode ierr; 755a46d3faSBarry Smith superlu_options_t options; 765a46d3faSBarry Smith 775a46d3faSBarry Smith PetscFunctionBegin; 785a46d3faSBarry Smith /* check if matrix is superlu_dist type */ 795a46d3faSBarry Smith if (A->ops->solve != MatSolve_SuperLU) PetscFunctionReturn(0); 805a46d3faSBarry Smith 815a46d3faSBarry Smith options = lu->options; 825a46d3faSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"SuperLU run parameters:\n");CHKERRQ(ierr); 835a46d3faSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Equil: %s\n",(options.Equil != NO) ? "YES": "NO");CHKERRQ(ierr); 845a46d3faSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," ColPerm: %D\n",options.ColPerm);CHKERRQ(ierr); 855a46d3faSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," IterRefine: %D\n",options.IterRefine);CHKERRQ(ierr); 865a46d3faSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," SymmetricMode: %s\n",(options.SymmetricMode != NO) ? "YES": "NO");CHKERRQ(ierr); 875a46d3faSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," DiagPivotThresh: %g\n",options.DiagPivotThresh);CHKERRQ(ierr); 885a46d3faSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," PivotGrowth: %s\n",(options.PivotGrowth != NO) ? "YES": "NO");CHKERRQ(ierr); 895a46d3faSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," ConditionNumber: %s\n",(options.ConditionNumber != NO) ? "YES": "NO");CHKERRQ(ierr); 905a46d3faSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," RowPerm: %D\n",options.RowPerm);CHKERRQ(ierr); 915a46d3faSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," ReplaceTinyPivot: %s\n",(options.ReplaceTinyPivot != NO) ? "YES": "NO");CHKERRQ(ierr); 925a46d3faSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," PrintStat: %s\n",(options.PrintStat != NO) ? "YES": "NO");CHKERRQ(ierr); 935a46d3faSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," lwork: %D\n",lu->lwork);CHKERRQ(ierr); 94d5f3da31SBarry Smith if (A->factortype == MAT_FACTOR_ILU){ 95cffbb591SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ILU_DropTol: %g\n",options.ILU_DropTol);CHKERRQ(ierr); 96cffbb591SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ILU_FillTol: %g\n",options.ILU_FillTol);CHKERRQ(ierr); 97cffbb591SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ILU_FillFactor: %g\n",options.ILU_FillFactor);CHKERRQ(ierr); 98cffbb591SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ILU_DropRule: %D\n",options.ILU_DropRule);CHKERRQ(ierr); 99cffbb591SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ILU_Norm: %D\n",options.ILU_Norm);CHKERRQ(ierr); 100cffbb591SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ILU_MILU: %D\n",options.ILU_MILU);CHKERRQ(ierr); 101cffbb591SHong Zhang } 1025a46d3faSBarry Smith PetscFunctionReturn(0); 1035a46d3faSBarry Smith } 1045a46d3faSBarry Smith 1055a46d3faSBarry Smith /* 1065a46d3faSBarry Smith These are the methods provided to REPLACE the corresponding methods of the 1075a46d3faSBarry Smith SeqAIJ matrix class. Other methods of SeqAIJ are not replaced 1085a46d3faSBarry Smith */ 1095a46d3faSBarry Smith #undef __FUNCT__ 1105a46d3faSBarry Smith #define __FUNCT__ "MatLUFactorNumeric_SuperLU" 1110481f469SBarry Smith PetscErrorCode MatLUFactorNumeric_SuperLU(Mat F,Mat A,const MatFactorInfo *info) 1125a46d3faSBarry Smith { 1131d5ca7f3SHong Zhang Mat_SuperLU *lu = (Mat_SuperLU*)F->spptr; 114cae5a23dSHong Zhang Mat_SeqAIJ *aa; 1155a46d3faSBarry Smith PetscErrorCode ierr; 1165a46d3faSBarry Smith PetscInt sinfo; 1175a46d3faSBarry Smith PetscReal ferr, berr; 1185a46d3faSBarry Smith NCformat *Ustore; 1195a46d3faSBarry Smith SCformat *Lstore; 1205a46d3faSBarry Smith 1215a46d3faSBarry Smith PetscFunctionBegin; 1225a46d3faSBarry Smith if (lu->flg == SAME_NONZERO_PATTERN){ /* successing numerical factorization */ 1235a46d3faSBarry Smith lu->options.Fact = SamePattern; 1245a46d3faSBarry Smith /* Ref: ~SuperLU_3.0/EXAMPLE/dlinsolx2.c */ 1255a46d3faSBarry Smith Destroy_SuperMatrix_Store(&lu->A); 126cae5a23dSHong Zhang if (lu->options.Equil){ 127cae5a23dSHong Zhang ierr = MatCopy_SeqAIJ(A,lu->A_dup,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 128cae5a23dSHong Zhang } 1295a46d3faSBarry Smith if ( lu->lwork >= 0 ) { 1305a46d3faSBarry Smith Destroy_SuperNode_Matrix(&lu->L); 1315a46d3faSBarry Smith Destroy_CompCol_Matrix(&lu->U); 1325a46d3faSBarry Smith lu->options.Fact = SamePattern; 1335a46d3faSBarry Smith } 1345a46d3faSBarry Smith } 1355a46d3faSBarry Smith 1365a46d3faSBarry Smith /* Create the SuperMatrix for lu->A=A^T: 1375a46d3faSBarry Smith Since SuperLU likes column-oriented matrices,we pass it the transpose, 1385a46d3faSBarry Smith and then solve A^T X = B in MatSolve(). */ 139cae5a23dSHong Zhang if (lu->options.Equil){ 140cae5a23dSHong Zhang aa = (Mat_SeqAIJ*)(lu->A_dup)->data; 141cae5a23dSHong Zhang } else { 142cae5a23dSHong Zhang aa = (Mat_SeqAIJ*)(A)->data; 143cae5a23dSHong Zhang } 1445a46d3faSBarry Smith #if defined(PETSC_USE_COMPLEX) 145d0f46423SBarry Smith zCreate_CompCol_Matrix(&lu->A,A->cmap->n,A->rmap->n,aa->nz,(doublecomplex*)aa->a,aa->j,aa->i, 1465a46d3faSBarry Smith SLU_NC,SLU_Z,SLU_GE); 1475a46d3faSBarry Smith #else 148d0f46423SBarry Smith dCreate_CompCol_Matrix(&lu->A,A->cmap->n,A->rmap->n,aa->nz,aa->a,aa->j,aa->i, 1495a46d3faSBarry Smith SLU_NC,SLU_D,SLU_GE); 1505a46d3faSBarry Smith #endif 1515a46d3faSBarry Smith 1525a46d3faSBarry Smith /* Numerical factorization */ 1535a46d3faSBarry Smith lu->B.ncol = 0; /* Indicate not to solve the system */ 154d5f3da31SBarry Smith if (F->factortype == MAT_FACTOR_LU){ 1555a46d3faSBarry Smith #if defined(PETSC_USE_COMPLEX) 1565a46d3faSBarry Smith zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 1575a46d3faSBarry Smith &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 1585d8b2efaSHong Zhang &lu->mem_usage, &lu->stat, &sinfo); 1595a46d3faSBarry Smith #else 1605a46d3faSBarry Smith dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 1615a46d3faSBarry Smith &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 1625d8b2efaSHong Zhang &lu->mem_usage, &lu->stat, &sinfo); 1635a46d3faSBarry Smith #endif 164d5f3da31SBarry Smith } else if (F->factortype == MAT_FACTOR_ILU){ 165cffbb591SHong Zhang /* Compute the incomplete factorization, condition number and pivot growth */ 166cffbb591SHong Zhang #if defined(PETSC_USE_COMPLEX) 167cffbb591SHong Zhang zgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r,lu->etree, lu->equed, lu->R, lu->C, 168cffbb591SHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, 1695d8b2efaSHong Zhang &lu->mem_usage, &lu->stat, &sinfo); 170cffbb591SHong Zhang #else 171cffbb591SHong Zhang dgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 172cffbb591SHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, 1735d8b2efaSHong Zhang &lu->mem_usage, &lu->stat, &sinfo); 174cffbb591SHong Zhang #endif 175cffbb591SHong Zhang } else { 176e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported"); 177cffbb591SHong Zhang } 1785a46d3faSBarry Smith if ( !sinfo || sinfo == lu->A.ncol+1 ) { 1795a46d3faSBarry Smith if ( lu->options.PivotGrowth ) 1805a46d3faSBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," Recip. pivot growth = %e\n", lu->rpg); 1815a46d3faSBarry Smith if ( lu->options.ConditionNumber ) 1825a46d3faSBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," Recip. condition number = %e\n", lu->rcond); 1835a46d3faSBarry Smith } else if ( sinfo > 0 ){ 1845a46d3faSBarry Smith if ( lu->lwork == -1 ) { 1855a46d3faSBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," ** Estimated memory: %D bytes\n", sinfo - lu->A.ncol); 1869308c137SBarry Smith } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot in row %D",sinfo); 1875a46d3faSBarry Smith } else { /* sinfo < 0 */ 188e32f2f54SBarry Smith SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_LIB, "info = %D, the %D-th argument in gssvx() had an illegal value", sinfo,-sinfo); 1895a46d3faSBarry Smith } 1905a46d3faSBarry Smith 1915a46d3faSBarry Smith if ( lu->options.PrintStat ) { 1925a46d3faSBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF,"MatLUFactorNumeric_SuperLU():\n"); 1935d8b2efaSHong Zhang StatPrint(&lu->stat); 1945a46d3faSBarry Smith Lstore = (SCformat *) lu->L.Store; 1955a46d3faSBarry Smith Ustore = (NCformat *) lu->U.Store; 1965a46d3faSBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," No of nonzeros in factor L = %D\n", Lstore->nnz); 1975a46d3faSBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," No of nonzeros in factor U = %D\n", Ustore->nnz); 1985a46d3faSBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," No of nonzeros in L+U = %D\n", Lstore->nnz + Ustore->nnz - lu->A.ncol); 1996da386baSHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," L\\U MB %.3f\ttotal MB needed %.3f\n", 2006da386baSHong Zhang lu->mem_usage.for_lu/1e6, lu->mem_usage.total_needed/1e6); 2015a46d3faSBarry Smith } 2025a46d3faSBarry Smith 2035a46d3faSBarry Smith lu->flg = SAME_NONZERO_PATTERN; 2041d5ca7f3SHong Zhang F->ops->solve = MatSolve_SuperLU; 2051d5ca7f3SHong Zhang F->ops->solvetranspose = MatSolveTranspose_SuperLU; 206*e0b74bf9SHong Zhang F->ops->matsolve = MatMatSolve_SuperLU; 2075a46d3faSBarry Smith PetscFunctionReturn(0); 2085a46d3faSBarry Smith } 2095a46d3faSBarry Smith 21014b396bbSKris Buschelman #undef __FUNCT__ 211f0c56d0fSKris Buschelman #define __FUNCT__ "MatDestroy_SuperLU" 212dfbe8321SBarry Smith PetscErrorCode MatDestroy_SuperLU(Mat A) 21314b396bbSKris Buschelman { 214dfbe8321SBarry Smith PetscErrorCode ierr; 215f0c56d0fSKris Buschelman Mat_SuperLU *lu=(Mat_SuperLU*)A->spptr; 21614b396bbSKris Buschelman 21714b396bbSKris Buschelman PetscFunctionBegin; 21875af56d4SHong Zhang if (lu->CleanUpSuperLU) { /* Free the SuperLU datastructures */ 21975af56d4SHong Zhang Destroy_SuperMatrix_Store(&lu->A); 22075af56d4SHong Zhang Destroy_SuperMatrix_Store(&lu->B); 22175af56d4SHong Zhang Destroy_SuperMatrix_Store(&lu->X); 2225d8b2efaSHong Zhang StatFree(&lu->stat); 2230e742b69SHong Zhang if ( lu->lwork >= 0 ) { 2240e742b69SHong Zhang Destroy_SuperNode_Matrix(&lu->L); 2250e742b69SHong Zhang Destroy_CompCol_Matrix(&lu->U); 2260e742b69SHong Zhang } 2270e742b69SHong Zhang } 2289ce50919SHong Zhang 2299ce50919SHong Zhang ierr = PetscFree(lu->etree);CHKERRQ(ierr); 2309ce50919SHong Zhang ierr = PetscFree(lu->perm_r);CHKERRQ(ierr); 2319ce50919SHong Zhang ierr = PetscFree(lu->perm_c);CHKERRQ(ierr); 2329ce50919SHong Zhang ierr = PetscFree(lu->R);CHKERRQ(ierr); 2339ce50919SHong Zhang ierr = PetscFree(lu->C);CHKERRQ(ierr); 2340e742b69SHong Zhang 235d954db57SHong Zhang /* clear composed functions */ 236d954db57SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatFactorGetSolverPackage_C","",PETSC_NULL);CHKERRQ(ierr); 237790a96dfSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatSuperluSetILUDropTol_C","",PETSC_NULL);CHKERRQ(ierr); 238d954db57SHong Zhang 239b24902e0SBarry Smith ierr = MatDestroy_SeqAIJ(A);CHKERRQ(ierr); 240cae5a23dSHong Zhang if (lu->A_dup){ierr = MatDestroy(lu->A_dup);CHKERRQ(ierr);} 241c31cb41cSBarry Smith ierr = PetscFree(lu->rhs_dup);CHKERRQ(ierr); 24214b396bbSKris Buschelman PetscFunctionReturn(0); 24314b396bbSKris Buschelman } 24414b396bbSKris Buschelman 24514b396bbSKris Buschelman #undef __FUNCT__ 246f0c56d0fSKris Buschelman #define __FUNCT__ "MatView_SuperLU" 247dfbe8321SBarry Smith PetscErrorCode MatView_SuperLU(Mat A,PetscViewer viewer) 24814b396bbSKris Buschelman { 249dfbe8321SBarry Smith PetscErrorCode ierr; 250ace3abfcSBarry Smith PetscBool iascii; 25114b396bbSKris Buschelman PetscViewerFormat format; 25214b396bbSKris Buschelman 25314b396bbSKris Buschelman PetscFunctionBegin; 2542692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 25532077d6dSBarry Smith if (iascii) { 25614b396bbSKris Buschelman ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 2572f59403fSHong Zhang if (format == PETSC_VIEWER_ASCII_INFO) { 258f0c56d0fSKris Buschelman ierr = MatFactorInfo_SuperLU(A,viewer);CHKERRQ(ierr); 25914b396bbSKris Buschelman } 26014b396bbSKris Buschelman } 26114b396bbSKris Buschelman PetscFunctionReturn(0); 26214b396bbSKris Buschelman } 26314b396bbSKris Buschelman 26414b396bbSKris Buschelman 26514b396bbSKris Buschelman #undef __FUNCT__ 266c29e884eSHong Zhang #define __FUNCT__ "MatSolve_SuperLU_Private" 267c29e884eSHong Zhang PetscErrorCode MatSolve_SuperLU_Private(Mat A,Vec b,Vec x) 26814b396bbSKris Buschelman { 269f0c56d0fSKris Buschelman Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr; 27075af56d4SHong Zhang PetscScalar *barray,*xarray; 271dfbe8321SBarry Smith PetscErrorCode ierr; 272cae5a23dSHong Zhang PetscInt info,i,n=x->map->n; 273da7a1d00SHong Zhang PetscReal ferr,berr; 27414b396bbSKris Buschelman 27514b396bbSKris Buschelman PetscFunctionBegin; 276937a6b0eSHong Zhang if ( lu->lwork == -1 ) { 277937a6b0eSHong Zhang PetscFunctionReturn(0); 278937a6b0eSHong Zhang } 279cae5a23dSHong Zhang 28075af56d4SHong Zhang lu->B.ncol = 1; /* Set the number of right-hand side */ 281cae5a23dSHong Zhang if (lu->options.Equil && !lu->rhs_dup){ 282cae5a23dSHong Zhang /* superlu overwrites b when Equil is used, thus create rhs_dup to keep user's b unchanged */ 283cae5a23dSHong Zhang ierr = PetscMalloc(n*sizeof(PetscScalar),&lu->rhs_dup);CHKERRQ(ierr); 284cae5a23dSHong Zhang } 285cae5a23dSHong Zhang if (lu->options.Equil){ 286cae5a23dSHong Zhang /* Copy b into rsh_dup */ 28775af56d4SHong Zhang ierr = VecGetArray(b,&barray);CHKERRQ(ierr); 288cae5a23dSHong Zhang ierr = PetscMemcpy(lu->rhs_dup,barray,n*sizeof(PetscScalar));CHKERRQ(ierr); 289cae5a23dSHong Zhang ierr = VecRestoreArray(b,&barray);CHKERRQ(ierr); 290cae5a23dSHong Zhang barray = lu->rhs_dup; 291cae5a23dSHong Zhang } else { 292cae5a23dSHong Zhang ierr = VecGetArray(b,&barray);CHKERRQ(ierr); 293cae5a23dSHong Zhang } 29475af56d4SHong Zhang ierr = VecGetArray(x,&xarray);CHKERRQ(ierr); 2955fe6150dSHong Zhang 2965fe6150dSHong Zhang #if defined(PETSC_USE_COMPLEX) 2975fe6150dSHong Zhang ((DNformat*)lu->B.Store)->nzval = (doublecomplex*)barray; 2985fe6150dSHong Zhang ((DNformat*)lu->X.Store)->nzval = (doublecomplex*)xarray; 2995fe6150dSHong Zhang #else 3005fe6150dSHong Zhang ((DNformat*)lu->B.Store)->nzval = barray; 30175af56d4SHong Zhang ((DNformat*)lu->X.Store)->nzval = xarray; 3025fe6150dSHong Zhang #endif 30375af56d4SHong Zhang 30475af56d4SHong Zhang lu->options.Fact = FACTORED; /* Indicate the factored form of A is supplied. */ 305d5f3da31SBarry Smith if (A->factortype == MAT_FACTOR_LU){ 3068914a3f7SHong Zhang #if defined(PETSC_USE_COMPLEX) 3078914a3f7SHong Zhang zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 3088914a3f7SHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 3095d8b2efaSHong Zhang &lu->mem_usage, &lu->stat, &info); 3108914a3f7SHong Zhang #else 31175af56d4SHong Zhang dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 31275af56d4SHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 3135d8b2efaSHong Zhang &lu->mem_usage, &lu->stat, &info); 3148914a3f7SHong Zhang #endif 315d5f3da31SBarry Smith } else if (A->factortype == MAT_FACTOR_ILU){ 316cffbb591SHong Zhang #if defined(PETSC_USE_COMPLEX) 317cffbb591SHong Zhang zgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 318cffbb591SHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, 3195d8b2efaSHong Zhang &lu->mem_usage, &lu->stat, &info); 320cffbb591SHong Zhang #else 321cffbb591SHong Zhang dgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 322cffbb591SHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, 3235d8b2efaSHong Zhang &lu->mem_usage, &lu->stat, &info); 324cffbb591SHong Zhang #endif 325cffbb591SHong Zhang } else { 326e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported"); 327cffbb591SHong Zhang } 328cae5a23dSHong Zhang if (!lu->options.Equil){ 32975af56d4SHong Zhang ierr = VecRestoreArray(b,&barray);CHKERRQ(ierr); 330cae5a23dSHong Zhang } 33175af56d4SHong Zhang ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr); 33275af56d4SHong Zhang 333958c9bccSBarry Smith if ( !info || info == lu->A.ncol+1 ) { 33475af56d4SHong Zhang if ( lu->options.IterRefine ) { 3358914a3f7SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF,"Iterative Refinement:\n"); 3368914a3f7SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," %8s%8s%16s%16s\n", "rhs", "Steps", "FERR", "BERR"); 33775af56d4SHong Zhang for (i = 0; i < 1; ++i) 3385d8b2efaSHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," %8d%8d%16e%16e\n", i+1, lu->stat.RefineSteps, ferr, berr); 33975af56d4SHong Zhang } 3408914a3f7SHong Zhang } else if ( info > 0 ){ 3418914a3f7SHong Zhang if ( lu->lwork == -1 ) { 34277431f27SBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," ** Estimated memory: %D bytes\n", info - lu->A.ncol); 3438914a3f7SHong Zhang } else { 34477431f27SBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," Warning: gssvx() returns info %D\n",info); 3458914a3f7SHong Zhang } 3468914a3f7SHong Zhang } else if (info < 0){ 347e32f2f54SBarry Smith SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_LIB, "info = %D, the %D-th argument in gssvx() had an illegal value", info,-info); 34814b396bbSKris Buschelman } 34914b396bbSKris Buschelman 3508914a3f7SHong Zhang if ( lu->options.PrintStat ) { 3518914a3f7SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF,"MatSolve__SuperLU():\n"); 3525d8b2efaSHong Zhang StatPrint(&lu->stat); 3538914a3f7SHong Zhang } 35475af56d4SHong Zhang PetscFunctionReturn(0); 35575af56d4SHong Zhang } 35614b396bbSKris Buschelman 35714b396bbSKris Buschelman #undef __FUNCT__ 358c29e884eSHong Zhang #define __FUNCT__ "MatSolve_SuperLU" 359c29e884eSHong Zhang PetscErrorCode MatSolve_SuperLU(Mat A,Vec b,Vec x) 360c29e884eSHong Zhang { 361c29e884eSHong Zhang Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr; 362c29e884eSHong Zhang PetscErrorCode ierr; 363c29e884eSHong Zhang 364c29e884eSHong Zhang PetscFunctionBegin; 365c29e884eSHong Zhang lu->options.Trans = TRANS; 366c29e884eSHong Zhang ierr = MatSolve_SuperLU_Private(A,b,x);CHKERRQ(ierr); 367c29e884eSHong Zhang PetscFunctionReturn(0); 368c29e884eSHong Zhang } 369c29e884eSHong Zhang 370c29e884eSHong Zhang #undef __FUNCT__ 371c7c1fe80SHong Zhang #define __FUNCT__ "MatSolveTranspose_SuperLU" 372c7c1fe80SHong Zhang PetscErrorCode MatSolveTranspose_SuperLU(Mat A,Vec b,Vec x) 373c7c1fe80SHong Zhang { 374c7c1fe80SHong Zhang Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr; 375c7c1fe80SHong Zhang PetscErrorCode ierr; 376c7c1fe80SHong Zhang 377c7c1fe80SHong Zhang PetscFunctionBegin; 378c7c1fe80SHong Zhang lu->options.Trans = NOTRANS; 379c29e884eSHong Zhang ierr = MatSolve_SuperLU_Private(A,b,x);CHKERRQ(ierr); 380c7c1fe80SHong Zhang PetscFunctionReturn(0); 381c7c1fe80SHong Zhang } 382c7c1fe80SHong Zhang 383*e0b74bf9SHong Zhang #undef __FUNCT__ 384*e0b74bf9SHong Zhang #define __FUNCT__ "MatMatSolve_SuperLU" 385*e0b74bf9SHong Zhang PetscErrorCode MatMatSolve_SuperLU(Mat A,Mat B,Mat X) 386*e0b74bf9SHong Zhang { 387*e0b74bf9SHong Zhang Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr; 388*e0b74bf9SHong Zhang 389*e0b74bf9SHong Zhang PetscFunctionBegin; 390*e0b74bf9SHong Zhang lu->options.Trans = TRANS; 391*e0b74bf9SHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatMatSolve_SuperLU() is not implemented yet"); 392*e0b74bf9SHong Zhang PetscFunctionReturn(0); 393*e0b74bf9SHong Zhang } 394*e0b74bf9SHong Zhang 39514b396bbSKris Buschelman /* 39614b396bbSKris Buschelman Note the r permutation is ignored 39714b396bbSKris Buschelman */ 39814b396bbSKris Buschelman #undef __FUNCT__ 399f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_SuperLU" 4000481f469SBarry Smith PetscErrorCode MatLUFactorSymbolic_SuperLU(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info) 40114b396bbSKris Buschelman { 4021d5ca7f3SHong Zhang Mat_SuperLU *lu = (Mat_SuperLU*)(F->spptr); 403b24902e0SBarry Smith 404b24902e0SBarry Smith PetscFunctionBegin; 405b24902e0SBarry Smith lu->flg = DIFFERENT_NONZERO_PATTERN; 406b24902e0SBarry Smith lu->CleanUpSuperLU = PETSC_TRUE; 4071d5ca7f3SHong Zhang F->ops->lufactornumeric = MatLUFactorNumeric_SuperLU; 408b24902e0SBarry Smith PetscFunctionReturn(0); 409b24902e0SBarry Smith } 410b24902e0SBarry Smith 41135bd34faSBarry Smith EXTERN_C_BEGIN 41235bd34faSBarry Smith #undef __FUNCT__ 4135ccb76cbSHong Zhang #define __FUNCT__ "MatSuperluSetILUDropTol_SuperLU" 4145ccb76cbSHong Zhang PetscErrorCode MatSuperluSetILUDropTol_SuperLU(Mat F,PetscReal dtol) 4155ccb76cbSHong Zhang { 4165ccb76cbSHong Zhang Mat_SuperLU *lu= (Mat_SuperLU*)F->spptr; 4175ccb76cbSHong Zhang 4185ccb76cbSHong Zhang PetscFunctionBegin; 4195ccb76cbSHong Zhang lu->options.ILU_DropTol = dtol; 4205ccb76cbSHong Zhang PetscFunctionReturn(0); 4215ccb76cbSHong Zhang } 4225ccb76cbSHong Zhang EXTERN_C_END 4235ccb76cbSHong Zhang 4245ccb76cbSHong Zhang #undef __FUNCT__ 4255ccb76cbSHong Zhang #define __FUNCT__ "MatSuperluSetILUDropTol" 4265ccb76cbSHong Zhang /*@ 4275ccb76cbSHong Zhang MatSuperluSetILUDropTol - Set SuperLU ILU drop tolerance 4285ccb76cbSHong Zhang Logically Collective on Mat 4295ccb76cbSHong Zhang 4305ccb76cbSHong Zhang Input Parameters: 4315ccb76cbSHong Zhang + F - the factored matrix obtained by calling MatGetFactor() from PETSc-SuperLU interface 4325ccb76cbSHong Zhang - dtol - drop tolerance 4335ccb76cbSHong Zhang 4345ccb76cbSHong Zhang Options Database: 4355ccb76cbSHong Zhang . -mat_superlu_ilu_droptol <dtol> 4365ccb76cbSHong Zhang 4375ccb76cbSHong Zhang Level: beginner 4385ccb76cbSHong Zhang 4395ccb76cbSHong Zhang References: SuperLU Users' Guide 4405ccb76cbSHong Zhang 4415ccb76cbSHong Zhang .seealso: MatGetFactor() 4425ccb76cbSHong Zhang @*/ 4435ccb76cbSHong Zhang PetscErrorCode MatSuperluSetILUDropTol(Mat F,PetscReal dtol) 4445ccb76cbSHong Zhang { 4455ccb76cbSHong Zhang PetscErrorCode ierr; 4465ccb76cbSHong Zhang 4475ccb76cbSHong Zhang PetscFunctionBegin; 4485ccb76cbSHong Zhang PetscValidHeaderSpecific(F,MAT_CLASSID,1); 4495ccb76cbSHong Zhang PetscValidLogicalCollectiveInt(F,dtol,2); 4505ccb76cbSHong Zhang ierr = PetscTryMethod(F,"MatSuperluSetILUDropTol_C",(Mat,PetscReal),(F,dtol));CHKERRQ(ierr); 4515ccb76cbSHong Zhang PetscFunctionReturn(0); 4525ccb76cbSHong Zhang } 4535ccb76cbSHong Zhang 4545ccb76cbSHong Zhang EXTERN_C_BEGIN 4555ccb76cbSHong Zhang #undef __FUNCT__ 45635bd34faSBarry Smith #define __FUNCT__ "MatFactorGetSolverPackage_seqaij_superlu" 45735bd34faSBarry Smith PetscErrorCode MatFactorGetSolverPackage_seqaij_superlu(Mat A,const MatSolverPackage *type) 45835bd34faSBarry Smith { 45935bd34faSBarry Smith PetscFunctionBegin; 4602692d6eeSBarry Smith *type = MATSOLVERSUPERLU; 46135bd34faSBarry Smith PetscFunctionReturn(0); 46235bd34faSBarry Smith } 46335bd34faSBarry Smith EXTERN_C_END 46435bd34faSBarry Smith 465b24902e0SBarry Smith 466b24902e0SBarry Smith /*MC 467ba992d64SSatish Balay MATSOLVERSUPERLU = "superlu" - A solver package providing solvers LU and ILU for sequential matrices 468b24902e0SBarry Smith via the external package SuperLU. 469b24902e0SBarry Smith 470e2e64c6bSBarry Smith Use ./configure --download-superlu to have PETSc installed with SuperLU 471b24902e0SBarry Smith 472b24902e0SBarry Smith Options Database Keys: 4732de22e2eSHong Zhang + -mat_superlu_equil: <FALSE> Equil (None) 4742de22e2eSHong Zhang . -mat_superlu_colperm <COLAMD> (choose one of) NATURAL MMD_ATA MMD_AT_PLUS_A COLAMD 4752de22e2eSHong Zhang . -mat_superlu_iterrefine <NOREFINE> (choose one of) NOREFINE SINGLE DOUBLE EXTRA 4762de22e2eSHong Zhang . -mat_superlu_symmetricmode: <FALSE> SymmetricMode (None) 4772de22e2eSHong Zhang . -mat_superlu_diagpivotthresh <1>: DiagPivotThresh (None) 4782de22e2eSHong Zhang . -mat_superlu_pivotgrowth: <FALSE> PivotGrowth (None) 4792de22e2eSHong Zhang . -mat_superlu_conditionnumber: <FALSE> ConditionNumber (None) 4802de22e2eSHong Zhang . -mat_superlu_rowperm <NOROWPERM> (choose one of) NOROWPERM LargeDiag 4812de22e2eSHong Zhang . -mat_superlu_replacetinypivot: <FALSE> ReplaceTinyPivot (None) 4822de22e2eSHong Zhang . -mat_superlu_printstat: <FALSE> PrintStat (None) 4832de22e2eSHong Zhang . -mat_superlu_lwork <0>: size of work array in bytes used by factorization (None) 4842de22e2eSHong Zhang . -mat_superlu_ilu_droptol <0>: ILU_DropTol (None) 4852de22e2eSHong Zhang . -mat_superlu_ilu_filltol <0>: ILU_FillTol (None) 4862de22e2eSHong Zhang . -mat_superlu_ilu_fillfactor <0>: ILU_FillFactor (None) 4872de22e2eSHong Zhang . -mat_superlu_ilu_droprull <0>: ILU_DropRule (None) 4882de22e2eSHong Zhang . -mat_superlu_ilu_norm <0>: ILU_Norm (None) 4892de22e2eSHong Zhang - -mat_superlu_ilu_milu <0>: ILU_MILU (None) 490b24902e0SBarry Smith 4912692d6eeSBarry Smith Notes: Do not confuse this with MATSOLVERSUPERLU_DIST which is for parallel sparse solves 4925c9eb25fSBarry Smith 493b24902e0SBarry Smith Level: beginner 494b24902e0SBarry Smith 495ba992d64SSatish Balay .seealso: PCLU, PCILU, MATSOLVERSUPERLU_DIST, MATSOLVERMUMPS, MATSOLVERSPOOLES, PCFactorSetMatSolverPackage(), MatSolverPackage 496b24902e0SBarry Smith M*/ 497b24902e0SBarry Smith 4985c9eb25fSBarry Smith EXTERN_C_BEGIN 499b24902e0SBarry Smith #undef __FUNCT__ 500b24902e0SBarry Smith #define __FUNCT__ "MatGetFactor_seqaij_superlu" 5015c9eb25fSBarry Smith PetscErrorCode MatGetFactor_seqaij_superlu(Mat A,MatFactorType ftype,Mat *F) 502b24902e0SBarry Smith { 50314b396bbSKris Buschelman Mat B; 504f0c56d0fSKris Buschelman Mat_SuperLU *lu; 5056849ba73SBarry Smith PetscErrorCode ierr; 5065d8b2efaSHong Zhang PetscInt indx,m=A->rmap->n,n=A->cmap->n; 507ace3abfcSBarry Smith PetscBool flg; 5085ca28756SHong Zhang const char *colperm[]={"NATURAL","MMD_ATA","MMD_AT_PLUS_A","COLAMD"}; /* MY_PERMC - not supported by the petsc interface yet */ 5095ca28756SHong Zhang const char *iterrefine[]={"NOREFINE", "SINGLE", "DOUBLE", "EXTRA"}; 5105ca28756SHong Zhang const char *rowperm[]={"NOROWPERM", "LargeDiag"}; /* MY_PERMC - not supported by the petsc interface yet */ 51114b396bbSKris Buschelman 51214b396bbSKris Buschelman PetscFunctionBegin; 5137adad957SLisandro Dalcin ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr); 514d0f46423SBarry Smith ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 5157adad957SLisandro Dalcin ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 516899d7b4fSKris Buschelman ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr); 517f0c56d0fSKris Buschelman 518cffbb591SHong Zhang if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU){ 519b24902e0SBarry Smith B->ops->lufactorsymbolic = MatLUFactorSymbolic_SuperLU; 520cffbb591SHong Zhang B->ops->ilufactorsymbolic = MatLUFactorSymbolic_SuperLU; 521b3a44c85SBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported"); 522cffbb591SHong Zhang 523b24902e0SBarry Smith B->ops->destroy = MatDestroy_SuperLU; 5243519fcdcSHong Zhang B->ops->view = MatView_SuperLU; 525d5f3da31SBarry Smith B->factortype = ftype; 52694b7f48cSBarry Smith B->assembled = PETSC_TRUE; /* required by -ksp_view */ 5275c9eb25fSBarry Smith B->preallocated = PETSC_TRUE; 52814b396bbSKris Buschelman 529b24902e0SBarry Smith ierr = PetscNewLog(B,Mat_SuperLU,&lu);CHKERRQ(ierr); 530cae5a23dSHong Zhang 531cffbb591SHong Zhang if (ftype == MAT_FACTOR_LU){ 5329ce50919SHong Zhang set_default_options(&lu->options); 5333d6581fbSHong Zhang /* Comments from SuperLU_4.0/SRC/dgssvx.c: 5343d6581fbSHong Zhang "Whether or not the system will be equilibrated depends on the 5353d6581fbSHong Zhang scaling of the matrix A, but if equilibration is used, A is 5363d6581fbSHong Zhang overwritten by diag(R)*A*diag(C) and B by diag(R)*B 5373d6581fbSHong Zhang (if options->Trans=NOTRANS) or diag(C)*B (if options->Trans = TRANS or CONJ)." 5383d6581fbSHong Zhang We set 'options.Equil = NO' as default because additional space is needed for it. 5393d6581fbSHong Zhang */ 5403d6581fbSHong Zhang lu->options.Equil = NO; 541cffbb591SHong Zhang } else if (ftype == MAT_FACTOR_ILU){ 542cffbb591SHong Zhang /* Set the default input options of ilu: 543cffbb591SHong Zhang options.Fact = DOFACT; 5443d6581fbSHong Zhang options.Equil = YES; // must be YES for ilu - don't know why 545cffbb591SHong Zhang options.ColPerm = COLAMD; 546cffbb591SHong Zhang options.DiagPivotThresh = 0.1; //different from complete LU 547cffbb591SHong Zhang options.Trans = NOTRANS; 548cffbb591SHong Zhang options.IterRefine = NOREFINE; 549cffbb591SHong Zhang options.SymmetricMode = NO; 550cffbb591SHong Zhang options.PivotGrowth = NO; 551cffbb591SHong Zhang options.ConditionNumber = NO; 552cffbb591SHong Zhang options.PrintStat = YES; 553cffbb591SHong Zhang options.RowPerm = LargeDiag; 554cffbb591SHong Zhang options.ILU_DropTol = 1e-4; 555cffbb591SHong Zhang options.ILU_FillTol = 1e-2; 556cffbb591SHong Zhang options.ILU_FillFactor = 10.0; 557cffbb591SHong Zhang options.ILU_DropRule = DROP_BASIC | DROP_AREA; 558cffbb591SHong Zhang options.ILU_Norm = INF_NORM; 559cffbb591SHong Zhang options.ILU_MILU = SMILU_2; 560cffbb591SHong Zhang */ 561cffbb591SHong Zhang ilu_set_default_options(&lu->options); 562cffbb591SHong Zhang } 5639ce50919SHong Zhang lu->options.PrintStat = NO; 5641a2e2f44SHong Zhang 5655d8b2efaSHong Zhang /* Initialize the statistics variables. */ 5665d8b2efaSHong Zhang StatInit(&lu->stat); 5678914a3f7SHong Zhang lu->lwork = 0; /* allocate space internally by system malloc */ 5689ce50919SHong Zhang 5697adad957SLisandro Dalcin ierr = PetscOptionsBegin(((PetscObject)A)->comm,((PetscObject)A)->prefix,"SuperLU Options","Mat");CHKERRQ(ierr); 570acfcf0e5SJed Brown ierr = PetscOptionsBool("-mat_superlu_equil","Equil","None",(PetscBool)lu->options.Equil,(PetscBool*)&lu->options.Equil,0);CHKERRQ(ierr); 5718914a3f7SHong Zhang ierr = PetscOptionsEList("-mat_superlu_colperm","ColPerm","None",colperm,4,colperm[3],&indx,&flg);CHKERRQ(ierr); 5729ce50919SHong Zhang if (flg) {lu->options.ColPerm = (colperm_t)indx;} 5738914a3f7SHong Zhang ierr = PetscOptionsEList("-mat_superlu_iterrefine","IterRefine","None",iterrefine,4,iterrefine[0],&indx,&flg);CHKERRQ(ierr); 5749ce50919SHong Zhang if (flg) { lu->options.IterRefine = (IterRefine_t)indx;} 575acfcf0e5SJed Brown ierr = PetscOptionsBool("-mat_superlu_symmetricmode","SymmetricMode","None",(PetscBool)lu->options.SymmetricMode,&flg,0);CHKERRQ(ierr); 5769ce50919SHong Zhang if (flg) lu->options.SymmetricMode = YES; 5778914a3f7SHong Zhang ierr = PetscOptionsReal("-mat_superlu_diagpivotthresh","DiagPivotThresh","None",lu->options.DiagPivotThresh,&lu->options.DiagPivotThresh,PETSC_NULL);CHKERRQ(ierr); 578acfcf0e5SJed Brown ierr = PetscOptionsBool("-mat_superlu_pivotgrowth","PivotGrowth","None",(PetscBool)lu->options.PivotGrowth,&flg,0);CHKERRQ(ierr); 5799ce50919SHong Zhang if (flg) lu->options.PivotGrowth = YES; 580acfcf0e5SJed Brown ierr = PetscOptionsBool("-mat_superlu_conditionnumber","ConditionNumber","None",(PetscBool)lu->options.ConditionNumber,&flg,0);CHKERRQ(ierr); 5819ce50919SHong Zhang if (flg) lu->options.ConditionNumber = YES; 582d7ebd59bSHong Zhang ierr = PetscOptionsEList("-mat_superlu_rowperm","rowperm","None",rowperm,2,rowperm[lu->options.RowPerm],&indx,&flg);CHKERRQ(ierr); 5839ce50919SHong Zhang if (flg) {lu->options.RowPerm = (rowperm_t)indx;} 584acfcf0e5SJed Brown ierr = PetscOptionsBool("-mat_superlu_replacetinypivot","ReplaceTinyPivot","None",(PetscBool)lu->options.ReplaceTinyPivot,&flg,0);CHKERRQ(ierr); 5859ce50919SHong Zhang if (flg) lu->options.ReplaceTinyPivot = YES; 586acfcf0e5SJed Brown ierr = PetscOptionsBool("-mat_superlu_printstat","PrintStat","None",(PetscBool)lu->options.PrintStat,&flg,0);CHKERRQ(ierr); 5879ce50919SHong Zhang if (flg) lu->options.PrintStat = YES; 5888914a3f7SHong Zhang ierr = PetscOptionsInt("-mat_superlu_lwork","size of work array in bytes used by factorization","None",lu->lwork,&lu->lwork,PETSC_NULL);CHKERRQ(ierr); 5895fe6150dSHong Zhang if (lu->lwork > 0 ){ 5905fe6150dSHong Zhang ierr = PetscMalloc(lu->lwork,&lu->work);CHKERRQ(ierr); 5915fe6150dSHong Zhang } else if (lu->lwork != 0 && lu->lwork != -1){ 59277431f27SBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," Warning: lwork %D is not supported by SUPERLU. The default lwork=0 is used.\n",lu->lwork); 5938914a3f7SHong Zhang lu->lwork = 0; 5948914a3f7SHong Zhang } 595cffbb591SHong Zhang /* ilu options */ 596cffbb591SHong Zhang ierr = PetscOptionsReal("-mat_superlu_ilu_droptol","ILU_DropTol","None",lu->options.ILU_DropTol,&lu->options.ILU_DropTol,PETSC_NULL);CHKERRQ(ierr); 597cffbb591SHong Zhang ierr = PetscOptionsReal("-mat_superlu_ilu_filltol","ILU_FillTol","None",lu->options.ILU_FillTol,&lu->options.ILU_FillTol,PETSC_NULL);CHKERRQ(ierr); 598cffbb591SHong Zhang ierr = PetscOptionsReal("-mat_superlu_ilu_fillfactor","ILU_FillFactor","None",lu->options.ILU_FillFactor,&lu->options.ILU_FillFactor,PETSC_NULL);CHKERRQ(ierr); 599cffbb591SHong Zhang ierr = PetscOptionsInt("-mat_superlu_ilu_droprull","ILU_DropRule","None",lu->options.ILU_DropRule,&lu->options.ILU_DropRule,PETSC_NULL);CHKERRQ(ierr); 600cffbb591SHong Zhang ierr = PetscOptionsInt("-mat_superlu_ilu_norm","ILU_Norm","None",lu->options.ILU_Norm,&indx,&flg);CHKERRQ(ierr); 601cffbb591SHong Zhang if (flg){ 602cffbb591SHong Zhang lu->options.ILU_Norm = (norm_t)indx; 603cffbb591SHong Zhang } 604cffbb591SHong Zhang ierr = PetscOptionsInt("-mat_superlu_ilu_milu","ILU_MILU","None",lu->options.ILU_MILU,&indx,&flg);CHKERRQ(ierr); 605cffbb591SHong Zhang if (flg){ 606cffbb591SHong Zhang lu->options.ILU_MILU = (milu_t)indx; 607cffbb591SHong Zhang } 6089ce50919SHong Zhang PetscOptionsEnd(); 60938abfdaeSHong Zhang if (lu->options.Equil == YES) { 61038abfdaeSHong Zhang /* superlu overwrites input matrix and rhs when Equil is used, thus create A_dup to keep user's A unchanged */ 61138abfdaeSHong Zhang ierr = MatDuplicate_SeqAIJ(A,MAT_COPY_VALUES,&lu->A_dup);CHKERRQ(ierr); 61238abfdaeSHong Zhang } 6139ce50919SHong Zhang 6145d8b2efaSHong Zhang /* Allocate spaces (notice sizes are for the transpose) */ 6155d8b2efaSHong Zhang ierr = PetscMalloc(m*sizeof(PetscInt),&lu->etree);CHKERRQ(ierr); 6165d8b2efaSHong Zhang ierr = PetscMalloc(n*sizeof(PetscInt),&lu->perm_r);CHKERRQ(ierr); 6175d8b2efaSHong Zhang ierr = PetscMalloc(m*sizeof(PetscInt),&lu->perm_c);CHKERRQ(ierr); 6185d8b2efaSHong Zhang ierr = PetscMalloc(n*sizeof(PetscScalar),&lu->R);CHKERRQ(ierr); 6195d8b2efaSHong Zhang ierr = PetscMalloc(m*sizeof(PetscScalar),&lu->C);CHKERRQ(ierr); 6205d8b2efaSHong Zhang 6215d8b2efaSHong Zhang /* create rhs and solution x without allocate space for .Store */ 6225d8b2efaSHong Zhang #if defined(PETSC_USE_COMPLEX) 6235d8b2efaSHong Zhang zCreate_Dense_Matrix(&lu->B, m, 1, PETSC_NULL, m, SLU_DN, SLU_Z, SLU_GE); 6245d8b2efaSHong Zhang zCreate_Dense_Matrix(&lu->X, m, 1, PETSC_NULL, m, SLU_DN, SLU_Z, SLU_GE); 6255d8b2efaSHong Zhang #else 6265d8b2efaSHong Zhang dCreate_Dense_Matrix(&lu->B, m, 1, PETSC_NULL, m, SLU_DN, SLU_D, SLU_GE); 6275d8b2efaSHong Zhang dCreate_Dense_Matrix(&lu->X, m, 1, PETSC_NULL, m, SLU_DN, SLU_D, SLU_GE); 6285d8b2efaSHong Zhang #endif 6295d8b2efaSHong Zhang 63075af56d4SHong Zhang #ifdef SUPERLU2 6315c9eb25fSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatCreateNull","MatCreateNull_SuperLU",(void(*)(void))MatCreateNull_SuperLU);CHKERRQ(ierr); 63275af56d4SHong Zhang #endif 63335bd34faSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_seqaij_superlu",MatFactorGetSolverPackage_seqaij_superlu);CHKERRQ(ierr); 6345ccb76cbSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSuperluSetILUDropTol_C","MatSuperluSetILUDropTol_SuperLU",MatSuperluSetILUDropTol_SuperLU);CHKERRQ(ierr); 6355c9eb25fSBarry Smith B->spptr = lu; 636899d7b4fSKris Buschelman *F = B; 63714b396bbSKris Buschelman PetscFunctionReturn(0); 63814b396bbSKris Buschelman } 6395c9eb25fSBarry Smith EXTERN_C_END 640d954db57SHong Zhang 641