Actual source code: aijcholmod.c
petsc-3.14.6 2021-03-30
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;
13: PetscErrorCode ierr;
16: MatMarkDiagonal_SeqAIJ(A);
17: adiag = aij->diag;
18: for (i=0,nz=0; i<m; i++) nz += ai[i+1] - adiag[i];
19: PetscMalloc2(m+1,&ci,nz,&cj);
20: if (values) {
21: vain = PETSC_TRUE;
22: PetscMalloc1(nz,&ca);
23: MatSeqAIJGetArrayRead(A,&aa);
24: }
25: for (i=0,k=0; i<m; i++) {
26: ci[i] = k;
27: for (j=adiag[i]; j<ai[i+1]; j++,k++) {
28: cj[k] = aj[j];
29: if (values) ca[k] = PetscConj(aa[j]);
30: }
31: }
32: ci[i] = k;
33: *aijalloc = PETSC_TRUE;
34: *valloc = vain;
35: if (values) {
36: MatSeqAIJRestoreArrayRead(A,&aa);
37: }
39: PetscMemzero(C,sizeof(*C));
41: C->nrow = (size_t)A->cmap->n;
42: C->ncol = (size_t)A->rmap->n;
43: C->nzmax = (size_t)nz;
44: C->p = ci;
45: C->i = cj;
46: C->x = values ? ca : 0;
47: C->stype = -1;
48: C->itype = CHOLMOD_INT_TYPE;
49: C->xtype = values ? CHOLMOD_SCALAR_TYPE : CHOLMOD_PATTERN;
50: C->dtype = CHOLMOD_DOUBLE;
51: C->sorted = 1;
52: C->packed = 1;
53: return(0);
54: }
56: static PetscErrorCode MatFactorGetSolverType_seqaij_cholmod(Mat A,MatSolverType *type)
57: {
59: *type = MATSOLVERCHOLMOD;
60: return(0);
61: }
63: /* Almost a copy of MatGetFactor_seqsbaij_cholmod, yuck */
64: PETSC_INTERN PetscErrorCode MatGetFactor_seqaij_cholmod(Mat A,MatFactorType ftype,Mat *F)
65: {
66: Mat B;
67: Mat_CHOLMOD *chol;
69: PetscInt m=A->rmap->n,n=A->cmap->n;
70: const char *prefix;
73: #if defined(PETSC_USE_COMPLEX)
74: if (!A->hermitian) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Only for hermitian matrices");
75: #endif
76: /* Create the factorization matrix F */
77: MatCreate(PetscObjectComm((PetscObject)A),&B);
78: MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,m,n);
79: PetscStrallocpy("cholmod",&((PetscObject)B)->type_name);
80: MatGetOptionsPrefix(A,&prefix);
81: MatSetOptionsPrefix(B,prefix);
82: MatSetUp(B);
83: PetscNewLog(B,&chol);
85: chol->Wrap = MatWrapCholmod_seqaij;
86: B->data = chol;
88: B->ops->getinfo = MatGetInfo_CHOLMOD;
89: B->ops->view = MatView_CHOLMOD;
90: B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_CHOLMOD;
91: B->ops->destroy = MatDestroy_CHOLMOD;
93: PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverType_C",MatFactorGetSolverType_seqaij_cholmod);
95: B->factortype = MAT_FACTOR_CHOLESKY;
96: B->assembled = PETSC_TRUE;
97: B->preallocated = PETSC_TRUE;
99: PetscFree(B->solvertype);
100: PetscStrallocpy(MATSOLVERCHOLMOD,&B->solvertype);
101: B->useordering = PETSC_TRUE;
102: CholmodStart(B);
103: *F = B;
104: return(0);
105: }