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 PetscErrorCode ierr; 14586621ddSJed Brown 15586621ddSJed Brown PetscFunctionBegin; 16586621ddSJed Brown ierr = MatMarkDiagonal_SeqAIJ(A);CHKERRQ(ierr); 17586621ddSJed Brown adiag = aij->diag; 18586621ddSJed Brown for (i=0,nz=0; i<m; i++) nz += ai[i+1] - adiag[i]; 19218db3c1SStefano Zampini ierr = PetscMalloc2(m+1,&ci,nz,&cj);CHKERRQ(ierr); 20218db3c1SStefano Zampini if (values) { 21218db3c1SStefano Zampini vain = PETSC_TRUE; 22218db3c1SStefano Zampini ierr = PetscMalloc1(nz,&ca);CHKERRQ(ierr); 23218db3c1SStefano Zampini ierr = MatSeqAIJGetArrayRead(A,&aa);CHKERRQ(ierr); 24218db3c1SStefano Zampini } 25586621ddSJed Brown for (i=0,k=0; i<m; i++) { 26586621ddSJed Brown ci[i] = k; 27586621ddSJed Brown for (j=adiag[i]; j<ai[i+1]; j++,k++) { 28586621ddSJed Brown cj[k] = aj[j]; 29218db3c1SStefano Zampini if (values) ca[k] = PetscConj(aa[j]); 30586621ddSJed Brown } 31586621ddSJed Brown } 32586621ddSJed Brown ci[i] = k; 33586621ddSJed Brown *aijalloc = PETSC_TRUE; 34218db3c1SStefano Zampini *valloc = vain; 35218db3c1SStefano Zampini if (values) { 36218db3c1SStefano Zampini ierr = MatSeqAIJRestoreArrayRead(A,&aa);CHKERRQ(ierr); 37218db3c1SStefano Zampini } 38586621ddSJed Brown 39586621ddSJed Brown ierr = PetscMemzero(C,sizeof(*C));CHKERRQ(ierr); 402205254eSKarl Rupp 41586621ddSJed Brown C->nrow = (size_t)A->cmap->n; 42586621ddSJed Brown C->ncol = (size_t)A->rmap->n; 43586621ddSJed Brown C->nzmax = (size_t)nz; 44586621ddSJed Brown C->p = ci; 45586621ddSJed Brown C->i = cj; 46586621ddSJed Brown C->x = values ? ca : 0; 47586621ddSJed Brown C->stype = -1; 48586621ddSJed Brown C->itype = CHOLMOD_INT_TYPE; 49586621ddSJed Brown C->xtype = values ? CHOLMOD_SCALAR_TYPE : CHOLMOD_PATTERN; 50586621ddSJed Brown C->dtype = CHOLMOD_DOUBLE; 51586621ddSJed Brown C->sorted = 1; 52586621ddSJed Brown C->packed = 1; 53586621ddSJed Brown PetscFunctionReturn(0); 54586621ddSJed Brown } 55586621ddSJed Brown 56ea799195SBarry Smith static PetscErrorCode MatFactorGetSolverType_seqaij_cholmod(Mat A,MatSolverType *type) 57586621ddSJed Brown { 58586621ddSJed Brown PetscFunctionBegin; 59586621ddSJed Brown *type = MATSOLVERCHOLMOD; 60586621ddSJed Brown PetscFunctionReturn(0); 61586621ddSJed Brown } 62586621ddSJed Brown 63586621ddSJed Brown /* Almost a copy of MatGetFactor_seqsbaij_cholmod, yuck */ 64db87b0f2SBarry Smith PETSC_INTERN PetscErrorCode MatGetFactor_seqaij_cholmod(Mat A,MatFactorType ftype,Mat *F) 65586621ddSJed Brown { 66586621ddSJed Brown Mat B; 67586621ddSJed Brown Mat_CHOLMOD *chol; 68586621ddSJed Brown PetscErrorCode ierr; 69586621ddSJed Brown PetscInt m=A->rmap->n,n=A->cmap->n; 70218db3c1SStefano Zampini const char *prefix; 71586621ddSJed Brown 72586621ddSJed Brown PetscFunctionBegin; 73218db3c1SStefano Zampini #if defined(PETSC_USE_COMPLEX) 74218db3c1SStefano Zampini if (!A->hermitian) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Only for hermitian matrices"); 75218db3c1SStefano Zampini #endif 76586621ddSJed Brown /* Create the factorization matrix F */ 77ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr); 78586621ddSJed Brown ierr = MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,m,n);CHKERRQ(ierr); 79ecddaf3cSBarry Smith ierr = PetscStrallocpy("cholmod",&((PetscObject)B)->type_name);CHKERRQ(ierr); 80218db3c1SStefano Zampini ierr = MatGetOptionsPrefix(A,&prefix);CHKERRQ(ierr); 81218db3c1SStefano Zampini ierr = MatSetOptionsPrefix(B,prefix);CHKERRQ(ierr); 82ecddaf3cSBarry Smith ierr = MatSetUp(B);CHKERRQ(ierr); 83b00a9115SJed Brown ierr = PetscNewLog(B,&chol);CHKERRQ(ierr); 842205254eSKarl Rupp 85586621ddSJed Brown chol->Wrap = MatWrapCholmod_seqaij; 866b8f6f9dSBarry Smith B->data = chol; 87586621ddSJed Brown 88218db3c1SStefano Zampini B->ops->getinfo = MatGetInfo_CHOLMOD; 89586621ddSJed Brown B->ops->view = MatView_CHOLMOD; 90586621ddSJed Brown B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_CHOLMOD; 91586621ddSJed Brown B->ops->destroy = MatDestroy_CHOLMOD; 922205254eSKarl Rupp 933ca39a21SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverType_C",MatFactorGetSolverType_seqaij_cholmod);CHKERRQ(ierr); 942205254eSKarl Rupp 95586621ddSJed Brown B->factortype = MAT_FACTOR_CHOLESKY; 96218db3c1SStefano Zampini B->assembled = PETSC_TRUE; 97586621ddSJed Brown B->preallocated = PETSC_TRUE; 98586621ddSJed Brown 9900c67f3bSHong Zhang ierr = PetscFree(B->solvertype);CHKERRQ(ierr); 10000c67f3bSHong Zhang ierr = PetscStrallocpy(MATSOLVERCHOLMOD,&B->solvertype);CHKERRQ(ierr); 101*f73b0415SBarry Smith B->canuseordering = PETSC_TRUE; 102586621ddSJed Brown ierr = CholmodStart(B);CHKERRQ(ierr); 103586621ddSJed Brown *F = B; 104586621ddSJed Brown PetscFunctionReturn(0); 105586621ddSJed Brown } 106