1 2 #include <../src/mat/impls/aij/seq/aij.h> 3 #include <../src/mat/impls/sbaij/seq/cholmod/cholmodimpl.h> 4 5 static PetscErrorCode MatWrapCholmod_seqaij(Mat A, PetscBool values, cholmod_sparse *C, PetscBool *aijalloc, PetscBool *valloc) { 6 Mat_SeqAIJ *aij = (Mat_SeqAIJ *)A->data; 7 const PetscScalar *aa; 8 PetscScalar *ca; 9 const PetscInt *ai = aij->i, *aj = aij->j, *adiag; 10 PetscInt m = A->rmap->n, i, j, k, nz, *ci, *cj; 11 PetscBool vain = PETSC_FALSE; 12 13 PetscFunctionBegin; 14 PetscCall(MatMarkDiagonal_SeqAIJ(A)); 15 adiag = aij->diag; 16 for (i = 0, nz = 0; i < m; i++) nz += ai[i + 1] - adiag[i]; 17 PetscCall(PetscMalloc2(m + 1, &ci, nz, &cj)); 18 if (values) { 19 vain = PETSC_TRUE; 20 PetscCall(PetscMalloc1(nz, &ca)); 21 PetscCall(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) { PetscCall(MatSeqAIJRestoreArrayRead(A, &aa)); } 34 35 PetscCall(PetscMemzero(C, sizeof(*C))); 36 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 PetscFunctionReturn(0); 50 } 51 52 static PetscErrorCode MatFactorGetSolverType_seqaij_cholmod(Mat A, MatSolverType *type) { 53 PetscFunctionBegin; 54 *type = MATSOLVERCHOLMOD; 55 PetscFunctionReturn(0); 56 } 57 58 /* Almost a copy of MatGetFactor_seqsbaij_cholmod, yuck */ 59 PETSC_INTERN PetscErrorCode MatGetFactor_seqaij_cholmod(Mat A, MatFactorType ftype, Mat *F) { 60 Mat B; 61 Mat_CHOLMOD *chol; 62 PetscInt m = A->rmap->n, n = A->cmap->n; 63 64 PetscFunctionBegin; 65 #if defined(PETSC_USE_COMPLEX) 66 PetscCheck(A->hermitian == PETSC_BOOL3_TRUE, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Only for Hermitian matrices"); 67 #endif 68 /* Create the factorization matrix F */ 69 PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B)); 70 PetscCall(MatSetSizes(B, PETSC_DECIDE, PETSC_DECIDE, m, n)); 71 PetscCall(PetscStrallocpy("cholmod", &((PetscObject)B)->type_name)); 72 PetscCall(MatSetUp(B)); 73 PetscCall(PetscNewLog(B, &chol)); 74 75 chol->Wrap = MatWrapCholmod_seqaij; 76 B->data = chol; 77 78 B->ops->getinfo = MatGetInfo_CHOLMOD; 79 B->ops->view = MatView_CHOLMOD; 80 B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_CHOLMOD; 81 B->ops->destroy = MatDestroy_CHOLMOD; 82 83 PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatFactorGetSolverType_C", MatFactorGetSolverType_seqaij_cholmod)); 84 85 B->factortype = MAT_FACTOR_CHOLESKY; 86 B->assembled = PETSC_TRUE; 87 B->preallocated = PETSC_TRUE; 88 89 PetscCall(PetscFree(B->solvertype)); 90 PetscCall(PetscStrallocpy(MATSOLVERCHOLMOD, &B->solvertype)); 91 B->canuseordering = PETSC_TRUE; 92 PetscCall(PetscStrallocpy(MATORDERINGEXTERNAL, (char **)&B->preferredordering[MAT_FACTOR_CHOLESKY])); 93 PetscCall(CholmodStart(B)); 94 *F = B; 95 PetscFunctionReturn(0); 96 } 97