MatICCFactorSymbolic#
Performs symbolic incomplete Cholesky factorization for a symmetric matrix. Use MatCholeskyFactorNumeric() to complete the factorization.
Synopsis#
#include "petscmat.h"
PetscErrorCode MatICCFactorSymbolic(Mat fact, Mat mat, IS perm, const MatFactorInfo *info)
Collective
Input Parameters#
fact - the factorized matrix obtained with
MatGetFactor()mat - the matrix to be factored
perm - row and column permutation
info - structure containing
levels - number of levels of fill.
expected fill - as ratio of original fill.
Notes#
Most users should employ the KSP interface for linear solvers
instead of working directly with matrix algebra routines such as this.
See, e.g., KSPCreate().
This uses the definition of level of fill as in Y. Saad [Saa03]
Fortran Note#
A valid (non-null) info argument must be provided
References#
Yousef Saad. Iterative Methods for Sparse Linear Systems. SIAM, 2nd edition, 2003. doi:10.1016/S1570-579X(01)80025-2.
See Also#
Matrices, Mat, MatGetFactor(), MatCholeskyFactorNumeric(), MatCholeskyFactor(), MatFactorInfo
Level#
developer
Location#
Implementations#
MatICCFactorSymbolic_SeqAIJ() in src/mat/impls/aij/seq/aijfact.c
MatICCFactorSymbolic_SeqAIJKokkos() in src/mat/impls/aij/seq/kokkos/aijkok.kokkos.cxx
MatICCFactorSymbolic_SeqAIJCUSPARSE() in src/mat/impls/aij/seq/seqcusparse/aijcusparse.cu
MatICCFactorSymbolic_SeqAIJHIPSPARSE() in src/mat/impls/aij/seq/seqhipsparse/aijhipsparse.hip.c
MatICCFactorSymbolic_SeqBAIJ() in src/mat/impls/baij/seq/baijfact.c
MatICCFactorSymbolic_SeqSBAIJ() in src/mat/impls/sbaij/seq/sbaijfact2.c
Index of all Mat routines
Table of Contents for all manual pages
Index of all manual pages