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) {
34: MatSeqAIJRestoreArrayRead(A,&aa);
35: }
37: PetscMemzero(C,sizeof(*C));
39: C->nrow = (size_t)A->cmap->n;
40: C->ncol = (size_t)A->rmap->n;
41: C->nzmax = (size_t)nz;
42: C->p = ci;
43: C->i = cj;
44: C->x = values ? ca : 0;
45: C->stype = -1;
46: C->itype = CHOLMOD_INT_TYPE;
47: C->xtype = values ? CHOLMOD_SCALAR_TYPE : CHOLMOD_PATTERN;
48: C->dtype = CHOLMOD_DOUBLE;
49: C->sorted = 1;
50: C->packed = 1;
51: return 0;
52: }
54: static PetscErrorCode MatFactorGetSolverType_seqaij_cholmod(Mat A,MatSolverType *type)
55: {
56: *type = MATSOLVERCHOLMOD;
57: return 0;
58: }
60: /* Almost a copy of MatGetFactor_seqsbaij_cholmod, yuck */
61: PETSC_INTERN PetscErrorCode MatGetFactor_seqaij_cholmod(Mat A,MatFactorType ftype,Mat *F)
62: {
63: Mat B;
64: Mat_CHOLMOD *chol;
65: PetscInt m=A->rmap->n,n=A->cmap->n;
66: const char *prefix;
68: #if defined(PETSC_USE_COMPLEX)
70: #endif
71: /* Create the factorization matrix F */
72: MatCreate(PetscObjectComm((PetscObject)A),&B);
73: MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,m,n);
74: PetscStrallocpy("cholmod",&((PetscObject)B)->type_name);
75: MatGetOptionsPrefix(A,&prefix);
76: MatSetOptionsPrefix(B,prefix);
77: MatSetUp(B);
78: PetscNewLog(B,&chol);
80: chol->Wrap = MatWrapCholmod_seqaij;
81: B->data = chol;
83: B->ops->getinfo = MatGetInfo_CHOLMOD;
84: B->ops->view = MatView_CHOLMOD;
85: B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_CHOLMOD;
86: B->ops->destroy = MatDestroy_CHOLMOD;
88: PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverType_C",MatFactorGetSolverType_seqaij_cholmod);
90: B->factortype = MAT_FACTOR_CHOLESKY;
91: B->assembled = PETSC_TRUE;
92: B->preallocated = PETSC_TRUE;
94: PetscFree(B->solvertype);
95: PetscStrallocpy(MATSOLVERCHOLMOD,&B->solvertype);
96: B->canuseordering = PETSC_TRUE;
97: PetscStrallocpy(MATORDERINGEXTERNAL,(char**)&B->preferredordering[MAT_FACTOR_CHOLESKY]);
98: CholmodStart(B);
99: *F = B;
100: return 0;
101: }