xref: /petsc/src/mat/impls/aij/seq/cholmod/aijcholmod.c (revision b94d7ded0a05f1bbd5e48daa6f92b28259c75b44)
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 
5218db3c1SStefano Zampini static PetscErrorCode MatWrapCholmod_seqaij(Mat A,PetscBool values,cholmod_sparse *C,PetscBool *aijalloc,PetscBool *valloc)
6586621ddSJed Brown {
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;
34218db3c1SStefano Zampini   if (values) {
359566063dSJacob Faibussowitsch     PetscCall(MatSeqAIJRestoreArrayRead(A,&aa));
36218db3c1SStefano Zampini   }
37586621ddSJed Brown 
389566063dSJacob Faibussowitsch   PetscCall(PetscMemzero(C,sizeof(*C)));
392205254eSKarl Rupp 
40586621ddSJed Brown   C->nrow   = (size_t)A->cmap->n;
41586621ddSJed Brown   C->ncol   = (size_t)A->rmap->n;
42586621ddSJed Brown   C->nzmax  = (size_t)nz;
43586621ddSJed Brown   C->p      = ci;
44586621ddSJed Brown   C->i      = cj;
45586621ddSJed Brown   C->x      = values ? ca : 0;
46586621ddSJed Brown   C->stype  = -1;
47586621ddSJed Brown   C->itype  = CHOLMOD_INT_TYPE;
48586621ddSJed Brown   C->xtype  = values ? CHOLMOD_SCALAR_TYPE : CHOLMOD_PATTERN;
49586621ddSJed Brown   C->dtype  = CHOLMOD_DOUBLE;
50586621ddSJed Brown   C->sorted = 1;
51586621ddSJed Brown   C->packed = 1;
52586621ddSJed Brown   PetscFunctionReturn(0);
53586621ddSJed Brown }
54586621ddSJed Brown 
55ea799195SBarry Smith static PetscErrorCode MatFactorGetSolverType_seqaij_cholmod(Mat A,MatSolverType *type)
56586621ddSJed Brown {
57586621ddSJed Brown   PetscFunctionBegin;
58586621ddSJed Brown   *type = MATSOLVERCHOLMOD;
59586621ddSJed Brown   PetscFunctionReturn(0);
60586621ddSJed Brown }
61586621ddSJed Brown 
62586621ddSJed Brown /* Almost a copy of MatGetFactor_seqsbaij_cholmod, yuck */
63db87b0f2SBarry Smith PETSC_INTERN PetscErrorCode MatGetFactor_seqaij_cholmod(Mat A,MatFactorType ftype,Mat *F)
64586621ddSJed Brown {
65586621ddSJed Brown   Mat            B;
66586621ddSJed Brown   Mat_CHOLMOD    *chol;
67586621ddSJed Brown   PetscInt       m=A->rmap->n,n=A->cmap->n;
68586621ddSJed Brown 
69586621ddSJed Brown   PetscFunctionBegin;
70218db3c1SStefano Zampini #if defined(PETSC_USE_COMPLEX)
71*b94d7dedSBarry Smith   PetscCheck(A->hermitian == PETSC_BOOL3_TRUE,PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Only for Hermitian matrices");
72218db3c1SStefano Zampini #endif
73586621ddSJed Brown   /* Create the factorization matrix F */
749566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)A),&B));
759566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,m,n));
769566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy("cholmod",&((PetscObject)B)->type_name));
779566063dSJacob Faibussowitsch   PetscCall(MatSetUp(B));
789566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(B,&chol));
792205254eSKarl Rupp 
80586621ddSJed Brown   chol->Wrap = MatWrapCholmod_seqaij;
816b8f6f9dSBarry Smith   B->data    = chol;
82586621ddSJed Brown 
83218db3c1SStefano Zampini   B->ops->getinfo                = MatGetInfo_CHOLMOD;
84586621ddSJed Brown   B->ops->view                   = MatView_CHOLMOD;
85586621ddSJed Brown   B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_CHOLMOD;
86586621ddSJed Brown   B->ops->destroy                = MatDestroy_CHOLMOD;
872205254eSKarl Rupp 
889566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverType_C",MatFactorGetSolverType_seqaij_cholmod));
892205254eSKarl Rupp 
90586621ddSJed Brown   B->factortype   = MAT_FACTOR_CHOLESKY;
91218db3c1SStefano Zampini   B->assembled    = PETSC_TRUE;
92586621ddSJed Brown   B->preallocated = PETSC_TRUE;
93586621ddSJed Brown 
949566063dSJacob Faibussowitsch   PetscCall(PetscFree(B->solvertype));
959566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(MATSOLVERCHOLMOD,&B->solvertype));
96f73b0415SBarry Smith   B->canuseordering = PETSC_TRUE;
979566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(MATORDERINGEXTERNAL,(char**)&B->preferredordering[MAT_FACTOR_CHOLESKY]));
989566063dSJacob Faibussowitsch   PetscCall(CholmodStart(B));
99586621ddSJed Brown   *F   = B;
100586621ddSJed Brown   PetscFunctionReturn(0);
101586621ddSJed Brown }
102