169e15a41SShri Abhyankar 269e15a41SShri Abhyankar /* 369e15a41SShri Abhyankar Provides an interface to the KLUv1.2 sparse solver 469e15a41SShri Abhyankar 569e15a41SShri Abhyankar When build with PETSC_USE_64BIT_INDICES this will use SuiteSparse_long as the 669e15a41SShri Abhyankar integer type in KLU, otherwise it will use int. This means 769e15a41SShri Abhyankar all integers in this file are simply declared as PetscInt. Also it means 869e15a41SShri Abhyankar that KLU SuiteSparse_long version MUST be built with 64 bit integers when used. 969e15a41SShri Abhyankar 1069e15a41SShri Abhyankar */ 1169e15a41SShri Abhyankar #include <../src/mat/impls/aij/seq/aij.h> 1269e15a41SShri Abhyankar 1369e15a41SShri Abhyankar #if defined(PETSC_USE_64BIT_INDICES) 1469e15a41SShri Abhyankar #define klu_K_defaults klu_l_defaults 1569e15a41SShri Abhyankar #define klu_K_analyze klu_l_analyze 1669e15a41SShri Abhyankar #define klu_K_analyze_given klu_l_analyze_given 1769e15a41SShri Abhyankar #define klu_K_free_symbolic klu_l_free_symbolic 1869e15a41SShri Abhyankar #define klu_K_free_numeric klu_l_free_numeric 1969e15a41SShri Abhyankar #define klu_K_common klu_l_common 2069e15a41SShri Abhyankar #define klu_K_symbolic klu_l_symbolic 2169e15a41SShri Abhyankar #define klu_K_numeric klu_l_numeric 2269e15a41SShri Abhyankar #if defined(PETSC_USE_COMPLEX) 2369e15a41SShri Abhyankar #define klu_K_factor klu_zl_factor 2469e15a41SShri Abhyankar #define klu_K_solve klu_zl_solve 2569e15a41SShri Abhyankar #define klu_K_tsolve klu_zl_tsolve 2669e15a41SShri Abhyankar #define klu_K_refactor klu_zl_refactor 2769e15a41SShri Abhyankar #define klu_K_sort klu_zl_sort 2869e15a41SShri Abhyankar #define klu_K_flops klu_zl_flops 2969e15a41SShri Abhyankar #define klu_K_rgrowth klu_zl_rgrowth 3069e15a41SShri Abhyankar #define klu_K_condest klu_zl_condest 3169e15a41SShri Abhyankar #define klu_K_rcond klu_zl_rcond 3269e15a41SShri Abhyankar #define klu_K_scale klu_zl_scale 3369e15a41SShri Abhyankar #else 3469e15a41SShri Abhyankar #define klu_K_factor klu_l_factor 3569e15a41SShri Abhyankar #define klu_K_solve klu_l_solve 3669e15a41SShri Abhyankar #define klu_K_tsolve klu_l_tsolve 3769e15a41SShri Abhyankar #define klu_K_refactor klu_l_refactor 3869e15a41SShri Abhyankar #define klu_K_sort klu_l_sort 3969e15a41SShri Abhyankar #define klu_K_flops klu_l_flops 4069e15a41SShri Abhyankar #define klu_K_rgrowth klu_l_rgrowth 4169e15a41SShri Abhyankar #define klu_K_condest klu_l_condest 4269e15a41SShri Abhyankar #define klu_K_rcond klu_l_rcond 4369e15a41SShri Abhyankar #define klu_K_scale klu_l_scale 4469e15a41SShri Abhyankar #endif 4569e15a41SShri Abhyankar #else 4669e15a41SShri Abhyankar #define klu_K_defaults klu_defaults 4769e15a41SShri Abhyankar #define klu_K_analyze klu_analyze 4869e15a41SShri Abhyankar #define klu_K_analyze_given klu_analyze_given 4969e15a41SShri Abhyankar #define klu_K_free_symbolic klu_free_symbolic 5069e15a41SShri Abhyankar #define klu_K_free_numeric klu_free_numeric 5169e15a41SShri Abhyankar #define klu_K_common klu_common 5269e15a41SShri Abhyankar #define klu_K_symbolic klu_symbolic 5369e15a41SShri Abhyankar #define klu_K_numeric klu_numeric 5469e15a41SShri Abhyankar #if defined(PETSC_USE_COMPLEX) 5569e15a41SShri Abhyankar #define klu_K_factor klu_z_factor 5669e15a41SShri Abhyankar #define klu_K_solve klu_z_solve 5769e15a41SShri Abhyankar #define klu_K_tsolve klu_z_tsolve 5869e15a41SShri Abhyankar #define klu_K_refactor klu_z_refactor 5969e15a41SShri Abhyankar #define klu_K_sort klu_z_sort 6069e15a41SShri Abhyankar #define klu_K_flops klu_z_flops 6169e15a41SShri Abhyankar #define klu_K_rgrowth klu_z_rgrowth 6269e15a41SShri Abhyankar #define klu_K_condest klu_z_condest 6369e15a41SShri Abhyankar #define klu_K_rcond klu_z_rcond 6469e15a41SShri Abhyankar #define klu_K_scale klu_z_scale 6569e15a41SShri Abhyankar #else 6669e15a41SShri Abhyankar #define klu_K_factor klu_factor 6769e15a41SShri Abhyankar #define klu_K_solve klu_solve 6869e15a41SShri Abhyankar #define klu_K_tsolve klu_tsolve 6969e15a41SShri Abhyankar #define klu_K_refactor klu_refactor 7069e15a41SShri Abhyankar #define klu_K_sort klu_sort 7169e15a41SShri Abhyankar #define klu_K_flops klu_flops 7269e15a41SShri Abhyankar #define klu_K_rgrowth klu_rgrowth 7369e15a41SShri Abhyankar #define klu_K_condest klu_condest 7469e15a41SShri Abhyankar #define klu_K_rcond klu_rcond 7569e15a41SShri Abhyankar #define klu_K_scale klu_scale 7669e15a41SShri Abhyankar #endif 7769e15a41SShri Abhyankar #endif 7869e15a41SShri Abhyankar 7969e15a41SShri Abhyankar #define SuiteSparse_long long long 8069e15a41SShri Abhyankar #define SuiteSparse_long_max LONG_LONG_MAX 8169e15a41SShri Abhyankar #define SuiteSparse_long_id "%lld" 8269e15a41SShri Abhyankar 8369e15a41SShri Abhyankar EXTERN_C_BEGIN 8469e15a41SShri Abhyankar #include <klu.h> 8569e15a41SShri Abhyankar EXTERN_C_END 8669e15a41SShri Abhyankar 8769e15a41SShri Abhyankar static const char *KluOrderingTypes[] = {"AMD","COLAMD","PETSC"}; 8869e15a41SShri Abhyankar static const char *scale[] ={"NONE","SUM","MAX"}; 8969e15a41SShri Abhyankar 9069e15a41SShri Abhyankar typedef struct { 9169e15a41SShri Abhyankar klu_K_common Common; 9269e15a41SShri Abhyankar klu_K_symbolic *Symbolic; 9369e15a41SShri Abhyankar klu_K_numeric *Numeric; 9469e15a41SShri Abhyankar PetscInt *perm_c,*perm_r; 9569e15a41SShri Abhyankar MatStructure flg; 9669e15a41SShri Abhyankar PetscBool PetscMatOrdering; 9769e15a41SShri Abhyankar PetscBool CleanUpKLU; 9869e15a41SShri Abhyankar } Mat_KLU; 9969e15a41SShri Abhyankar 10069e15a41SShri Abhyankar #undef __FUNCT__ 10169e15a41SShri Abhyankar #define __FUNCT__ "MatDestroy_KLU" 10269e15a41SShri Abhyankar static PetscErrorCode MatDestroy_KLU(Mat A) 10369e15a41SShri Abhyankar { 10469e15a41SShri Abhyankar PetscErrorCode ierr; 105*569c4379SBarry Smith Mat_KLU *lu=(Mat_KLU*)A->data; 10669e15a41SShri Abhyankar 10769e15a41SShri Abhyankar PetscFunctionBegin; 108*569c4379SBarry Smith if (lu->CleanUpKLU) { 10969e15a41SShri Abhyankar klu_K_free_symbolic(&lu->Symbolic,&lu->Common); 11069e15a41SShri Abhyankar klu_K_free_numeric(&lu->Numeric,&lu->Common); 111828fd075SShri Abhyankar ierr = PetscFree2(lu->perm_r,lu->perm_c);CHKERRQ(ierr); 11269e15a41SShri Abhyankar } 113*569c4379SBarry Smith ierr = PetscFree(A->data);CHKERRQ(ierr); 11469e15a41SShri Abhyankar PetscFunctionReturn(0); 11569e15a41SShri Abhyankar } 11669e15a41SShri Abhyankar 11769e15a41SShri Abhyankar #undef __FUNCT__ 11869e15a41SShri Abhyankar #define __FUNCT__ "MatSolveTranspose_KLU" 11969e15a41SShri Abhyankar static PetscErrorCode MatSolveTranspose_KLU(Mat A,Vec b,Vec x) 12069e15a41SShri Abhyankar { 121*569c4379SBarry Smith Mat_KLU *lu = (Mat_KLU*)A->data; 12269e15a41SShri Abhyankar PetscScalar *xa; 12369e15a41SShri Abhyankar PetscErrorCode ierr; 12469e15a41SShri Abhyankar PetscInt status; 12569e15a41SShri Abhyankar 12669e15a41SShri Abhyankar PetscFunctionBegin; 12769e15a41SShri Abhyankar /* KLU uses a column major format, solve Ax = b by klu_*_solve */ 12869e15a41SShri Abhyankar /* ----------------------------------*/ 12969e15a41SShri Abhyankar ierr = VecCopy(b,x); /* klu_solve stores the solution in rhs */ 13069e15a41SShri Abhyankar ierr = VecGetArray(x,&xa); 13169e15a41SShri Abhyankar status = klu_K_solve(lu->Symbolic,lu->Numeric,A->rmap->n,1,(PetscReal*)xa,&lu->Common); 132d9ca1df4SBarry Smith if (status != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"KLU Solve failed"); 13369e15a41SShri Abhyankar ierr = VecRestoreArray(x,&xa);CHKERRQ(ierr); 13469e15a41SShri Abhyankar PetscFunctionReturn(0); 13569e15a41SShri Abhyankar } 13669e15a41SShri Abhyankar 13769e15a41SShri Abhyankar #undef __FUNCT__ 13869e15a41SShri Abhyankar #define __FUNCT__ "MatSolve_KLU" 13969e15a41SShri Abhyankar static PetscErrorCode MatSolve_KLU(Mat A,Vec b,Vec x) 14069e15a41SShri Abhyankar { 141*569c4379SBarry Smith Mat_KLU *lu = (Mat_KLU*)A->data; 14269e15a41SShri Abhyankar PetscScalar *xa; 14369e15a41SShri Abhyankar PetscErrorCode ierr; 14469e15a41SShri Abhyankar PetscInt status; 14569e15a41SShri Abhyankar 14669e15a41SShri Abhyankar PetscFunctionBegin; 14769e15a41SShri Abhyankar /* KLU uses a column major format, solve A^Tx = b by klu_*_tsolve */ 14869e15a41SShri Abhyankar /* ----------------------------------*/ 14969e15a41SShri Abhyankar ierr = VecCopy(b,x); /* klu_solve stores the solution in rhs */ 15069e15a41SShri Abhyankar ierr = VecGetArray(x,&xa); 15169e15a41SShri Abhyankar #if defined(PETSC_USE_COMPLEX) 15269e15a41SShri Abhyankar PetscInt conj_solve=1; 15369e15a41SShri Abhyankar status = klu_K_tsolve(lu->Symbolic,lu->Numeric,A->rmap->n,1,(PetscReal*)xa,conj_solve,&lu->Common); /* conjugate solve */ 15469e15a41SShri Abhyankar #else 15569e15a41SShri Abhyankar status = klu_K_tsolve(lu->Symbolic,lu->Numeric,A->rmap->n,1,xa,&lu->Common); 15669e15a41SShri Abhyankar #endif 157d9ca1df4SBarry Smith if (status != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"KLU Solve failed"); 15869e15a41SShri Abhyankar ierr = VecRestoreArray(x,&xa);CHKERRQ(ierr); 15969e15a41SShri Abhyankar PetscFunctionReturn(0); 16069e15a41SShri Abhyankar } 16169e15a41SShri Abhyankar 16269e15a41SShri Abhyankar #undef __FUNCT__ 16369e15a41SShri Abhyankar #define __FUNCT__ "MatLUFactorNumeric_KLU" 16469e15a41SShri Abhyankar static PetscErrorCode MatLUFactorNumeric_KLU(Mat F,Mat A,const MatFactorInfo *info) 16569e15a41SShri Abhyankar { 166*569c4379SBarry Smith Mat_KLU *lu = (Mat_KLU*)(F)->data; 16769e15a41SShri Abhyankar Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 16869e15a41SShri Abhyankar PetscInt *ai = a->i,*aj=a->j; 16969e15a41SShri Abhyankar PetscScalar *av = a->a; 17069e15a41SShri Abhyankar 17169e15a41SShri Abhyankar PetscFunctionBegin; 17269e15a41SShri Abhyankar /* numeric factorization of A' */ 17369e15a41SShri Abhyankar /* ----------------------------*/ 17469e15a41SShri Abhyankar 17569e15a41SShri Abhyankar if (lu->flg == SAME_NONZERO_PATTERN && lu->Numeric) { 17669e15a41SShri Abhyankar klu_K_free_numeric(&lu->Numeric,&lu->Common); 17769e15a41SShri Abhyankar } 17869e15a41SShri Abhyankar lu->Numeric = klu_K_factor(ai,aj,(PetscReal*)av,lu->Symbolic,&lu->Common); 1796c4ed002SBarry Smith if(!lu->Numeric) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"KLU Numeric factorization failed"); 18069e15a41SShri Abhyankar 18169e15a41SShri Abhyankar lu->flg = SAME_NONZERO_PATTERN; 18269e15a41SShri Abhyankar lu->CleanUpKLU = PETSC_TRUE; 18369e15a41SShri Abhyankar F->ops->solve = MatSolve_KLU; 18469e15a41SShri Abhyankar F->ops->solvetranspose = MatSolveTranspose_KLU; 18569e15a41SShri Abhyankar PetscFunctionReturn(0); 18669e15a41SShri Abhyankar } 18769e15a41SShri Abhyankar 18869e15a41SShri Abhyankar #undef __FUNCT__ 18969e15a41SShri Abhyankar #define __FUNCT__ "MatLUFactorSymbolic_KLU" 19069e15a41SShri Abhyankar static PetscErrorCode MatLUFactorSymbolic_KLU(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info) 19169e15a41SShri Abhyankar { 19269e15a41SShri Abhyankar Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 193*569c4379SBarry Smith Mat_KLU *lu = (Mat_KLU*)(F->data); 19469e15a41SShri Abhyankar PetscErrorCode ierr; 19569e15a41SShri Abhyankar PetscInt i,*ai = a->i,*aj = a->j,m=A->rmap->n,n=A->cmap->n; 19669e15a41SShri Abhyankar const PetscInt *ra,*ca; 19769e15a41SShri Abhyankar 19869e15a41SShri Abhyankar PetscFunctionBegin; 19969e15a41SShri Abhyankar if (lu->PetscMatOrdering) { 20069e15a41SShri Abhyankar ierr = ISGetIndices(r,&ra);CHKERRQ(ierr); 20169e15a41SShri Abhyankar ierr = ISGetIndices(c,&ca);CHKERRQ(ierr); 20234c5654cSShri Abhyankar ierr = PetscMalloc2(m,&lu->perm_r,n,&lu->perm_c);CHKERRQ(ierr); 20369e15a41SShri Abhyankar /* we cannot simply memcpy on 64 bit archs */ 20469e15a41SShri Abhyankar for (i = 0; i < m; i++) lu->perm_r[i] = ra[i]; 20569e15a41SShri Abhyankar for (i = 0; i < n; i++) lu->perm_c[i] = ca[i]; 20669e15a41SShri Abhyankar ierr = ISRestoreIndices(r,&ra);CHKERRQ(ierr); 20769e15a41SShri Abhyankar ierr = ISRestoreIndices(c,&ca);CHKERRQ(ierr); 20869e15a41SShri Abhyankar } 20969e15a41SShri Abhyankar 21069e15a41SShri Abhyankar /* symbolic factorization of A' */ 21169e15a41SShri Abhyankar /* ---------------------------------------------------------------------- */ 21269e15a41SShri Abhyankar if (lu->PetscMatOrdering) { /* use Petsc ordering */ 21369e15a41SShri Abhyankar lu->Symbolic = klu_K_analyze_given(n,ai,aj,lu->perm_c,lu->perm_r,&lu->Common); 21469e15a41SShri Abhyankar } else { /* use klu internal ordering */ 21569e15a41SShri Abhyankar lu->Symbolic = klu_K_analyze(n,ai,aj,&lu->Common); 21669e15a41SShri Abhyankar } 2176c4ed002SBarry Smith if (!lu->Symbolic) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"KLU Symbolic Factorization failed"); 21869e15a41SShri Abhyankar 21969e15a41SShri Abhyankar lu->flg = DIFFERENT_NONZERO_PATTERN; 22069e15a41SShri Abhyankar lu->CleanUpKLU = PETSC_TRUE; 22169e15a41SShri Abhyankar (F)->ops->lufactornumeric = MatLUFactorNumeric_KLU; 22269e15a41SShri Abhyankar PetscFunctionReturn(0); 22369e15a41SShri Abhyankar } 22469e15a41SShri Abhyankar 22569e15a41SShri Abhyankar #undef __FUNCT__ 22669e15a41SShri Abhyankar #define __FUNCT__ "MatFactorInfo_KLU" 22769e15a41SShri Abhyankar static PetscErrorCode MatFactorInfo_KLU(Mat A,PetscViewer viewer) 22869e15a41SShri Abhyankar { 229*569c4379SBarry Smith Mat_KLU *lu= (Mat_KLU*)A->data; 23069e15a41SShri Abhyankar klu_K_numeric *Numeric=(klu_K_numeric*)lu->Numeric; 23169e15a41SShri Abhyankar PetscErrorCode ierr; 23269e15a41SShri Abhyankar 23369e15a41SShri Abhyankar PetscFunctionBegin; 23469e15a41SShri Abhyankar /* check if matrix is KLU type */ 23569e15a41SShri Abhyankar if (A->ops->solve != MatSolve_KLU) PetscFunctionReturn(0); 23669e15a41SShri Abhyankar 23769e15a41SShri Abhyankar ierr = PetscViewerASCIIPrintf(viewer,"KLU stats:\n");CHKERRQ(ierr); 23869e15a41SShri Abhyankar ierr = PetscViewerASCIIPrintf(viewer," Number of diagonal blocks: %d\n",Numeric->nblocks); 23969e15a41SShri Abhyankar ierr = PetscViewerASCIIPrintf(viewer," Total nonzeros=%d\n",Numeric->lnz+Numeric->unz);CHKERRQ(ierr); 24069e15a41SShri Abhyankar 24169e15a41SShri Abhyankar ierr = PetscViewerASCIIPrintf(viewer,"KLU runtime parameters:\n");CHKERRQ(ierr); 24269e15a41SShri Abhyankar 24369e15a41SShri Abhyankar /* Control parameters used by numeric factorization */ 24469e15a41SShri Abhyankar ierr = PetscViewerASCIIPrintf(viewer," Partial pivoting tolerance: %g\n",lu->Common.tol);CHKERRQ(ierr); 24569e15a41SShri Abhyankar /* BTF preordering */ 24669e15a41SShri Abhyankar ierr = PetscViewerASCIIPrintf(viewer," BTF preordering enabled: %d\n",lu->Common.btf);CHKERRQ(ierr); 24769e15a41SShri Abhyankar /* mat ordering */ 24869e15a41SShri Abhyankar if (!lu->PetscMatOrdering) { 24969e15a41SShri Abhyankar ierr = PetscViewerASCIIPrintf(viewer," Ordering: %s (not using the PETSc ordering)\n",KluOrderingTypes[(int)lu->Common.ordering]);CHKERRQ(ierr); 25069e15a41SShri Abhyankar } else { 25169e15a41SShri Abhyankar ierr = PetscViewerASCIIPrintf(viewer," Using PETSc ordering\n");CHKERRQ(ierr); 25269e15a41SShri Abhyankar } 25369e15a41SShri Abhyankar /* matrix row scaling */ 25469e15a41SShri Abhyankar ierr = PetscViewerASCIIPrintf(viewer, " Matrix row scaling: %s\n",scale[(int)lu->Common.scale]);CHKERRQ(ierr); 25569e15a41SShri Abhyankar PetscFunctionReturn(0); 25669e15a41SShri Abhyankar } 25769e15a41SShri Abhyankar 25869e15a41SShri Abhyankar #undef __FUNCT__ 25969e15a41SShri Abhyankar #define __FUNCT__ "MatView_KLU" 26069e15a41SShri Abhyankar static PetscErrorCode MatView_KLU(Mat A,PetscViewer viewer) 26169e15a41SShri Abhyankar { 26269e15a41SShri Abhyankar PetscErrorCode ierr; 26369e15a41SShri Abhyankar PetscBool iascii; 26469e15a41SShri Abhyankar PetscViewerFormat format; 26569e15a41SShri Abhyankar 26669e15a41SShri Abhyankar PetscFunctionBegin; 26769e15a41SShri Abhyankar ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 26869e15a41SShri Abhyankar if (iascii) { 26969e15a41SShri Abhyankar ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 27069e15a41SShri Abhyankar if (format == PETSC_VIEWER_ASCII_INFO) { 27169e15a41SShri Abhyankar ierr = MatFactorInfo_KLU(A,viewer);CHKERRQ(ierr); 27269e15a41SShri Abhyankar } 27369e15a41SShri Abhyankar } 27469e15a41SShri Abhyankar PetscFunctionReturn(0); 27569e15a41SShri Abhyankar } 27669e15a41SShri Abhyankar 27769e15a41SShri Abhyankar #undef __FUNCT__ 27869e15a41SShri Abhyankar #define __FUNCT__ "MatFactorGetSolverPackage_seqaij_klu" 27969e15a41SShri Abhyankar PetscErrorCode MatFactorGetSolverPackage_seqaij_klu(Mat A,const MatSolverPackage *type) 28069e15a41SShri Abhyankar { 28169e15a41SShri Abhyankar PetscFunctionBegin; 28269e15a41SShri Abhyankar *type = MATSOLVERKLU; 28369e15a41SShri Abhyankar PetscFunctionReturn(0); 28469e15a41SShri Abhyankar } 28569e15a41SShri Abhyankar 28669e15a41SShri Abhyankar 28769e15a41SShri Abhyankar /*MC 28869e15a41SShri Abhyankar MATSOLVERKLU = "klu" - A matrix type providing direct solvers (LU) for sequential matrices 28969e15a41SShri Abhyankar via the external package KLU. 29069e15a41SShri Abhyankar 291057b0833SShri Abhyankar ./configure --download-suitesparse to install PETSc to use KLU 29269e15a41SShri Abhyankar 293c2b89b5dSBarry Smith Use -pc_type lu -pc_factor_mat_solver_package klu to us this direct solver 294c2b89b5dSBarry Smith 29569e15a41SShri Abhyankar Consult KLU documentation for more information on the options database keys below. 29669e15a41SShri Abhyankar 29769e15a41SShri Abhyankar Options Database Keys: 29869e15a41SShri Abhyankar + -mat_klu_pivot_tol <0.001> - Partial pivoting tolerance 29969e15a41SShri Abhyankar . -mat_klu_use_btf <1> - Use BTF preordering 30069e15a41SShri Abhyankar . -mat_klu_ordering <AMD> - KLU reordering scheme to reduce fill-in (choose one of) AMD COLAMD PETSC 30169e15a41SShri Abhyankar - -mat_klu_row_scale <NONE> - Matrix row scaling (choose one of) NONE SUM MAX 30269e15a41SShri Abhyankar 303a364b7d2SBarry Smith Note: KLU is part of SuiteSparse http://faculty.cse.tamu.edu/davis/suitesparse.html 304a364b7d2SBarry Smith 30569e15a41SShri Abhyankar Level: beginner 30669e15a41SShri Abhyankar 30769e15a41SShri Abhyankar .seealso: PCLU, MATSOLVERUMFPACK, MATSOLVERCHOLMOD, PCFactorSetMatSolverPackage(), MatSolverPackage 30869e15a41SShri Abhyankar M*/ 30969e15a41SShri Abhyankar 31069e15a41SShri Abhyankar #undef __FUNCT__ 31169e15a41SShri Abhyankar #define __FUNCT__ "MatGetFactor_seqaij_klu" 312db87b0f2SBarry Smith PETSC_INTERN PetscErrorCode MatGetFactor_seqaij_klu(Mat A,MatFactorType ftype,Mat *F) 31369e15a41SShri Abhyankar { 31469e15a41SShri Abhyankar Mat B; 31569e15a41SShri Abhyankar Mat_KLU *lu; 31669e15a41SShri Abhyankar PetscErrorCode ierr; 31769e15a41SShri Abhyankar PetscInt m=A->rmap->n,n=A->cmap->n,idx,status; 31869e15a41SShri Abhyankar PetscBool flg; 31969e15a41SShri Abhyankar 32069e15a41SShri Abhyankar PetscFunctionBegin; 32169e15a41SShri Abhyankar /* Create the factorization matrix F */ 32269e15a41SShri Abhyankar ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr); 32369e15a41SShri Abhyankar ierr = MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,m,n);CHKERRQ(ierr); 324*569c4379SBarry Smith ierr = PetscStrallocpy("klu",&((PetscObject)B)->type_name);CHKERRQ(ierr); 325*569c4379SBarry Smith ierr = MatSetUp(B);CHKERRQ(ierr); 326*569c4379SBarry Smith 32769e15a41SShri Abhyankar ierr = PetscNewLog(B,&lu);CHKERRQ(ierr); 32869e15a41SShri Abhyankar 329*569c4379SBarry Smith B->data = lu; 330*569c4379SBarry Smith B->ops->getinfo = MatGetInfo_External; 33169e15a41SShri Abhyankar B->ops->lufactorsymbolic = MatLUFactorSymbolic_KLU; 33269e15a41SShri Abhyankar B->ops->destroy = MatDestroy_KLU; 33369e15a41SShri Abhyankar B->ops->view = MatView_KLU; 33469e15a41SShri Abhyankar 33569e15a41SShri Abhyankar ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_seqaij_klu);CHKERRQ(ierr); 33669e15a41SShri Abhyankar 33769e15a41SShri Abhyankar B->factortype = MAT_FACTOR_LU; 33869e15a41SShri Abhyankar B->assembled = PETSC_TRUE; /* required by -ksp_view */ 33969e15a41SShri Abhyankar B->preallocated = PETSC_TRUE; 34069e15a41SShri Abhyankar 34100c67f3bSHong Zhang ierr = PetscFree(B->solvertype);CHKERRQ(ierr); 34200c67f3bSHong Zhang ierr = PetscStrallocpy(MATSOLVERKLU,&B->solvertype);CHKERRQ(ierr); 34300c67f3bSHong Zhang 34469e15a41SShri Abhyankar /* initializations */ 34569e15a41SShri Abhyankar /* ------------------------------------------------*/ 34669e15a41SShri Abhyankar /* get the default control parameters */ 34769e15a41SShri Abhyankar status = klu_K_defaults(&lu->Common); 3486c4ed002SBarry Smith if(status <= 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"KLU Initialization failed"); 3496c4ed002SBarry Smith 35069e15a41SShri Abhyankar lu->Common.scale = 0; /* No row scaling */ 35169e15a41SShri Abhyankar 35269e15a41SShri Abhyankar ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"KLU Options","Mat");CHKERRQ(ierr); 35369e15a41SShri Abhyankar /* Partial pivoting tolerance */ 35469e15a41SShri Abhyankar ierr = PetscOptionsReal("-mat_klu_pivot_tol","Partial pivoting tolerance","None",lu->Common.tol,&lu->Common.tol,NULL);CHKERRQ(ierr); 35569e15a41SShri Abhyankar /* BTF pre-ordering */ 35669e15a41SShri Abhyankar ierr = PetscOptionsInt("-mat_klu_use_btf","Enable BTF preordering","None",lu->Common.btf,&lu->Common.btf,NULL);CHKERRQ(ierr); 35769e15a41SShri Abhyankar /* Matrix reordering */ 35869e15a41SShri Abhyankar ierr = PetscOptionsEList("-mat_klu_ordering","Internal ordering method","None",KluOrderingTypes,sizeof(KluOrderingTypes)/sizeof(KluOrderingTypes[0]),KluOrderingTypes[0],&idx,&flg);CHKERRQ(ierr); 35969e15a41SShri Abhyankar if (flg) { 36069e15a41SShri Abhyankar if ((int)idx == 2) lu->PetscMatOrdering = PETSC_TRUE; /* use Petsc mat ordering (note: size is for the transpose, and PETSc r = Klu perm_c) */ 36169e15a41SShri Abhyankar else lu->Common.ordering = (int)idx; 36269e15a41SShri Abhyankar } 36369e15a41SShri Abhyankar /* Matrix row scaling */ 36469e15a41SShri Abhyankar ierr = PetscOptionsEList("-mat_klu_row_scale","Matrix row scaling","None",scale,3,scale[0],&idx,&flg);CHKERRQ(ierr); 36569e15a41SShri Abhyankar PetscOptionsEnd(); 36669e15a41SShri Abhyankar *F = B; 36769e15a41SShri Abhyankar PetscFunctionReturn(0); 36869e15a41SShri Abhyankar } 369