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 59e475b0dSSatish Balay When build 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 { 23ec55ff42SBarry Smith PetscErrorCode ierr; 24641875f9SMatthew G Knepley 25641875f9SMatthew G Knepley PetscFunctionBegin; 26641875f9SMatthew G Knepley if (status > CHOLMOD_OK) { 27e49b6e0cSMatthew G. Knepley ierr = PetscInfo4(static_F,"CHOLMOD warning %d at %s:%d: %s\n",status,file,line,message);CHKERRV(ierr); 28641875f9SMatthew G Knepley } else if (status == CHOLMOD_OK) { /* Documentation says this can happen, but why? */ 29e49b6e0cSMatthew G. Knepley ierr = PetscInfo3(static_F,"CHOLMOD OK at %s:%d: %s\n",file,line,message);CHKERRV(ierr); 30641875f9SMatthew G Knepley } else { 31e49b6e0cSMatthew G. Knepley ierr = PetscErrorPrintf("CHOLMOD error %d at %s:%d: %s\n",status,file,line,message);CHKERRV(ierr); 32641875f9SMatthew G Knepley } 33641875f9SMatthew G Knepley PetscFunctionReturnVoid(); 34641875f9SMatthew G Knepley } 35641875f9SMatthew G Knepley 367087cfbeSBarry Smith PetscErrorCode CholmodStart(Mat F) 37641875f9SMatthew G Knepley { 38641875f9SMatthew G Knepley PetscErrorCode ierr; 396b8f6f9dSBarry Smith Mat_CHOLMOD *chol=(Mat_CHOLMOD*)F->data; 40641875f9SMatthew G Knepley cholmod_common *c; 41ace3abfcSBarry Smith PetscBool flg; 42641875f9SMatthew G Knepley 43641875f9SMatthew G Knepley PetscFunctionBegin; 44641875f9SMatthew G Knepley if (chol->common) PetscFunctionReturn(0); 45854ce69bSBarry Smith ierr = PetscMalloc1(1,&chol->common);CHKERRQ(ierr); 46641875f9SMatthew G Knepley ierr = !cholmod_X_start(chol->common);CHKERRQ(ierr); 4726fbe8dcSKarl Rupp 48641875f9SMatthew G Knepley c = chol->common; 49641875f9SMatthew G Knepley c->error_handler = CholmodErrorHandler; 50641875f9SMatthew G Knepley 51641875f9SMatthew G Knepley #define CHOLMOD_OPTION_DOUBLE(name,help) do { \ 52641875f9SMatthew G Knepley PetscReal tmp = (PetscReal)c->name; \ 538afaa268SBarry Smith ierr = PetscOptionsReal("-mat_cholmod_" #name,help,"None",tmp,&tmp,NULL);CHKERRQ(ierr); \ 54641875f9SMatthew G Knepley c->name = (double)tmp; \ 55641875f9SMatthew G Knepley } while (0) 5626fbe8dcSKarl Rupp 57641875f9SMatthew G Knepley #define CHOLMOD_OPTION_INT(name,help) do { \ 58641875f9SMatthew G Knepley PetscInt tmp = (PetscInt)c->name; \ 598afaa268SBarry Smith ierr = PetscOptionsInt("-mat_cholmod_" #name,help,"None",tmp,&tmp,NULL);CHKERRQ(ierr); \ 60641875f9SMatthew G Knepley c->name = (int)tmp; \ 61641875f9SMatthew G Knepley } while (0) 6226fbe8dcSKarl Rupp 63641875f9SMatthew G Knepley #define CHOLMOD_OPTION_SIZE_T(name,help) do { \ 64641875f9SMatthew G Knepley PetscInt tmp = (PetscInt)c->name; \ 658afaa268SBarry Smith ierr = PetscOptionsInt("-mat_cholmod_" #name,help,"None",tmp,&tmp,NULL);CHKERRQ(ierr); \ 66ce94432eSBarry Smith if (tmp < 0) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_OUTOFRANGE,"value must be positive"); \ 67641875f9SMatthew G Knepley c->name = (size_t)tmp; \ 68641875f9SMatthew G Knepley } while (0) 6926fbe8dcSKarl Rupp 70b9eaa5e8SBarry Smith #define CHOLMOD_OPTION_BOOL(name,help) do { \ 71ace3abfcSBarry Smith PetscBool tmp = (PetscBool) !!c->name; \ 728afaa268SBarry Smith ierr = PetscOptionsBool("-mat_cholmod_" #name,help,"None",tmp,&tmp,NULL);CHKERRQ(ierr); \ 73641875f9SMatthew G Knepley c->name = (int)tmp; \ 74641875f9SMatthew G Knepley } while (0) 75641875f9SMatthew G Knepley 76ce94432eSBarry Smith ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)F),((PetscObject)F)->prefix,"CHOLMOD Options","Mat");CHKERRQ(ierr); 77641875f9SMatthew G Knepley /* CHOLMOD handles first-time packing and refactor-packing separately, but we usually want them to be the same. */ 78ace3abfcSBarry Smith chol->pack = (PetscBool)c->final_pack; 7926fbe8dcSKarl Rupp 80b9eaa5e8SBarry Smith #if defined(PETSC_USE_SUITESPARSE_GPU) 81b9eaa5e8SBarry Smith c->useGPU = 1; 82b9eaa5e8SBarry Smith CHOLMOD_OPTION_INT(useGPU,"Use GPU for BLAS 1, otherwise 0"); 83b9eaa5e8SBarry Smith #endif 84b9eaa5e8SBarry Smith 858afaa268SBarry Smith ierr = PetscOptionsBool("-mat_cholmod_pack","Pack factors after factorization [disable for frequent repeat factorization]","None",chol->pack,&chol->pack,NULL);CHKERRQ(ierr); 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}; 958afaa268SBarry Smith ierr = PetscOptionsEnum("-mat_cholmod_factor","Factorization method","None",list,(PetscEnum)c->supernodal,(PetscEnum*)&c->supernodal,NULL);CHKERRQ(ierr); 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; 109641875f9SMatthew G Knepley ierr = PetscOptionsRealArray("-mat_cholmod_zrelax","3 real supernodal relaxed amalgamation parameters","None",tmp,&n,&flg);CHKERRQ(ierr); 110ce94432eSBarry Smith if (flg && n != 3) SETERRQ(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]}; 115641875f9SMatthew G Knepley ierr = PetscOptionsIntArray("-mat_cholmod_nrelax","3 size_t supernodal relaxed amalgamation parameters","None",tmp,&n,&flg);CHKERRQ(ierr); 116ce94432eSBarry Smith if (flg && n != 3) SETERRQ(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"); 122641875f9SMatthew G Knepley ierr = PetscOptionsEnd();CHKERRQ(ierr); 123641875f9SMatthew G Knepley PetscFunctionReturn(0); 124641875f9SMatthew G Knepley } 125641875f9SMatthew G Knepley 126ace3abfcSBarry Smith static PetscErrorCode MatWrapCholmod_seqsbaij(Mat A,PetscBool values,cholmod_sparse *C,PetscBool *aijalloc) 127641875f9SMatthew G Knepley { 128641875f9SMatthew G Knepley Mat_SeqSBAIJ *sbaij = (Mat_SeqSBAIJ*)A->data; 129641875f9SMatthew G Knepley PetscErrorCode ierr; 130641875f9SMatthew G Knepley 131641875f9SMatthew G Knepley PetscFunctionBegin; 132641875f9SMatthew G Knepley ierr = PetscMemzero(C,sizeof(*C));CHKERRQ(ierr); 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; 139641875f9SMatthew G Knepley C->x = sbaij->a; 140641875f9SMatthew G Knepley C->stype = -1; 141641875f9SMatthew G Knepley C->itype = CHOLMOD_INT_TYPE; 142641875f9SMatthew G Knepley C->xtype = CHOLMOD_SCALAR_TYPE; 143641875f9SMatthew G Knepley C->dtype = CHOLMOD_DOUBLE; 144641875f9SMatthew G Knepley C->sorted = 1; 145641875f9SMatthew G Knepley C->packed = 1; 146641875f9SMatthew G Knepley *aijalloc = PETSC_FALSE; 147641875f9SMatthew G Knepley PetscFunctionReturn(0); 148641875f9SMatthew G Knepley } 149641875f9SMatthew G Knepley 150d9ca1df4SBarry Smith static PetscErrorCode VecWrapCholmodRead(Vec X,cholmod_dense *Y) 151641875f9SMatthew G Knepley { 152641875f9SMatthew G Knepley PetscErrorCode ierr; 153d9ca1df4SBarry Smith const PetscScalar *x; 154641875f9SMatthew G Knepley PetscInt n; 155641875f9SMatthew G Knepley 156641875f9SMatthew G Knepley PetscFunctionBegin; 157641875f9SMatthew G Knepley ierr = PetscMemzero(Y,sizeof(*Y));CHKERRQ(ierr); 158d9ca1df4SBarry Smith ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 159641875f9SMatthew G Knepley ierr = VecGetSize(X,&n);CHKERRQ(ierr); 16026fbe8dcSKarl Rupp 161641875f9SMatthew G Knepley Y->x = (double*)x; 162641875f9SMatthew G Knepley Y->nrow = n; 163641875f9SMatthew G Knepley Y->ncol = 1; 164641875f9SMatthew G Knepley Y->nzmax = n; 165641875f9SMatthew G Knepley Y->d = n; 166641875f9SMatthew G Knepley Y->x = (double*)x; 167641875f9SMatthew G Knepley Y->xtype = CHOLMOD_SCALAR_TYPE; 168641875f9SMatthew G Knepley Y->dtype = CHOLMOD_DOUBLE; 169641875f9SMatthew G Knepley PetscFunctionReturn(0); 170641875f9SMatthew G Knepley } 171641875f9SMatthew G Knepley 172d9ca1df4SBarry Smith static PetscErrorCode VecUnWrapCholmodRead(Vec X,cholmod_dense *Y) 173d9ca1df4SBarry Smith { 174d9ca1df4SBarry Smith PetscErrorCode ierr; 175d9ca1df4SBarry Smith 176d9ca1df4SBarry Smith PetscFunctionBegin; 177b7be58f8SBarry Smith ierr = VecRestoreArrayRead(X,NULL);CHKERRQ(ierr); 178d9ca1df4SBarry Smith PetscFunctionReturn(0); 179d9ca1df4SBarry Smith } 180d9ca1df4SBarry Smith 181eb9872f6SBarry Smith PETSC_INTERN PetscErrorCode MatDestroy_CHOLMOD(Mat F) 182641875f9SMatthew G Knepley { 183641875f9SMatthew G Knepley PetscErrorCode ierr; 1846b8f6f9dSBarry Smith Mat_CHOLMOD *chol=(Mat_CHOLMOD*)F->data; 185641875f9SMatthew G Knepley 186641875f9SMatthew G Knepley PetscFunctionBegin; 187641875f9SMatthew G Knepley ierr = !cholmod_X_free_factor(&chol->factor,chol->common);CHKERRQ(ierr); 188641875f9SMatthew G Knepley ierr = !cholmod_X_finish(chol->common);CHKERRQ(ierr); 189641875f9SMatthew G Knepley ierr = PetscFree(chol->common);CHKERRQ(ierr); 190641875f9SMatthew G Knepley ierr = PetscFree(chol->matrix);CHKERRQ(ierr); 1916b8f6f9dSBarry Smith ierr = PetscFree(F->data);CHKERRQ(ierr); 192641875f9SMatthew G Knepley PetscFunctionReturn(0); 193641875f9SMatthew G Knepley } 194641875f9SMatthew G Knepley 195641875f9SMatthew G Knepley static PetscErrorCode MatSolve_CHOLMOD(Mat,Vec,Vec); 196641875f9SMatthew G Knepley 197fcd503bcSBarry Smith /*static const char *const CholmodOrderingMethods[] = {"User","AMD","METIS","NESDIS(default)","Natural","NESDIS(small=20000)","NESDIS(small=4,no constrained)","NESDIS()"};*/ 198641875f9SMatthew G Knepley 199*860c79edSBarry Smith static PetscErrorCode MatView_Info_CHOLMOD(Mat F,PetscViewer viewer) 200641875f9SMatthew G Knepley { 2016b8f6f9dSBarry Smith Mat_CHOLMOD *chol = (Mat_CHOLMOD*)F->data; 202641875f9SMatthew G Knepley const cholmod_common *c = chol->common; 203641875f9SMatthew G Knepley PetscErrorCode ierr; 204641875f9SMatthew G Knepley PetscInt i; 205641875f9SMatthew G Knepley 206641875f9SMatthew G Knepley PetscFunctionBegin; 207641875f9SMatthew G Knepley if (F->ops->solve != MatSolve_CHOLMOD) PetscFunctionReturn(0); 208641875f9SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer,"CHOLMOD run parameters:\n");CHKERRQ(ierr); 209641875f9SMatthew G Knepley ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 210641875f9SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer,"Pack factors after symbolic factorization: %s\n",chol->pack ? "TRUE" : "FALSE");CHKERRQ(ierr); 211641875f9SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer,"Common.dbound %g (Smallest absolute value of diagonal entries of D)\n",c->dbound);CHKERRQ(ierr); 212641875f9SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer,"Common.grow0 %g\n",c->grow0);CHKERRQ(ierr); 213641875f9SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer,"Common.grow1 %g\n",c->grow1);CHKERRQ(ierr); 214641875f9SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer,"Common.grow2 %u\n",(unsigned)c->grow2);CHKERRQ(ierr); 215641875f9SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer,"Common.maxrank %u\n",(unsigned)c->maxrank);CHKERRQ(ierr); 216641875f9SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer,"Common.supernodal_switch %g\n",c->supernodal_switch);CHKERRQ(ierr); 217641875f9SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer,"Common.supernodal %d\n",c->supernodal);CHKERRQ(ierr); 218641875f9SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer,"Common.final_asis %d\n",c->final_asis);CHKERRQ(ierr); 219641875f9SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer,"Common.final_super %d\n",c->final_super);CHKERRQ(ierr); 220641875f9SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer,"Common.final_ll %d\n",c->final_ll);CHKERRQ(ierr); 221641875f9SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer,"Common.final_pack %d\n",c->final_pack);CHKERRQ(ierr); 222641875f9SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer,"Common.final_monotonic %d\n",c->final_monotonic);CHKERRQ(ierr); 223641875f9SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer,"Common.final_resymbol %d\n",c->final_resymbol);CHKERRQ(ierr); 224641875f9SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer,"Common.zrelax [%g,%g,%g]\n",c->zrelax[0],c->zrelax[1],c->zrelax[2]);CHKERRQ(ierr); 225641875f9SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer,"Common.nrelax [%u,%u,%u]\n",(unsigned)c->nrelax[0],(unsigned)c->nrelax[1],(unsigned)c->nrelax[2]);CHKERRQ(ierr); 226641875f9SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer,"Common.prefer_upper %d\n",c->prefer_upper);CHKERRQ(ierr); 227641875f9SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer,"Common.print %d\n",c->print);CHKERRQ(ierr); 228641875f9SMatthew G Knepley for (i=0; i<c->nmethods; i++) { 229641875f9SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer,"Ordering method %D%s:\n",i,i==c->selected ? " [SELECTED]" : "");CHKERRQ(ierr); 230641875f9SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer," lnz %g, fl %g, prune_dense %g, prune_dense2 %g\n", 231641875f9SMatthew G Knepley c->method[i].lnz,c->method[i].fl,c->method[i].prune_dense,c->method[i].prune_dense2);CHKERRQ(ierr); 232641875f9SMatthew G Knepley } 233641875f9SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer,"Common.postorder %d\n",c->postorder);CHKERRQ(ierr); 234641875f9SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer,"Common.default_nesdis %d (use NESDIS instead of METIS for nested dissection)\n",c->default_nesdis);CHKERRQ(ierr); 235641875f9SMatthew G Knepley /* Statistics */ 236641875f9SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer,"Common.fl %g (flop count from most recent analysis)\n",c->fl);CHKERRQ(ierr); 237641875f9SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer,"Common.lnz %g (fundamental nz in L)\n",c->lnz);CHKERRQ(ierr); 238641875f9SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer,"Common.anz %g\n",c->anz);CHKERRQ(ierr); 239641875f9SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer,"Common.modfl %g (flop count from most recent update)\n",c->modfl);CHKERRQ(ierr); 240641875f9SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer,"Common.malloc_count %g (number of live objects)\n",(double)c->malloc_count);CHKERRQ(ierr);CHKERRQ(ierr); 241641875f9SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer,"Common.memory_usage %g (peak memory usage in bytes)\n",(double)c->memory_usage);CHKERRQ(ierr);CHKERRQ(ierr); 242641875f9SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer,"Common.memory_inuse %g (current memory usage in bytes)\n",(double)c->memory_inuse);CHKERRQ(ierr);CHKERRQ(ierr); 243641875f9SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer,"Common.nrealloc_col %g (number of column reallocations)\n",c->nrealloc_col);CHKERRQ(ierr);CHKERRQ(ierr); 244641875f9SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer,"Common.nrealloc_factor %g (number of factor reallocations due to column reallocations)\n",c->nrealloc_factor);CHKERRQ(ierr);CHKERRQ(ierr); 245641875f9SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer,"Common.ndbounds_hit %g (number of times diagonal was modified by dbound)\n",c->ndbounds_hit);CHKERRQ(ierr);CHKERRQ(ierr); 246641875f9SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer,"Common.rowfacfl %g (number of flops in last call to cholmod_rowfac)\n",c->rowfacfl);CHKERRQ(ierr);CHKERRQ(ierr); 247641875f9SMatthew G Knepley ierr = PetscViewerASCIIPrintf(viewer,"Common.aatfl %g (number of flops to compute A(:,f)*A(:,f)')\n",c->aatfl);CHKERRQ(ierr);CHKERRQ(ierr); 248b9eaa5e8SBarry Smith #if defined(PETSC_USE_SUITESPARSE_GPU) 249b9eaa5e8SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Common.useGPU %d\n",c->useGPU);CHKERRQ(ierr);CHKERRQ(ierr); 250b9eaa5e8SBarry Smith #endif 251641875f9SMatthew G Knepley ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 252641875f9SMatthew G Knepley PetscFunctionReturn(0); 253641875f9SMatthew G Knepley } 254641875f9SMatthew G Knepley 255eb9872f6SBarry Smith PETSC_INTERN PetscErrorCode MatView_CHOLMOD(Mat F,PetscViewer viewer) 256641875f9SMatthew G Knepley { 257641875f9SMatthew G Knepley PetscErrorCode ierr; 258ace3abfcSBarry Smith PetscBool iascii; 259641875f9SMatthew G Knepley PetscViewerFormat format; 260641875f9SMatthew G Knepley 261641875f9SMatthew G Knepley PetscFunctionBegin; 262251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 263641875f9SMatthew G Knepley if (iascii) { 264641875f9SMatthew G Knepley ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 265641875f9SMatthew G Knepley if (format == PETSC_VIEWER_ASCII_INFO) { 266*860c79edSBarry Smith ierr = MatView_Info_CHOLMOD(F,viewer);CHKERRQ(ierr); 267641875f9SMatthew G Knepley } 268641875f9SMatthew G Knepley } 269641875f9SMatthew G Knepley PetscFunctionReturn(0); 270641875f9SMatthew G Knepley } 271641875f9SMatthew G Knepley 272641875f9SMatthew G Knepley static PetscErrorCode MatSolve_CHOLMOD(Mat F,Vec B,Vec X) 273641875f9SMatthew G Knepley { 2746b8f6f9dSBarry Smith Mat_CHOLMOD *chol = (Mat_CHOLMOD*)F->data; 275641875f9SMatthew G Knepley cholmod_dense cholB,*cholX; 276641875f9SMatthew G Knepley PetscScalar *x; 277641875f9SMatthew G Knepley PetscErrorCode ierr; 278641875f9SMatthew G Knepley 279641875f9SMatthew G Knepley PetscFunctionBegin; 280d9ca1df4SBarry Smith ierr = VecWrapCholmodRead(B,&cholB);CHKERRQ(ierr); 281641875f9SMatthew G Knepley static_F = F; 282641875f9SMatthew G Knepley cholX = cholmod_X_solve(CHOLMOD_A,chol->factor,&cholB,chol->common); 283641875f9SMatthew G Knepley if (!cholX) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CHOLMOD failed"); 284d9ca1df4SBarry Smith ierr = VecUnWrapCholmodRead(B,&cholB);CHKERRQ(ierr); 285641875f9SMatthew G Knepley ierr = VecGetArray(X,&x);CHKERRQ(ierr); 286641875f9SMatthew G Knepley ierr = PetscMemcpy(x,cholX->x,cholX->nrow*sizeof(*x));CHKERRQ(ierr); 287641875f9SMatthew G Knepley ierr = !cholmod_X_free_dense(&cholX,chol->common);CHKERRQ(ierr); 288641875f9SMatthew G Knepley ierr = VecRestoreArray(X,&x);CHKERRQ(ierr); 289641875f9SMatthew G Knepley PetscFunctionReturn(0); 290641875f9SMatthew G Knepley } 291641875f9SMatthew G Knepley 292641875f9SMatthew G Knepley static PetscErrorCode MatCholeskyFactorNumeric_CHOLMOD(Mat F,Mat A,const MatFactorInfo *info) 293641875f9SMatthew G Knepley { 2946b8f6f9dSBarry Smith Mat_CHOLMOD *chol = (Mat_CHOLMOD*)F->data; 295641875f9SMatthew G Knepley cholmod_sparse cholA; 296ace3abfcSBarry Smith PetscBool aijalloc; 297641875f9SMatthew G Knepley PetscErrorCode ierr; 298641875f9SMatthew G Knepley 299641875f9SMatthew G Knepley PetscFunctionBegin; 300641875f9SMatthew G Knepley ierr = (*chol->Wrap)(A,PETSC_TRUE,&cholA,&aijalloc);CHKERRQ(ierr); 301641875f9SMatthew G Knepley static_F = F; 302641875f9SMatthew G Knepley ierr = !cholmod_X_factorize(&cholA,chol->factor,chol->common); 303ce94432eSBarry Smith if (ierr) SETERRQ1(PetscObjectComm((PetscObject)F),PETSC_ERR_LIB,"CHOLMOD factorization failed with status %d",chol->common->status); 304ce94432eSBarry Smith if (chol->common->status == CHOLMOD_NOT_POSDEF) SETERRQ1(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); 305641875f9SMatthew G Knepley 306641875f9SMatthew G Knepley if (aijalloc) {ierr = PetscFree3(cholA.p,cholA.i,cholA.x);CHKERRQ(ierr);} 307641875f9SMatthew G Knepley 308641875f9SMatthew G Knepley F->ops->solve = MatSolve_CHOLMOD; 309641875f9SMatthew G Knepley F->ops->solvetranspose = MatSolve_CHOLMOD; 310641875f9SMatthew G Knepley PetscFunctionReturn(0); 311641875f9SMatthew G Knepley } 312641875f9SMatthew G Knepley 313eb9872f6SBarry Smith PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_CHOLMOD(Mat F,Mat A,IS perm,const MatFactorInfo *info) 314641875f9SMatthew G Knepley { 3156b8f6f9dSBarry Smith Mat_CHOLMOD *chol = (Mat_CHOLMOD*)F->data; 316641875f9SMatthew G Knepley PetscErrorCode ierr; 317641875f9SMatthew G Knepley cholmod_sparse cholA; 318ace3abfcSBarry Smith PetscBool aijalloc; 319641875f9SMatthew G Knepley PetscInt *fset = 0; 320641875f9SMatthew G Knepley size_t fsize = 0; 321641875f9SMatthew G Knepley 322641875f9SMatthew G Knepley PetscFunctionBegin; 323641875f9SMatthew G Knepley ierr = (*chol->Wrap)(A,PETSC_FALSE,&cholA,&aijalloc);CHKERRQ(ierr); 324641875f9SMatthew G Knepley static_F = F; 325641875f9SMatthew G Knepley if (chol->factor) { 326641875f9SMatthew G Knepley ierr = !cholmod_X_resymbol(&cholA,fset,fsize,(int)chol->pack,chol->factor,chol->common); 327ce94432eSBarry Smith if (ierr) SETERRQ1(PetscObjectComm((PetscObject)F),PETSC_ERR_LIB,"CHOLMOD analysis failed with status %d",chol->common->status); 328641875f9SMatthew G Knepley } else if (perm) { 329641875f9SMatthew G Knepley const PetscInt *ip; 330641875f9SMatthew G Knepley ierr = ISGetIndices(perm,&ip);CHKERRQ(ierr); 331641875f9SMatthew G Knepley chol->factor = cholmod_X_analyze_p(&cholA,(PetscInt*)ip,fset,fsize,chol->common); 332ce94432eSBarry Smith if (!chol->factor) SETERRQ1(PetscObjectComm((PetscObject)F),PETSC_ERR_LIB,"CHOLMOD analysis failed with status %d",chol->common->status); 333641875f9SMatthew G Knepley ierr = ISRestoreIndices(perm,&ip);CHKERRQ(ierr); 334641875f9SMatthew G Knepley } else { 335641875f9SMatthew G Knepley chol->factor = cholmod_X_analyze(&cholA,chol->common); 336ce94432eSBarry Smith if (!chol->factor) SETERRQ1(PetscObjectComm((PetscObject)F),PETSC_ERR_LIB,"CHOLMOD analysis failed with status %d",chol->common->status); 337641875f9SMatthew G Knepley } 338641875f9SMatthew G Knepley 339641875f9SMatthew G Knepley if (aijalloc) {ierr = PetscFree3(cholA.p,cholA.i,cholA.x);CHKERRQ(ierr);} 340641875f9SMatthew G Knepley 341641875f9SMatthew G Knepley F->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_CHOLMOD; 342641875f9SMatthew G Knepley PetscFunctionReturn(0); 343641875f9SMatthew G Knepley } 344641875f9SMatthew G Knepley 345db87b0f2SBarry Smith static PetscErrorCode MatFactorGetSolverPackage_seqsbaij_cholmod(Mat A,const MatSolverPackage *type) 346641875f9SMatthew G Knepley { 347641875f9SMatthew G Knepley PetscFunctionBegin; 348641875f9SMatthew G Knepley *type = MATSOLVERCHOLMOD; 349641875f9SMatthew G Knepley PetscFunctionReturn(0); 350641875f9SMatthew G Knepley } 351641875f9SMatthew G Knepley 352641875f9SMatthew G Knepley /*MC 353641875f9SMatthew G Knepley MATSOLVERCHOLMOD = "cholmod" - A matrix type providing direct solvers (Cholesky) for sequential matrices 354641875f9SMatthew G Knepley via the external package CHOLMOD. 355641875f9SMatthew G Knepley 356c2b89b5dSBarry Smith Use ./configure --download-suitesparse to install PETSc to use CHOLMOD 357c2b89b5dSBarry Smith 358c2b89b5dSBarry Smith Use -pc_type lu -pc_factor_mat_solver_package cholmod to use this direct solver 359641875f9SMatthew G Knepley 360641875f9SMatthew G Knepley Consult CHOLMOD documentation for more information about the Common parameters 361641875f9SMatthew G Knepley which correspond to the options database keys below. 362641875f9SMatthew G Knepley 363641875f9SMatthew G Knepley Options Database Keys: 364e08999f5SMatthew G Knepley + -mat_cholmod_dbound <0> - Minimum absolute value of diagonal entries of D (None) 365e08999f5SMatthew G Knepley . -mat_cholmod_grow0 <1.2> - Global growth ratio when factors are modified (None) 366e08999f5SMatthew G Knepley . -mat_cholmod_grow1 <1.2> - Column growth ratio when factors are modified (None) 367e08999f5SMatthew G Knepley . -mat_cholmod_grow2 <5> - Affine column growth constant when factors are modified (None) 368e08999f5SMatthew G Knepley . -mat_cholmod_maxrank <8> - Max rank of update, larger values are faster but use more memory [2,4,8] (None) 369e08999f5SMatthew G Knepley . -mat_cholmod_factor <AUTO> - (choose one of) SIMPLICIAL AUTO SUPERNODAL 370e08999f5SMatthew G Knepley . -mat_cholmod_supernodal_switch <40> - flop/nnz_L threshold for switching to supernodal factorization (None) 371e08999f5SMatthew G Knepley . -mat_cholmod_final_asis <TRUE> - Leave factors "as is" (None) 372e08999f5SMatthew G Knepley . -mat_cholmod_final_pack <TRUE> - Pack the columns when finished (use FALSE if the factors will be updated later) (None) 373e08999f5SMatthew G Knepley . -mat_cholmod_zrelax <0.8> - 3 real supernodal relaxed amalgamation parameters (None) 374e08999f5SMatthew G Knepley . -mat_cholmod_nrelax <4> - 3 size_t supernodal relaxed amalgamation parameters (None) 375e08999f5SMatthew G Knepley . -mat_cholmod_prefer_upper <TRUE> - Work with upper triangular form (faster when using fill-reducing ordering, slower in natural ordering) (None) 376e08999f5SMatthew G Knepley - -mat_cholmod_print <3> - Verbosity level (None) 377641875f9SMatthew G Knepley 378641875f9SMatthew G Knepley Level: beginner 379641875f9SMatthew G Knepley 380a364b7d2SBarry Smith Note: CHOLMOD is part of SuiteSparse http://faculty.cse.tamu.edu/davis/suitesparse.html 381a364b7d2SBarry Smith 382641875f9SMatthew G Knepley .seealso: PCCHOLESKY, PCFactorSetMatSolverPackage(), MatSolverPackage 383641875f9SMatthew G Knepley M*/ 384b2573a8aSBarry Smith 385db87b0f2SBarry Smith PETSC_INTERN PetscErrorCode MatGetFactor_seqsbaij_cholmod(Mat A,MatFactorType ftype,Mat *F) 386641875f9SMatthew G Knepley { 387641875f9SMatthew G Knepley Mat B; 388641875f9SMatthew G Knepley Mat_CHOLMOD *chol; 389641875f9SMatthew G Knepley PetscErrorCode ierr; 390641875f9SMatthew G Knepley PetscInt m=A->rmap->n,n=A->cmap->n,bs; 391641875f9SMatthew G Knepley 392641875f9SMatthew G Knepley PetscFunctionBegin; 393641875f9SMatthew G Knepley if (ftype != MAT_FACTOR_CHOLESKY) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"CHOLMOD cannot do %s factorization with SBAIJ, only %s", 394641875f9SMatthew G Knepley MatFactorTypes[ftype],MatFactorTypes[MAT_FACTOR_CHOLESKY]); 395641875f9SMatthew G Knepley ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr); 396ce94432eSBarry Smith if (bs != 1) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"CHOLMOD only supports block size=1, given %D",bs); 397641875f9SMatthew G Knepley /* Create the factorization matrix F */ 398ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr); 399641875f9SMatthew G Knepley ierr = MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,m,n);CHKERRQ(ierr); 4006b8f6f9dSBarry Smith ierr = PetscStrallocpy("cholmod",&((PetscObject)B)->type_name);CHKERRQ(ierr); 4016b8f6f9dSBarry Smith ierr = MatSetUp(B);CHKERRQ(ierr); 402b00a9115SJed Brown ierr = PetscNewLog(B,&chol);CHKERRQ(ierr); 40326fbe8dcSKarl Rupp 404641875f9SMatthew G Knepley chol->Wrap = MatWrapCholmod_seqsbaij; 4056b8f6f9dSBarry Smith B->data = chol; 406641875f9SMatthew G Knepley 4076b8f6f9dSBarry Smith B->ops->getinfo = MatGetInfo_External; 408641875f9SMatthew G Knepley B->ops->view = MatView_CHOLMOD; 409641875f9SMatthew G Knepley B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_CHOLMOD; 410641875f9SMatthew G Knepley B->ops->destroy = MatDestroy_CHOLMOD; 411bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_seqsbaij_cholmod);CHKERRQ(ierr); 412641875f9SMatthew G Knepley B->factortype = MAT_FACTOR_CHOLESKY; 413641875f9SMatthew G Knepley B->assembled = PETSC_TRUE; /* required by -ksp_view */ 414641875f9SMatthew G Knepley B->preallocated = PETSC_TRUE; 415641875f9SMatthew G Knepley 416641875f9SMatthew G Knepley ierr = CholmodStart(B);CHKERRQ(ierr); 41700c67f3bSHong Zhang 41800c67f3bSHong Zhang ierr = PetscFree(B->solvertype);CHKERRQ(ierr); 41900c67f3bSHong Zhang ierr = PetscStrallocpy(MATSOLVERCHOLMOD,&B->solvertype);CHKERRQ(ierr); 42000c67f3bSHong Zhang 421641875f9SMatthew G Knepley *F = B; 422641875f9SMatthew G Knepley PetscFunctionReturn(0); 423641875f9SMatthew G Knepley } 424