Actual source code: aijcholmod.c
2: #include <../src/mat/impls/aij/seq/aij.h>
3: #include <../src/mat/impls/sbaij/seq/cholmod/cholmodimpl.h>
7: static PetscErrorCode MatWrapCholmod_seqaij(Mat A,PetscBool values,cholmod_sparse *C,PetscBool *aijalloc)
8: {
9: Mat_SeqAIJ *aij = (Mat_SeqAIJ*)A->data;
10: const PetscInt *ai = aij->i,*aj = aij->j,*adiag;
11: const MatScalar *aa = aij->a;
12: PetscInt m = A->rmap->n,i,j,k,nz,*ci,*cj;
13: PetscScalar *ca;
14: PetscErrorCode ierr;
17: MatMarkDiagonal_SeqAIJ(A);
18: adiag = aij->diag;
19: for (i=0,nz=0; i<m; i++) nz += ai[i+1] - adiag[i];
20: PetscMalloc3(m+1,PetscInt,&ci,nz,PetscInt,&cj,values?nz:0,PetscScalar,&ca);
21: for (i=0,k=0; i<m; i++) {
22: ci[i] = k;
23: for (j=adiag[i]; j<ai[i+1]; j++,k++) {
24: cj[k] = aj[j];
25: if (values) ca[k] = aa[j];
26: }
27: }
28: ci[i] = k;
29: *aijalloc = PETSC_TRUE;
31: PetscMemzero(C,sizeof(*C));
32: C->nrow = (size_t)A->cmap->n;
33: C->ncol = (size_t)A->rmap->n;
34: C->nzmax = (size_t)nz;
35: C->p = ci;
36: C->i = cj;
37: C->x = values ? ca : 0;
38: C->stype = -1;
39: C->itype = CHOLMOD_INT_TYPE;
40: C->xtype = values ? CHOLMOD_SCALAR_TYPE : CHOLMOD_PATTERN;
41: C->dtype = CHOLMOD_DOUBLE;
42: C->sorted = 1;
43: C->packed = 1;
44: return(0);
45: }
50: PetscErrorCode MatFactorGetSolverPackage_seqaij_cholmod(Mat A,const MatSolverPackage *type)
51: {
53: *type = MATSOLVERCHOLMOD;
54: return(0);
55: }
61: /* Almost a copy of MatGetFactor_seqsbaij_cholmod, yuck */
62: PetscErrorCode MatGetFactor_seqaij_cholmod(Mat A,MatFactorType ftype,Mat *F)
63: {
64: Mat B;
65: Mat_CHOLMOD *chol;
67: PetscInt m=A->rmap->n,n=A->cmap->n;
70: if (ftype != MAT_FACTOR_CHOLESKY) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"CHOLMOD cannot do %s factorization with AIJ, only %s",
71: MatFactorTypes[ftype],MatFactorTypes[MAT_FACTOR_CHOLESKY]);
72: /* Create the factorization matrix F */
73: MatCreate(((PetscObject)A)->comm,&B);
74: MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,m,n);
75: MatSetType(B,((PetscObject)A)->type_name);
76: MatSeqAIJSetPreallocation(B,0,PETSC_NULL);
77: PetscNewLog(B,Mat_CHOLMOD,&chol);
78: chol->Wrap = MatWrapCholmod_seqaij;
79: chol->Destroy = MatDestroy_SeqAIJ;
80: B->spptr = chol;
82: B->ops->view = MatView_CHOLMOD;
83: B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_CHOLMOD;
84: B->ops->destroy = MatDestroy_CHOLMOD;
85: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_seqaij_cholmod",MatFactorGetSolverPackage_seqaij_cholmod);
86: B->factortype = MAT_FACTOR_CHOLESKY;
87: B->assembled = PETSC_TRUE; /* required by -ksp_view */
88: B->preallocated = PETSC_TRUE;
90: CholmodStart(B);
91: *F = B;
92: return(0);
93: }