Actual source code: aijcholmod.c
2: #include <../src/mat/impls/aij/seq/aij.h>
3: #include <../src/mat/impls/sbaij/seq/cholmod/cholmodimpl.h>
5: static PetscErrorCode MatWrapCholmod_seqaij(Mat A, PetscBool values, cholmod_sparse *C, PetscBool *aijalloc, PetscBool *valloc)
6: {
7: Mat_SeqAIJ *aij = (Mat_SeqAIJ *)A->data;
8: const PetscScalar *aa;
9: PetscScalar *ca;
10: const PetscInt *ai = aij->i, *aj = aij->j, *adiag;
11: PetscInt m = A->rmap->n, i, j, k, nz, *ci, *cj;
12: PetscBool vain = PETSC_FALSE;
14: MatMarkDiagonal_SeqAIJ(A);
15: adiag = aij->diag;
16: for (i = 0, nz = 0; i < m; i++) nz += ai[i + 1] - adiag[i];
17: PetscMalloc2(m + 1, &ci, nz, &cj);
18: if (values) {
19: vain = PETSC_TRUE;
20: PetscMalloc1(nz, &ca);
21: MatSeqAIJGetArrayRead(A, &aa);
22: }
23: for (i = 0, k = 0; i < m; i++) {
24: ci[i] = k;
25: for (j = adiag[i]; j < ai[i + 1]; j++, k++) {
26: cj[k] = aj[j];
27: if (values) ca[k] = PetscConj(aa[j]);
28: }
29: }
30: ci[i] = k;
31: *aijalloc = PETSC_TRUE;
32: *valloc = vain;
33: if (values) MatSeqAIJRestoreArrayRead(A, &aa);
35: PetscMemzero(C, sizeof(*C));
37: C->nrow = (size_t)A->cmap->n;
38: C->ncol = (size_t)A->rmap->n;
39: C->nzmax = (size_t)nz;
40: C->p = ci;
41: C->i = cj;
42: C->x = values ? ca : 0;
43: C->stype = -1;
44: C->itype = CHOLMOD_INT_TYPE;
45: C->xtype = values ? CHOLMOD_SCALAR_TYPE : CHOLMOD_PATTERN;
46: C->dtype = CHOLMOD_DOUBLE;
47: C->sorted = 1;
48: C->packed = 1;
49: return 0;
50: }
52: static PetscErrorCode MatFactorGetSolverType_seqaij_cholmod(Mat A, MatSolverType *type)
53: {
54: *type = MATSOLVERCHOLMOD;
55: return 0;
56: }
58: /* Almost a copy of MatGetFactor_seqsbaij_cholmod, yuck */
59: PETSC_INTERN PetscErrorCode MatGetFactor_seqaij_cholmod(Mat A, MatFactorType ftype, Mat *F)
60: {
61: Mat B;
62: Mat_CHOLMOD *chol;
63: PetscInt m = A->rmap->n, n = A->cmap->n;
65: #if defined(PETSC_USE_COMPLEX)
67: #endif
68: /* Create the factorization matrix F */
69: MatCreate(PetscObjectComm((PetscObject)A), &B);
70: MatSetSizes(B, PETSC_DECIDE, PETSC_DECIDE, m, n);
71: PetscStrallocpy("cholmod", &((PetscObject)B)->type_name);
72: MatSetUp(B);
73: PetscNew(&chol);
75: chol->Wrap = MatWrapCholmod_seqaij;
76: B->data = chol;
78: B->ops->getinfo = MatGetInfo_CHOLMOD;
79: B->ops->view = MatView_CHOLMOD;
80: B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_CHOLMOD;
81: B->ops->destroy = MatDestroy_CHOLMOD;
83: PetscObjectComposeFunction((PetscObject)B, "MatFactorGetSolverType_C", MatFactorGetSolverType_seqaij_cholmod);
85: B->factortype = MAT_FACTOR_CHOLESKY;
86: B->assembled = PETSC_TRUE;
87: B->preallocated = PETSC_TRUE;
89: PetscFree(B->solvertype);
90: PetscStrallocpy(MATSOLVERCHOLMOD, &B->solvertype);
91: B->canuseordering = PETSC_TRUE;
92: PetscStrallocpy(MATORDERINGEXTERNAL, (char **)&B->preferredordering[MAT_FACTOR_CHOLESKY]);
93: CholmodStart(B);
94: *F = B;
95: return 0;
96: }