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 34641875f9SMatthew G Knepley #define CHOLMOD_OPTION_DOUBLE(name,help) do { \ 35641875f9SMatthew G Knepley PetscReal tmp = (PetscReal)c->name; \ 369566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-mat_cholmod_" #name,help,"None",tmp,&tmp,NULL)); \ 37641875f9SMatthew G Knepley c->name = (double)tmp; \ 38641875f9SMatthew G Knepley } while (0) 3926fbe8dcSKarl Rupp 40641875f9SMatthew G Knepley #define CHOLMOD_OPTION_INT(name,help) do { \ 41641875f9SMatthew G Knepley PetscInt tmp = (PetscInt)c->name; \ 429566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-mat_cholmod_" #name,help,"None",tmp,&tmp,NULL)); \ 43641875f9SMatthew G Knepley c->name = (int)tmp; \ 44641875f9SMatthew G Knepley } while (0) 4526fbe8dcSKarl Rupp 46641875f9SMatthew G Knepley #define CHOLMOD_OPTION_SIZE_T(name,help) do { \ 4754b3d318SStefano Zampini PetscReal tmp = (PetscInt)c->name; \ 489566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-mat_cholmod_" #name,help,"None",tmp,&tmp,NULL)); \ 4908401ef6SPierre Jolivet PetscCheck(tmp >= 0,PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_OUTOFRANGE,"value must be positive"); \ 50641875f9SMatthew G Knepley c->name = (size_t)tmp; \ 51641875f9SMatthew G Knepley } while (0) 5226fbe8dcSKarl Rupp 53b9eaa5e8SBarry Smith #define CHOLMOD_OPTION_BOOL(name,help) do { \ 54ace3abfcSBarry Smith PetscBool tmp = (PetscBool) !!c->name; \ 559566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-mat_cholmod_" #name,help,"None",tmp,&tmp,NULL)); \ 56641875f9SMatthew G Knepley c->name = (int)tmp; \ 57641875f9SMatthew G Knepley } while (0) 58641875f9SMatthew G Knepley 5926cc229bSBarry Smith static PetscErrorCode CholmodSetOptions(Mat F) 6026cc229bSBarry Smith { 6126cc229bSBarry Smith Mat_CHOLMOD *chol=(Mat_CHOLMOD*)F->data; 6226cc229bSBarry Smith cholmod_common *c=chol->common; 6326cc229bSBarry Smith PetscBool flg; 6426cc229bSBarry Smith 6526cc229bSBarry Smith PetscFunctionBegin; 66d0609cedSBarry Smith PetscOptionsBegin(PetscObjectComm((PetscObject)F),((PetscObject)F)->prefix,"CHOLMOD Options","Mat"); 6754b3d318SStefano Zampini CHOLMOD_OPTION_INT(nmethods,"Number of different ordering methods to try"); 6826fbe8dcSKarl Rupp 69b9eaa5e8SBarry Smith #if defined(PETSC_USE_SUITESPARSE_GPU) 70b9eaa5e8SBarry Smith c->useGPU = 1; 71b9eaa5e8SBarry Smith CHOLMOD_OPTION_INT(useGPU,"Use GPU for BLAS 1, otherwise 0"); 7254b3d318SStefano Zampini CHOLMOD_OPTION_SIZE_T(maxGpuMemBytes,"Maximum memory to allocate on the GPU"); 7354b3d318SStefano Zampini CHOLMOD_OPTION_DOUBLE(maxGpuMemFraction,"Fraction of available GPU memory to allocate"); 74b9eaa5e8SBarry Smith #endif 75b9eaa5e8SBarry Smith 7654b3d318SStefano Zampini /* CHOLMOD handles first-time packing and refactor-packing separately, but we usually want them to be the same. */ 7754b3d318SStefano Zampini chol->pack = (PetscBool)c->final_pack; 789566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-mat_cholmod_pack","Pack factors after factorization [disable for frequent repeat factorization]","None",chol->pack,&chol->pack,NULL)); 79641875f9SMatthew G Knepley c->final_pack = (int)chol->pack; 80641875f9SMatthew G Knepley 81641875f9SMatthew G Knepley CHOLMOD_OPTION_DOUBLE(dbound,"Minimum absolute value of diagonal entries of D"); 82641875f9SMatthew G Knepley CHOLMOD_OPTION_DOUBLE(grow0,"Global growth ratio when factors are modified"); 83641875f9SMatthew G Knepley CHOLMOD_OPTION_DOUBLE(grow1,"Column growth ratio when factors are modified"); 84641875f9SMatthew G Knepley CHOLMOD_OPTION_SIZE_T(grow2,"Affine column growth constant when factors are modified"); 85641875f9SMatthew G Knepley CHOLMOD_OPTION_SIZE_T(maxrank,"Max rank of update, larger values are faster but use more memory [2,4,8]"); 86641875f9SMatthew G Knepley { 87641875f9SMatthew G Knepley static const char *const list[] = {"SIMPLICIAL","AUTO","SUPERNODAL","MatCholmodFactorType","MAT_CHOLMOD_FACTOR_",0}; 889566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-mat_cholmod_factor","Factorization method","None",list,(PetscEnum)c->supernodal,(PetscEnum*)&c->supernodal,NULL)); 89641875f9SMatthew G Knepley } 90641875f9SMatthew G Knepley if (c->supernodal) CHOLMOD_OPTION_DOUBLE(supernodal_switch,"flop/nnz_L threshold for switching to supernodal factorization"); 91b9eaa5e8SBarry Smith CHOLMOD_OPTION_BOOL(final_asis,"Leave factors \"as is\""); 92b9eaa5e8SBarry Smith CHOLMOD_OPTION_BOOL(final_pack,"Pack the columns when finished (use FALSE if the factors will be updated later)"); 93641875f9SMatthew G Knepley if (!c->final_asis) { 94b9eaa5e8SBarry Smith CHOLMOD_OPTION_BOOL(final_super,"Leave supernodal factors instead of converting to simplicial"); 95b9eaa5e8SBarry Smith CHOLMOD_OPTION_BOOL(final_ll,"Turn LDL' factorization into LL'"); 96b9eaa5e8SBarry Smith CHOLMOD_OPTION_BOOL(final_monotonic,"Ensure columns are monotonic when done"); 97b9eaa5e8SBarry Smith CHOLMOD_OPTION_BOOL(final_resymbol,"Remove numerically zero values resulting from relaxed supernodal amalgamation"); 98641875f9SMatthew G Knepley } 99641875f9SMatthew G Knepley { 100641875f9SMatthew G Knepley PetscReal tmp[] = {(PetscReal)c->zrelax[0],(PetscReal)c->zrelax[1],(PetscReal)c->zrelax[2]}; 101641875f9SMatthew G Knepley PetscInt n = 3; 1029566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-mat_cholmod_zrelax","3 real supernodal relaxed amalgamation parameters","None",tmp,&n,&flg)); 103aed4548fSBarry Smith PetscCheck(!flg || n == 3,PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_OUTOFRANGE,"must provide exactly 3 parameters to -mat_cholmod_zrelax"); 104641875f9SMatthew G Knepley if (flg) while (n--) c->zrelax[n] = (double)tmp[n]; 105641875f9SMatthew G Knepley } 106641875f9SMatthew G Knepley { 107641875f9SMatthew G Knepley PetscInt n,tmp[] = {(PetscInt)c->nrelax[0],(PetscInt)c->nrelax[1],(PetscInt)c->nrelax[2]}; 1089566063dSJacob Faibussowitsch PetscCall(PetscOptionsIntArray("-mat_cholmod_nrelax","3 size_t supernodal relaxed amalgamation parameters","None",tmp,&n,&flg)); 109aed4548fSBarry Smith PetscCheck(!flg || n == 3,PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_OUTOFRANGE,"must provide exactly 3 parameters to -mat_cholmod_nrelax"); 110641875f9SMatthew G Knepley if (flg) while (n--) c->nrelax[n] = (size_t)tmp[n]; 111641875f9SMatthew G Knepley } 112b9eaa5e8SBarry Smith CHOLMOD_OPTION_BOOL(prefer_upper,"Work with upper triangular form [faster when using fill-reducing ordering, slower in natural ordering]"); 113b9eaa5e8SBarry Smith CHOLMOD_OPTION_BOOL(default_nesdis,"Use NESDIS instead of METIS for nested dissection"); 114641875f9SMatthew G Knepley CHOLMOD_OPTION_INT(print,"Verbosity level"); 115d0609cedSBarry Smith PetscOptionsEnd(); 116641875f9SMatthew G Knepley PetscFunctionReturn(0); 117641875f9SMatthew G Knepley } 118641875f9SMatthew G Knepley 11926cc229bSBarry Smith PetscErrorCode CholmodStart(Mat F) 12026cc229bSBarry Smith { 12126cc229bSBarry Smith Mat_CHOLMOD *chol=(Mat_CHOLMOD*)F->data; 12226cc229bSBarry Smith cholmod_common *c; 12326cc229bSBarry Smith 12426cc229bSBarry Smith PetscFunctionBegin; 12526cc229bSBarry Smith if (chol->common) PetscFunctionReturn(0); 12626cc229bSBarry Smith PetscCall(PetscMalloc1(1,&chol->common)); 12726cc229bSBarry Smith PetscCall(!cholmod_X_start(chol->common)); 12826cc229bSBarry Smith 12926cc229bSBarry Smith c = chol->common; 13026cc229bSBarry Smith c->error_handler = CholmodErrorHandler; 13126cc229bSBarry Smith PetscFunctionReturn(0); 13226cc229bSBarry Smith } 13326cc229bSBarry Smith 134218db3c1SStefano Zampini static PetscErrorCode MatWrapCholmod_seqsbaij(Mat A,PetscBool values,cholmod_sparse *C,PetscBool *aijalloc,PetscBool *valloc) 135641875f9SMatthew G Knepley { 136641875f9SMatthew G Knepley Mat_SeqSBAIJ *sbaij = (Mat_SeqSBAIJ*)A->data; 137218db3c1SStefano Zampini PetscBool vallocin = PETSC_FALSE; 138641875f9SMatthew G Knepley 139641875f9SMatthew G Knepley PetscFunctionBegin; 1409566063dSJacob Faibussowitsch PetscCall(PetscMemzero(C,sizeof(*C))); 141641875f9SMatthew 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 */ 142641875f9SMatthew G Knepley C->nrow = (size_t)A->cmap->n; 143641875f9SMatthew G Knepley C->ncol = (size_t)A->rmap->n; 144641875f9SMatthew G Knepley C->nzmax = (size_t)sbaij->maxnz; 145641875f9SMatthew G Knepley C->p = sbaij->i; 146641875f9SMatthew G Knepley C->i = sbaij->j; 147218db3c1SStefano Zampini if (values) { 148218db3c1SStefano Zampini #if defined(PETSC_USE_COMPLEX) 149218db3c1SStefano Zampini /* we need to pass CHOLMOD the conjugate matrix */ 150218db3c1SStefano Zampini PetscScalar *v; 151218db3c1SStefano Zampini PetscInt i; 152218db3c1SStefano Zampini 1539566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(sbaij->maxnz,&v)); 154218db3c1SStefano Zampini for (i = 0; i < sbaij->maxnz; i++) v[i] = PetscConj(sbaij->a[i]); 155218db3c1SStefano Zampini C->x = v; 156218db3c1SStefano Zampini vallocin = PETSC_TRUE; 157218db3c1SStefano Zampini #else 158641875f9SMatthew G Knepley C->x = sbaij->a; 159218db3c1SStefano Zampini #endif 160218db3c1SStefano Zampini } 161641875f9SMatthew G Knepley C->stype = -1; 162641875f9SMatthew G Knepley C->itype = CHOLMOD_INT_TYPE; 163218db3c1SStefano Zampini C->xtype = values ? CHOLMOD_SCALAR_TYPE : CHOLMOD_PATTERN; 164641875f9SMatthew G Knepley C->dtype = CHOLMOD_DOUBLE; 165641875f9SMatthew G Knepley C->sorted = 1; 166641875f9SMatthew G Knepley C->packed = 1; 167641875f9SMatthew G Knepley *aijalloc = PETSC_FALSE; 168218db3c1SStefano Zampini *valloc = vallocin; 169641875f9SMatthew G Knepley PetscFunctionReturn(0); 170641875f9SMatthew G Knepley } 171641875f9SMatthew G Knepley 172218db3c1SStefano Zampini #define GET_ARRAY_READ 0 173218db3c1SStefano Zampini #define GET_ARRAY_WRITE 1 174218db3c1SStefano Zampini 175a2fc1e05SToby Isaac PetscErrorCode VecWrapCholmod(Vec X,PetscInt rw,cholmod_dense *Y) 176641875f9SMatthew G Knepley { 177218db3c1SStefano Zampini PetscScalar *x; 178641875f9SMatthew G Knepley PetscInt n; 179641875f9SMatthew G Knepley 180641875f9SMatthew G Knepley PetscFunctionBegin; 1819566063dSJacob Faibussowitsch PetscCall(PetscMemzero(Y,sizeof(*Y))); 182218db3c1SStefano Zampini switch (rw) { 183218db3c1SStefano Zampini case GET_ARRAY_READ: 1849566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X,(const PetscScalar**)&x)); 185218db3c1SStefano Zampini break; 186218db3c1SStefano Zampini case GET_ARRAY_WRITE: 1879566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(X,&x)); 188218db3c1SStefano Zampini break; 189218db3c1SStefano Zampini default: 19098921bdaSJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Case %" PetscInt_FMT " not handled",rw); 191218db3c1SStefano Zampini break; 192218db3c1SStefano Zampini } 1939566063dSJacob Faibussowitsch PetscCall(VecGetSize(X,&n)); 19426fbe8dcSKarl Rupp 195218db3c1SStefano Zampini Y->x = x; 196641875f9SMatthew G Knepley Y->nrow = n; 197641875f9SMatthew G Knepley Y->ncol = 1; 198641875f9SMatthew G Knepley Y->nzmax = n; 199641875f9SMatthew G Knepley Y->d = n; 200641875f9SMatthew G Knepley Y->xtype = CHOLMOD_SCALAR_TYPE; 201641875f9SMatthew G Knepley Y->dtype = CHOLMOD_DOUBLE; 202641875f9SMatthew G Knepley PetscFunctionReturn(0); 203641875f9SMatthew G Knepley } 204641875f9SMatthew G Knepley 205a2fc1e05SToby Isaac PetscErrorCode VecUnWrapCholmod(Vec X,PetscInt rw,cholmod_dense *Y) 206d9ca1df4SBarry Smith { 207d9ca1df4SBarry Smith PetscFunctionBegin; 208218db3c1SStefano Zampini switch (rw) { 209218db3c1SStefano Zampini case GET_ARRAY_READ: 2109566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X,(const PetscScalar**)&Y->x)); 211218db3c1SStefano Zampini break; 212218db3c1SStefano Zampini case GET_ARRAY_WRITE: 2139566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(X,(PetscScalar**)&Y->x)); 214218db3c1SStefano Zampini break; 215218db3c1SStefano Zampini default: 21698921bdaSJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Case %" PetscInt_FMT " not handled",rw); 217218db3c1SStefano Zampini break; 218218db3c1SStefano Zampini } 219218db3c1SStefano Zampini PetscFunctionReturn(0); 220218db3c1SStefano Zampini } 221218db3c1SStefano Zampini 222a2fc1e05SToby Isaac PetscErrorCode MatDenseWrapCholmod(Mat X,PetscInt rw,cholmod_dense *Y) 223218db3c1SStefano Zampini { 224218db3c1SStefano Zampini PetscScalar *x; 225218db3c1SStefano Zampini PetscInt m,n,lda; 226218db3c1SStefano Zampini 227218db3c1SStefano Zampini PetscFunctionBegin; 2289566063dSJacob Faibussowitsch PetscCall(PetscMemzero(Y,sizeof(*Y))); 229218db3c1SStefano Zampini switch (rw) { 230218db3c1SStefano Zampini case GET_ARRAY_READ: 2319566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayRead(X,(const PetscScalar**)&x)); 232218db3c1SStefano Zampini break; 233218db3c1SStefano Zampini case GET_ARRAY_WRITE: 2349566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayWrite(X,&x)); 235218db3c1SStefano Zampini break; 236218db3c1SStefano Zampini default: 23798921bdaSJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Case %" PetscInt_FMT " not handled",rw); 238218db3c1SStefano Zampini break; 239218db3c1SStefano Zampini } 2409566063dSJacob Faibussowitsch PetscCall(MatDenseGetLDA(X,&lda)); 2419566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(X,&m,&n)); 242218db3c1SStefano Zampini 243218db3c1SStefano Zampini Y->x = x; 244218db3c1SStefano Zampini Y->nrow = m; 245218db3c1SStefano Zampini Y->ncol = n; 246218db3c1SStefano Zampini Y->nzmax = lda*n; 247218db3c1SStefano Zampini Y->d = lda; 248218db3c1SStefano Zampini Y->xtype = CHOLMOD_SCALAR_TYPE; 249218db3c1SStefano Zampini Y->dtype = CHOLMOD_DOUBLE; 250218db3c1SStefano Zampini PetscFunctionReturn(0); 251218db3c1SStefano Zampini } 252218db3c1SStefano Zampini 253a2fc1e05SToby Isaac PetscErrorCode MatDenseUnWrapCholmod(Mat X,PetscInt rw,cholmod_dense *Y) 254218db3c1SStefano Zampini { 255218db3c1SStefano Zampini PetscFunctionBegin; 256218db3c1SStefano Zampini switch (rw) { 257218db3c1SStefano Zampini case GET_ARRAY_READ: 2589566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayRead(X,(const PetscScalar**)&Y->x)); 259218db3c1SStefano Zampini break; 260218db3c1SStefano Zampini case GET_ARRAY_WRITE: 261218db3c1SStefano Zampini /* we don't have MatDenseRestoreArrayWrite */ 2629566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArray(X,(PetscScalar**)&Y->x)); 263218db3c1SStefano Zampini break; 264218db3c1SStefano Zampini default: 26598921bdaSJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Case %" PetscInt_FMT " not handled",rw); 266218db3c1SStefano Zampini break; 267218db3c1SStefano Zampini } 268d9ca1df4SBarry Smith PetscFunctionReturn(0); 269d9ca1df4SBarry Smith } 270d9ca1df4SBarry Smith 271eb9872f6SBarry Smith PETSC_INTERN PetscErrorCode MatDestroy_CHOLMOD(Mat F) 272641875f9SMatthew G Knepley { 2736b8f6f9dSBarry Smith Mat_CHOLMOD *chol=(Mat_CHOLMOD*)F->data; 274641875f9SMatthew G Knepley 275641875f9SMatthew G Knepley PetscFunctionBegin; 276a2fc1e05SToby Isaac if (chol->spqrfact) { 2779566063dSJacob Faibussowitsch PetscCall(!SuiteSparseQR_C_free(&chol->spqrfact, chol->common)); 278a2fc1e05SToby Isaac } 279a2fc1e05SToby Isaac if (chol->factor) { 2809566063dSJacob Faibussowitsch PetscCall(!cholmod_X_free_factor(&chol->factor,chol->common)); 281a2fc1e05SToby Isaac } 282a2fc1e05SToby Isaac if (chol->common->itype == CHOLMOD_INT) { 2839566063dSJacob Faibussowitsch PetscCall(!cholmod_finish(chol->common)); 284a2fc1e05SToby Isaac } else { 2859566063dSJacob Faibussowitsch PetscCall(!cholmod_l_finish(chol->common)); 286a2fc1e05SToby Isaac } 2879566063dSJacob Faibussowitsch PetscCall(PetscFree(chol->common)); 2889566063dSJacob Faibussowitsch PetscCall(PetscFree(chol->matrix)); 2899566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)F,"MatFactorGetSolverType_C",NULL)); 2902e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)F,"MatQRFactorSymbolic_C",NULL)); 2912e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)F,"MatQRFactorNumeric_C",NULL)); 2929566063dSJacob Faibussowitsch PetscCall(PetscFree(F->data)); 293641875f9SMatthew G Knepley PetscFunctionReturn(0); 294641875f9SMatthew G Knepley } 295641875f9SMatthew G Knepley 296641875f9SMatthew G Knepley static PetscErrorCode MatSolve_CHOLMOD(Mat,Vec,Vec); 297218db3c1SStefano Zampini static PetscErrorCode MatMatSolve_CHOLMOD(Mat,Mat,Mat); 298641875f9SMatthew G Knepley 299fcd503bcSBarry Smith /*static const char *const CholmodOrderingMethods[] = {"User","AMD","METIS","NESDIS(default)","Natural","NESDIS(small=20000)","NESDIS(small=4,no constrained)","NESDIS()"};*/ 300641875f9SMatthew G Knepley 301860c79edSBarry Smith static PetscErrorCode MatView_Info_CHOLMOD(Mat F,PetscViewer viewer) 302641875f9SMatthew G Knepley { 3036b8f6f9dSBarry Smith Mat_CHOLMOD *chol = (Mat_CHOLMOD*)F->data; 304641875f9SMatthew G Knepley const cholmod_common *c = chol->common; 305641875f9SMatthew G Knepley PetscInt i; 306641875f9SMatthew G Knepley 307641875f9SMatthew G Knepley PetscFunctionBegin; 308641875f9SMatthew G Knepley if (F->ops->solve != MatSolve_CHOLMOD) PetscFunctionReturn(0); 3099566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"CHOLMOD run parameters:\n")); 3109566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 3119566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Pack factors after symbolic factorization: %s\n",chol->pack ? "TRUE" : "FALSE")); 3129566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Common.dbound %g (Smallest absolute value of diagonal entries of D)\n",c->dbound)); 3139566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Common.grow0 %g\n",c->grow0)); 3149566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Common.grow1 %g\n",c->grow1)); 3159566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Common.grow2 %u\n",(unsigned)c->grow2)); 3169566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Common.maxrank %u\n",(unsigned)c->maxrank)); 3179566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Common.supernodal_switch %g\n",c->supernodal_switch)); 3189566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Common.supernodal %d\n",c->supernodal)); 3199566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Common.final_asis %d\n",c->final_asis)); 3209566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Common.final_super %d\n",c->final_super)); 3219566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Common.final_ll %d\n",c->final_ll)); 3229566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Common.final_pack %d\n",c->final_pack)); 3239566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Common.final_monotonic %d\n",c->final_monotonic)); 3249566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Common.final_resymbol %d\n",c->final_resymbol)); 3259566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Common.zrelax [%g,%g,%g]\n",c->zrelax[0],c->zrelax[1],c->zrelax[2])); 3269566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Common.nrelax [%u,%u,%u]\n",(unsigned)c->nrelax[0],(unsigned)c->nrelax[1],(unsigned)c->nrelax[2])); 3279566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Common.prefer_upper %d\n",c->prefer_upper)); 3289566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Common.print %d\n",c->print)); 329641875f9SMatthew G Knepley for (i=0; i<c->nmethods; i++) { 3309566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Ordering method %" PetscInt_FMT "%s:\n",i,i==c->selected ? " [SELECTED]" : "")); 331d0609cedSBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer," lnz %g, fl %g, prune_dense %g, prune_dense2 %g\n", 332d0609cedSBarry Smith c->method[i].lnz,c->method[i].fl,c->method[i].prune_dense,c->method[i].prune_dense2)); 333641875f9SMatthew G Knepley } 3349566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Common.postorder %d\n",c->postorder)); 3359566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Common.default_nesdis %d (use NESDIS instead of METIS for nested dissection)\n",c->default_nesdis)); 336641875f9SMatthew G Knepley /* Statistics */ 3379566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Common.fl %g (flop count from most recent analysis)\n",c->fl)); 3389566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Common.lnz %g (fundamental nz in L)\n",c->lnz)); 3399566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Common.anz %g\n",c->anz)); 3409566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Common.modfl %g (flop count from most recent update)\n",c->modfl)); 3419566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Common.malloc_count %g (number of live objects)\n",(double)c->malloc_count)); 3429566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Common.memory_usage %g (peak memory usage in bytes)\n",(double)c->memory_usage)); 3439566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Common.memory_inuse %g (current memory usage in bytes)\n",(double)c->memory_inuse)); 3449566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Common.nrealloc_col %g (number of column reallocations)\n",c->nrealloc_col)); 3459566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Common.nrealloc_factor %g (number of factor reallocations due to column reallocations)\n",c->nrealloc_factor)); 3469566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Common.ndbounds_hit %g (number of times diagonal was modified by dbound)\n",c->ndbounds_hit)); 3479566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Common.rowfacfl %g (number of flops in last call to cholmod_rowfac)\n",c->rowfacfl)); 3489566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Common.aatfl %g (number of flops to compute A(:,f)*A(:,f)')\n",c->aatfl)); 349b9eaa5e8SBarry Smith #if defined(PETSC_USE_SUITESPARSE_GPU) 3509566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Common.useGPU %d\n",c->useGPU)); 351b9eaa5e8SBarry Smith #endif 3529566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 353641875f9SMatthew G Knepley PetscFunctionReturn(0); 354641875f9SMatthew G Knepley } 355641875f9SMatthew G Knepley 356eb9872f6SBarry Smith PETSC_INTERN PetscErrorCode MatView_CHOLMOD(Mat F,PetscViewer viewer) 357641875f9SMatthew G Knepley { 358ace3abfcSBarry Smith PetscBool iascii; 359641875f9SMatthew G Knepley PetscViewerFormat format; 360641875f9SMatthew G Knepley 361641875f9SMatthew G Knepley PetscFunctionBegin; 3629566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii)); 363641875f9SMatthew G Knepley if (iascii) { 3649566063dSJacob Faibussowitsch PetscCall(PetscViewerGetFormat(viewer,&format)); 365641875f9SMatthew G Knepley if (format == PETSC_VIEWER_ASCII_INFO) { 3669566063dSJacob Faibussowitsch PetscCall(MatView_Info_CHOLMOD(F,viewer)); 367641875f9SMatthew G Knepley } 368641875f9SMatthew G Knepley } 369641875f9SMatthew G Knepley PetscFunctionReturn(0); 370641875f9SMatthew G Knepley } 371641875f9SMatthew G Knepley 372641875f9SMatthew G Knepley static PetscErrorCode MatSolve_CHOLMOD(Mat F,Vec B,Vec X) 373641875f9SMatthew G Knepley { 3746b8f6f9dSBarry Smith Mat_CHOLMOD *chol = (Mat_CHOLMOD*)F->data; 375218db3c1SStefano Zampini cholmod_dense cholB,cholX,*X_handle,*Y_handle = NULL,*E_handle = NULL; 376641875f9SMatthew G Knepley 377641875f9SMatthew G Knepley PetscFunctionBegin; 378641875f9SMatthew G Knepley static_F = F; 3799566063dSJacob Faibussowitsch PetscCall(VecWrapCholmod(B,GET_ARRAY_READ,&cholB)); 3809566063dSJacob Faibussowitsch PetscCall(VecWrapCholmod(X,GET_ARRAY_WRITE,&cholX)); 381218db3c1SStefano Zampini X_handle = &cholX; 3829566063dSJacob Faibussowitsch PetscCall(!cholmod_X_solve2(CHOLMOD_A,chol->factor,&cholB,NULL,&X_handle,NULL,&Y_handle,&E_handle,chol->common)); 3839566063dSJacob Faibussowitsch PetscCall(!cholmod_X_free_dense(&Y_handle,chol->common)); 3849566063dSJacob Faibussowitsch PetscCall(!cholmod_X_free_dense(&E_handle,chol->common)); 3859566063dSJacob Faibussowitsch PetscCall(VecUnWrapCholmod(B,GET_ARRAY_READ,&cholB)); 3869566063dSJacob Faibussowitsch PetscCall(VecUnWrapCholmod(X,GET_ARRAY_WRITE,&cholX)); 3879566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(4.0*chol->common->lnz)); 388218db3c1SStefano Zampini PetscFunctionReturn(0); 389218db3c1SStefano Zampini } 390218db3c1SStefano Zampini 391218db3c1SStefano Zampini static PetscErrorCode MatMatSolve_CHOLMOD(Mat F,Mat B,Mat X) 392218db3c1SStefano Zampini { 393218db3c1SStefano Zampini Mat_CHOLMOD *chol = (Mat_CHOLMOD*)F->data; 394218db3c1SStefano Zampini cholmod_dense cholB,cholX,*X_handle,*Y_handle = NULL,*E_handle = NULL; 395218db3c1SStefano Zampini 396218db3c1SStefano Zampini PetscFunctionBegin; 397218db3c1SStefano Zampini static_F = F; 3989566063dSJacob Faibussowitsch PetscCall(MatDenseWrapCholmod(B,GET_ARRAY_READ,&cholB)); 3999566063dSJacob Faibussowitsch PetscCall(MatDenseWrapCholmod(X,GET_ARRAY_WRITE,&cholX)); 400218db3c1SStefano Zampini X_handle = &cholX; 4019566063dSJacob Faibussowitsch PetscCall(!cholmod_X_solve2(CHOLMOD_A,chol->factor,&cholB,NULL,&X_handle,NULL,&Y_handle,&E_handle,chol->common)); 4029566063dSJacob Faibussowitsch PetscCall(!cholmod_X_free_dense(&Y_handle,chol->common)); 4039566063dSJacob Faibussowitsch PetscCall(!cholmod_X_free_dense(&E_handle,chol->common)); 4049566063dSJacob Faibussowitsch PetscCall(MatDenseUnWrapCholmod(B,GET_ARRAY_READ,&cholB)); 4059566063dSJacob Faibussowitsch PetscCall(MatDenseUnWrapCholmod(X,GET_ARRAY_WRITE,&cholX)); 4069566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(4.0*B->cmap->n*chol->common->lnz)); 407641875f9SMatthew G Knepley PetscFunctionReturn(0); 408641875f9SMatthew G Knepley } 409641875f9SMatthew G Knepley 410641875f9SMatthew G Knepley static PetscErrorCode MatCholeskyFactorNumeric_CHOLMOD(Mat F,Mat A,const MatFactorInfo *info) 411641875f9SMatthew G Knepley { 4126b8f6f9dSBarry Smith Mat_CHOLMOD *chol = (Mat_CHOLMOD*)F->data; 413641875f9SMatthew G Knepley cholmod_sparse cholA; 414218db3c1SStefano Zampini PetscBool aijalloc,valloc; 415d0609cedSBarry Smith int err; 416641875f9SMatthew G Knepley 417641875f9SMatthew G Knepley PetscFunctionBegin; 4189566063dSJacob Faibussowitsch PetscCall((*chol->Wrap)(A,PETSC_TRUE,&cholA,&aijalloc,&valloc)); 419641875f9SMatthew G Knepley static_F = F; 420d0609cedSBarry Smith err = !cholmod_X_factorize(&cholA,chol->factor,chol->common); 421d0609cedSBarry Smith PetscCheck(!err,PetscObjectComm((PetscObject)F),PETSC_ERR_LIB,"CHOLMOD factorization failed with status %d",chol->common->status); 42208401ef6SPierre 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); 423641875f9SMatthew G Knepley 4249566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(chol->common->fl)); 4259566063dSJacob Faibussowitsch if (aijalloc) PetscCall(PetscFree2(cholA.p,cholA.i)); 4269566063dSJacob Faibussowitsch if (valloc) PetscCall(PetscFree(cholA.x)); 427218db3c1SStefano Zampini #if defined(PETSC_USE_SUITESPARSE_GPU) 4289566063dSJacob 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)); 429218db3c1SStefano Zampini #endif 430641875f9SMatthew G Knepley 431641875f9SMatthew G Knepley F->ops->solve = MatSolve_CHOLMOD; 432641875f9SMatthew G Knepley F->ops->solvetranspose = MatSolve_CHOLMOD; 433218db3c1SStefano Zampini F->ops->matsolve = MatMatSolve_CHOLMOD; 434218db3c1SStefano Zampini F->ops->matsolvetranspose = MatMatSolve_CHOLMOD; 435641875f9SMatthew G Knepley PetscFunctionReturn(0); 436641875f9SMatthew G Knepley } 437641875f9SMatthew G Knepley 438eb9872f6SBarry Smith PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_CHOLMOD(Mat F,Mat A,IS perm,const MatFactorInfo *info) 439641875f9SMatthew G Knepley { 4406b8f6f9dSBarry Smith Mat_CHOLMOD *chol = (Mat_CHOLMOD*)F->data; 441d0609cedSBarry Smith int err; 442641875f9SMatthew G Knepley cholmod_sparse cholA; 443218db3c1SStefano Zampini PetscBool aijalloc,valloc; 444641875f9SMatthew G Knepley PetscInt *fset = 0; 445641875f9SMatthew G Knepley size_t fsize = 0; 446641875f9SMatthew G Knepley 447641875f9SMatthew G Knepley PetscFunctionBegin; 44826cc229bSBarry Smith /* Set options to F */ 44926cc229bSBarry Smith PetscCall(CholmodSetOptions(F)); 45026cc229bSBarry Smith 4519566063dSJacob Faibussowitsch PetscCall((*chol->Wrap)(A,PETSC_FALSE,&cholA,&aijalloc,&valloc)); 452641875f9SMatthew G Knepley static_F = F; 453641875f9SMatthew G Knepley if (chol->factor) { 454d0609cedSBarry Smith err = !cholmod_X_resymbol(&cholA,fset,fsize,(int)chol->pack,chol->factor,chol->common); 455d0609cedSBarry Smith PetscCheck(!err,PetscObjectComm((PetscObject)F),PETSC_ERR_LIB,"CHOLMOD analysis failed with status %d",chol->common->status); 456641875f9SMatthew G Knepley } else if (perm) { 457641875f9SMatthew G Knepley const PetscInt *ip; 4589566063dSJacob Faibussowitsch PetscCall(ISGetIndices(perm,&ip)); 459641875f9SMatthew G Knepley chol->factor = cholmod_X_analyze_p(&cholA,(PetscInt*)ip,fset,fsize,chol->common); 46028b400f6SJacob Faibussowitsch PetscCheck(chol->factor,PetscObjectComm((PetscObject)F),PETSC_ERR_LIB,"CHOLMOD analysis failed using PETSc ordering with status %d",chol->common->status); 4619566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(perm,&ip)); 462641875f9SMatthew G Knepley } else { 463641875f9SMatthew G Knepley chol->factor = cholmod_X_analyze(&cholA,chol->common); 46428b400f6SJacob Faibussowitsch PetscCheck(chol->factor,PetscObjectComm((PetscObject)F),PETSC_ERR_LIB,"CHOLMOD analysis failed using internal ordering with status %d",chol->common->status); 465641875f9SMatthew G Knepley } 466641875f9SMatthew G Knepley 4679566063dSJacob Faibussowitsch if (aijalloc) PetscCall(PetscFree2(cholA.p,cholA.i)); 4689566063dSJacob Faibussowitsch if (valloc) PetscCall(PetscFree(cholA.x)); 469641875f9SMatthew G Knepley 470641875f9SMatthew G Knepley F->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_CHOLMOD; 471641875f9SMatthew G Knepley PetscFunctionReturn(0); 472641875f9SMatthew G Knepley } 473641875f9SMatthew G Knepley 474ea799195SBarry Smith static PetscErrorCode MatFactorGetSolverType_seqsbaij_cholmod(Mat A,MatSolverType *type) 475641875f9SMatthew G Knepley { 476641875f9SMatthew G Knepley PetscFunctionBegin; 477641875f9SMatthew G Knepley *type = MATSOLVERCHOLMOD; 478641875f9SMatthew G Knepley PetscFunctionReturn(0); 479641875f9SMatthew G Knepley } 480641875f9SMatthew G Knepley 481218db3c1SStefano Zampini PETSC_INTERN PetscErrorCode MatGetInfo_CHOLMOD(Mat F,MatInfoType flag,MatInfo *info) 482218db3c1SStefano Zampini { 483218db3c1SStefano Zampini Mat_CHOLMOD *chol = (Mat_CHOLMOD*)F->data; 484218db3c1SStefano Zampini 485218db3c1SStefano Zampini PetscFunctionBegin; 486218db3c1SStefano Zampini info->block_size = 1.0; 487218db3c1SStefano Zampini info->nz_allocated = chol->common->lnz; 488218db3c1SStefano Zampini info->nz_used = chol->common->lnz; 489218db3c1SStefano Zampini info->nz_unneeded = 0.0; 490218db3c1SStefano Zampini info->assemblies = 0.0; 491218db3c1SStefano Zampini info->mallocs = 0.0; 492218db3c1SStefano Zampini info->memory = chol->common->memory_inuse; 493218db3c1SStefano Zampini info->fill_ratio_given = 0; 494218db3c1SStefano Zampini info->fill_ratio_needed = 0; 495218db3c1SStefano Zampini info->factor_mallocs = chol->common->malloc_count; 496218db3c1SStefano Zampini PetscFunctionReturn(0); 497218db3c1SStefano Zampini } 498218db3c1SStefano Zampini 499641875f9SMatthew G Knepley /*MC 500ee300463SSatish Balay MATSOLVERCHOLMOD 501ee300463SSatish Balay 502ee300463SSatish Balay A matrix type providing direct solvers (Cholesky) for sequential matrices 503641875f9SMatthew G Knepley via the external package CHOLMOD. 504641875f9SMatthew G Knepley 505c2b89b5dSBarry Smith Use ./configure --download-suitesparse to install PETSc to use CHOLMOD 506c2b89b5dSBarry Smith 507f29f8b16SBarry Smith Use -pc_type cholesky -pc_factor_mat_solver_type cholmod to use this direct solver 508641875f9SMatthew G Knepley 509641875f9SMatthew G Knepley Consult CHOLMOD documentation for more information about the Common parameters 510641875f9SMatthew G Knepley which correspond to the options database keys below. 511641875f9SMatthew G Knepley 512641875f9SMatthew G Knepley Options Database Keys: 513e08999f5SMatthew G Knepley + -mat_cholmod_dbound <0> - Minimum absolute value of diagonal entries of D (None) 514e08999f5SMatthew G Knepley . -mat_cholmod_grow0 <1.2> - Global growth ratio when factors are modified (None) 515e08999f5SMatthew G Knepley . -mat_cholmod_grow1 <1.2> - Column growth ratio when factors are modified (None) 516e08999f5SMatthew G Knepley . -mat_cholmod_grow2 <5> - Affine column growth constant when factors are modified (None) 517e08999f5SMatthew G Knepley . -mat_cholmod_maxrank <8> - Max rank of update, larger values are faster but use more memory [2,4,8] (None) 518e08999f5SMatthew G Knepley . -mat_cholmod_factor <AUTO> - (choose one of) SIMPLICIAL AUTO SUPERNODAL 519e08999f5SMatthew G Knepley . -mat_cholmod_supernodal_switch <40> - flop/nnz_L threshold for switching to supernodal factorization (None) 520e08999f5SMatthew G Knepley . -mat_cholmod_final_asis <TRUE> - Leave factors "as is" (None) 521e08999f5SMatthew G Knepley . -mat_cholmod_final_pack <TRUE> - Pack the columns when finished (use FALSE if the factors will be updated later) (None) 522e08999f5SMatthew G Knepley . -mat_cholmod_zrelax <0.8> - 3 real supernodal relaxed amalgamation parameters (None) 523e08999f5SMatthew G Knepley . -mat_cholmod_nrelax <4> - 3 size_t supernodal relaxed amalgamation parameters (None) 524e08999f5SMatthew G Knepley . -mat_cholmod_prefer_upper <TRUE> - Work with upper triangular form (faster when using fill-reducing ordering, slower in natural ordering) (None) 5252c7c0729SBarry Smith . -mat_cholmod_print <3> - Verbosity level (None) 526ee300463SSatish Balay - -mat_ordering_type internal - Use the ordering provided by Cholmod 527641875f9SMatthew G Knepley 528641875f9SMatthew G Knepley Level: beginner 529641875f9SMatthew G Knepley 530a364b7d2SBarry Smith Note: CHOLMOD is part of SuiteSparse http://faculty.cse.tamu.edu/davis/suitesparse.html 531a364b7d2SBarry Smith 532db781477SPatrick Sanan .seealso: `PCCHOLESKY`, `PCFactorSetMatSolverType()`, `MatSolverType` 533641875f9SMatthew G Knepley M*/ 534b2573a8aSBarry Smith 535db87b0f2SBarry Smith PETSC_INTERN PetscErrorCode MatGetFactor_seqsbaij_cholmod(Mat A,MatFactorType ftype,Mat *F) 536641875f9SMatthew G Knepley { 537641875f9SMatthew G Knepley Mat B; 538641875f9SMatthew G Knepley Mat_CHOLMOD *chol; 539641875f9SMatthew G Knepley PetscInt m=A->rmap->n,n=A->cmap->n,bs; 540641875f9SMatthew G Knepley 541641875f9SMatthew G Knepley PetscFunctionBegin; 5429566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(A,&bs)); 54308401ef6SPierre Jolivet PetscCheck(bs == 1,PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"CHOLMOD only supports block size=1, given %" PetscInt_FMT,bs); 544218db3c1SStefano Zampini #if defined(PETSC_USE_COMPLEX) 545*b94d7dedSBarry Smith PetscCheck(A->hermitian == PETSC_BOOL3_TRUE,PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Only for Hermitian matrices"); 546218db3c1SStefano Zampini #endif 547641875f9SMatthew G Knepley /* Create the factorization matrix F */ 5489566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A),&B)); 5499566063dSJacob Faibussowitsch PetscCall(MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,m,n)); 5509566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy("cholmod",&((PetscObject)B)->type_name)); 5519566063dSJacob Faibussowitsch PetscCall(MatSetUp(B)); 5529566063dSJacob Faibussowitsch PetscCall(PetscNewLog(B,&chol)); 55326fbe8dcSKarl Rupp 554641875f9SMatthew G Knepley chol->Wrap = MatWrapCholmod_seqsbaij; 5556b8f6f9dSBarry Smith B->data = chol; 556641875f9SMatthew G Knepley 557218db3c1SStefano Zampini B->ops->getinfo = MatGetInfo_CHOLMOD; 558641875f9SMatthew G Knepley B->ops->view = MatView_CHOLMOD; 559641875f9SMatthew G Knepley B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_CHOLMOD; 560641875f9SMatthew G Knepley B->ops->destroy = MatDestroy_CHOLMOD; 5619566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverType_C",MatFactorGetSolverType_seqsbaij_cholmod)); 562641875f9SMatthew G Knepley B->factortype = MAT_FACTOR_CHOLESKY; 563218db3c1SStefano Zampini B->assembled = PETSC_TRUE; 564641875f9SMatthew G Knepley B->preallocated = PETSC_TRUE; 565641875f9SMatthew G Knepley 5669566063dSJacob Faibussowitsch PetscCall(CholmodStart(B)); 56700c67f3bSHong Zhang 5689566063dSJacob Faibussowitsch PetscCall(PetscFree(B->solvertype)); 5699566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(MATSOLVERCHOLMOD,&B->solvertype)); 570f73b0415SBarry Smith B->canuseordering = PETSC_TRUE; 5719566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(MATORDERINGEXTERNAL,(char**)&B->preferredordering[MAT_FACTOR_CHOLESKY])); 572641875f9SMatthew G Knepley *F = B; 573641875f9SMatthew G Knepley PetscFunctionReturn(0); 574641875f9SMatthew G Knepley } 575