xref: /petsc/src/ksp/pc/impls/factor/cholesky/cholesky.c (revision 2e956fe4fc852fabc23b437482e1fb7b77fddb0d)
19b54502bSHong Zhang 
29b54502bSHong Zhang /*
39b54502bSHong Zhang    Defines a direct factorization preconditioner for any Mat implementation
49b54502bSHong Zhang    Note: this need not be consided a preconditioner since it supplies
59b54502bSHong Zhang          a direct solver.
69b54502bSHong Zhang */
7c6db04a5SJed Brown #include <../src/ksp/pc/impls/factor/factor.h>         /*I "petscpc.h" I*/
89b54502bSHong Zhang 
99b54502bSHong Zhang typedef struct {
10075768bcSBarry Smith   PC_Factor hdr;
119b54502bSHong Zhang   IS        row,col;                 /* index sets used for reordering */
129b54502bSHong Zhang } PC_Cholesky;
139b54502bSHong Zhang 
144416b707SBarry Smith static PetscErrorCode PCSetFromOptions_Cholesky(PetscOptionItems *PetscOptionsObject,PC pc)
159b54502bSHong Zhang {
169b54502bSHong Zhang   PetscFunctionBegin;
17d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject,"Cholesky options");
189566063dSJacob Faibussowitsch   PetscCall(PCSetFromOptions_Factor(PetscOptionsObject,pc));
19d0609cedSBarry Smith   PetscOptionsHeadEnd();
209b54502bSHong Zhang   PetscFunctionReturn(0);
219b54502bSHong Zhang }
229b54502bSHong Zhang 
239b54502bSHong Zhang static PetscErrorCode PCSetUp_Cholesky(PC pc)
249b54502bSHong Zhang {
25ace3abfcSBarry Smith   PetscBool              flg;
269b54502bSHong Zhang   PC_Cholesky            *dir = (PC_Cholesky*)pc->data;
27ea799195SBarry Smith   MatSolverType          stype;
2800e125f8SBarry Smith   MatFactorError         err;
29f023e1d5SPierre Jolivet   const char             *prefix;
309b54502bSHong Zhang 
319b54502bSHong Zhang   PetscFunctionBegin;
32a7f29c02SPierre Jolivet   if (!((PetscObject)pc->pmat)->prefix) {
33f023e1d5SPierre Jolivet     PetscCall(PCGetOptionsPrefix(pc,&prefix));
34f023e1d5SPierre Jolivet     PetscCall(MatSetOptionsPrefix(pc->pmat,prefix));
35a7f29c02SPierre Jolivet   }
36c6e4fdc6SHong Zhang   pc->failedreason = PC_NOERROR;
373d1c1ea0SBarry Smith   if (dir->hdr.reusefill && pc->setupcalled) ((PC_Factor*)dir)->info.fill = dir->hdr.actualfill;
389b54502bSHong Zhang 
399566063dSJacob Faibussowitsch   PetscCall(MatSetErrorIfFailure(pc->pmat,pc->erroriffailure));
403d1c1ea0SBarry Smith   if (dir->hdr.inplace) {
419b54502bSHong Zhang     if (dir->row && dir->col && (dir->row != dir->col)) {
429566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&dir->row));
439b54502bSHong Zhang     }
449566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&dir->col));
452c7c0729SBarry Smith     /* should only get reordering if the factor matrix uses it but cannot determine because MatGetFactor() not called */
469566063dSJacob Faibussowitsch     PetscCall(PCFactorSetDefaultOrdering_Factor(pc));
479566063dSJacob Faibussowitsch     PetscCall(MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col));
489b54502bSHong Zhang     if (dir->col && (dir->row != dir->col)) {  /* only use row ordering for SBAIJ */
499566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&dir->col));
509b54502bSHong Zhang     }
519566063dSJacob Faibussowitsch     if (dir->row) PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row));
529566063dSJacob Faibussowitsch     PetscCall(MatCholeskyFactor(pc->pmat,dir->row,&((PC_Factor*)dir)->info));
539566063dSJacob Faibussowitsch     PetscCall(MatFactorGetError(pc->pmat,&err));
5400e125f8SBarry Smith     if (err) { /* Factor() fails */
5500e125f8SBarry Smith       pc->failedreason = (PCFailedReason)err;
566baea169SHong Zhang       PetscFunctionReturn(0);
576baea169SHong Zhang     }
582fa5cd67SKarl Rupp 
59075768bcSBarry Smith     ((PC_Factor*)dir)->fact = pc->pmat;
609b54502bSHong Zhang   } else {
619b54502bSHong Zhang     MatInfo info;
6200e125f8SBarry Smith 
639b54502bSHong Zhang     if (!pc->setupcalled) {
64f73b0415SBarry Smith       PetscBool canuseordering;
652c7c0729SBarry Smith       if (!((PC_Factor*)dir)->fact) {
669566063dSJacob Faibussowitsch         PetscCall(MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_CHOLESKY,&((PC_Factor*)dir)->fact));
679566063dSJacob Faibussowitsch         PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact));
682c7c0729SBarry Smith       }
699566063dSJacob Faibussowitsch       PetscCall(MatFactorGetCanUseOrdering(((PC_Factor*)dir)->fact,&canuseordering));
70f73b0415SBarry Smith       if (canuseordering) {
719566063dSJacob Faibussowitsch         PetscCall(PCFactorSetDefaultOrdering_Factor(pc));
729566063dSJacob Faibussowitsch         PetscCall(MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col));
739bfd6278SHong Zhang         /* check if dir->row == dir->col */
744ac6704cSBarry Smith         if (dir->row) {
759566063dSJacob Faibussowitsch           PetscCall(ISEqual(dir->row,dir->col,&flg));
7628b400f6SJacob Faibussowitsch           PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"row and column permutations must be equal");
774ac6704cSBarry Smith         }
789566063dSJacob Faibussowitsch         PetscCall(ISDestroy(&dir->col)); /* only pass one ordering into CholeskyFactor */
799bfd6278SHong Zhang 
8090d69ab7SBarry Smith         flg  = PETSC_FALSE;
819566063dSJacob Faibussowitsch         PetscCall(PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&flg,NULL));
829b54502bSHong Zhang         if (flg) {
839b54502bSHong Zhang           PetscReal tol = 1.e-10;
849566063dSJacob Faibussowitsch           PetscCall(PetscOptionsGetReal(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&tol,NULL));
859566063dSJacob Faibussowitsch           PetscCall(MatReorderForNonzeroDiagonal(pc->pmat,tol,dir->row,dir->row));
869b54502bSHong Zhang         }
879566063dSJacob Faibussowitsch         PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row));
88a1f19f5aSHong Zhang       }
899566063dSJacob Faibussowitsch       PetscCall(MatCholeskyFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,&((PC_Factor*)dir)->info));
909566063dSJacob Faibussowitsch       PetscCall(MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info));
913d1c1ea0SBarry Smith       dir->hdr.actualfill = info.fill_ratio_needed;
929b54502bSHong Zhang     } else if (pc->flag != SAME_NONZERO_PATTERN) {
933d1c1ea0SBarry Smith       if (!dir->hdr.reuseordering) {
94f73b0415SBarry Smith         PetscBool canuseordering;
959566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&((PC_Factor*)dir)->fact));
969566063dSJacob Faibussowitsch         PetscCall(MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_CHOLESKY,&((PC_Factor*)dir)->fact));
979566063dSJacob Faibussowitsch         PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact));
989566063dSJacob Faibussowitsch         PetscCall(MatFactorGetCanUseOrdering(((PC_Factor*)dir)->fact,&canuseordering));
99f73b0415SBarry Smith         if (canuseordering) {
1009566063dSJacob Faibussowitsch           PetscCall(ISDestroy(&dir->row));
1019566063dSJacob Faibussowitsch           PetscCall(PCFactorSetDefaultOrdering_Factor(pc));
1029566063dSJacob Faibussowitsch           PetscCall(MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col));
1039566063dSJacob Faibussowitsch           PetscCall(ISDestroy(&dir->col)); /* only use dir->row ordering in CholeskyFactor */
1049d61c4f5SHong Zhang 
10590d69ab7SBarry Smith           flg  = PETSC_FALSE;
1069566063dSJacob Faibussowitsch           PetscCall(PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&flg,NULL));
1079b54502bSHong Zhang           if (flg) {
1089b54502bSHong Zhang             PetscReal tol = 1.e-10;
1099566063dSJacob Faibussowitsch             PetscCall(PetscOptionsGetReal(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&tol,NULL));
1109566063dSJacob Faibussowitsch             PetscCall(MatReorderForNonzeroDiagonal(pc->pmat,tol,dir->row,dir->row));
1119b54502bSHong Zhang           }
1129566063dSJacob Faibussowitsch           PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row));
1139b54502bSHong Zhang         }
1142c7c0729SBarry Smith       }
1159566063dSJacob Faibussowitsch       PetscCall(MatCholeskyFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,&((PC_Factor*)dir)->info));
1169566063dSJacob Faibussowitsch       PetscCall(MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info));
1173d1c1ea0SBarry Smith       dir->hdr.actualfill = info.fill_ratio_needed;
1189566063dSJacob Faibussowitsch       PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact));
11904545d6dSBarry Smith     } else {
1209566063dSJacob Faibussowitsch       PetscCall(MatFactorGetError(((PC_Factor*)dir)->fact,&err));
121160a8794SBarry Smith       if (err == MAT_FACTOR_NUMERIC_ZEROPIVOT) {
1229566063dSJacob Faibussowitsch         PetscCall(MatFactorClearError(((PC_Factor*)dir)->fact));
123b8b68cfdSBarry Smith         pc->failedreason = PC_NOERROR;
12404545d6dSBarry Smith       }
1259b54502bSHong Zhang     }
1269566063dSJacob Faibussowitsch     PetscCall(MatFactorGetError(((PC_Factor*)dir)->fact,&err));
12700e125f8SBarry Smith     if (err) { /* FactorSymbolic() fails */
12800e125f8SBarry Smith       pc->failedreason = (PCFailedReason)err;
1296baea169SHong Zhang       PetscFunctionReturn(0);
1306baea169SHong Zhang     }
1316baea169SHong Zhang 
1329566063dSJacob Faibussowitsch     PetscCall(MatCholeskyFactorNumeric(((PC_Factor*)dir)->fact,pc->pmat,&((PC_Factor*)dir)->info));
1339566063dSJacob Faibussowitsch     PetscCall(MatFactorGetError(((PC_Factor*)dir)->fact,&err));
13400e125f8SBarry Smith     if (err) { /* FactorNumeric() fails */
13500e125f8SBarry Smith       pc->failedreason = (PCFailedReason)err;
1366baea169SHong Zhang     }
1379b54502bSHong Zhang   }
13800c67f3bSHong Zhang 
1399566063dSJacob Faibussowitsch   PetscCall(PCFactorGetMatSolverType(pc,&stype));
14000c67f3bSHong Zhang   if (!stype) {
141ea799195SBarry Smith     MatSolverType solverpackage;
1429566063dSJacob Faibussowitsch     PetscCall(MatFactorGetSolverType(((PC_Factor*)dir)->fact,&solverpackage));
1439566063dSJacob Faibussowitsch     PetscCall(PCFactorSetMatSolverType(pc,solverpackage));
14400c67f3bSHong Zhang   }
1459b54502bSHong Zhang   PetscFunctionReturn(0);
1469b54502bSHong Zhang }
1479b54502bSHong Zhang 
14869d2c0f9SBarry Smith static PetscErrorCode PCReset_Cholesky(PC pc)
1499b54502bSHong Zhang {
1509b54502bSHong Zhang   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;
1519b54502bSHong Zhang 
1529b54502bSHong Zhang   PetscFunctionBegin;
1539566063dSJacob Faibussowitsch   if (!dir->hdr.inplace && ((PC_Factor*)dir)->fact) PetscCall(MatDestroy(&((PC_Factor*)dir)->fact));
1549566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&dir->row));
1559566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&dir->col));
15669d2c0f9SBarry Smith   PetscFunctionReturn(0);
15769d2c0f9SBarry Smith }
15869d2c0f9SBarry Smith 
15969d2c0f9SBarry Smith static PetscErrorCode PCDestroy_Cholesky(PC pc)
16069d2c0f9SBarry Smith {
16169d2c0f9SBarry Smith   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;
16269d2c0f9SBarry Smith 
16369d2c0f9SBarry Smith   PetscFunctionBegin;
1649566063dSJacob Faibussowitsch   PetscCall(PCReset_Cholesky(pc));
1659566063dSJacob Faibussowitsch   PetscCall(PetscFree(((PC_Factor*)dir)->ordering));
1669566063dSJacob Faibussowitsch   PetscCall(PetscFree(((PC_Factor*)dir)->solvertype));
167*2e956fe4SStefano Zampini   PetscCall(PCFactorClearComposedFunctions(pc));
1689566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc->data));
1699b54502bSHong Zhang   PetscFunctionReturn(0);
1709b54502bSHong Zhang }
1719b54502bSHong Zhang 
1729b54502bSHong Zhang static PetscErrorCode PCApply_Cholesky(PC pc,Vec x,Vec y)
1739b54502bSHong Zhang {
1749b54502bSHong Zhang   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;
1759b54502bSHong Zhang 
1769b54502bSHong Zhang   PetscFunctionBegin;
1773d1c1ea0SBarry Smith   if (dir->hdr.inplace) {
1789566063dSJacob Faibussowitsch     PetscCall(MatSolve(pc->pmat,x,y));
1792fa5cd67SKarl Rupp   } else {
1809566063dSJacob Faibussowitsch     PetscCall(MatSolve(((PC_Factor*)dir)->fact,x,y));
1812fa5cd67SKarl Rupp   }
1829b54502bSHong Zhang   PetscFunctionReturn(0);
1839b54502bSHong Zhang }
1849b54502bSHong Zhang 
1857b6e2003SPierre Jolivet static PetscErrorCode PCMatApply_Cholesky(PC pc,Mat X,Mat Y)
1867b6e2003SPierre Jolivet {
1877b6e2003SPierre Jolivet   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;
1887b6e2003SPierre Jolivet 
1897b6e2003SPierre Jolivet   PetscFunctionBegin;
1907b6e2003SPierre Jolivet   if (dir->hdr.inplace) {
1919566063dSJacob Faibussowitsch     PetscCall(MatMatSolve(pc->pmat,X,Y));
1927b6e2003SPierre Jolivet   } else {
1939566063dSJacob Faibussowitsch     PetscCall(MatMatSolve(((PC_Factor*)dir)->fact,X,Y));
1947b6e2003SPierre Jolivet   }
1957b6e2003SPierre Jolivet   PetscFunctionReturn(0);
1967b6e2003SPierre Jolivet }
1977b6e2003SPierre Jolivet 
1983d72fe1eSJed Brown static PetscErrorCode PCApplySymmetricLeft_Cholesky(PC pc,Vec x,Vec y)
1993d72fe1eSJed Brown {
2003d72fe1eSJed Brown   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;
2013d72fe1eSJed Brown 
2023d72fe1eSJed Brown   PetscFunctionBegin;
2033d72fe1eSJed Brown   if (dir->hdr.inplace) {
2049566063dSJacob Faibussowitsch     PetscCall(MatForwardSolve(pc->pmat,x,y));
2053d72fe1eSJed Brown   } else {
2069566063dSJacob Faibussowitsch     PetscCall(MatForwardSolve(((PC_Factor*)dir)->fact,x,y));
2073d72fe1eSJed Brown   }
2083d72fe1eSJed Brown   PetscFunctionReturn(0);
2093d72fe1eSJed Brown }
2103d72fe1eSJed Brown 
2113d72fe1eSJed Brown static PetscErrorCode PCApplySymmetricRight_Cholesky(PC pc,Vec x,Vec y)
2123d72fe1eSJed Brown {
2133d72fe1eSJed Brown   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;
2143d72fe1eSJed Brown 
2153d72fe1eSJed Brown   PetscFunctionBegin;
2163d72fe1eSJed Brown   if (dir->hdr.inplace) {
2179566063dSJacob Faibussowitsch     PetscCall(MatBackwardSolve(pc->pmat,x,y));
2183d72fe1eSJed Brown   } else {
2199566063dSJacob Faibussowitsch     PetscCall(MatBackwardSolve(((PC_Factor*)dir)->fact,x,y));
2203d72fe1eSJed Brown   }
2213d72fe1eSJed Brown   PetscFunctionReturn(0);
2223d72fe1eSJed Brown }
2233d72fe1eSJed Brown 
2249b54502bSHong Zhang static PetscErrorCode PCApplyTranspose_Cholesky(PC pc,Vec x,Vec y)
2259b54502bSHong Zhang {
2269b54502bSHong Zhang   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;
2279b54502bSHong Zhang 
2289b54502bSHong Zhang   PetscFunctionBegin;
2293d1c1ea0SBarry Smith   if (dir->hdr.inplace) {
2309566063dSJacob Faibussowitsch     PetscCall(MatSolveTranspose(pc->pmat,x,y));
2312fa5cd67SKarl Rupp   } else {
2329566063dSJacob Faibussowitsch     PetscCall(MatSolveTranspose(((PC_Factor*)dir)->fact,x,y));
2332fa5cd67SKarl Rupp   }
2349b54502bSHong Zhang   PetscFunctionReturn(0);
2359b54502bSHong Zhang }
2369b54502bSHong Zhang 
2379b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/
2389b54502bSHong Zhang 
2399b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/
2409b54502bSHong Zhang 
2419b54502bSHong Zhang /*@
2422401956bSBarry Smith    PCFactorSetReuseOrdering - When similar matrices are factored, this
2439b54502bSHong Zhang    causes the ordering computed in the first factor to be used for all
2449b54502bSHong Zhang    following factors.
2459b54502bSHong Zhang 
246ad4df100SBarry Smith    Logically Collective on PC
2479b54502bSHong Zhang 
2489b54502bSHong Zhang    Input Parameters:
2499b54502bSHong Zhang +  pc - the preconditioner context
2509b54502bSHong Zhang -  flag - PETSC_TRUE to reuse else PETSC_FALSE
2519b54502bSHong Zhang 
2529b54502bSHong Zhang    Options Database Key:
2532401956bSBarry Smith .  -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering()
2549b54502bSHong Zhang 
2559b54502bSHong Zhang    Level: intermediate
2569b54502bSHong Zhang 
257db781477SPatrick Sanan .seealso: `PCFactorSetReuseFill()`
2589b54502bSHong Zhang @*/
2597087cfbeSBarry Smith PetscErrorCode  PCFactorSetReuseOrdering(PC pc,PetscBool flag)
2609b54502bSHong Zhang {
2619b54502bSHong Zhang   PetscFunctionBegin;
2620700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
263acfcf0e5SJed Brown   PetscValidLogicalCollectiveBool(pc,flag,2);
264cac4c232SBarry Smith   PetscTryMethod(pc,"PCFactorSetReuseOrdering_C",(PC,PetscBool),(pc,flag));
2659b54502bSHong Zhang   PetscFunctionReturn(0);
2669b54502bSHong Zhang }
2679b54502bSHong Zhang 
2689b54502bSHong Zhang /*MC
26996fc60bcSBarry Smith    PCCHOLESKY - Uses a direct solver, based on Cholesky factorization, as a preconditioner
2709b54502bSHong Zhang 
2719b54502bSHong Zhang    Options Database Keys:
2722401956bSBarry Smith +  -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering()
2733ca39a21SBarry Smith .  -pc_factor_mat_solver_type - Actives PCFactorSetMatSolverType() to choose the direct solver, like superlu
2742401956bSBarry Smith .  -pc_factor_reuse_fill - Activates PCFactorSetReuseFill()
27555ba2a51SBarry Smith .  -pc_factor_fill <fill> - Sets fill amount
2762401956bSBarry Smith .  -pc_factor_in_place - Activates in-place factorization
277145b9266SHong Zhang -  -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine
2789b54502bSHong Zhang 
27995452b02SPatrick Sanan    Notes:
28095452b02SPatrick Sanan     Not all options work for all matrix formats
2819b54502bSHong Zhang 
2829b54502bSHong Zhang    Level: beginner
2839b54502bSHong Zhang 
28495452b02SPatrick Sanan    Notes:
28595452b02SPatrick Sanan     Usually this will compute an "exact" solution in one iteration and does
2869b54502bSHong Zhang           not need a Krylov method (i.e. you can use -ksp_type preonly, or
2879b54502bSHong Zhang           KSPSetType(ksp,KSPPREONLY) for the Krylov method
2889b54502bSHong Zhang 
289db781477SPatrick Sanan .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`,
290db781477SPatrick Sanan           `PCILU`, `PCLU`, `PCICC`, `PCFactorSetReuseOrdering()`, `PCFactorSetReuseFill()`, `PCFactorGetMatrix()`,
291db781477SPatrick Sanan           `PCFactorSetFill()`, `PCFactorSetShiftNonzero()`, `PCFactorSetShiftType()`, `PCFactorSetShiftAmount()`
292db781477SPatrick Sanan           `PCFactorSetUseInPlace()`, `PCFactorGetUseInPlace()`, `PCFactorSetMatOrderingType()`
2939b54502bSHong Zhang 
2949b54502bSHong Zhang M*/
2959b54502bSHong Zhang 
2968cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_Cholesky(PC pc)
2979b54502bSHong Zhang {
2989b54502bSHong Zhang   PC_Cholesky    *dir;
2999b54502bSHong Zhang 
3009b54502bSHong Zhang   PetscFunctionBegin;
3019566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(pc,&dir));
3023d1c1ea0SBarry Smith   pc->data = (void*)dir;
3039566063dSJacob Faibussowitsch   PetscCall(PCFactorInitialize(pc,MAT_FACTOR_CHOLESKY));
3042fa5cd67SKarl Rupp 
305075768bcSBarry Smith   ((PC_Factor*)dir)->info.fill  = 5.0;
3062fa5cd67SKarl Rupp 
3079b54502bSHong Zhang   pc->ops->destroy             = PCDestroy_Cholesky;
30869d2c0f9SBarry Smith   pc->ops->reset               = PCReset_Cholesky;
3099b54502bSHong Zhang   pc->ops->apply               = PCApply_Cholesky;
3107b6e2003SPierre Jolivet   pc->ops->matapply            = PCMatApply_Cholesky;
3113d72fe1eSJed Brown   pc->ops->applysymmetricleft  = PCApplySymmetricLeft_Cholesky;
3123d72fe1eSJed Brown   pc->ops->applysymmetricright = PCApplySymmetricRight_Cholesky;
3139b54502bSHong Zhang   pc->ops->applytranspose      = PCApplyTranspose_Cholesky;
3149b54502bSHong Zhang   pc->ops->setup               = PCSetUp_Cholesky;
3159b54502bSHong Zhang   pc->ops->setfromoptions      = PCSetFromOptions_Cholesky;
31692e08861SBarry Smith   pc->ops->view                = PCView_Factor;
3170a545947SLisandro Dalcin   pc->ops->applyrichardson     = NULL;
3189b54502bSHong Zhang   PetscFunctionReturn(0);
3199b54502bSHong Zhang }
320