xref: /petsc/src/ksp/pc/impls/factor/cholesky/cholesky.c (revision dbbe0bcd3f3a8fbab5a45420dc06f8387e5764c6)
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 
14*dbbe0bcdSBarry Smith static PetscErrorCode PCSetFromOptions_Cholesky(PC pc,PetscOptionItems *PetscOptionsObject)
159b54502bSHong Zhang {
169b54502bSHong Zhang   PetscFunctionBegin;
17d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject,"Cholesky options");
18*dbbe0bcdSBarry Smith   PetscCall(PCSetFromOptions_Factor(pc,PetscOptionsObject));
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;
32c6e4fdc6SHong Zhang   pc->failedreason = PC_NOERROR;
333d1c1ea0SBarry Smith   if (dir->hdr.reusefill && pc->setupcalled) ((PC_Factor*)dir)->info.fill = dir->hdr.actualfill;
349b54502bSHong Zhang 
3526cc229bSBarry Smith   PetscCall(PCGetOptionsPrefix(pc,&prefix));
3626cc229bSBarry Smith   PetscCall(MatSetOptionsPrefixFactor(pc->pmat,prefix));
3726cc229bSBarry Smith 
389566063dSJacob Faibussowitsch   PetscCall(MatSetErrorIfFailure(pc->pmat,pc->erroriffailure));
393d1c1ea0SBarry Smith   if (dir->hdr.inplace) {
409b54502bSHong Zhang     if (dir->row && dir->col && (dir->row != dir->col)) {
419566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&dir->row));
429b54502bSHong Zhang     }
439566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&dir->col));
442c7c0729SBarry Smith     /* should only get reordering if the factor matrix uses it but cannot determine because MatGetFactor() not called */
459566063dSJacob Faibussowitsch     PetscCall(PCFactorSetDefaultOrdering_Factor(pc));
469566063dSJacob Faibussowitsch     PetscCall(MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col));
479b54502bSHong Zhang     if (dir->col && (dir->row != dir->col)) {  /* only use row ordering for SBAIJ */
489566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&dir->col));
499b54502bSHong Zhang     }
509566063dSJacob Faibussowitsch     if (dir->row) PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row));
519566063dSJacob Faibussowitsch     PetscCall(MatCholeskyFactor(pc->pmat,dir->row,&((PC_Factor*)dir)->info));
529566063dSJacob Faibussowitsch     PetscCall(MatFactorGetError(pc->pmat,&err));
5300e125f8SBarry Smith     if (err) { /* Factor() fails */
5400e125f8SBarry Smith       pc->failedreason = (PCFailedReason)err;
556baea169SHong Zhang       PetscFunctionReturn(0);
566baea169SHong Zhang     }
572fa5cd67SKarl Rupp 
58075768bcSBarry Smith     ((PC_Factor*)dir)->fact = pc->pmat;
599b54502bSHong Zhang   } else {
609b54502bSHong Zhang     MatInfo info;
6100e125f8SBarry Smith 
629b54502bSHong Zhang     if (!pc->setupcalled) {
63f73b0415SBarry Smith       PetscBool canuseordering;
642c7c0729SBarry Smith       if (!((PC_Factor*)dir)->fact) {
659566063dSJacob Faibussowitsch         PetscCall(MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_CHOLESKY,&((PC_Factor*)dir)->fact));
669566063dSJacob Faibussowitsch         PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact));
672c7c0729SBarry Smith       }
689566063dSJacob Faibussowitsch       PetscCall(MatFactorGetCanUseOrdering(((PC_Factor*)dir)->fact,&canuseordering));
69f73b0415SBarry Smith       if (canuseordering) {
709566063dSJacob Faibussowitsch         PetscCall(PCFactorSetDefaultOrdering_Factor(pc));
719566063dSJacob Faibussowitsch         PetscCall(MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col));
729bfd6278SHong Zhang         /* check if dir->row == dir->col */
734ac6704cSBarry Smith         if (dir->row) {
749566063dSJacob Faibussowitsch           PetscCall(ISEqual(dir->row,dir->col,&flg));
7528b400f6SJacob Faibussowitsch           PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"row and column permutations must be equal");
764ac6704cSBarry Smith         }
779566063dSJacob Faibussowitsch         PetscCall(ISDestroy(&dir->col)); /* only pass one ordering into CholeskyFactor */
789bfd6278SHong Zhang 
7990d69ab7SBarry Smith         flg  = PETSC_FALSE;
809566063dSJacob Faibussowitsch         PetscCall(PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&flg,NULL));
819b54502bSHong Zhang         if (flg) {
829b54502bSHong Zhang           PetscReal tol = 1.e-10;
839566063dSJacob Faibussowitsch           PetscCall(PetscOptionsGetReal(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&tol,NULL));
849566063dSJacob Faibussowitsch           PetscCall(MatReorderForNonzeroDiagonal(pc->pmat,tol,dir->row,dir->row));
859b54502bSHong Zhang         }
869566063dSJacob Faibussowitsch         PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row));
87a1f19f5aSHong Zhang       }
889566063dSJacob Faibussowitsch       PetscCall(MatCholeskyFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,&((PC_Factor*)dir)->info));
899566063dSJacob Faibussowitsch       PetscCall(MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info));
903d1c1ea0SBarry Smith       dir->hdr.actualfill = info.fill_ratio_needed;
919b54502bSHong Zhang     } else if (pc->flag != SAME_NONZERO_PATTERN) {
923d1c1ea0SBarry Smith       if (!dir->hdr.reuseordering) {
93f73b0415SBarry Smith         PetscBool canuseordering;
949566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&((PC_Factor*)dir)->fact));
959566063dSJacob Faibussowitsch         PetscCall(MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_CHOLESKY,&((PC_Factor*)dir)->fact));
969566063dSJacob Faibussowitsch         PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact));
979566063dSJacob Faibussowitsch         PetscCall(MatFactorGetCanUseOrdering(((PC_Factor*)dir)->fact,&canuseordering));
98f73b0415SBarry Smith         if (canuseordering) {
999566063dSJacob Faibussowitsch           PetscCall(ISDestroy(&dir->row));
1009566063dSJacob Faibussowitsch           PetscCall(PCFactorSetDefaultOrdering_Factor(pc));
1019566063dSJacob Faibussowitsch           PetscCall(MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col));
1029566063dSJacob Faibussowitsch           PetscCall(ISDestroy(&dir->col)); /* only use dir->row ordering in CholeskyFactor */
1039d61c4f5SHong Zhang 
10490d69ab7SBarry Smith           flg  = PETSC_FALSE;
1059566063dSJacob Faibussowitsch           PetscCall(PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&flg,NULL));
1069b54502bSHong Zhang           if (flg) {
1079b54502bSHong Zhang             PetscReal tol = 1.e-10;
1089566063dSJacob Faibussowitsch             PetscCall(PetscOptionsGetReal(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&tol,NULL));
1099566063dSJacob Faibussowitsch             PetscCall(MatReorderForNonzeroDiagonal(pc->pmat,tol,dir->row,dir->row));
1109b54502bSHong Zhang           }
1119566063dSJacob Faibussowitsch           PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row));
1129b54502bSHong Zhang         }
1132c7c0729SBarry Smith       }
1149566063dSJacob Faibussowitsch       PetscCall(MatCholeskyFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,&((PC_Factor*)dir)->info));
1159566063dSJacob Faibussowitsch       PetscCall(MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info));
1163d1c1ea0SBarry Smith       dir->hdr.actualfill = info.fill_ratio_needed;
1179566063dSJacob Faibussowitsch       PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact));
11804545d6dSBarry Smith     } else {
1199566063dSJacob Faibussowitsch       PetscCall(MatFactorGetError(((PC_Factor*)dir)->fact,&err));
120160a8794SBarry Smith       if (err == MAT_FACTOR_NUMERIC_ZEROPIVOT) {
1219566063dSJacob Faibussowitsch         PetscCall(MatFactorClearError(((PC_Factor*)dir)->fact));
122b8b68cfdSBarry Smith         pc->failedreason = PC_NOERROR;
12304545d6dSBarry Smith       }
1249b54502bSHong Zhang     }
1259566063dSJacob Faibussowitsch     PetscCall(MatFactorGetError(((PC_Factor*)dir)->fact,&err));
12600e125f8SBarry Smith     if (err) { /* FactorSymbolic() fails */
12700e125f8SBarry Smith       pc->failedreason = (PCFailedReason)err;
1286baea169SHong Zhang       PetscFunctionReturn(0);
1296baea169SHong Zhang     }
1306baea169SHong Zhang 
1319566063dSJacob Faibussowitsch     PetscCall(MatCholeskyFactorNumeric(((PC_Factor*)dir)->fact,pc->pmat,&((PC_Factor*)dir)->info));
1329566063dSJacob Faibussowitsch     PetscCall(MatFactorGetError(((PC_Factor*)dir)->fact,&err));
13300e125f8SBarry Smith     if (err) { /* FactorNumeric() fails */
13400e125f8SBarry Smith       pc->failedreason = (PCFailedReason)err;
1356baea169SHong Zhang     }
1369b54502bSHong Zhang   }
13700c67f3bSHong Zhang 
1389566063dSJacob Faibussowitsch   PetscCall(PCFactorGetMatSolverType(pc,&stype));
13900c67f3bSHong Zhang   if (!stype) {
140ea799195SBarry Smith     MatSolverType solverpackage;
1419566063dSJacob Faibussowitsch     PetscCall(MatFactorGetSolverType(((PC_Factor*)dir)->fact,&solverpackage));
1429566063dSJacob Faibussowitsch     PetscCall(PCFactorSetMatSolverType(pc,solverpackage));
14300c67f3bSHong Zhang   }
1449b54502bSHong Zhang   PetscFunctionReturn(0);
1459b54502bSHong Zhang }
1469b54502bSHong Zhang 
14769d2c0f9SBarry Smith static PetscErrorCode PCReset_Cholesky(PC pc)
1489b54502bSHong Zhang {
1499b54502bSHong Zhang   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;
1509b54502bSHong Zhang 
1519b54502bSHong Zhang   PetscFunctionBegin;
1529566063dSJacob Faibussowitsch   if (!dir->hdr.inplace && ((PC_Factor*)dir)->fact) PetscCall(MatDestroy(&((PC_Factor*)dir)->fact));
1539566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&dir->row));
1549566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&dir->col));
15569d2c0f9SBarry Smith   PetscFunctionReturn(0);
15669d2c0f9SBarry Smith }
15769d2c0f9SBarry Smith 
15869d2c0f9SBarry Smith static PetscErrorCode PCDestroy_Cholesky(PC pc)
15969d2c0f9SBarry Smith {
16069d2c0f9SBarry Smith   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;
16169d2c0f9SBarry Smith 
16269d2c0f9SBarry Smith   PetscFunctionBegin;
1639566063dSJacob Faibussowitsch   PetscCall(PCReset_Cholesky(pc));
1649566063dSJacob Faibussowitsch   PetscCall(PetscFree(((PC_Factor*)dir)->ordering));
1659566063dSJacob Faibussowitsch   PetscCall(PetscFree(((PC_Factor*)dir)->solvertype));
1662e956fe4SStefano Zampini   PetscCall(PCFactorClearComposedFunctions(pc));
1679566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc->data));
1689b54502bSHong Zhang   PetscFunctionReturn(0);
1699b54502bSHong Zhang }
1709b54502bSHong Zhang 
1719b54502bSHong Zhang static PetscErrorCode PCApply_Cholesky(PC pc,Vec x,Vec y)
1729b54502bSHong Zhang {
1739b54502bSHong Zhang   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;
1749b54502bSHong Zhang 
1759b54502bSHong Zhang   PetscFunctionBegin;
1763d1c1ea0SBarry Smith   if (dir->hdr.inplace) {
1779566063dSJacob Faibussowitsch     PetscCall(MatSolve(pc->pmat,x,y));
1782fa5cd67SKarl Rupp   } else {
1799566063dSJacob Faibussowitsch     PetscCall(MatSolve(((PC_Factor*)dir)->fact,x,y));
1802fa5cd67SKarl Rupp   }
1819b54502bSHong Zhang   PetscFunctionReturn(0);
1829b54502bSHong Zhang }
1839b54502bSHong Zhang 
1847b6e2003SPierre Jolivet static PetscErrorCode PCMatApply_Cholesky(PC pc,Mat X,Mat Y)
1857b6e2003SPierre Jolivet {
1867b6e2003SPierre Jolivet   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;
1877b6e2003SPierre Jolivet 
1887b6e2003SPierre Jolivet   PetscFunctionBegin;
1897b6e2003SPierre Jolivet   if (dir->hdr.inplace) {
1909566063dSJacob Faibussowitsch     PetscCall(MatMatSolve(pc->pmat,X,Y));
1917b6e2003SPierre Jolivet   } else {
1929566063dSJacob Faibussowitsch     PetscCall(MatMatSolve(((PC_Factor*)dir)->fact,X,Y));
1937b6e2003SPierre Jolivet   }
1947b6e2003SPierre Jolivet   PetscFunctionReturn(0);
1957b6e2003SPierre Jolivet }
1967b6e2003SPierre Jolivet 
1973d72fe1eSJed Brown static PetscErrorCode PCApplySymmetricLeft_Cholesky(PC pc,Vec x,Vec y)
1983d72fe1eSJed Brown {
1993d72fe1eSJed Brown   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;
2003d72fe1eSJed Brown 
2013d72fe1eSJed Brown   PetscFunctionBegin;
2023d72fe1eSJed Brown   if (dir->hdr.inplace) {
2039566063dSJacob Faibussowitsch     PetscCall(MatForwardSolve(pc->pmat,x,y));
2043d72fe1eSJed Brown   } else {
2059566063dSJacob Faibussowitsch     PetscCall(MatForwardSolve(((PC_Factor*)dir)->fact,x,y));
2063d72fe1eSJed Brown   }
2073d72fe1eSJed Brown   PetscFunctionReturn(0);
2083d72fe1eSJed Brown }
2093d72fe1eSJed Brown 
2103d72fe1eSJed Brown static PetscErrorCode PCApplySymmetricRight_Cholesky(PC pc,Vec x,Vec y)
2113d72fe1eSJed Brown {
2123d72fe1eSJed Brown   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;
2133d72fe1eSJed Brown 
2143d72fe1eSJed Brown   PetscFunctionBegin;
2153d72fe1eSJed Brown   if (dir->hdr.inplace) {
2169566063dSJacob Faibussowitsch     PetscCall(MatBackwardSolve(pc->pmat,x,y));
2173d72fe1eSJed Brown   } else {
2189566063dSJacob Faibussowitsch     PetscCall(MatBackwardSolve(((PC_Factor*)dir)->fact,x,y));
2193d72fe1eSJed Brown   }
2203d72fe1eSJed Brown   PetscFunctionReturn(0);
2213d72fe1eSJed Brown }
2223d72fe1eSJed Brown 
2239b54502bSHong Zhang static PetscErrorCode PCApplyTranspose_Cholesky(PC pc,Vec x,Vec y)
2249b54502bSHong Zhang {
2259b54502bSHong Zhang   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;
2269b54502bSHong Zhang 
2279b54502bSHong Zhang   PetscFunctionBegin;
2283d1c1ea0SBarry Smith   if (dir->hdr.inplace) {
2299566063dSJacob Faibussowitsch     PetscCall(MatSolveTranspose(pc->pmat,x,y));
2302fa5cd67SKarl Rupp   } else {
2319566063dSJacob Faibussowitsch     PetscCall(MatSolveTranspose(((PC_Factor*)dir)->fact,x,y));
2322fa5cd67SKarl Rupp   }
2339b54502bSHong Zhang   PetscFunctionReturn(0);
2349b54502bSHong Zhang }
2359b54502bSHong Zhang 
2369b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/
2379b54502bSHong Zhang 
2389b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/
2399b54502bSHong Zhang 
2409b54502bSHong Zhang /*@
2412401956bSBarry Smith    PCFactorSetReuseOrdering - When similar matrices are factored, this
2429b54502bSHong Zhang    causes the ordering computed in the first factor to be used for all
2439b54502bSHong Zhang    following factors.
2449b54502bSHong Zhang 
245ad4df100SBarry Smith    Logically Collective on PC
2469b54502bSHong Zhang 
2479b54502bSHong Zhang    Input Parameters:
2489b54502bSHong Zhang +  pc - the preconditioner context
2499b54502bSHong Zhang -  flag - PETSC_TRUE to reuse else PETSC_FALSE
2509b54502bSHong Zhang 
2519b54502bSHong Zhang    Options Database Key:
2522401956bSBarry Smith .  -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering()
2539b54502bSHong Zhang 
2549b54502bSHong Zhang    Level: intermediate
2559b54502bSHong Zhang 
256db781477SPatrick Sanan .seealso: `PCFactorSetReuseFill()`
2579b54502bSHong Zhang @*/
2587087cfbeSBarry Smith PetscErrorCode  PCFactorSetReuseOrdering(PC pc,PetscBool flag)
2599b54502bSHong Zhang {
2609b54502bSHong Zhang   PetscFunctionBegin;
2610700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
262acfcf0e5SJed Brown   PetscValidLogicalCollectiveBool(pc,flag,2);
263cac4c232SBarry Smith   PetscTryMethod(pc,"PCFactorSetReuseOrdering_C",(PC,PetscBool),(pc,flag));
2649b54502bSHong Zhang   PetscFunctionReturn(0);
2659b54502bSHong Zhang }
2669b54502bSHong Zhang 
2679b54502bSHong Zhang /*MC
26896fc60bcSBarry Smith    PCCHOLESKY - Uses a direct solver, based on Cholesky factorization, as a preconditioner
2699b54502bSHong Zhang 
2709b54502bSHong Zhang    Options Database Keys:
2712401956bSBarry Smith +  -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering()
2723ca39a21SBarry Smith .  -pc_factor_mat_solver_type - Actives PCFactorSetMatSolverType() to choose the direct solver, like superlu
2732401956bSBarry Smith .  -pc_factor_reuse_fill - Activates PCFactorSetReuseFill()
27455ba2a51SBarry Smith .  -pc_factor_fill <fill> - Sets fill amount
2752401956bSBarry Smith .  -pc_factor_in_place - Activates in-place factorization
276145b9266SHong Zhang -  -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine
2779b54502bSHong Zhang 
27895452b02SPatrick Sanan    Notes:
27995452b02SPatrick Sanan     Not all options work for all matrix formats
2809b54502bSHong Zhang 
2819b54502bSHong Zhang    Level: beginner
2829b54502bSHong Zhang 
28395452b02SPatrick Sanan    Notes:
28495452b02SPatrick Sanan     Usually this will compute an "exact" solution in one iteration and does
2859b54502bSHong Zhang           not need a Krylov method (i.e. you can use -ksp_type preonly, or
2869b54502bSHong Zhang           KSPSetType(ksp,KSPPREONLY) for the Krylov method
2879b54502bSHong Zhang 
288db781477SPatrick Sanan .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`,
289db781477SPatrick Sanan           `PCILU`, `PCLU`, `PCICC`, `PCFactorSetReuseOrdering()`, `PCFactorSetReuseFill()`, `PCFactorGetMatrix()`,
290db781477SPatrick Sanan           `PCFactorSetFill()`, `PCFactorSetShiftNonzero()`, `PCFactorSetShiftType()`, `PCFactorSetShiftAmount()`
291db781477SPatrick Sanan           `PCFactorSetUseInPlace()`, `PCFactorGetUseInPlace()`, `PCFactorSetMatOrderingType()`
2929b54502bSHong Zhang 
2939b54502bSHong Zhang M*/
2949b54502bSHong Zhang 
2958cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_Cholesky(PC pc)
2969b54502bSHong Zhang {
2979b54502bSHong Zhang   PC_Cholesky    *dir;
2989b54502bSHong Zhang 
2999b54502bSHong Zhang   PetscFunctionBegin;
3009566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(pc,&dir));
3013d1c1ea0SBarry Smith   pc->data = (void*)dir;
3029566063dSJacob Faibussowitsch   PetscCall(PCFactorInitialize(pc,MAT_FACTOR_CHOLESKY));
3032fa5cd67SKarl Rupp 
304075768bcSBarry Smith   ((PC_Factor*)dir)->info.fill  = 5.0;
3052fa5cd67SKarl Rupp 
3069b54502bSHong Zhang   pc->ops->destroy             = PCDestroy_Cholesky;
30769d2c0f9SBarry Smith   pc->ops->reset               = PCReset_Cholesky;
3089b54502bSHong Zhang   pc->ops->apply               = PCApply_Cholesky;
3097b6e2003SPierre Jolivet   pc->ops->matapply            = PCMatApply_Cholesky;
3103d72fe1eSJed Brown   pc->ops->applysymmetricleft  = PCApplySymmetricLeft_Cholesky;
3113d72fe1eSJed Brown   pc->ops->applysymmetricright = PCApplySymmetricRight_Cholesky;
3129b54502bSHong Zhang   pc->ops->applytranspose      = PCApplyTranspose_Cholesky;
3139b54502bSHong Zhang   pc->ops->setup               = PCSetUp_Cholesky;
3149b54502bSHong Zhang   pc->ops->setfromoptions      = PCSetFromOptions_Cholesky;
31592e08861SBarry Smith   pc->ops->view                = PCView_Factor;
3160a545947SLisandro Dalcin   pc->ops->applyrichardson     = NULL;
3179b54502bSHong Zhang   PetscFunctionReturn(0);
3189b54502bSHong Zhang }
319