1586621ddSJed Brown 2c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h> 3c6db04a5SJed Brown #include <../src/mat/impls/sbaij/seq/cholmod/cholmodimpl.h> 4586621ddSJed Brown 59371c9d4SSatish Balay static PetscErrorCode MatWrapCholmod_seqaij(Mat A, PetscBool values, cholmod_sparse *C, PetscBool *aijalloc, PetscBool *valloc) { 6586621ddSJed Brown Mat_SeqAIJ *aij = (Mat_SeqAIJ *)A->data; 7218db3c1SStefano Zampini const PetscScalar *aa; 8586621ddSJed Brown PetscScalar *ca; 9218db3c1SStefano Zampini const PetscInt *ai = aij->i, *aj = aij->j, *adiag; 10218db3c1SStefano Zampini PetscInt m = A->rmap->n, i, j, k, nz, *ci, *cj; 11218db3c1SStefano Zampini PetscBool vain = PETSC_FALSE; 12586621ddSJed Brown 13586621ddSJed Brown PetscFunctionBegin; 149566063dSJacob Faibussowitsch PetscCall(MatMarkDiagonal_SeqAIJ(A)); 15586621ddSJed Brown adiag = aij->diag; 16586621ddSJed Brown for (i = 0, nz = 0; i < m; i++) nz += ai[i + 1] - adiag[i]; 179566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(m + 1, &ci, nz, &cj)); 18218db3c1SStefano Zampini if (values) { 19218db3c1SStefano Zampini vain = PETSC_TRUE; 209566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nz, &ca)); 219566063dSJacob Faibussowitsch PetscCall(MatSeqAIJGetArrayRead(A, &aa)); 22218db3c1SStefano Zampini } 23586621ddSJed Brown for (i = 0, k = 0; i < m; i++) { 24586621ddSJed Brown ci[i] = k; 25586621ddSJed Brown for (j = adiag[i]; j < ai[i + 1]; j++, k++) { 26586621ddSJed Brown cj[k] = aj[j]; 27218db3c1SStefano Zampini if (values) ca[k] = PetscConj(aa[j]); 28586621ddSJed Brown } 29586621ddSJed Brown } 30586621ddSJed Brown ci[i] = k; 31586621ddSJed Brown *aijalloc = PETSC_TRUE; 32218db3c1SStefano Zampini *valloc = vain; 33*48a46eb9SPierre Jolivet if (values) PetscCall(MatSeqAIJRestoreArrayRead(A, &aa)); 34586621ddSJed Brown 359566063dSJacob Faibussowitsch PetscCall(PetscMemzero(C, sizeof(*C))); 362205254eSKarl Rupp 37586621ddSJed Brown C->nrow = (size_t)A->cmap->n; 38586621ddSJed Brown C->ncol = (size_t)A->rmap->n; 39586621ddSJed Brown C->nzmax = (size_t)nz; 40586621ddSJed Brown C->p = ci; 41586621ddSJed Brown C->i = cj; 42586621ddSJed Brown C->x = values ? ca : 0; 43586621ddSJed Brown C->stype = -1; 44586621ddSJed Brown C->itype = CHOLMOD_INT_TYPE; 45586621ddSJed Brown C->xtype = values ? CHOLMOD_SCALAR_TYPE : CHOLMOD_PATTERN; 46586621ddSJed Brown C->dtype = CHOLMOD_DOUBLE; 47586621ddSJed Brown C->sorted = 1; 48586621ddSJed Brown C->packed = 1; 49586621ddSJed Brown PetscFunctionReturn(0); 50586621ddSJed Brown } 51586621ddSJed Brown 529371c9d4SSatish Balay static PetscErrorCode MatFactorGetSolverType_seqaij_cholmod(Mat A, MatSolverType *type) { 53586621ddSJed Brown PetscFunctionBegin; 54586621ddSJed Brown *type = MATSOLVERCHOLMOD; 55586621ddSJed Brown PetscFunctionReturn(0); 56586621ddSJed Brown } 57586621ddSJed Brown 58586621ddSJed Brown /* Almost a copy of MatGetFactor_seqsbaij_cholmod, yuck */ 599371c9d4SSatish Balay PETSC_INTERN PetscErrorCode MatGetFactor_seqaij_cholmod(Mat A, MatFactorType ftype, Mat *F) { 60586621ddSJed Brown Mat B; 61586621ddSJed Brown Mat_CHOLMOD *chol; 62586621ddSJed Brown PetscInt m = A->rmap->n, n = A->cmap->n; 63586621ddSJed Brown 64586621ddSJed Brown PetscFunctionBegin; 65218db3c1SStefano Zampini #if defined(PETSC_USE_COMPLEX) 66b94d7dedSBarry Smith PetscCheck(A->hermitian == PETSC_BOOL3_TRUE, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Only for Hermitian matrices"); 67218db3c1SStefano Zampini #endif 68586621ddSJed Brown /* Create the factorization matrix F */ 699566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B)); 709566063dSJacob Faibussowitsch PetscCall(MatSetSizes(B, PETSC_DECIDE, PETSC_DECIDE, m, n)); 719566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy("cholmod", &((PetscObject)B)->type_name)); 729566063dSJacob Faibussowitsch PetscCall(MatSetUp(B)); 739566063dSJacob Faibussowitsch PetscCall(PetscNewLog(B, &chol)); 742205254eSKarl Rupp 75586621ddSJed Brown chol->Wrap = MatWrapCholmod_seqaij; 766b8f6f9dSBarry Smith B->data = chol; 77586621ddSJed Brown 78218db3c1SStefano Zampini B->ops->getinfo = MatGetInfo_CHOLMOD; 79586621ddSJed Brown B->ops->view = MatView_CHOLMOD; 80586621ddSJed Brown B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_CHOLMOD; 81586621ddSJed Brown B->ops->destroy = MatDestroy_CHOLMOD; 822205254eSKarl Rupp 839566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatFactorGetSolverType_C", MatFactorGetSolverType_seqaij_cholmod)); 842205254eSKarl Rupp 85586621ddSJed Brown B->factortype = MAT_FACTOR_CHOLESKY; 86218db3c1SStefano Zampini B->assembled = PETSC_TRUE; 87586621ddSJed Brown B->preallocated = PETSC_TRUE; 88586621ddSJed Brown 899566063dSJacob Faibussowitsch PetscCall(PetscFree(B->solvertype)); 909566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(MATSOLVERCHOLMOD, &B->solvertype)); 91f73b0415SBarry Smith B->canuseordering = PETSC_TRUE; 929566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(MATORDERINGEXTERNAL, (char **)&B->preferredordering[MAT_FACTOR_CHOLESKY])); 939566063dSJacob Faibussowitsch PetscCall(CholmodStart(B)); 94586621ddSJed Brown *F = B; 95586621ddSJed Brown PetscFunctionReturn(0); 96586621ddSJed Brown } 97