xref: /petsc/src/mat/impls/aij/seq/cholmod/aijcholmod.c (revision 03e5aca47aa9bd9d5a48c628e4311d43ab2119ff)
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 
5d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatWrapCholmod_seqaij(Mat A, PetscBool values, cholmod_sparse *C, PetscBool *aijalloc, PetscBool *valloc)
6d71ae5a4SJacob Faibussowitsch {
7586621ddSJed Brown   Mat_SeqAIJ        *aij = (Mat_SeqAIJ *)A->data;
8218db3c1SStefano Zampini   const PetscScalar *aa;
9586621ddSJed Brown   PetscScalar       *ca;
10218db3c1SStefano Zampini   const PetscInt    *ai = aij->i, *aj = aij->j, *adiag;
11218db3c1SStefano Zampini   PetscInt           m    = A->rmap->n, i, j, k, nz, *ci, *cj;
12218db3c1SStefano Zampini   PetscBool          vain = PETSC_FALSE;
13586621ddSJed Brown 
14586621ddSJed Brown   PetscFunctionBegin;
159566063dSJacob Faibussowitsch   PetscCall(MatMarkDiagonal_SeqAIJ(A));
16586621ddSJed Brown   adiag = aij->diag;
17586621ddSJed Brown   for (i = 0, nz = 0; i < m; i++) nz += ai[i + 1] - adiag[i];
189566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(m + 1, &ci, nz, &cj));
19218db3c1SStefano Zampini   if (values) {
20218db3c1SStefano Zampini     vain = PETSC_TRUE;
219566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nz, &ca));
229566063dSJacob Faibussowitsch     PetscCall(MatSeqAIJGetArrayRead(A, &aa));
23218db3c1SStefano Zampini   }
24586621ddSJed Brown   for (i = 0, k = 0; i < m; i++) {
25586621ddSJed Brown     ci[i] = k;
26586621ddSJed Brown     for (j = adiag[i]; j < ai[i + 1]; j++, k++) {
27586621ddSJed Brown       cj[k] = aj[j];
28218db3c1SStefano Zampini       if (values) ca[k] = PetscConj(aa[j]);
29586621ddSJed Brown     }
30586621ddSJed Brown   }
31586621ddSJed Brown   ci[i]     = k;
32586621ddSJed Brown   *aijalloc = PETSC_TRUE;
33218db3c1SStefano Zampini   *valloc   = vain;
3448a46eb9SPierre Jolivet   if (values) PetscCall(MatSeqAIJRestoreArrayRead(A, &aa));
35586621ddSJed Brown 
369566063dSJacob Faibussowitsch   PetscCall(PetscMemzero(C, sizeof(*C)));
372205254eSKarl Rupp 
38586621ddSJed Brown   C->nrow   = (size_t)A->cmap->n;
39586621ddSJed Brown   C->ncol   = (size_t)A->rmap->n;
40586621ddSJed Brown   C->nzmax  = (size_t)nz;
41586621ddSJed Brown   C->p      = ci;
42586621ddSJed Brown   C->i      = cj;
43586621ddSJed Brown   C->x      = values ? ca : 0;
44586621ddSJed Brown   C->stype  = -1;
45586621ddSJed Brown   C->itype  = CHOLMOD_INT_TYPE;
46586621ddSJed Brown   C->xtype  = values ? CHOLMOD_SCALAR_TYPE : CHOLMOD_PATTERN;
47586621ddSJed Brown   C->dtype  = CHOLMOD_DOUBLE;
48586621ddSJed Brown   C->sorted = 1;
49586621ddSJed Brown   C->packed = 1;
503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
51586621ddSJed Brown }
52586621ddSJed Brown 
53d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatFactorGetSolverType_seqaij_cholmod(Mat A, MatSolverType *type)
54d71ae5a4SJacob Faibussowitsch {
55586621ddSJed Brown   PetscFunctionBegin;
56586621ddSJed Brown   *type = MATSOLVERCHOLMOD;
573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
58586621ddSJed Brown }
59586621ddSJed Brown 
60586621ddSJed Brown /* Almost a copy of MatGetFactor_seqsbaij_cholmod, yuck */
61d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode MatGetFactor_seqaij_cholmod(Mat A, MatFactorType ftype, Mat *F)
62d71ae5a4SJacob Faibussowitsch {
63586621ddSJed Brown   Mat          B;
64586621ddSJed Brown   Mat_CHOLMOD *chol;
65586621ddSJed Brown   PetscInt     m = A->rmap->n, n = A->cmap->n;
66586621ddSJed Brown 
67586621ddSJed Brown   PetscFunctionBegin;
68218db3c1SStefano Zampini #if defined(PETSC_USE_COMPLEX)
69*03e5aca4SStefano Zampini   if (A->hermitian != PETSC_BOOL3_TRUE) {
70*03e5aca4SStefano Zampini     PetscCall(PetscInfo(A, "Only for Hermitian matrices.\n"));
71*03e5aca4SStefano Zampini     *F = NULL;
72*03e5aca4SStefano Zampini     PetscFunctionReturn(PETSC_SUCCESS);
73*03e5aca4SStefano Zampini   }
74218db3c1SStefano Zampini #endif
75586621ddSJed Brown   /* Create the factorization matrix F */
769566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B));
779566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(B, PETSC_DECIDE, PETSC_DECIDE, m, n));
789566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy("cholmod", &((PetscObject)B)->type_name));
799566063dSJacob Faibussowitsch   PetscCall(MatSetUp(B));
804dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&chol));
812205254eSKarl Rupp 
82586621ddSJed Brown   chol->Wrap = MatWrapCholmod_seqaij;
836b8f6f9dSBarry Smith   B->data    = chol;
84586621ddSJed Brown 
85218db3c1SStefano Zampini   B->ops->getinfo                = MatGetInfo_CHOLMOD;
86586621ddSJed Brown   B->ops->view                   = MatView_CHOLMOD;
87586621ddSJed Brown   B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_CHOLMOD;
88586621ddSJed Brown   B->ops->destroy                = MatDestroy_CHOLMOD;
892205254eSKarl Rupp 
909566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatFactorGetSolverType_C", MatFactorGetSolverType_seqaij_cholmod));
912205254eSKarl Rupp 
92586621ddSJed Brown   B->factortype   = MAT_FACTOR_CHOLESKY;
93218db3c1SStefano Zampini   B->assembled    = PETSC_TRUE;
94586621ddSJed Brown   B->preallocated = PETSC_TRUE;
95586621ddSJed Brown 
969566063dSJacob Faibussowitsch   PetscCall(PetscFree(B->solvertype));
979566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(MATSOLVERCHOLMOD, &B->solvertype));
98f73b0415SBarry Smith   B->canuseordering = PETSC_TRUE;
999566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(MATORDERINGEXTERNAL, (char **)&B->preferredordering[MAT_FACTOR_CHOLESKY]));
1009566063dSJacob Faibussowitsch   PetscCall(CholmodStart(B));
101586621ddSJed Brown   *F = B;
1023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
103586621ddSJed Brown }
104