xref: /petsc/src/mat/impls/aij/seq/cholmod/aijcholmod.c (revision ae1ee55146a7ad071171b861759b85940c7e5c67)
1c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h>
2c6db04a5SJed Brown #include <../src/mat/impls/sbaij/seq/cholmod/cholmodimpl.h>
3586621ddSJed Brown 
MatWrapCholmod_seqaij(Mat A,PetscBool values,cholmod_sparse * C,PetscBool * aijalloc,PetscBool * valloc)4d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatWrapCholmod_seqaij(Mat A, PetscBool values, cholmod_sparse *C, PetscBool *aijalloc, PetscBool *valloc)
5d71ae5a4SJacob Faibussowitsch {
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;
14421480d9SBarry Smith   PetscCall(MatGetDiagonalMarkers_SeqAIJ(A, &adiag, NULL));
15586621ddSJed Brown   for (i = 0, nz = 0; i < m; i++) nz += ai[i + 1] - adiag[i];
169566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(m + 1, &ci, nz, &cj));
17218db3c1SStefano Zampini   if (values) {
18218db3c1SStefano Zampini     vain = PETSC_TRUE;
199566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nz, &ca));
209566063dSJacob Faibussowitsch     PetscCall(MatSeqAIJGetArrayRead(A, &aa));
21218db3c1SStefano Zampini   }
22586621ddSJed Brown   for (i = 0, k = 0; i < m; i++) {
23586621ddSJed Brown     ci[i] = k;
24586621ddSJed Brown     for (j = adiag[i]; j < ai[i + 1]; j++, k++) {
25586621ddSJed Brown       cj[k] = aj[j];
26218db3c1SStefano Zampini       if (values) ca[k] = PetscConj(aa[j]);
27586621ddSJed Brown     }
28586621ddSJed Brown   }
29586621ddSJed Brown   ci[i]     = k;
30586621ddSJed Brown   *aijalloc = PETSC_TRUE;
31218db3c1SStefano Zampini   *valloc   = vain;
3248a46eb9SPierre Jolivet   if (values) PetscCall(MatSeqAIJRestoreArrayRead(A, &aa));
33586621ddSJed Brown 
349566063dSJacob Faibussowitsch   PetscCall(PetscMemzero(C, sizeof(*C)));
352205254eSKarl Rupp 
36586621ddSJed Brown   C->nrow   = (size_t)A->cmap->n;
37586621ddSJed Brown   C->ncol   = (size_t)A->rmap->n;
38586621ddSJed Brown   C->nzmax  = (size_t)nz;
39586621ddSJed Brown   C->p      = ci;
40586621ddSJed Brown   C->i      = cj;
41586621ddSJed Brown   C->x      = values ? ca : 0;
42586621ddSJed Brown   C->stype  = -1;
43586621ddSJed Brown   C->itype  = CHOLMOD_INT_TYPE;
44586621ddSJed Brown   C->xtype  = values ? CHOLMOD_SCALAR_TYPE : CHOLMOD_PATTERN;
45586621ddSJed Brown   C->dtype  = CHOLMOD_DOUBLE;
46586621ddSJed Brown   C->sorted = 1;
47586621ddSJed Brown   C->packed = 1;
483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
49586621ddSJed Brown }
50586621ddSJed Brown 
MatFactorGetSolverType_seqaij_cholmod(Mat A,MatSolverType * type)51d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatFactorGetSolverType_seqaij_cholmod(Mat A, MatSolverType *type)
52d71ae5a4SJacob Faibussowitsch {
53586621ddSJed Brown   PetscFunctionBegin;
54586621ddSJed Brown   *type = MATSOLVERCHOLMOD;
553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
56586621ddSJed Brown }
57586621ddSJed Brown 
58586621ddSJed Brown /* Almost a copy of MatGetFactor_seqsbaij_cholmod, yuck */
MatGetFactor_seqaij_cholmod(Mat A,MatFactorType ftype,Mat * F)59d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode MatGetFactor_seqaij_cholmod(Mat A, MatFactorType ftype, Mat *F)
60d71ae5a4SJacob Faibussowitsch {
61586621ddSJed Brown   Mat          B;
62586621ddSJed Brown   Mat_CHOLMOD *chol;
63586621ddSJed Brown   PetscInt     m = A->rmap->n, n = A->cmap->n;
64586621ddSJed Brown 
65586621ddSJed Brown   PetscFunctionBegin;
66*fc2fb351SPierre Jolivet   if (PetscDefined(USE_COMPLEX) && A->hermitian != PETSC_BOOL3_TRUE) {
6703e5aca4SStefano Zampini     PetscCall(PetscInfo(A, "Only for Hermitian matrices.\n"));
6803e5aca4SStefano Zampini     *F = NULL;
6903e5aca4SStefano Zampini     PetscFunctionReturn(PETSC_SUCCESS);
7003e5aca4SStefano Zampini   }
71586621ddSJed Brown   /* Create the factorization matrix F */
729566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B));
739566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(B, PETSC_DECIDE, PETSC_DECIDE, m, n));
749566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy("cholmod", &((PetscObject)B)->type_name));
759566063dSJacob Faibussowitsch   PetscCall(MatSetUp(B));
764dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&chol));
772205254eSKarl Rupp 
78586621ddSJed Brown   chol->Wrap = MatWrapCholmod_seqaij;
796b8f6f9dSBarry Smith   B->data    = chol;
80586621ddSJed Brown 
81218db3c1SStefano Zampini   B->ops->getinfo                = MatGetInfo_CHOLMOD;
82586621ddSJed Brown   B->ops->view                   = MatView_CHOLMOD;
83586621ddSJed Brown   B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_CHOLMOD;
84586621ddSJed Brown   B->ops->destroy                = MatDestroy_CHOLMOD;
852205254eSKarl Rupp 
869566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatFactorGetSolverType_C", MatFactorGetSolverType_seqaij_cholmod));
872205254eSKarl Rupp 
88586621ddSJed Brown   B->factortype   = MAT_FACTOR_CHOLESKY;
89218db3c1SStefano Zampini   B->assembled    = PETSC_TRUE;
90586621ddSJed Brown   B->preallocated = PETSC_TRUE;
91586621ddSJed Brown 
929566063dSJacob Faibussowitsch   PetscCall(PetscFree(B->solvertype));
939566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(MATSOLVERCHOLMOD, &B->solvertype));
94f73b0415SBarry Smith   B->canuseordering = PETSC_TRUE;
959566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(MATORDERINGEXTERNAL, (char **)&B->preferredordering[MAT_FACTOR_CHOLESKY]));
969566063dSJacob Faibussowitsch   PetscCall(CholmodStart(B));
97586621ddSJed Brown   *F = B;
983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
99586621ddSJed Brown }
100