Actual source code: cholmodimpl.h
1: #pragma once
3: #include <petscsys.h>
5: #if defined(PETSC_USE_COMPLEX)
6: #define CHOLMOD_SCALAR_TYPE CHOLMOD_COMPLEX
7: #else
8: #define CHOLMOD_SCALAR_TYPE CHOLMOD_REAL
9: #endif
11: #if defined(PETSC_USE_64BIT_INDICES)
12: #define CHOLMOD_INT_TYPE CHOLMOD_LONG
13: #define cholmod_X_start cholmod_l_start
14: #define cholmod_X_analyze cholmod_l_analyze
15: /* the type casts are needed because PetscInt is long long while SuiteSparse_long is long and compilers warn even when they are identical */
16: #define cholmod_X_analyze_p(a, b, c, d, f) cholmod_l_analyze_p(a, (SuiteSparse_long *)b, (SuiteSparse_long *)c, d, f)
17: #define cholmod_X_copy cholmod_l_copy
18: #define cholmod_X_factorize cholmod_l_factorize
19: #define cholmod_X_finish cholmod_l_finish
20: #define cholmod_X_free_factor cholmod_l_free_factor
21: #define cholmod_X_free_dense cholmod_l_free_dense
22: #define cholmod_X_resymbol(a, b, c, d, f, e) cholmod_l_resymbol(a, (SuiteSparse_long *)b, c, d, f, e)
23: #define cholmod_X_solve cholmod_l_solve
24: #define cholmod_X_solve2 cholmod_l_solve2
25: #define cholmod_X_check_sparse cholmod_l_check_sparse
26: #else
27: #define CHOLMOD_INT_TYPE CHOLMOD_INT
28: #define cholmod_X_start cholmod_start
29: #define cholmod_X_analyze cholmod_analyze
30: #define cholmod_X_analyze_p cholmod_analyze_p
31: #define cholmod_X_copy cholmod_copy
32: #define cholmod_X_factorize cholmod_factorize
33: #define cholmod_X_finish cholmod_finish
34: #define cholmod_X_free_factor cholmod_free_factor
35: #define cholmod_X_free_dense cholmod_free_dense
36: #define cholmod_X_resymbol cholmod_resymbol
37: #define cholmod_X_solve cholmod_solve
38: #define cholmod_X_solve2 cholmod_solve2
39: #define cholmod_X_check_sparse cholmod_check_sparse
40: #endif
42: EXTERN_C_BEGIN
43: #include <cholmod.h>
44: #include <SuiteSparseQR_C.h>
45: EXTERN_C_END
47: typedef struct {
48: PetscErrorCode (*Wrap)(Mat, PetscBool, cholmod_sparse *, PetscBool *, PetscBool *);
49: cholmod_sparse *matrix;
50: cholmod_factor *factor;
51: cholmod_common *common;
52: SuiteSparseQR_C_factorization *spqrfact;
53: PetscScalar scale;
54: PetscBool pack;
55: PetscBool normal;
56: } Mat_CHOLMOD;
58: PETSC_INTERN PetscErrorCode CholmodStart(Mat);
59: PETSC_INTERN PetscErrorCode MatView_CHOLMOD(Mat, PetscViewer);
60: PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_CHOLMOD(Mat, Mat, IS, const MatFactorInfo *);
61: PETSC_INTERN PetscErrorCode MatGetInfo_CHOLMOD(Mat, MatInfoType, MatInfo *);
62: PETSC_INTERN PetscErrorCode MatDestroy_CHOLMOD(Mat);
64: PETSC_INTERN PetscErrorCode VecWrapCholmod(Vec, PetscInt, cholmod_dense *);
65: PETSC_INTERN PetscErrorCode VecUnWrapCholmod(Vec, PetscInt, cholmod_dense *);
66: PETSC_INTERN PetscErrorCode MatDenseWrapCholmod(Mat, PetscInt, cholmod_dense *);
67: PETSC_INTERN PetscErrorCode MatDenseUnWrapCholmod(Mat, PetscInt, cholmod_dense *);