1be1d678aSKris Buschelman #define PETSCMAT_DLL 214b396bbSKris Buschelman 35a46d3faSBarry Smith /* -------------------------------------------------------------------- 45a46d3faSBarry Smith 55a46d3faSBarry Smith This file implements a subclass of the SeqAIJ matrix class that uses 65a46d3faSBarry Smith the SuperLU 3.0 sparse solver. You can use this as a starting point for 75a46d3faSBarry Smith implementing your own subclass of a PETSc matrix class. 85a46d3faSBarry Smith 95a46d3faSBarry Smith This demonstrates a way to make an implementation inheritence of a PETSc 105a46d3faSBarry Smith matrix type. This means constructing a new matrix type (SuperLU) by changing some 115a46d3faSBarry Smith of the methods of the previous type (SeqAIJ), adding additional data, and possibly 125a46d3faSBarry Smith additional method. (See any book on object oriented programming). 1314b396bbSKris Buschelman */ 1414b396bbSKris Buschelman 155a46d3faSBarry Smith /* 165a46d3faSBarry Smith Defines the data structure for the base matrix type (SeqAIJ) 175a46d3faSBarry Smith */ 1814b396bbSKris Buschelman #include "src/mat/impls/aij/seq/aij.h" 1914b396bbSKris Buschelman 205a46d3faSBarry Smith /* 215a46d3faSBarry Smith SuperLU include files 225a46d3faSBarry Smith */ 2314b396bbSKris Buschelman EXTERN_C_BEGIN 2414b396bbSKris Buschelman #if defined(PETSC_USE_COMPLEX) 2566f57492SHong Zhang #include "slu_zdefs.h" 2614b396bbSKris Buschelman #else 2766f57492SHong Zhang #include "slu_ddefs.h" 2814b396bbSKris Buschelman #endif 2966f57492SHong Zhang #include "slu_util.h" 3014b396bbSKris Buschelman EXTERN_C_END 3114b396bbSKris Buschelman 325a46d3faSBarry Smith /* 335a46d3faSBarry Smith This is the data we are "ADDING" to the SeqAIJ matrix type to get the SuperLU matrix type 345a46d3faSBarry Smith */ 3514b396bbSKris Buschelman typedef struct { 3675af56d4SHong Zhang SuperMatrix A,L,U,B,X; 3775af56d4SHong Zhang superlu_options_t options; 38da7a1d00SHong Zhang PetscInt *perm_c; /* column permutation vector */ 39da7a1d00SHong Zhang PetscInt *perm_r; /* row permutations from partial pivoting */ 40da7a1d00SHong Zhang PetscInt *etree; 41da7a1d00SHong Zhang PetscReal *R, *C; 4275af56d4SHong Zhang char equed[1]; 43da7a1d00SHong Zhang PetscInt lwork; 4475af56d4SHong Zhang void *work; 45da7a1d00SHong Zhang PetscReal rpg, rcond; 4675af56d4SHong Zhang mem_usage_t mem_usage; 4775af56d4SHong Zhang MatStructure flg; 4814b396bbSKris Buschelman 495a46d3faSBarry Smith /* 505a46d3faSBarry Smith This is where the methods for the superclass (SeqAIJ) are kept while we 515a46d3faSBarry Smith reset the pointers in the function table to the new (SuperLU) versions 525a46d3faSBarry Smith */ 536849ba73SBarry Smith PetscErrorCode (*MatDuplicate)(Mat,MatDuplicateOption,Mat*); 546849ba73SBarry Smith PetscErrorCode (*MatView)(Mat,PetscViewer); 556849ba73SBarry Smith PetscErrorCode (*MatAssemblyEnd)(Mat,MatAssemblyType); 566849ba73SBarry Smith PetscErrorCode (*MatLUFactorSymbolic)(Mat,IS,IS,MatFactorInfo*,Mat*); 576849ba73SBarry Smith PetscErrorCode (*MatDestroy)(Mat); 5814b396bbSKris Buschelman 5914b396bbSKris Buschelman /* Flag to clean up (non-global) SuperLU objects during Destroy */ 6014b396bbSKris Buschelman PetscTruth CleanUpSuperLU; 61f0c56d0fSKris Buschelman } Mat_SuperLU; 6214b396bbSKris Buschelman 63e0e586b9SHong Zhang extern PetscErrorCode MatFactorInfo_SuperLU(Mat,PetscViewer); 64e0e586b9SHong Zhang extern PetscErrorCode MatLUFactorNumeric_SuperLU(Mat,MatFactorInfo *,Mat *); 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); 69e0e586b9SHong Zhang extern PetscErrorCode MatSolveTranspose_SuperLU(Mat,Vec,Vec); 70e0e586b9SHong Zhang extern PetscErrorCode MatLUFactorSymbolic_SuperLU(Mat,IS,IS,MatFactorInfo *,Mat *); 71e0e586b9SHong Zhang extern PetscErrorCode MatDuplicate_SuperLU(Mat, MatDuplicateOption, Mat *); 72e0e586b9SHong Zhang 735a46d3faSBarry Smith /* 745a46d3faSBarry Smith Takes a SuperLU matrix (that is a SeqAIJ matrix with the additional SuperLU data-structures 755a46d3faSBarry Smith and methods) and converts it back to a regular SeqAIJ matrix. 765a46d3faSBarry Smith */ 775a46d3faSBarry Smith EXTERN_C_BEGIN 785a46d3faSBarry Smith #undef __FUNCT__ 795a46d3faSBarry Smith #define __FUNCT__ "MatConvert_SuperLU_SeqAIJ" 805a46d3faSBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SuperLU_SeqAIJ(Mat A,MatType type,MatReuse reuse,Mat *newmat) 815a46d3faSBarry Smith { 825a46d3faSBarry Smith PetscErrorCode ierr; 835a46d3faSBarry Smith Mat B=*newmat; 845a46d3faSBarry Smith Mat_SuperLU *lu=(Mat_SuperLU *)A->spptr; 8514b396bbSKris Buschelman 865a46d3faSBarry Smith PetscFunctionBegin; 875a46d3faSBarry Smith if (reuse == MAT_INITIAL_MATRIX) { 885a46d3faSBarry Smith ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr); 895a46d3faSBarry Smith } 905a46d3faSBarry Smith /* Reset the original SeqAIJ function pointers */ 915a46d3faSBarry Smith B->ops->duplicate = lu->MatDuplicate; 925a46d3faSBarry Smith B->ops->view = lu->MatView; 935a46d3faSBarry Smith B->ops->assemblyend = lu->MatAssemblyEnd; 945a46d3faSBarry Smith B->ops->lufactorsymbolic = lu->MatLUFactorSymbolic; 955a46d3faSBarry Smith B->ops->destroy = lu->MatDestroy; 9615d897acSHong Zhang ierr = PetscFree(lu);CHKERRQ(ierr); 9715d897acSHong Zhang A->spptr = PETSC_NULL; 985a46d3faSBarry Smith 995a46d3faSBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaij_superlu_C","",PETSC_NULL);CHKERRQ(ierr); 1005a46d3faSBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_superlu_seqaij_C","",PETSC_NULL);CHKERRQ(ierr); 1015a46d3faSBarry Smith 1025a46d3faSBarry Smith /* change the type name back to its original value */ 1035a46d3faSBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJ);CHKERRQ(ierr); 1045a46d3faSBarry Smith *newmat = B; 1055a46d3faSBarry Smith PetscFunctionReturn(0); 1065a46d3faSBarry Smith } 1075a46d3faSBarry Smith EXTERN_C_END 108b0592601SKris Buschelman 109b0592601SKris Buschelman EXTERN_C_BEGIN 1105a46d3faSBarry Smith #undef __FUNCT__ 1115a46d3faSBarry Smith #define __FUNCT__ "MatConvert_SeqAIJ_SuperLU" 1125a46d3faSBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqAIJ_SuperLU(Mat A,MatType type,MatReuse reuse,Mat *newmat) 1135a46d3faSBarry Smith { 1145a46d3faSBarry Smith PetscErrorCode ierr; 1155a46d3faSBarry Smith Mat B=*newmat; 1165a46d3faSBarry Smith Mat_SuperLU *lu; 1175a46d3faSBarry Smith 1185a46d3faSBarry Smith PetscFunctionBegin; 1195a46d3faSBarry Smith if (reuse == MAT_INITIAL_MATRIX){ 1205a46d3faSBarry Smith ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr); 1215a46d3faSBarry Smith } 1225a46d3faSBarry Smith 12338f2d2fdSLisandro Dalcin ierr = PetscNewLog(B,Mat_SuperLU,&lu);CHKERRQ(ierr); 1245a46d3faSBarry Smith /* save the original SeqAIJ methods that we are changing */ 1255a46d3faSBarry Smith lu->MatDuplicate = A->ops->duplicate; 1265a46d3faSBarry Smith lu->MatView = A->ops->view; 1275a46d3faSBarry Smith lu->MatAssemblyEnd = A->ops->assemblyend; 1285a46d3faSBarry Smith lu->MatLUFactorSymbolic = A->ops->lufactorsymbolic; 1295a46d3faSBarry Smith lu->MatDestroy = A->ops->destroy; 1305a46d3faSBarry Smith lu->CleanUpSuperLU = PETSC_FALSE; 1315a46d3faSBarry Smith 1325a46d3faSBarry Smith /* add to the matrix the location for all the SuperLU data is to be stored */ 13315d897acSHong Zhang B->spptr = (void*)lu; /* attach Mat_SuperLU to B->spptr is a bad design! */ 1345a46d3faSBarry Smith 1355a46d3faSBarry Smith /* set the methods in the function table to the SuperLU versions */ 1365a46d3faSBarry Smith B->ops->duplicate = MatDuplicate_SuperLU; 1375a46d3faSBarry Smith B->ops->view = MatView_SuperLU; 1385a46d3faSBarry Smith B->ops->assemblyend = MatAssemblyEnd_SuperLU; 1395a46d3faSBarry Smith B->ops->lufactorsymbolic = MatLUFactorSymbolic_SuperLU; 1405a46d3faSBarry Smith B->ops->choleskyfactorsymbolic = 0; 1415a46d3faSBarry Smith B->ops->destroy = MatDestroy_SuperLU; 1425a46d3faSBarry Smith 1435a46d3faSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_superlu_C", 1445a46d3faSBarry Smith "MatConvert_SeqAIJ_SuperLU",MatConvert_SeqAIJ_SuperLU);CHKERRQ(ierr); 1455a46d3faSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_superlu_seqaij_C", 1465a46d3faSBarry Smith "MatConvert_SuperLU_SeqAIJ",MatConvert_SuperLU_SeqAIJ);CHKERRQ(ierr); 1475a46d3faSBarry Smith ierr = PetscInfo(0,"Using SuperLU for SeqAIJ LU factorization and solves.\n");CHKERRQ(ierr); 1485a46d3faSBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)B,MATSUPERLU);CHKERRQ(ierr); 1495a46d3faSBarry Smith *newmat = B; 1505a46d3faSBarry Smith PetscFunctionReturn(0); 1515a46d3faSBarry Smith } 152b0592601SKris Buschelman EXTERN_C_END 15314b396bbSKris Buschelman 1545a46d3faSBarry Smith /* 1555a46d3faSBarry Smith Utility function 1565a46d3faSBarry Smith */ 1575a46d3faSBarry Smith #undef __FUNCT__ 1585a46d3faSBarry Smith #define __FUNCT__ "MatFactorInfo_SuperLU" 1595a46d3faSBarry Smith PetscErrorCode MatFactorInfo_SuperLU(Mat A,PetscViewer viewer) 1605a46d3faSBarry Smith { 1615a46d3faSBarry Smith Mat_SuperLU *lu= (Mat_SuperLU*)A->spptr; 1625a46d3faSBarry Smith PetscErrorCode ierr; 1635a46d3faSBarry Smith superlu_options_t options; 1645a46d3faSBarry Smith 1655a46d3faSBarry Smith PetscFunctionBegin; 1665a46d3faSBarry Smith /* check if matrix is superlu_dist type */ 1675a46d3faSBarry Smith if (A->ops->solve != MatSolve_SuperLU) PetscFunctionReturn(0); 1685a46d3faSBarry Smith 1695a46d3faSBarry Smith options = lu->options; 1705a46d3faSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"SuperLU run parameters:\n");CHKERRQ(ierr); 1715a46d3faSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Equil: %s\n",(options.Equil != NO) ? "YES": "NO");CHKERRQ(ierr); 1725a46d3faSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," ColPerm: %D\n",options.ColPerm);CHKERRQ(ierr); 1735a46d3faSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," IterRefine: %D\n",options.IterRefine);CHKERRQ(ierr); 1745a46d3faSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," SymmetricMode: %s\n",(options.SymmetricMode != NO) ? "YES": "NO");CHKERRQ(ierr); 1755a46d3faSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," DiagPivotThresh: %g\n",options.DiagPivotThresh);CHKERRQ(ierr); 1765a46d3faSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," PivotGrowth: %s\n",(options.PivotGrowth != NO) ? "YES": "NO");CHKERRQ(ierr); 1775a46d3faSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," ConditionNumber: %s\n",(options.ConditionNumber != NO) ? "YES": "NO");CHKERRQ(ierr); 1785a46d3faSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," RowPerm: %D\n",options.RowPerm);CHKERRQ(ierr); 1795a46d3faSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," ReplaceTinyPivot: %s\n",(options.ReplaceTinyPivot != NO) ? "YES": "NO");CHKERRQ(ierr); 1805a46d3faSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," PrintStat: %s\n",(options.PrintStat != NO) ? "YES": "NO");CHKERRQ(ierr); 1815a46d3faSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," lwork: %D\n",lu->lwork);CHKERRQ(ierr); 1825a46d3faSBarry Smith 1835a46d3faSBarry Smith PetscFunctionReturn(0); 1845a46d3faSBarry Smith } 1855a46d3faSBarry Smith 1865a46d3faSBarry Smith /* 1875a46d3faSBarry Smith These are the methods provided to REPLACE the corresponding methods of the 1885a46d3faSBarry Smith SeqAIJ matrix class. Other methods of SeqAIJ are not replaced 1895a46d3faSBarry Smith */ 1905a46d3faSBarry Smith #undef __FUNCT__ 1915a46d3faSBarry Smith #define __FUNCT__ "MatLUFactorNumeric_SuperLU" 1925a46d3faSBarry Smith PetscErrorCode MatLUFactorNumeric_SuperLU(Mat A,MatFactorInfo *info,Mat *F) 1935a46d3faSBarry Smith { 1945a46d3faSBarry Smith Mat_SeqAIJ *aa = (Mat_SeqAIJ*)(A)->data; 1955a46d3faSBarry Smith Mat_SuperLU *lu = (Mat_SuperLU*)(*F)->spptr; 1965a46d3faSBarry Smith PetscErrorCode ierr; 1975a46d3faSBarry Smith PetscInt sinfo; 1985a46d3faSBarry Smith SuperLUStat_t stat; 1995a46d3faSBarry Smith PetscReal ferr, berr; 2005a46d3faSBarry Smith NCformat *Ustore; 2015a46d3faSBarry Smith SCformat *Lstore; 2025a46d3faSBarry Smith 2035a46d3faSBarry Smith PetscFunctionBegin; 2045a46d3faSBarry Smith if (lu->flg == SAME_NONZERO_PATTERN){ /* successing numerical factorization */ 2055a46d3faSBarry Smith lu->options.Fact = SamePattern; 2065a46d3faSBarry Smith /* Ref: ~SuperLU_3.0/EXAMPLE/dlinsolx2.c */ 2075a46d3faSBarry Smith Destroy_SuperMatrix_Store(&lu->A); 2085a46d3faSBarry Smith if ( lu->lwork >= 0 ) { 2095a46d3faSBarry Smith Destroy_SuperNode_Matrix(&lu->L); 2105a46d3faSBarry Smith Destroy_CompCol_Matrix(&lu->U); 2115a46d3faSBarry Smith lu->options.Fact = SamePattern; 2125a46d3faSBarry Smith } 2135a46d3faSBarry Smith } 2145a46d3faSBarry Smith 2155a46d3faSBarry Smith /* Create the SuperMatrix for lu->A=A^T: 2165a46d3faSBarry Smith Since SuperLU likes column-oriented matrices,we pass it the transpose, 2175a46d3faSBarry Smith and then solve A^T X = B in MatSolve(). */ 2185a46d3faSBarry Smith #if defined(PETSC_USE_COMPLEX) 2195a46d3faSBarry Smith zCreate_CompCol_Matrix(&lu->A,A->cmap.n,A->rmap.n,aa->nz,(doublecomplex*)aa->a,aa->j,aa->i, 2205a46d3faSBarry Smith SLU_NC,SLU_Z,SLU_GE); 2215a46d3faSBarry Smith #else 2225a46d3faSBarry Smith dCreate_CompCol_Matrix(&lu->A,A->cmap.n,A->rmap.n,aa->nz,aa->a,aa->j,aa->i, 2235a46d3faSBarry Smith SLU_NC,SLU_D,SLU_GE); 2245a46d3faSBarry Smith #endif 2255a46d3faSBarry Smith 2265a46d3faSBarry Smith /* Initialize the statistics variables. */ 2275a46d3faSBarry Smith StatInit(&stat); 2285a46d3faSBarry Smith 2295a46d3faSBarry Smith /* Numerical factorization */ 2305a46d3faSBarry Smith lu->B.ncol = 0; /* Indicate not to solve the system */ 2315a46d3faSBarry Smith #if defined(PETSC_USE_COMPLEX) 2325a46d3faSBarry Smith zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 2335a46d3faSBarry Smith &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 2345a46d3faSBarry Smith &lu->mem_usage, &stat, &sinfo); 2355a46d3faSBarry Smith #else 2365a46d3faSBarry Smith dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 2375a46d3faSBarry Smith &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 2385a46d3faSBarry Smith &lu->mem_usage, &stat, &sinfo); 2395a46d3faSBarry Smith #endif 2405a46d3faSBarry Smith if ( !sinfo || sinfo == lu->A.ncol+1 ) { 2415a46d3faSBarry Smith if ( lu->options.PivotGrowth ) 2425a46d3faSBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," Recip. pivot growth = %e\n", lu->rpg); 2435a46d3faSBarry Smith if ( lu->options.ConditionNumber ) 2445a46d3faSBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," Recip. condition number = %e\n", lu->rcond); 2455a46d3faSBarry Smith } else if ( sinfo > 0 ){ 2465a46d3faSBarry Smith if ( lu->lwork == -1 ) { 2475a46d3faSBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," ** Estimated memory: %D bytes\n", sinfo - lu->A.ncol); 2485a46d3faSBarry Smith } else { 2495a46d3faSBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," Warning: gssvx() returns info %D\n",sinfo); 2505a46d3faSBarry Smith } 2515a46d3faSBarry Smith } else { /* sinfo < 0 */ 2525a46d3faSBarry Smith SETERRQ2(PETSC_ERR_LIB, "info = %D, the %D-th argument in gssvx() had an illegal value", sinfo,-sinfo); 2535a46d3faSBarry Smith } 2545a46d3faSBarry Smith 2555a46d3faSBarry Smith if ( lu->options.PrintStat ) { 2565a46d3faSBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF,"MatLUFactorNumeric_SuperLU():\n"); 2575a46d3faSBarry Smith StatPrint(&stat); 2585a46d3faSBarry Smith Lstore = (SCformat *) lu->L.Store; 2595a46d3faSBarry Smith Ustore = (NCformat *) lu->U.Store; 2605a46d3faSBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," No of nonzeros in factor L = %D\n", Lstore->nnz); 2615a46d3faSBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," No of nonzeros in factor U = %D\n", Ustore->nnz); 2625a46d3faSBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," No of nonzeros in L+U = %D\n", Lstore->nnz + Ustore->nnz - lu->A.ncol); 2635a46d3faSBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," L\\U MB %.3f\ttotal MB needed %.3f\texpansions %D\n", 2645a46d3faSBarry Smith lu->mem_usage.for_lu/1e6, lu->mem_usage.total_needed/1e6, 2655a46d3faSBarry Smith lu->mem_usage.expansions); 2665a46d3faSBarry Smith } 2675a46d3faSBarry Smith StatFree(&stat); 2685a46d3faSBarry Smith 2695a46d3faSBarry Smith lu->flg = SAME_NONZERO_PATTERN; 2705a46d3faSBarry Smith PetscFunctionReturn(0); 2715a46d3faSBarry Smith } 2725a46d3faSBarry Smith 27314b396bbSKris Buschelman #undef __FUNCT__ 274f0c56d0fSKris Buschelman #define __FUNCT__ "MatDestroy_SuperLU" 275dfbe8321SBarry Smith PetscErrorCode MatDestroy_SuperLU(Mat A) 27614b396bbSKris Buschelman { 277dfbe8321SBarry Smith PetscErrorCode ierr; 278f0c56d0fSKris Buschelman Mat_SuperLU *lu=(Mat_SuperLU*)A->spptr; 27914b396bbSKris Buschelman 28014b396bbSKris Buschelman PetscFunctionBegin; 28175af56d4SHong Zhang if (lu->CleanUpSuperLU) { /* Free the SuperLU datastructures */ 28275af56d4SHong Zhang Destroy_SuperMatrix_Store(&lu->A); 28375af56d4SHong Zhang Destroy_SuperMatrix_Store(&lu->B); 28475af56d4SHong Zhang Destroy_SuperMatrix_Store(&lu->X); 2859ce50919SHong Zhang 2869ce50919SHong Zhang ierr = PetscFree(lu->etree);CHKERRQ(ierr); 2879ce50919SHong Zhang ierr = PetscFree(lu->perm_r);CHKERRQ(ierr); 2889ce50919SHong Zhang ierr = PetscFree(lu->perm_c);CHKERRQ(ierr); 2899ce50919SHong Zhang ierr = PetscFree(lu->R);CHKERRQ(ierr); 2909ce50919SHong Zhang ierr = PetscFree(lu->C);CHKERRQ(ierr); 29175af56d4SHong Zhang if ( lu->lwork >= 0 ) { 292fb3e25aaSKris Buschelman Destroy_SuperNode_Matrix(&lu->L); 293fb3e25aaSKris Buschelman Destroy_CompCol_Matrix(&lu->U); 29475af56d4SHong Zhang } 295fb3e25aaSKris Buschelman } 296ceb03754SKris Buschelman ierr = MatConvert_SuperLU_SeqAIJ(A,MATSEQAIJ,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr); 297fb3e25aaSKris Buschelman ierr = (*A->ops->destroy)(A);CHKERRQ(ierr); 29814b396bbSKris Buschelman PetscFunctionReturn(0); 29914b396bbSKris Buschelman } 30014b396bbSKris Buschelman 30114b396bbSKris Buschelman #undef __FUNCT__ 302f0c56d0fSKris Buschelman #define __FUNCT__ "MatView_SuperLU" 303dfbe8321SBarry Smith PetscErrorCode MatView_SuperLU(Mat A,PetscViewer viewer) 30414b396bbSKris Buschelman { 305dfbe8321SBarry Smith PetscErrorCode ierr; 30632077d6dSBarry Smith PetscTruth iascii; 30714b396bbSKris Buschelman PetscViewerFormat format; 308f0c56d0fSKris Buschelman Mat_SuperLU *lu=(Mat_SuperLU*)(A->spptr); 30914b396bbSKris Buschelman 31014b396bbSKris Buschelman PetscFunctionBegin; 31114b396bbSKris Buschelman ierr = (*lu->MatView)(A,viewer);CHKERRQ(ierr); 31214b396bbSKris Buschelman 31332077d6dSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 31432077d6dSBarry Smith if (iascii) { 31514b396bbSKris Buschelman ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 3162f59403fSHong Zhang if (format == PETSC_VIEWER_ASCII_INFO) { 317f0c56d0fSKris Buschelman ierr = MatFactorInfo_SuperLU(A,viewer);CHKERRQ(ierr); 31814b396bbSKris Buschelman } 31914b396bbSKris Buschelman } 32014b396bbSKris Buschelman PetscFunctionReturn(0); 32114b396bbSKris Buschelman } 32214b396bbSKris Buschelman 32314b396bbSKris Buschelman #undef __FUNCT__ 324f0c56d0fSKris Buschelman #define __FUNCT__ "MatAssemblyEnd_SuperLU" 325dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_SuperLU(Mat A,MatAssemblyType mode) { 326dfbe8321SBarry Smith PetscErrorCode ierr; 327f0c56d0fSKris Buschelman Mat_SuperLU *lu=(Mat_SuperLU*)(A->spptr); 32814b396bbSKris Buschelman 32914b396bbSKris Buschelman PetscFunctionBegin; 33014b396bbSKris Buschelman ierr = (*lu->MatAssemblyEnd)(A,mode);CHKERRQ(ierr); 331b0592601SKris Buschelman lu->MatLUFactorSymbolic = A->ops->lufactorsymbolic; 332f0c56d0fSKris Buschelman A->ops->lufactorsymbolic = MatLUFactorSymbolic_SuperLU; 33314b396bbSKris Buschelman PetscFunctionReturn(0); 33414b396bbSKris Buschelman } 33514b396bbSKris Buschelman 33614b396bbSKris Buschelman 33714b396bbSKris Buschelman #undef __FUNCT__ 338c29e884eSHong Zhang #define __FUNCT__ "MatSolve_SuperLU_Private" 339c29e884eSHong Zhang PetscErrorCode MatSolve_SuperLU_Private(Mat A,Vec b,Vec x) 34014b396bbSKris Buschelman { 341f0c56d0fSKris Buschelman Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr; 34275af56d4SHong Zhang PetscScalar *barray,*xarray; 343dfbe8321SBarry Smith PetscErrorCode ierr; 344da7a1d00SHong Zhang PetscInt info,i; 34575af56d4SHong Zhang SuperLUStat_t stat; 346da7a1d00SHong Zhang PetscReal ferr,berr; 34714b396bbSKris Buschelman 34814b396bbSKris Buschelman PetscFunctionBegin; 349937a6b0eSHong Zhang if ( lu->lwork == -1 ) { 350937a6b0eSHong Zhang PetscFunctionReturn(0); 351937a6b0eSHong Zhang } 35275af56d4SHong Zhang lu->B.ncol = 1; /* Set the number of right-hand side */ 35375af56d4SHong Zhang ierr = VecGetArray(b,&barray);CHKERRQ(ierr); 35475af56d4SHong Zhang ierr = VecGetArray(x,&xarray);CHKERRQ(ierr); 3555fe6150dSHong Zhang 3565fe6150dSHong Zhang #if defined(PETSC_USE_COMPLEX) 3575fe6150dSHong Zhang ((DNformat*)lu->B.Store)->nzval = (doublecomplex*)barray; 3585fe6150dSHong Zhang ((DNformat*)lu->X.Store)->nzval = (doublecomplex*)xarray; 3595fe6150dSHong Zhang #else 3605fe6150dSHong Zhang ((DNformat*)lu->B.Store)->nzval = barray; 36175af56d4SHong Zhang ((DNformat*)lu->X.Store)->nzval = xarray; 3625fe6150dSHong Zhang #endif 36375af56d4SHong Zhang 36475af56d4SHong Zhang /* Initialize the statistics variables. */ 36575af56d4SHong Zhang StatInit(&stat); 36675af56d4SHong Zhang 36775af56d4SHong Zhang lu->options.Fact = FACTORED; /* Indicate the factored form of A is supplied. */ 3688914a3f7SHong Zhang #if defined(PETSC_USE_COMPLEX) 3698914a3f7SHong Zhang zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 3708914a3f7SHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 3718914a3f7SHong Zhang &lu->mem_usage, &stat, &info); 3728914a3f7SHong Zhang #else 37375af56d4SHong Zhang dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 37475af56d4SHong Zhang &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 37575af56d4SHong Zhang &lu->mem_usage, &stat, &info); 3768914a3f7SHong Zhang #endif 37775af56d4SHong Zhang ierr = VecRestoreArray(b,&barray);CHKERRQ(ierr); 37875af56d4SHong Zhang ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr); 37975af56d4SHong Zhang 380958c9bccSBarry Smith if ( !info || info == lu->A.ncol+1 ) { 38175af56d4SHong Zhang if ( lu->options.IterRefine ) { 3828914a3f7SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF,"Iterative Refinement:\n"); 3838914a3f7SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," %8s%8s%16s%16s\n", "rhs", "Steps", "FERR", "BERR"); 38475af56d4SHong Zhang for (i = 0; i < 1; ++i) 3858914a3f7SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," %8d%8d%16e%16e\n", i+1, stat.RefineSteps, ferr, berr); 38675af56d4SHong Zhang } 3878914a3f7SHong Zhang } else if ( info > 0 ){ 3888914a3f7SHong Zhang if ( lu->lwork == -1 ) { 38977431f27SBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," ** Estimated memory: %D bytes\n", info - lu->A.ncol); 3908914a3f7SHong Zhang } else { 39177431f27SBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," Warning: gssvx() returns info %D\n",info); 3928914a3f7SHong Zhang } 3938914a3f7SHong Zhang } else if (info < 0){ 39477431f27SBarry Smith SETERRQ2(PETSC_ERR_LIB, "info = %D, the %D-th argument in gssvx() had an illegal value", info,-info); 39514b396bbSKris Buschelman } 39614b396bbSKris Buschelman 3978914a3f7SHong Zhang if ( lu->options.PrintStat ) { 3988914a3f7SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF,"MatSolve__SuperLU():\n"); 3998914a3f7SHong Zhang StatPrint(&stat); 4008914a3f7SHong Zhang } 40175af56d4SHong Zhang StatFree(&stat); 40275af56d4SHong Zhang PetscFunctionReturn(0); 40375af56d4SHong Zhang } 40414b396bbSKris Buschelman 40514b396bbSKris Buschelman #undef __FUNCT__ 406c29e884eSHong Zhang #define __FUNCT__ "MatSolve_SuperLU" 407c29e884eSHong Zhang PetscErrorCode MatSolve_SuperLU(Mat A,Vec b,Vec x) 408c29e884eSHong Zhang { 409c29e884eSHong Zhang Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr; 410c29e884eSHong Zhang PetscErrorCode ierr; 411c29e884eSHong Zhang 412c29e884eSHong Zhang PetscFunctionBegin; 413c29e884eSHong Zhang lu->options.Trans = TRANS; 414c29e884eSHong Zhang ierr = MatSolve_SuperLU_Private(A,b,x);CHKERRQ(ierr); 415c29e884eSHong Zhang PetscFunctionReturn(0); 416c29e884eSHong Zhang } 417c29e884eSHong Zhang 418c29e884eSHong Zhang #undef __FUNCT__ 419c7c1fe80SHong Zhang #define __FUNCT__ "MatSolveTranspose_SuperLU" 420c7c1fe80SHong Zhang PetscErrorCode MatSolveTranspose_SuperLU(Mat A,Vec b,Vec x) 421c7c1fe80SHong Zhang { 422c7c1fe80SHong Zhang Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr; 423c7c1fe80SHong Zhang PetscErrorCode ierr; 424c7c1fe80SHong Zhang 425c7c1fe80SHong Zhang PetscFunctionBegin; 426c7c1fe80SHong Zhang lu->options.Trans = NOTRANS; 427c29e884eSHong Zhang ierr = MatSolve_SuperLU_Private(A,b,x);CHKERRQ(ierr); 428c7c1fe80SHong Zhang PetscFunctionReturn(0); 429c7c1fe80SHong Zhang } 430c7c1fe80SHong Zhang 43114b396bbSKris Buschelman 43214b396bbSKris Buschelman /* 43314b396bbSKris Buschelman Note the r permutation is ignored 43414b396bbSKris Buschelman */ 43514b396bbSKris Buschelman #undef __FUNCT__ 436f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_SuperLU" 437dfbe8321SBarry Smith PetscErrorCode MatLUFactorSymbolic_SuperLU(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F) 43814b396bbSKris Buschelman { 43914b396bbSKris Buschelman Mat B; 440f0c56d0fSKris Buschelman Mat_SuperLU *lu; 4416849ba73SBarry Smith PetscErrorCode ierr; 4422a4c71feSBarry Smith PetscInt m=A->rmap.n,n=A->cmap.n,indx; 4439ce50919SHong Zhang PetscTruth flg; 4445ca28756SHong Zhang const char *colperm[]={"NATURAL","MMD_ATA","MMD_AT_PLUS_A","COLAMD"}; /* MY_PERMC - not supported by the petsc interface yet */ 4455ca28756SHong Zhang const char *iterrefine[]={"NOREFINE", "SINGLE", "DOUBLE", "EXTRA"}; 4465ca28756SHong Zhang const char *rowperm[]={"NOROWPERM", "LargeDiag"}; /* MY_PERMC - not supported by the petsc interface yet */ 44714b396bbSKris Buschelman 44814b396bbSKris Buschelman PetscFunctionBegin; 449*7adad957SLisandro Dalcin ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr); 4502a4c71feSBarry Smith ierr = MatSetSizes(B,A->rmap.n,A->cmap.n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 451*7adad957SLisandro Dalcin ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 452899d7b4fSKris Buschelman ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr); 453f0c56d0fSKris Buschelman 454f0c56d0fSKris Buschelman B->ops->lufactornumeric = MatLUFactorNumeric_SuperLU; 455f0c56d0fSKris Buschelman B->ops->solve = MatSolve_SuperLU; 456c7c1fe80SHong Zhang B->ops->solvetranspose = MatSolveTranspose_SuperLU; 45714b396bbSKris Buschelman B->factor = FACTOR_LU; 45894b7f48cSBarry Smith B->assembled = PETSC_TRUE; /* required by -ksp_view */ 45914b396bbSKris Buschelman 460f0c56d0fSKris Buschelman lu = (Mat_SuperLU*)(B->spptr); 4619ce50919SHong Zhang 4629ce50919SHong Zhang /* Set SuperLU options */ 4639ce50919SHong Zhang /* the default values for options argument: 4649ce50919SHong Zhang options.Fact = DOFACT; 4659ce50919SHong Zhang options.Equil = YES; 4669ce50919SHong Zhang options.ColPerm = COLAMD; 4679ce50919SHong Zhang options.DiagPivotThresh = 1.0; 4689ce50919SHong Zhang options.Trans = NOTRANS; 4699ce50919SHong Zhang options.IterRefine = NOREFINE; 4709ce50919SHong Zhang options.SymmetricMode = NO; 4719ce50919SHong Zhang options.PivotGrowth = NO; 4729ce50919SHong Zhang options.ConditionNumber = NO; 4739ce50919SHong Zhang options.PrintStat = YES; 4749ce50919SHong Zhang */ 4759ce50919SHong Zhang set_default_options(&lu->options); 4768914a3f7SHong Zhang /* equilibration causes error in solve(), thus not supported here. See dgssvx.c for possible reason. */ 4778914a3f7SHong Zhang lu->options.Equil = NO; 4789ce50919SHong Zhang lu->options.PrintStat = NO; 4798914a3f7SHong Zhang lu->lwork = 0; /* allocate space internally by system malloc */ 4809ce50919SHong Zhang 481*7adad957SLisandro Dalcin ierr = PetscOptionsBegin(((PetscObject)A)->comm,((PetscObject)A)->prefix,"SuperLU Options","Mat");CHKERRQ(ierr); 4828914a3f7SHong Zhang ierr = PetscOptionsEList("-mat_superlu_colperm","ColPerm","None",colperm,4,colperm[3],&indx,&flg);CHKERRQ(ierr); 4839ce50919SHong Zhang if (flg) {lu->options.ColPerm = (colperm_t)indx;} 4848914a3f7SHong Zhang ierr = PetscOptionsEList("-mat_superlu_iterrefine","IterRefine","None",iterrefine,4,iterrefine[0],&indx,&flg);CHKERRQ(ierr); 4859ce50919SHong Zhang if (flg) { lu->options.IterRefine = (IterRefine_t)indx;} 4864dc4c822SBarry Smith ierr = PetscOptionsTruth("-mat_superlu_symmetricmode","SymmetricMode","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr); 4879ce50919SHong Zhang if (flg) lu->options.SymmetricMode = YES; 4888914a3f7SHong Zhang ierr = PetscOptionsReal("-mat_superlu_diagpivotthresh","DiagPivotThresh","None",lu->options.DiagPivotThresh,&lu->options.DiagPivotThresh,PETSC_NULL);CHKERRQ(ierr); 4894dc4c822SBarry Smith ierr = PetscOptionsTruth("-mat_superlu_pivotgrowth","PivotGrowth","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr); 4909ce50919SHong Zhang if (flg) lu->options.PivotGrowth = YES; 4914dc4c822SBarry Smith ierr = PetscOptionsTruth("-mat_superlu_conditionnumber","ConditionNumber","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr); 4929ce50919SHong Zhang if (flg) lu->options.ConditionNumber = YES; 4938914a3f7SHong Zhang ierr = PetscOptionsEList("-mat_superlu_rowperm","rowperm","None",rowperm,2,rowperm[0],&indx,&flg);CHKERRQ(ierr); 4949ce50919SHong Zhang if (flg) {lu->options.RowPerm = (rowperm_t)indx;} 4954dc4c822SBarry Smith ierr = PetscOptionsTruth("-mat_superlu_replacetinypivot","ReplaceTinyPivot","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr); 4969ce50919SHong Zhang if (flg) lu->options.ReplaceTinyPivot = YES; 4974dc4c822SBarry Smith ierr = PetscOptionsTruth("-mat_superlu_printstat","PrintStat","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr); 4989ce50919SHong Zhang if (flg) lu->options.PrintStat = YES; 4998914a3f7SHong Zhang ierr = PetscOptionsInt("-mat_superlu_lwork","size of work array in bytes used by factorization","None",lu->lwork,&lu->lwork,PETSC_NULL);CHKERRQ(ierr); 5005fe6150dSHong Zhang if (lu->lwork > 0 ){ 5015fe6150dSHong Zhang ierr = PetscMalloc(lu->lwork,&lu->work);CHKERRQ(ierr); 5025fe6150dSHong Zhang } else if (lu->lwork != 0 && lu->lwork != -1){ 50377431f27SBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," Warning: lwork %D is not supported by SUPERLU. The default lwork=0 is used.\n",lu->lwork); 5048914a3f7SHong Zhang lu->lwork = 0; 5058914a3f7SHong Zhang } 5069ce50919SHong Zhang PetscOptionsEnd(); 5079ce50919SHong Zhang 50875af56d4SHong Zhang #ifdef SUPERLU2 5092ecb5974SKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatCreateNull","MatCreateNull_SuperLU", 510f0c56d0fSKris Buschelman (void(*)(void))MatCreateNull_SuperLU);CHKERRQ(ierr); 51175af56d4SHong Zhang #endif 51214b396bbSKris Buschelman 51375af56d4SHong Zhang /* Allocate spaces (notice sizes are for the transpose) */ 514da7a1d00SHong Zhang ierr = PetscMalloc(m*sizeof(PetscInt),&lu->etree);CHKERRQ(ierr); 515da7a1d00SHong Zhang ierr = PetscMalloc(n*sizeof(PetscInt),&lu->perm_r);CHKERRQ(ierr); 516da7a1d00SHong Zhang ierr = PetscMalloc(m*sizeof(PetscInt),&lu->perm_c);CHKERRQ(ierr); 517da7a1d00SHong Zhang ierr = PetscMalloc(n*sizeof(PetscInt),&lu->R);CHKERRQ(ierr); 518da7a1d00SHong Zhang ierr = PetscMalloc(m*sizeof(PetscInt),&lu->C);CHKERRQ(ierr); 51975af56d4SHong Zhang 5209ce50919SHong Zhang /* create rhs and solution x without allocate space for .Store */ 5215fe6150dSHong Zhang #if defined(PETSC_USE_COMPLEX) 522937a6b0eSHong Zhang zCreate_Dense_Matrix(&lu->B, m, 1, PETSC_NULL, m, SLU_DN, SLU_Z, SLU_GE); 523937a6b0eSHong Zhang zCreate_Dense_Matrix(&lu->X, m, 1, PETSC_NULL, m, SLU_DN, SLU_Z, SLU_GE); 5245fe6150dSHong Zhang #else 52575af56d4SHong Zhang dCreate_Dense_Matrix(&lu->B, m, 1, PETSC_NULL, m, SLU_DN, SLU_D, SLU_GE); 52675af56d4SHong Zhang dCreate_Dense_Matrix(&lu->X, m, 1, PETSC_NULL, m, SLU_DN, SLU_D, SLU_GE); 5275fe6150dSHong Zhang #endif 52814b396bbSKris Buschelman 52914b396bbSKris Buschelman lu->flg = DIFFERENT_NONZERO_PATTERN; 530e740cb95SKris Buschelman lu->CleanUpSuperLU = PETSC_TRUE; 531899d7b4fSKris Buschelman 532899d7b4fSKris Buschelman *F = B; 53338f2d2fdSLisandro Dalcin ierr = PetscLogObjectMemory(B,(A->rmap.n+A->cmap.n)*sizeof(PetscInt));CHKERRQ(ierr); 53414b396bbSKris Buschelman PetscFunctionReturn(0); 53514b396bbSKris Buschelman } 53614b396bbSKris Buschelman 537f0c56d0fSKris Buschelman 538f0c56d0fSKris Buschelman #undef __FUNCT__ 539f0c56d0fSKris Buschelman #define __FUNCT__ "MatDuplicate_SuperLU" 540dfbe8321SBarry Smith PetscErrorCode MatDuplicate_SuperLU(Mat A, MatDuplicateOption op, Mat *M) { 541da7a1d00SHong Zhang PetscErrorCode ierr; 5428f340917SKris Buschelman Mat_SuperLU *lu=(Mat_SuperLU *)A->spptr; 5438f340917SKris Buschelman 54414b396bbSKris Buschelman PetscFunctionBegin; 5458f340917SKris Buschelman ierr = (*lu->MatDuplicate)(A,op,M);CHKERRQ(ierr); 5463f953163SKris Buschelman ierr = PetscMemcpy((*M)->spptr,lu,sizeof(Mat_SuperLU));CHKERRQ(ierr); 54714b396bbSKris Buschelman PetscFunctionReturn(0); 54814b396bbSKris Buschelman } 54914b396bbSKris Buschelman 550b0592601SKris Buschelman 55124b6179bSKris Buschelman /*MC 552fafad747SKris Buschelman MATSUPERLU - MATSUPERLU = "superlu" - A matrix type providing direct solvers (LU) for sequential matrices 55324b6179bSKris Buschelman via the external package SuperLU. 55424b6179bSKris Buschelman 55524b6179bSKris Buschelman If SuperLU is installed (see the manual for 55624b6179bSKris Buschelman instructions on how to declare the existence of external packages), 55724b6179bSKris Buschelman a matrix type can be constructed which invokes SuperLU solvers. 55824b6179bSKris Buschelman After calling MatCreate(...,A), simply call MatSetType(A,MATSUPERLU). 55924b6179bSKris Buschelman 56007d81ca4SBarry Smith This matrix inherits from MATSEQAIJ. As a result, MatSeqAIJSetPreallocation() is 56128b08bd3SKris Buschelman supported for this matrix type. One can also call MatConvert for an inplace conversion to or from 56228b08bd3SKris Buschelman the MATSEQAIJ type without data copy. 56324b6179bSKris Buschelman 56424b6179bSKris Buschelman Options Database Keys: 5650bad9183SKris Buschelman + -mat_type superlu - sets the matrix type to "superlu" during a call to MatSetFromOptions() 56651678082SBarry Smith . -mat_superlu_ordering <0,1,2,3> - 0: natural ordering, 56724b6179bSKris Buschelman 1: MMD applied to A'*A, 56824b6179bSKris Buschelman 2: MMD applied to A'+A, 56924b6179bSKris Buschelman 3: COLAMD, approximate minimum degree column ordering 57051678082SBarry Smith . -mat_superlu_iterrefine - have SuperLU do iterative refinement after the triangular solve 57151678082SBarry Smith choices: NOREFINE, SINGLE, DOUBLE, EXTRA; default is NOREFINE 57251678082SBarry Smith - -mat_superlu_printstat - print SuperLU statistics about the factorization 57324b6179bSKris Buschelman 57424b6179bSKris Buschelman Level: beginner 57524b6179bSKris Buschelman 57624b6179bSKris Buschelman .seealso: PCLU 57724b6179bSKris Buschelman M*/ 57824b6179bSKris Buschelman 5795a46d3faSBarry Smith /* 5805a46d3faSBarry Smith Constructor for the new derived matrix class. It simply creates the base 5815a46d3faSBarry Smith matrix class and then adds the additional information/methods needed by SuperLU. 5825a46d3faSBarry Smith */ 583b0592601SKris Buschelman EXTERN_C_BEGIN 584b0592601SKris Buschelman #undef __FUNCT__ 585f0c56d0fSKris Buschelman #define __FUNCT__ "MatCreate_SuperLU" 586be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_SuperLU(Mat A) 587dfbe8321SBarry Smith { 588dfbe8321SBarry Smith PetscErrorCode ierr; 589b0592601SKris Buschelman 590b0592601SKris Buschelman PetscFunctionBegin; 591b0592601SKris Buschelman ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr); 592ceb03754SKris Buschelman ierr = MatConvert_SeqAIJ_SuperLU(A,MATSUPERLU,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr); 59314b396bbSKris Buschelman PetscFunctionReturn(0); 59414b396bbSKris Buschelman } 59514b396bbSKris Buschelman EXTERN_C_END 596