1641875f9SMatthew G Knepley 2641875f9SMatthew G Knepley /* 3d515b9b4SShri Abhyankar Provides an interface to the CHOLMOD sparse solver available through SuiteSparse version 4.2.1 4641875f9SMatthew G Knepley 58999bf53SRichard Mills When built with PETSC_USE_64BIT_INDICES this will use Suitesparse_long as the 6641875f9SMatthew G Knepley integer type in UMFPACK, otherwise it will use int. This means 7641875f9SMatthew G Knepley all integers in this file as simply declared as PetscInt. Also it means 89e475b0dSSatish Balay that one cannot use 64BIT_INDICES on 32bit machines [as Suitesparse_long is 32bit only] 9641875f9SMatthew G Knepley 10641875f9SMatthew G Knepley */ 11641875f9SMatthew G Knepley 12c6db04a5SJed Brown #include <../src/mat/impls/sbaij/seq/sbaij.h> 13c6db04a5SJed Brown #include <../src/mat/impls/sbaij/seq/cholmod/cholmodimpl.h> 14641875f9SMatthew G Knepley 15641875f9SMatthew G Knepley /* 16641875f9SMatthew G Knepley This is a terrible hack, but it allows the error handler to retain a context. 17641875f9SMatthew G Knepley Note that this hack really cannot be made both reentrant and concurrent. 18641875f9SMatthew G Knepley */ 19641875f9SMatthew G Knepley static Mat static_F; 20641875f9SMatthew G Knepley 21641875f9SMatthew G Knepley static void CholmodErrorHandler(int status,const char *file,int line,const char *message) 22641875f9SMatthew G Knepley { 23641875f9SMatthew G Knepley PetscFunctionBegin; 24641875f9SMatthew G Knepley if (status > CHOLMOD_OK) { 259566063dSJacob Faibussowitsch PetscCallVoid(PetscInfo(static_F,"CHOLMOD warning %d at %s:%d: %s\n",status,file,line,message)); 26641875f9SMatthew G Knepley } else if (status == CHOLMOD_OK) { /* Documentation says this can happen, but why? */ 279566063dSJacob Faibussowitsch PetscCallVoid(PetscInfo(static_F,"CHOLMOD OK at %s:%d: %s\n",file,line,message)); 28641875f9SMatthew G Knepley } else { 299566063dSJacob Faibussowitsch PetscCallVoid(PetscErrorPrintf("CHOLMOD error %d at %s:%d: %s\n",status,file,line,message)); 30641875f9SMatthew G Knepley } 31641875f9SMatthew G Knepley PetscFunctionReturnVoid(); 32641875f9SMatthew G Knepley } 33641875f9SMatthew G Knepley 347087cfbeSBarry Smith PetscErrorCode CholmodStart(Mat F) 35641875f9SMatthew G Knepley { 366b8f6f9dSBarry Smith Mat_CHOLMOD *chol=(Mat_CHOLMOD*)F->data; 37641875f9SMatthew G Knepley cholmod_common *c; 38ace3abfcSBarry Smith PetscBool flg; 39641875f9SMatthew G Knepley 40641875f9SMatthew G Knepley PetscFunctionBegin; 41641875f9SMatthew G Knepley if (chol->common) PetscFunctionReturn(0); 429566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(1,&chol->common)); 439566063dSJacob Faibussowitsch PetscCall(!cholmod_X_start(chol->common)); 4426fbe8dcSKarl Rupp 45641875f9SMatthew G Knepley c = chol->common; 46641875f9SMatthew G Knepley c->error_handler = CholmodErrorHandler; 47641875f9SMatthew G Knepley 48641875f9SMatthew G Knepley #define CHOLMOD_OPTION_DOUBLE(name,help) do { \ 49641875f9SMatthew G Knepley PetscReal tmp = (PetscReal)c->name; \ 509566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-mat_cholmod_" #name,help,"None",tmp,&tmp,NULL)); \ 51641875f9SMatthew G Knepley c->name = (double)tmp; \ 52641875f9SMatthew G Knepley } while (0) 5326fbe8dcSKarl Rupp 54641875f9SMatthew G Knepley #define CHOLMOD_OPTION_INT(name,help) do { \ 55641875f9SMatthew G Knepley PetscInt tmp = (PetscInt)c->name; \ 569566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-mat_cholmod_" #name,help,"None",tmp,&tmp,NULL)); \ 57641875f9SMatthew G Knepley c->name = (int)tmp; \ 58641875f9SMatthew G Knepley } while (0) 5926fbe8dcSKarl Rupp 60641875f9SMatthew G Knepley #define CHOLMOD_OPTION_SIZE_T(name,help) do { \ 6154b3d318SStefano Zampini PetscReal tmp = (PetscInt)c->name; \ 629566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-mat_cholmod_" #name,help,"None",tmp,&tmp,NULL)); \ 6308401ef6SPierre Jolivet PetscCheck(tmp >= 0,PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_OUTOFRANGE,"value must be positive"); \ 64641875f9SMatthew G Knepley c->name = (size_t)tmp; \ 65641875f9SMatthew G Knepley } while (0) 6626fbe8dcSKarl Rupp 67b9eaa5e8SBarry Smith #define CHOLMOD_OPTION_BOOL(name,help) do { \ 68ace3abfcSBarry Smith PetscBool tmp = (PetscBool) !!c->name; \ 699566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-mat_cholmod_" #name,help,"None",tmp,&tmp,NULL)); \ 70641875f9SMatthew G Knepley c->name = (int)tmp; \ 71641875f9SMatthew G Knepley } while (0) 72641875f9SMatthew G Knepley 73d0609cedSBarry Smith PetscOptionsBegin(PetscObjectComm((PetscObject)F),((PetscObject)F)->prefix,"CHOLMOD Options","Mat"); 7454b3d318SStefano Zampini CHOLMOD_OPTION_INT(nmethods,"Number of different ordering methods to try"); 7526fbe8dcSKarl Rupp 76b9eaa5e8SBarry Smith #if defined(PETSC_USE_SUITESPARSE_GPU) 77b9eaa5e8SBarry Smith c->useGPU = 1; 78b9eaa5e8SBarry Smith CHOLMOD_OPTION_INT(useGPU,"Use GPU for BLAS 1, otherwise 0"); 7954b3d318SStefano Zampini CHOLMOD_OPTION_SIZE_T(maxGpuMemBytes,"Maximum memory to allocate on the GPU"); 8054b3d318SStefano Zampini CHOLMOD_OPTION_DOUBLE(maxGpuMemFraction,"Fraction of available GPU memory to allocate"); 81b9eaa5e8SBarry Smith #endif 82b9eaa5e8SBarry Smith 8354b3d318SStefano Zampini /* CHOLMOD handles first-time packing and refactor-packing separately, but we usually want them to be the same. */ 8454b3d318SStefano Zampini chol->pack = (PetscBool)c->final_pack; 859566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-mat_cholmod_pack","Pack factors after factorization [disable for frequent repeat factorization]","None",chol->pack,&chol->pack,NULL)); 86641875f9SMatthew G Knepley c->final_pack = (int)chol->pack; 87641875f9SMatthew G Knepley 88641875f9SMatthew G Knepley CHOLMOD_OPTION_DOUBLE(dbound,"Minimum absolute value of diagonal entries of D"); 89641875f9SMatthew G Knepley CHOLMOD_OPTION_DOUBLE(grow0,"Global growth ratio when factors are modified"); 90641875f9SMatthew G Knepley CHOLMOD_OPTION_DOUBLE(grow1,"Column growth ratio when factors are modified"); 91641875f9SMatthew G Knepley CHOLMOD_OPTION_SIZE_T(grow2,"Affine column growth constant when factors are modified"); 92641875f9SMatthew G Knepley CHOLMOD_OPTION_SIZE_T(maxrank,"Max rank of update, larger values are faster but use more memory [2,4,8]"); 93641875f9SMatthew G Knepley { 94641875f9SMatthew G Knepley static const char *const list[] = {"SIMPLICIAL","AUTO","SUPERNODAL","MatCholmodFactorType","MAT_CHOLMOD_FACTOR_",0}; 959566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-mat_cholmod_factor","Factorization method","None",list,(PetscEnum)c->supernodal,(PetscEnum*)&c->supernodal,NULL)); 96641875f9SMatthew G Knepley } 97641875f9SMatthew G Knepley if (c->supernodal) CHOLMOD_OPTION_DOUBLE(supernodal_switch,"flop/nnz_L threshold for switching to supernodal factorization"); 98b9eaa5e8SBarry Smith CHOLMOD_OPTION_BOOL(final_asis,"Leave factors \"as is\""); 99b9eaa5e8SBarry Smith CHOLMOD_OPTION_BOOL(final_pack,"Pack the columns when finished (use FALSE if the factors will be updated later)"); 100641875f9SMatthew G Knepley if (!c->final_asis) { 101b9eaa5e8SBarry Smith CHOLMOD_OPTION_BOOL(final_super,"Leave supernodal factors instead of converting to simplicial"); 102b9eaa5e8SBarry Smith CHOLMOD_OPTION_BOOL(final_ll,"Turn LDL' factorization into LL'"); 103b9eaa5e8SBarry Smith CHOLMOD_OPTION_BOOL(final_monotonic,"Ensure columns are monotonic when done"); 104b9eaa5e8SBarry Smith CHOLMOD_OPTION_BOOL(final_resymbol,"Remove numerically zero values resulting from relaxed supernodal amalgamation"); 105641875f9SMatthew G Knepley } 106641875f9SMatthew G Knepley { 107641875f9SMatthew G Knepley PetscReal tmp[] = {(PetscReal)c->zrelax[0],(PetscReal)c->zrelax[1],(PetscReal)c->zrelax[2]}; 108641875f9SMatthew G Knepley PetscInt n = 3; 1099566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-mat_cholmod_zrelax","3 real supernodal relaxed amalgamation parameters","None",tmp,&n,&flg)); 110aed4548fSBarry Smith PetscCheck(!flg || n == 3,PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_OUTOFRANGE,"must provide exactly 3 parameters to -mat_cholmod_zrelax"); 111641875f9SMatthew G Knepley if (flg) while (n--) c->zrelax[n] = (double)tmp[n]; 112641875f9SMatthew G Knepley } 113641875f9SMatthew G Knepley { 114641875f9SMatthew G Knepley PetscInt n,tmp[] = {(PetscInt)c->nrelax[0],(PetscInt)c->nrelax[1],(PetscInt)c->nrelax[2]}; 1159566063dSJacob Faibussowitsch PetscCall(PetscOptionsIntArray("-mat_cholmod_nrelax","3 size_t supernodal relaxed amalgamation parameters","None",tmp,&n,&flg)); 116aed4548fSBarry Smith PetscCheck(!flg || n == 3,PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_OUTOFRANGE,"must provide exactly 3 parameters to -mat_cholmod_nrelax"); 117641875f9SMatthew G Knepley if (flg) while (n--) c->nrelax[n] = (size_t)tmp[n]; 118641875f9SMatthew G Knepley } 119b9eaa5e8SBarry Smith CHOLMOD_OPTION_BOOL(prefer_upper,"Work with upper triangular form [faster when using fill-reducing ordering, slower in natural ordering]"); 120b9eaa5e8SBarry Smith CHOLMOD_OPTION_BOOL(default_nesdis,"Use NESDIS instead of METIS for nested dissection"); 121641875f9SMatthew G Knepley CHOLMOD_OPTION_INT(print,"Verbosity level"); 122d0609cedSBarry Smith PetscOptionsEnd(); 123641875f9SMatthew G Knepley PetscFunctionReturn(0); 124641875f9SMatthew G Knepley } 125641875f9SMatthew G Knepley 126218db3c1SStefano Zampini static PetscErrorCode MatWrapCholmod_seqsbaij(Mat A,PetscBool values,cholmod_sparse *C,PetscBool *aijalloc,PetscBool *valloc) 127641875f9SMatthew G Knepley { 128641875f9SMatthew G Knepley Mat_SeqSBAIJ *sbaij = (Mat_SeqSBAIJ*)A->data; 129218db3c1SStefano Zampini PetscBool vallocin = PETSC_FALSE; 130641875f9SMatthew G Knepley 131641875f9SMatthew G Knepley PetscFunctionBegin; 1329566063dSJacob Faibussowitsch PetscCall(PetscMemzero(C,sizeof(*C))); 133641875f9SMatthew G Knepley /* CHOLMOD uses column alignment, SBAIJ stores the upper factor, so we pass it on as a lower factor, swapping the meaning of row and column */ 134641875f9SMatthew G Knepley C->nrow = (size_t)A->cmap->n; 135641875f9SMatthew G Knepley C->ncol = (size_t)A->rmap->n; 136641875f9SMatthew G Knepley C->nzmax = (size_t)sbaij->maxnz; 137641875f9SMatthew G Knepley C->p = sbaij->i; 138641875f9SMatthew G Knepley C->i = sbaij->j; 139218db3c1SStefano Zampini if (values) { 140218db3c1SStefano Zampini #if defined(PETSC_USE_COMPLEX) 141218db3c1SStefano Zampini /* we need to pass CHOLMOD the conjugate matrix */ 142218db3c1SStefano Zampini PetscScalar *v; 143218db3c1SStefano Zampini PetscInt i; 144218db3c1SStefano Zampini 1459566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(sbaij->maxnz,&v)); 146218db3c1SStefano Zampini for (i = 0; i < sbaij->maxnz; i++) v[i] = PetscConj(sbaij->a[i]); 147218db3c1SStefano Zampini C->x = v; 148218db3c1SStefano Zampini vallocin = PETSC_TRUE; 149218db3c1SStefano Zampini #else 150641875f9SMatthew G Knepley C->x = sbaij->a; 151218db3c1SStefano Zampini #endif 152218db3c1SStefano Zampini } 153641875f9SMatthew G Knepley C->stype = -1; 154641875f9SMatthew G Knepley C->itype = CHOLMOD_INT_TYPE; 155218db3c1SStefano Zampini C->xtype = values ? CHOLMOD_SCALAR_TYPE : CHOLMOD_PATTERN; 156641875f9SMatthew G Knepley C->dtype = CHOLMOD_DOUBLE; 157641875f9SMatthew G Knepley C->sorted = 1; 158641875f9SMatthew G Knepley C->packed = 1; 159641875f9SMatthew G Knepley *aijalloc = PETSC_FALSE; 160218db3c1SStefano Zampini *valloc = vallocin; 161641875f9SMatthew G Knepley PetscFunctionReturn(0); 162641875f9SMatthew G Knepley } 163641875f9SMatthew G Knepley 164218db3c1SStefano Zampini #define GET_ARRAY_READ 0 165218db3c1SStefano Zampini #define GET_ARRAY_WRITE 1 166218db3c1SStefano Zampini 167a2fc1e05SToby Isaac PetscErrorCode VecWrapCholmod(Vec X,PetscInt rw,cholmod_dense *Y) 168641875f9SMatthew G Knepley { 169218db3c1SStefano Zampini PetscScalar *x; 170641875f9SMatthew G Knepley PetscInt n; 171641875f9SMatthew G Knepley 172641875f9SMatthew G Knepley PetscFunctionBegin; 1739566063dSJacob Faibussowitsch PetscCall(PetscMemzero(Y,sizeof(*Y))); 174218db3c1SStefano Zampini switch (rw) { 175218db3c1SStefano Zampini case GET_ARRAY_READ: 1769566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X,(const PetscScalar**)&x)); 177218db3c1SStefano Zampini break; 178218db3c1SStefano Zampini case GET_ARRAY_WRITE: 1799566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(X,&x)); 180218db3c1SStefano Zampini break; 181218db3c1SStefano Zampini default: 18298921bdaSJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Case %" PetscInt_FMT " not handled",rw); 183218db3c1SStefano Zampini break; 184218db3c1SStefano Zampini } 1859566063dSJacob Faibussowitsch PetscCall(VecGetSize(X,&n)); 18626fbe8dcSKarl Rupp 187218db3c1SStefano Zampini Y->x = x; 188641875f9SMatthew G Knepley Y->nrow = n; 189641875f9SMatthew G Knepley Y->ncol = 1; 190641875f9SMatthew G Knepley Y->nzmax = n; 191641875f9SMatthew G Knepley Y->d = n; 192641875f9SMatthew G Knepley Y->xtype = CHOLMOD_SCALAR_TYPE; 193641875f9SMatthew G Knepley Y->dtype = CHOLMOD_DOUBLE; 194641875f9SMatthew G Knepley PetscFunctionReturn(0); 195641875f9SMatthew G Knepley } 196641875f9SMatthew G Knepley 197a2fc1e05SToby Isaac PetscErrorCode VecUnWrapCholmod(Vec X,PetscInt rw,cholmod_dense *Y) 198d9ca1df4SBarry Smith { 199d9ca1df4SBarry Smith PetscFunctionBegin; 200218db3c1SStefano Zampini switch (rw) { 201218db3c1SStefano Zampini case GET_ARRAY_READ: 2029566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X,(const PetscScalar**)&Y->x)); 203218db3c1SStefano Zampini break; 204218db3c1SStefano Zampini case GET_ARRAY_WRITE: 2059566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(X,(PetscScalar**)&Y->x)); 206218db3c1SStefano Zampini break; 207218db3c1SStefano Zampini default: 20898921bdaSJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Case %" PetscInt_FMT " not handled",rw); 209218db3c1SStefano Zampini break; 210218db3c1SStefano Zampini } 211218db3c1SStefano Zampini PetscFunctionReturn(0); 212218db3c1SStefano Zampini } 213218db3c1SStefano Zampini 214a2fc1e05SToby Isaac PetscErrorCode MatDenseWrapCholmod(Mat X,PetscInt rw,cholmod_dense *Y) 215218db3c1SStefano Zampini { 216218db3c1SStefano Zampini PetscScalar *x; 217218db3c1SStefano Zampini PetscInt m,n,lda; 218218db3c1SStefano Zampini 219218db3c1SStefano Zampini PetscFunctionBegin; 2209566063dSJacob Faibussowitsch PetscCall(PetscMemzero(Y,sizeof(*Y))); 221218db3c1SStefano Zampini switch (rw) { 222218db3c1SStefano Zampini case GET_ARRAY_READ: 2239566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayRead(X,(const PetscScalar**)&x)); 224218db3c1SStefano Zampini break; 225218db3c1SStefano Zampini case GET_ARRAY_WRITE: 2269566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayWrite(X,&x)); 227218db3c1SStefano Zampini break; 228218db3c1SStefano Zampini default: 22998921bdaSJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Case %" PetscInt_FMT " not handled",rw); 230218db3c1SStefano Zampini break; 231218db3c1SStefano Zampini } 2329566063dSJacob Faibussowitsch PetscCall(MatDenseGetLDA(X,&lda)); 2339566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(X,&m,&n)); 234218db3c1SStefano Zampini 235218db3c1SStefano Zampini Y->x = x; 236218db3c1SStefano Zampini Y->nrow = m; 237218db3c1SStefano Zampini Y->ncol = n; 238218db3c1SStefano Zampini Y->nzmax = lda*n; 239218db3c1SStefano Zampini Y->d = lda; 240218db3c1SStefano Zampini Y->xtype = CHOLMOD_SCALAR_TYPE; 241218db3c1SStefano Zampini Y->dtype = CHOLMOD_DOUBLE; 242218db3c1SStefano Zampini PetscFunctionReturn(0); 243218db3c1SStefano Zampini } 244218db3c1SStefano Zampini 245a2fc1e05SToby Isaac PetscErrorCode MatDenseUnWrapCholmod(Mat X,PetscInt rw,cholmod_dense *Y) 246218db3c1SStefano Zampini { 247218db3c1SStefano Zampini PetscFunctionBegin; 248218db3c1SStefano Zampini switch (rw) { 249218db3c1SStefano Zampini case GET_ARRAY_READ: 2509566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayRead(X,(const PetscScalar**)&Y->x)); 251218db3c1SStefano Zampini break; 252218db3c1SStefano Zampini case GET_ARRAY_WRITE: 253218db3c1SStefano Zampini /* we don't have MatDenseRestoreArrayWrite */ 2549566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArray(X,(PetscScalar**)&Y->x)); 255218db3c1SStefano Zampini break; 256218db3c1SStefano Zampini default: 25798921bdaSJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Case %" PetscInt_FMT " not handled",rw); 258218db3c1SStefano Zampini break; 259218db3c1SStefano Zampini } 260d9ca1df4SBarry Smith PetscFunctionReturn(0); 261d9ca1df4SBarry Smith } 262d9ca1df4SBarry Smith 263eb9872f6SBarry Smith PETSC_INTERN PetscErrorCode MatDestroy_CHOLMOD(Mat F) 264641875f9SMatthew G Knepley { 2656b8f6f9dSBarry Smith Mat_CHOLMOD *chol=(Mat_CHOLMOD*)F->data; 266641875f9SMatthew G Knepley 267641875f9SMatthew G Knepley PetscFunctionBegin; 268a2fc1e05SToby Isaac if (chol->spqrfact) { 2699566063dSJacob Faibussowitsch PetscCall(!SuiteSparseQR_C_free(&chol->spqrfact, chol->common)); 270a2fc1e05SToby Isaac } 271a2fc1e05SToby Isaac if (chol->factor) { 2729566063dSJacob Faibussowitsch PetscCall(!cholmod_X_free_factor(&chol->factor,chol->common)); 273a2fc1e05SToby Isaac } 274a2fc1e05SToby Isaac if (chol->common->itype == CHOLMOD_INT) { 2759566063dSJacob Faibussowitsch PetscCall(!cholmod_finish(chol->common)); 276a2fc1e05SToby Isaac } else { 2779566063dSJacob Faibussowitsch PetscCall(!cholmod_l_finish(chol->common)); 278a2fc1e05SToby Isaac } 2799566063dSJacob Faibussowitsch PetscCall(PetscFree(chol->common)); 2809566063dSJacob Faibussowitsch PetscCall(PetscFree(chol->matrix)); 2819566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)F,"MatFactorGetSolverType_C",NULL)); 282*2e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)F,"MatQRFactorSymbolic_C",NULL)); 283*2e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)F,"MatQRFactorNumeric_C",NULL)); 2849566063dSJacob Faibussowitsch PetscCall(PetscFree(F->data)); 285641875f9SMatthew G Knepley PetscFunctionReturn(0); 286641875f9SMatthew G Knepley } 287641875f9SMatthew G Knepley 288641875f9SMatthew G Knepley static PetscErrorCode MatSolve_CHOLMOD(Mat,Vec,Vec); 289218db3c1SStefano Zampini static PetscErrorCode MatMatSolve_CHOLMOD(Mat,Mat,Mat); 290641875f9SMatthew G Knepley 291fcd503bcSBarry Smith /*static const char *const CholmodOrderingMethods[] = {"User","AMD","METIS","NESDIS(default)","Natural","NESDIS(small=20000)","NESDIS(small=4,no constrained)","NESDIS()"};*/ 292641875f9SMatthew G Knepley 293860c79edSBarry Smith static PetscErrorCode MatView_Info_CHOLMOD(Mat F,PetscViewer viewer) 294641875f9SMatthew G Knepley { 2956b8f6f9dSBarry Smith Mat_CHOLMOD *chol = (Mat_CHOLMOD*)F->data; 296641875f9SMatthew G Knepley const cholmod_common *c = chol->common; 297641875f9SMatthew G Knepley PetscInt i; 298641875f9SMatthew G Knepley 299641875f9SMatthew G Knepley PetscFunctionBegin; 300641875f9SMatthew G Knepley if (F->ops->solve != MatSolve_CHOLMOD) PetscFunctionReturn(0); 3019566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"CHOLMOD run parameters:\n")); 3029566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 3039566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Pack factors after symbolic factorization: %s\n",chol->pack ? "TRUE" : "FALSE")); 3049566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Common.dbound %g (Smallest absolute value of diagonal entries of D)\n",c->dbound)); 3059566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Common.grow0 %g\n",c->grow0)); 3069566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Common.grow1 %g\n",c->grow1)); 3079566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Common.grow2 %u\n",(unsigned)c->grow2)); 3089566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Common.maxrank %u\n",(unsigned)c->maxrank)); 3099566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Common.supernodal_switch %g\n",c->supernodal_switch)); 3109566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Common.supernodal %d\n",c->supernodal)); 3119566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Common.final_asis %d\n",c->final_asis)); 3129566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Common.final_super %d\n",c->final_super)); 3139566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Common.final_ll %d\n",c->final_ll)); 3149566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Common.final_pack %d\n",c->final_pack)); 3159566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Common.final_monotonic %d\n",c->final_monotonic)); 3169566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Common.final_resymbol %d\n",c->final_resymbol)); 3179566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Common.zrelax [%g,%g,%g]\n",c->zrelax[0],c->zrelax[1],c->zrelax[2])); 3189566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Common.nrelax [%u,%u,%u]\n",(unsigned)c->nrelax[0],(unsigned)c->nrelax[1],(unsigned)c->nrelax[2])); 3199566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Common.prefer_upper %d\n",c->prefer_upper)); 3209566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Common.print %d\n",c->print)); 321641875f9SMatthew G Knepley for (i=0; i<c->nmethods; i++) { 3229566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Ordering method %" PetscInt_FMT "%s:\n",i,i==c->selected ? " [SELECTED]" : "")); 323d0609cedSBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer," lnz %g, fl %g, prune_dense %g, prune_dense2 %g\n", 324d0609cedSBarry Smith c->method[i].lnz,c->method[i].fl,c->method[i].prune_dense,c->method[i].prune_dense2)); 325641875f9SMatthew G Knepley } 3269566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Common.postorder %d\n",c->postorder)); 3279566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Common.default_nesdis %d (use NESDIS instead of METIS for nested dissection)\n",c->default_nesdis)); 328641875f9SMatthew G Knepley /* Statistics */ 3299566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Common.fl %g (flop count from most recent analysis)\n",c->fl)); 3309566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Common.lnz %g (fundamental nz in L)\n",c->lnz)); 3319566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Common.anz %g\n",c->anz)); 3329566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Common.modfl %g (flop count from most recent update)\n",c->modfl)); 3339566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Common.malloc_count %g (number of live objects)\n",(double)c->malloc_count)); 3349566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Common.memory_usage %g (peak memory usage in bytes)\n",(double)c->memory_usage)); 3359566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Common.memory_inuse %g (current memory usage in bytes)\n",(double)c->memory_inuse)); 3369566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Common.nrealloc_col %g (number of column reallocations)\n",c->nrealloc_col)); 3379566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Common.nrealloc_factor %g (number of factor reallocations due to column reallocations)\n",c->nrealloc_factor)); 3389566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Common.ndbounds_hit %g (number of times diagonal was modified by dbound)\n",c->ndbounds_hit)); 3399566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Common.rowfacfl %g (number of flops in last call to cholmod_rowfac)\n",c->rowfacfl)); 3409566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Common.aatfl %g (number of flops to compute A(:,f)*A(:,f)')\n",c->aatfl)); 341b9eaa5e8SBarry Smith #if defined(PETSC_USE_SUITESPARSE_GPU) 3429566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Common.useGPU %d\n",c->useGPU)); 343b9eaa5e8SBarry Smith #endif 3449566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 345641875f9SMatthew G Knepley PetscFunctionReturn(0); 346641875f9SMatthew G Knepley } 347641875f9SMatthew G Knepley 348eb9872f6SBarry Smith PETSC_INTERN PetscErrorCode MatView_CHOLMOD(Mat F,PetscViewer viewer) 349641875f9SMatthew G Knepley { 350ace3abfcSBarry Smith PetscBool iascii; 351641875f9SMatthew G Knepley PetscViewerFormat format; 352641875f9SMatthew G Knepley 353641875f9SMatthew G Knepley PetscFunctionBegin; 3549566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii)); 355641875f9SMatthew G Knepley if (iascii) { 3569566063dSJacob Faibussowitsch PetscCall(PetscViewerGetFormat(viewer,&format)); 357641875f9SMatthew G Knepley if (format == PETSC_VIEWER_ASCII_INFO) { 3589566063dSJacob Faibussowitsch PetscCall(MatView_Info_CHOLMOD(F,viewer)); 359641875f9SMatthew G Knepley } 360641875f9SMatthew G Knepley } 361641875f9SMatthew G Knepley PetscFunctionReturn(0); 362641875f9SMatthew G Knepley } 363641875f9SMatthew G Knepley 364641875f9SMatthew G Knepley static PetscErrorCode MatSolve_CHOLMOD(Mat F,Vec B,Vec X) 365641875f9SMatthew G Knepley { 3666b8f6f9dSBarry Smith Mat_CHOLMOD *chol = (Mat_CHOLMOD*)F->data; 367218db3c1SStefano Zampini cholmod_dense cholB,cholX,*X_handle,*Y_handle = NULL,*E_handle = NULL; 368641875f9SMatthew G Knepley 369641875f9SMatthew G Knepley PetscFunctionBegin; 370641875f9SMatthew G Knepley static_F = F; 3719566063dSJacob Faibussowitsch PetscCall(VecWrapCholmod(B,GET_ARRAY_READ,&cholB)); 3729566063dSJacob Faibussowitsch PetscCall(VecWrapCholmod(X,GET_ARRAY_WRITE,&cholX)); 373218db3c1SStefano Zampini X_handle = &cholX; 3749566063dSJacob Faibussowitsch PetscCall(!cholmod_X_solve2(CHOLMOD_A,chol->factor,&cholB,NULL,&X_handle,NULL,&Y_handle,&E_handle,chol->common)); 3759566063dSJacob Faibussowitsch PetscCall(!cholmod_X_free_dense(&Y_handle,chol->common)); 3769566063dSJacob Faibussowitsch PetscCall(!cholmod_X_free_dense(&E_handle,chol->common)); 3779566063dSJacob Faibussowitsch PetscCall(VecUnWrapCholmod(B,GET_ARRAY_READ,&cholB)); 3789566063dSJacob Faibussowitsch PetscCall(VecUnWrapCholmod(X,GET_ARRAY_WRITE,&cholX)); 3799566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(4.0*chol->common->lnz)); 380218db3c1SStefano Zampini PetscFunctionReturn(0); 381218db3c1SStefano Zampini } 382218db3c1SStefano Zampini 383218db3c1SStefano Zampini static PetscErrorCode MatMatSolve_CHOLMOD(Mat F,Mat B,Mat X) 384218db3c1SStefano Zampini { 385218db3c1SStefano Zampini Mat_CHOLMOD *chol = (Mat_CHOLMOD*)F->data; 386218db3c1SStefano Zampini cholmod_dense cholB,cholX,*X_handle,*Y_handle = NULL,*E_handle = NULL; 387218db3c1SStefano Zampini 388218db3c1SStefano Zampini PetscFunctionBegin; 389218db3c1SStefano Zampini static_F = F; 3909566063dSJacob Faibussowitsch PetscCall(MatDenseWrapCholmod(B,GET_ARRAY_READ,&cholB)); 3919566063dSJacob Faibussowitsch PetscCall(MatDenseWrapCholmod(X,GET_ARRAY_WRITE,&cholX)); 392218db3c1SStefano Zampini X_handle = &cholX; 3939566063dSJacob Faibussowitsch PetscCall(!cholmod_X_solve2(CHOLMOD_A,chol->factor,&cholB,NULL,&X_handle,NULL,&Y_handle,&E_handle,chol->common)); 3949566063dSJacob Faibussowitsch PetscCall(!cholmod_X_free_dense(&Y_handle,chol->common)); 3959566063dSJacob Faibussowitsch PetscCall(!cholmod_X_free_dense(&E_handle,chol->common)); 3969566063dSJacob Faibussowitsch PetscCall(MatDenseUnWrapCholmod(B,GET_ARRAY_READ,&cholB)); 3979566063dSJacob Faibussowitsch PetscCall(MatDenseUnWrapCholmod(X,GET_ARRAY_WRITE,&cholX)); 3989566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(4.0*B->cmap->n*chol->common->lnz)); 399641875f9SMatthew G Knepley PetscFunctionReturn(0); 400641875f9SMatthew G Knepley } 401641875f9SMatthew G Knepley 402641875f9SMatthew G Knepley static PetscErrorCode MatCholeskyFactorNumeric_CHOLMOD(Mat F,Mat A,const MatFactorInfo *info) 403641875f9SMatthew G Knepley { 4046b8f6f9dSBarry Smith Mat_CHOLMOD *chol = (Mat_CHOLMOD*)F->data; 405641875f9SMatthew G Knepley cholmod_sparse cholA; 406218db3c1SStefano Zampini PetscBool aijalloc,valloc; 407d0609cedSBarry Smith int err; 408641875f9SMatthew G Knepley 409641875f9SMatthew G Knepley PetscFunctionBegin; 4109566063dSJacob Faibussowitsch PetscCall((*chol->Wrap)(A,PETSC_TRUE,&cholA,&aijalloc,&valloc)); 411641875f9SMatthew G Knepley static_F = F; 412d0609cedSBarry Smith err = !cholmod_X_factorize(&cholA,chol->factor,chol->common); 413d0609cedSBarry Smith PetscCheck(!err,PetscObjectComm((PetscObject)F),PETSC_ERR_LIB,"CHOLMOD factorization failed with status %d",chol->common->status); 41408401ef6SPierre Jolivet PetscCheck(chol->common->status != CHOLMOD_NOT_POSDEF,PetscObjectComm((PetscObject)F),PETSC_ERR_MAT_CH_ZRPVT,"CHOLMOD detected that the matrix is not positive definite, failure at column %u",(unsigned)chol->factor->minor); 415641875f9SMatthew G Knepley 4169566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(chol->common->fl)); 4179566063dSJacob Faibussowitsch if (aijalloc) PetscCall(PetscFree2(cholA.p,cholA.i)); 4189566063dSJacob Faibussowitsch if (valloc) PetscCall(PetscFree(cholA.x)); 419218db3c1SStefano Zampini #if defined(PETSC_USE_SUITESPARSE_GPU) 4209566063dSJacob Faibussowitsch PetscCall(PetscLogGpuTimeAdd(chol->common->CHOLMOD_GPU_GEMM_TIME + chol->common->CHOLMOD_GPU_SYRK_TIME + chol->common->CHOLMOD_GPU_TRSM_TIME + chol->common->CHOLMOD_GPU_POTRF_TIME)); 421218db3c1SStefano Zampini #endif 422641875f9SMatthew G Knepley 423641875f9SMatthew G Knepley F->ops->solve = MatSolve_CHOLMOD; 424641875f9SMatthew G Knepley F->ops->solvetranspose = MatSolve_CHOLMOD; 425218db3c1SStefano Zampini F->ops->matsolve = MatMatSolve_CHOLMOD; 426218db3c1SStefano Zampini F->ops->matsolvetranspose = MatMatSolve_CHOLMOD; 427641875f9SMatthew G Knepley PetscFunctionReturn(0); 428641875f9SMatthew G Knepley } 429641875f9SMatthew G Knepley 430eb9872f6SBarry Smith PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_CHOLMOD(Mat F,Mat A,IS perm,const MatFactorInfo *info) 431641875f9SMatthew G Knepley { 4326b8f6f9dSBarry Smith Mat_CHOLMOD *chol = (Mat_CHOLMOD*)F->data; 433d0609cedSBarry Smith int err; 434641875f9SMatthew G Knepley cholmod_sparse cholA; 435218db3c1SStefano Zampini PetscBool aijalloc,valloc; 436641875f9SMatthew G Knepley PetscInt *fset = 0; 437641875f9SMatthew G Knepley size_t fsize = 0; 438641875f9SMatthew G Knepley 439641875f9SMatthew G Knepley PetscFunctionBegin; 4409566063dSJacob Faibussowitsch PetscCall((*chol->Wrap)(A,PETSC_FALSE,&cholA,&aijalloc,&valloc)); 441641875f9SMatthew G Knepley static_F = F; 442641875f9SMatthew G Knepley if (chol->factor) { 443d0609cedSBarry Smith err = !cholmod_X_resymbol(&cholA,fset,fsize,(int)chol->pack,chol->factor,chol->common); 444d0609cedSBarry Smith PetscCheck(!err,PetscObjectComm((PetscObject)F),PETSC_ERR_LIB,"CHOLMOD analysis failed with status %d",chol->common->status); 445641875f9SMatthew G Knepley } else if (perm) { 446641875f9SMatthew G Knepley const PetscInt *ip; 4479566063dSJacob Faibussowitsch PetscCall(ISGetIndices(perm,&ip)); 448641875f9SMatthew G Knepley chol->factor = cholmod_X_analyze_p(&cholA,(PetscInt*)ip,fset,fsize,chol->common); 44928b400f6SJacob Faibussowitsch PetscCheck(chol->factor,PetscObjectComm((PetscObject)F),PETSC_ERR_LIB,"CHOLMOD analysis failed using PETSc ordering with status %d",chol->common->status); 4509566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(perm,&ip)); 451641875f9SMatthew G Knepley } else { 452641875f9SMatthew G Knepley chol->factor = cholmod_X_analyze(&cholA,chol->common); 45328b400f6SJacob Faibussowitsch PetscCheck(chol->factor,PetscObjectComm((PetscObject)F),PETSC_ERR_LIB,"CHOLMOD analysis failed using internal ordering with status %d",chol->common->status); 454641875f9SMatthew G Knepley } 455641875f9SMatthew G Knepley 4569566063dSJacob Faibussowitsch if (aijalloc) PetscCall(PetscFree2(cholA.p,cholA.i)); 4579566063dSJacob Faibussowitsch if (valloc) PetscCall(PetscFree(cholA.x)); 458641875f9SMatthew G Knepley 459641875f9SMatthew G Knepley F->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_CHOLMOD; 460641875f9SMatthew G Knepley PetscFunctionReturn(0); 461641875f9SMatthew G Knepley } 462641875f9SMatthew G Knepley 463ea799195SBarry Smith static PetscErrorCode MatFactorGetSolverType_seqsbaij_cholmod(Mat A,MatSolverType *type) 464641875f9SMatthew G Knepley { 465641875f9SMatthew G Knepley PetscFunctionBegin; 466641875f9SMatthew G Knepley *type = MATSOLVERCHOLMOD; 467641875f9SMatthew G Knepley PetscFunctionReturn(0); 468641875f9SMatthew G Knepley } 469641875f9SMatthew G Knepley 470218db3c1SStefano Zampini PETSC_INTERN PetscErrorCode MatGetInfo_CHOLMOD(Mat F,MatInfoType flag,MatInfo *info) 471218db3c1SStefano Zampini { 472218db3c1SStefano Zampini Mat_CHOLMOD *chol = (Mat_CHOLMOD*)F->data; 473218db3c1SStefano Zampini 474218db3c1SStefano Zampini PetscFunctionBegin; 475218db3c1SStefano Zampini info->block_size = 1.0; 476218db3c1SStefano Zampini info->nz_allocated = chol->common->lnz; 477218db3c1SStefano Zampini info->nz_used = chol->common->lnz; 478218db3c1SStefano Zampini info->nz_unneeded = 0.0; 479218db3c1SStefano Zampini info->assemblies = 0.0; 480218db3c1SStefano Zampini info->mallocs = 0.0; 481218db3c1SStefano Zampini info->memory = chol->common->memory_inuse; 482218db3c1SStefano Zampini info->fill_ratio_given = 0; 483218db3c1SStefano Zampini info->fill_ratio_needed = 0; 484218db3c1SStefano Zampini info->factor_mallocs = chol->common->malloc_count; 485218db3c1SStefano Zampini PetscFunctionReturn(0); 486218db3c1SStefano Zampini } 487218db3c1SStefano Zampini 488641875f9SMatthew G Knepley /*MC 489ee300463SSatish Balay MATSOLVERCHOLMOD 490ee300463SSatish Balay 491ee300463SSatish Balay A matrix type providing direct solvers (Cholesky) for sequential matrices 492641875f9SMatthew G Knepley via the external package CHOLMOD. 493641875f9SMatthew G Knepley 494c2b89b5dSBarry Smith Use ./configure --download-suitesparse to install PETSc to use CHOLMOD 495c2b89b5dSBarry Smith 496f29f8b16SBarry Smith Use -pc_type cholesky -pc_factor_mat_solver_type cholmod to use this direct solver 497641875f9SMatthew G Knepley 498641875f9SMatthew G Knepley Consult CHOLMOD documentation for more information about the Common parameters 499641875f9SMatthew G Knepley which correspond to the options database keys below. 500641875f9SMatthew G Knepley 501641875f9SMatthew G Knepley Options Database Keys: 502e08999f5SMatthew G Knepley + -mat_cholmod_dbound <0> - Minimum absolute value of diagonal entries of D (None) 503e08999f5SMatthew G Knepley . -mat_cholmod_grow0 <1.2> - Global growth ratio when factors are modified (None) 504e08999f5SMatthew G Knepley . -mat_cholmod_grow1 <1.2> - Column growth ratio when factors are modified (None) 505e08999f5SMatthew G Knepley . -mat_cholmod_grow2 <5> - Affine column growth constant when factors are modified (None) 506e08999f5SMatthew G Knepley . -mat_cholmod_maxrank <8> - Max rank of update, larger values are faster but use more memory [2,4,8] (None) 507e08999f5SMatthew G Knepley . -mat_cholmod_factor <AUTO> - (choose one of) SIMPLICIAL AUTO SUPERNODAL 508e08999f5SMatthew G Knepley . -mat_cholmod_supernodal_switch <40> - flop/nnz_L threshold for switching to supernodal factorization (None) 509e08999f5SMatthew G Knepley . -mat_cholmod_final_asis <TRUE> - Leave factors "as is" (None) 510e08999f5SMatthew G Knepley . -mat_cholmod_final_pack <TRUE> - Pack the columns when finished (use FALSE if the factors will be updated later) (None) 511e08999f5SMatthew G Knepley . -mat_cholmod_zrelax <0.8> - 3 real supernodal relaxed amalgamation parameters (None) 512e08999f5SMatthew G Knepley . -mat_cholmod_nrelax <4> - 3 size_t supernodal relaxed amalgamation parameters (None) 513e08999f5SMatthew G Knepley . -mat_cholmod_prefer_upper <TRUE> - Work with upper triangular form (faster when using fill-reducing ordering, slower in natural ordering) (None) 5142c7c0729SBarry Smith . -mat_cholmod_print <3> - Verbosity level (None) 515ee300463SSatish Balay - -mat_ordering_type internal - Use the ordering provided by Cholmod 516641875f9SMatthew G Knepley 517641875f9SMatthew G Knepley Level: beginner 518641875f9SMatthew G Knepley 519a364b7d2SBarry Smith Note: CHOLMOD is part of SuiteSparse http://faculty.cse.tamu.edu/davis/suitesparse.html 520a364b7d2SBarry Smith 521db781477SPatrick Sanan .seealso: `PCCHOLESKY`, `PCFactorSetMatSolverType()`, `MatSolverType` 522641875f9SMatthew G Knepley M*/ 523b2573a8aSBarry Smith 524db87b0f2SBarry Smith PETSC_INTERN PetscErrorCode MatGetFactor_seqsbaij_cholmod(Mat A,MatFactorType ftype,Mat *F) 525641875f9SMatthew G Knepley { 526641875f9SMatthew G Knepley Mat B; 527641875f9SMatthew G Knepley Mat_CHOLMOD *chol; 528641875f9SMatthew G Knepley PetscInt m=A->rmap->n,n=A->cmap->n,bs; 529218db3c1SStefano Zampini const char *prefix; 530641875f9SMatthew G Knepley 531641875f9SMatthew G Knepley PetscFunctionBegin; 5329566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(A,&bs)); 53308401ef6SPierre Jolivet PetscCheck(bs == 1,PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"CHOLMOD only supports block size=1, given %" PetscInt_FMT,bs); 534218db3c1SStefano Zampini #if defined(PETSC_USE_COMPLEX) 53528b400f6SJacob Faibussowitsch PetscCheck(A->hermitian,PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Only for Hermitian matrices"); 536218db3c1SStefano Zampini #endif 537641875f9SMatthew G Knepley /* Create the factorization matrix F */ 5389566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A),&B)); 5399566063dSJacob Faibussowitsch PetscCall(MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,m,n)); 5409566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy("cholmod",&((PetscObject)B)->type_name)); 5419566063dSJacob Faibussowitsch PetscCall(MatGetOptionsPrefix(A,&prefix)); 5429566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(B,prefix)); 5439566063dSJacob Faibussowitsch PetscCall(MatSetUp(B)); 5449566063dSJacob Faibussowitsch PetscCall(PetscNewLog(B,&chol)); 54526fbe8dcSKarl Rupp 546641875f9SMatthew G Knepley chol->Wrap = MatWrapCholmod_seqsbaij; 5476b8f6f9dSBarry Smith B->data = chol; 548641875f9SMatthew G Knepley 549218db3c1SStefano Zampini B->ops->getinfo = MatGetInfo_CHOLMOD; 550641875f9SMatthew G Knepley B->ops->view = MatView_CHOLMOD; 551641875f9SMatthew G Knepley B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_CHOLMOD; 552641875f9SMatthew G Knepley B->ops->destroy = MatDestroy_CHOLMOD; 5539566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverType_C",MatFactorGetSolverType_seqsbaij_cholmod)); 554641875f9SMatthew G Knepley B->factortype = MAT_FACTOR_CHOLESKY; 555218db3c1SStefano Zampini B->assembled = PETSC_TRUE; 556641875f9SMatthew G Knepley B->preallocated = PETSC_TRUE; 557641875f9SMatthew G Knepley 5589566063dSJacob Faibussowitsch PetscCall(CholmodStart(B)); 55900c67f3bSHong Zhang 5609566063dSJacob Faibussowitsch PetscCall(PetscFree(B->solvertype)); 5619566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(MATSOLVERCHOLMOD,&B->solvertype)); 562f73b0415SBarry Smith B->canuseordering = PETSC_TRUE; 5639566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(MATORDERINGEXTERNAL,(char**)&B->preferredordering[MAT_FACTOR_CHOLESKY])); 564641875f9SMatthew G Knepley *F = B; 565641875f9SMatthew G Knepley PetscFunctionReturn(0); 566641875f9SMatthew G Knepley } 567