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 5ace3abfcSBarry Smith static PetscErrorCode MatWrapCholmod_seqaij(Mat A,PetscBool values,cholmod_sparse *C,PetscBool *aijalloc) 6586621ddSJed Brown { 7586621ddSJed Brown Mat_SeqAIJ *aij = (Mat_SeqAIJ*)A->data; 8586621ddSJed Brown const PetscInt *ai = aij->i,*aj = aij->j,*adiag; 9586621ddSJed Brown const MatScalar *aa = aij->a; 10586621ddSJed Brown PetscInt m = A->rmap->n,i,j,k,nz,*ci,*cj; 11586621ddSJed Brown PetscScalar *ca; 12586621ddSJed Brown PetscErrorCode ierr; 13586621ddSJed Brown 14586621ddSJed Brown PetscFunctionBegin; 15586621ddSJed Brown ierr = MatMarkDiagonal_SeqAIJ(A);CHKERRQ(ierr); 16586621ddSJed Brown adiag = aij->diag; 17586621ddSJed Brown for (i=0,nz=0; i<m; i++) nz += ai[i+1] - adiag[i]; 18dcca6d9dSJed Brown ierr = PetscMalloc3(m+1,&ci,nz,&cj,values ? nz : 0,&ca);CHKERRQ(ierr); 19586621ddSJed Brown for (i=0,k=0; i<m; i++) { 20586621ddSJed Brown ci[i] = k; 21586621ddSJed Brown for (j=adiag[i]; j<ai[i+1]; j++,k++) { 22586621ddSJed Brown cj[k] = aj[j]; 23586621ddSJed Brown if (values) ca[k] = aa[j]; 24586621ddSJed Brown } 25586621ddSJed Brown } 26586621ddSJed Brown ci[i] = k; 27586621ddSJed Brown *aijalloc = PETSC_TRUE; 28586621ddSJed Brown 29586621ddSJed Brown ierr = PetscMemzero(C,sizeof(*C));CHKERRQ(ierr); 302205254eSKarl Rupp 31586621ddSJed Brown C->nrow = (size_t)A->cmap->n; 32586621ddSJed Brown C->ncol = (size_t)A->rmap->n; 33586621ddSJed Brown C->nzmax = (size_t)nz; 34586621ddSJed Brown C->p = ci; 35586621ddSJed Brown C->i = cj; 36586621ddSJed Brown C->x = values ? ca : 0; 37586621ddSJed Brown C->stype = -1; 38586621ddSJed Brown C->itype = CHOLMOD_INT_TYPE; 39586621ddSJed Brown C->xtype = values ? CHOLMOD_SCALAR_TYPE : CHOLMOD_PATTERN; 40586621ddSJed Brown C->dtype = CHOLMOD_DOUBLE; 41586621ddSJed Brown C->sorted = 1; 42586621ddSJed Brown C->packed = 1; 43586621ddSJed Brown PetscFunctionReturn(0); 44586621ddSJed Brown } 45586621ddSJed Brown 46*ea799195SBarry Smith static PetscErrorCode MatFactorGetSolverType_seqaij_cholmod(Mat A,MatSolverType *type) 47586621ddSJed Brown { 48586621ddSJed Brown PetscFunctionBegin; 49586621ddSJed Brown *type = MATSOLVERCHOLMOD; 50586621ddSJed Brown PetscFunctionReturn(0); 51586621ddSJed Brown } 52586621ddSJed Brown 53586621ddSJed Brown /* Almost a copy of MatGetFactor_seqsbaij_cholmod, yuck */ 54db87b0f2SBarry Smith PETSC_INTERN PetscErrorCode MatGetFactor_seqaij_cholmod(Mat A,MatFactorType ftype,Mat *F) 55586621ddSJed Brown { 56586621ddSJed Brown Mat B; 57586621ddSJed Brown Mat_CHOLMOD *chol; 58586621ddSJed Brown PetscErrorCode ierr; 59586621ddSJed Brown PetscInt m=A->rmap->n,n=A->cmap->n; 60586621ddSJed Brown 61586621ddSJed Brown PetscFunctionBegin; 62ecddaf3cSBarry Smith if (ftype != MAT_FACTOR_CHOLESKY) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"CHOLMOD cannot do %s factorization with AIJ, only %s",MatFactorTypes[ftype],MatFactorTypes[MAT_FACTOR_CHOLESKY]); 63586621ddSJed Brown /* Create the factorization matrix F */ 64ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr); 65586621ddSJed Brown ierr = MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,m,n);CHKERRQ(ierr); 66ecddaf3cSBarry Smith ierr = PetscStrallocpy("cholmod",&((PetscObject)B)->type_name);CHKERRQ(ierr); 67ecddaf3cSBarry Smith ierr = MatSetUp(B);CHKERRQ(ierr); 68b00a9115SJed Brown ierr = PetscNewLog(B,&chol);CHKERRQ(ierr); 692205254eSKarl Rupp 70586621ddSJed Brown chol->Wrap = MatWrapCholmod_seqaij; 716b8f6f9dSBarry Smith B->data = chol; 72586621ddSJed Brown 736b8f6f9dSBarry Smith B->ops->getinfo = MatGetInfo_External; 741aef8b4cSStefano Zampini B->ops->matsolve = NULL; 75586621ddSJed Brown B->ops->view = MatView_CHOLMOD; 76586621ddSJed Brown B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_CHOLMOD; 77586621ddSJed Brown B->ops->destroy = MatDestroy_CHOLMOD; 786b8f6f9dSBarry Smith B->ops->getinfo = MatGetInfo_External; 792205254eSKarl Rupp 803ca39a21SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverType_C",MatFactorGetSolverType_seqaij_cholmod);CHKERRQ(ierr); 812205254eSKarl Rupp 82586621ddSJed Brown B->factortype = MAT_FACTOR_CHOLESKY; 83586621ddSJed Brown B->assembled = PETSC_TRUE; /* required by -ksp_view */ 84586621ddSJed Brown B->preallocated = PETSC_TRUE; 85586621ddSJed Brown 8600c67f3bSHong Zhang ierr = PetscFree(B->solvertype);CHKERRQ(ierr); 8700c67f3bSHong Zhang ierr = PetscStrallocpy(MATSOLVERCHOLMOD,&B->solvertype);CHKERRQ(ierr); 8800c67f3bSHong Zhang 89586621ddSJed Brown ierr = CholmodStart(B);CHKERRQ(ierr); 90586621ddSJed Brown *F = B; 91586621ddSJed Brown PetscFunctionReturn(0); 92586621ddSJed Brown } 93