xref: /petsc/src/ksp/pc/impls/factor/cholesky/cholesky.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
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;
17*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHead(PetscOptionsObject,"Cholesky options"));
18*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PCSetFromOptions_Factor(PetscOptionsObject,pc));
19*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsTail());
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;
299b54502bSHong Zhang 
309b54502bSHong Zhang   PetscFunctionBegin;
31c6e4fdc6SHong Zhang   pc->failedreason = PC_NOERROR;
323d1c1ea0SBarry Smith   if (dir->hdr.reusefill && pc->setupcalled) ((PC_Factor*)dir)->info.fill = dir->hdr.actualfill;
339b54502bSHong Zhang 
34*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetErrorIfFailure(pc->pmat,pc->erroriffailure));
353d1c1ea0SBarry Smith   if (dir->hdr.inplace) {
369b54502bSHong Zhang     if (dir->row && dir->col && (dir->row != dir->col)) {
37*5f80ce2aSJacob Faibussowitsch       CHKERRQ(ISDestroy(&dir->row));
389b54502bSHong Zhang     }
39*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISDestroy(&dir->col));
402c7c0729SBarry Smith     /* should only get reordering if the factor matrix uses it but cannot determine because MatGetFactor() not called */
41*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PCFactorSetDefaultOrdering_Factor(pc));
42*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col));
439b54502bSHong Zhang     if (dir->col && (dir->row != dir->col)) {  /* only use row ordering for SBAIJ */
44*5f80ce2aSJacob Faibussowitsch       CHKERRQ(ISDestroy(&dir->col));
459b54502bSHong Zhang     }
46*5f80ce2aSJacob Faibussowitsch     if (dir->row) CHKERRQ(PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row));
47*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCholeskyFactor(pc->pmat,dir->row,&((PC_Factor*)dir)->info));
48*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatFactorGetError(pc->pmat,&err));
4900e125f8SBarry Smith     if (err) { /* Factor() fails */
5000e125f8SBarry Smith       pc->failedreason = (PCFailedReason)err;
516baea169SHong Zhang       PetscFunctionReturn(0);
526baea169SHong Zhang     }
532fa5cd67SKarl Rupp 
54075768bcSBarry Smith     ((PC_Factor*)dir)->fact = pc->pmat;
559b54502bSHong Zhang   } else {
569b54502bSHong Zhang     MatInfo info;
5700e125f8SBarry Smith 
589b54502bSHong Zhang     if (!pc->setupcalled) {
59f73b0415SBarry Smith       PetscBool canuseordering;
602c7c0729SBarry Smith       if (!((PC_Factor*)dir)->fact) {
61*5f80ce2aSJacob Faibussowitsch         CHKERRQ(MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_CHOLESKY,&((PC_Factor*)dir)->fact));
62*5f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact));
632c7c0729SBarry Smith       }
64*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatFactorGetCanUseOrdering(((PC_Factor*)dir)->fact,&canuseordering));
65f73b0415SBarry Smith       if (canuseordering) {
66*5f80ce2aSJacob Faibussowitsch         CHKERRQ(PCFactorSetDefaultOrdering_Factor(pc));
67*5f80ce2aSJacob Faibussowitsch         CHKERRQ(MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col));
689bfd6278SHong Zhang         /* check if dir->row == dir->col */
694ac6704cSBarry Smith         if (dir->row) {
70*5f80ce2aSJacob Faibussowitsch           CHKERRQ(ISEqual(dir->row,dir->col,&flg));
712c71b3e2SJacob Faibussowitsch           PetscCheckFalse(!flg,PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"row and column permutations must be equal");
724ac6704cSBarry Smith         }
73*5f80ce2aSJacob Faibussowitsch         CHKERRQ(ISDestroy(&dir->col)); /* only pass one ordering into CholeskyFactor */
749bfd6278SHong Zhang 
7590d69ab7SBarry Smith         flg  = PETSC_FALSE;
76*5f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&flg,NULL));
779b54502bSHong Zhang         if (flg) {
789b54502bSHong Zhang           PetscReal tol = 1.e-10;
79*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscOptionsGetReal(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&tol,NULL));
80*5f80ce2aSJacob Faibussowitsch           CHKERRQ(MatReorderForNonzeroDiagonal(pc->pmat,tol,dir->row,dir->row));
819b54502bSHong Zhang         }
82*5f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row));
83a1f19f5aSHong Zhang       }
84*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatCholeskyFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,&((PC_Factor*)dir)->info));
85*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info));
863d1c1ea0SBarry Smith       dir->hdr.actualfill = info.fill_ratio_needed;
879b54502bSHong Zhang     } else if (pc->flag != SAME_NONZERO_PATTERN) {
883d1c1ea0SBarry Smith       if (!dir->hdr.reuseordering) {
89f73b0415SBarry Smith         PetscBool canuseordering;
90*5f80ce2aSJacob Faibussowitsch         CHKERRQ(MatDestroy(&((PC_Factor*)dir)->fact));
91*5f80ce2aSJacob Faibussowitsch         CHKERRQ(MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_CHOLESKY,&((PC_Factor*)dir)->fact));
92*5f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact));
93*5f80ce2aSJacob Faibussowitsch         CHKERRQ(MatFactorGetCanUseOrdering(((PC_Factor*)dir)->fact,&canuseordering));
94f73b0415SBarry Smith         if (canuseordering) {
95*5f80ce2aSJacob Faibussowitsch           CHKERRQ(ISDestroy(&dir->row));
96*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PCFactorSetDefaultOrdering_Factor(pc));
97*5f80ce2aSJacob Faibussowitsch           CHKERRQ(MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col));
98*5f80ce2aSJacob Faibussowitsch           CHKERRQ(ISDestroy(&dir->col)); /* only use dir->row ordering in CholeskyFactor */
999d61c4f5SHong Zhang 
10090d69ab7SBarry Smith           flg  = PETSC_FALSE;
101*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&flg,NULL));
1029b54502bSHong Zhang           if (flg) {
1039b54502bSHong Zhang             PetscReal tol = 1.e-10;
104*5f80ce2aSJacob Faibussowitsch             CHKERRQ(PetscOptionsGetReal(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&tol,NULL));
105*5f80ce2aSJacob Faibussowitsch             CHKERRQ(MatReorderForNonzeroDiagonal(pc->pmat,tol,dir->row,dir->row));
1069b54502bSHong Zhang           }
107*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row));
1089b54502bSHong Zhang         }
1092c7c0729SBarry Smith       }
110*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatCholeskyFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,&((PC_Factor*)dir)->info));
111*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info));
1123d1c1ea0SBarry Smith       dir->hdr.actualfill = info.fill_ratio_needed;
113*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact));
11404545d6dSBarry Smith     } else {
115*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatFactorGetError(((PC_Factor*)dir)->fact,&err));
116160a8794SBarry Smith       if (err == MAT_FACTOR_NUMERIC_ZEROPIVOT) {
117*5f80ce2aSJacob Faibussowitsch         CHKERRQ(MatFactorClearError(((PC_Factor*)dir)->fact));
118b8b68cfdSBarry Smith         pc->failedreason = PC_NOERROR;
11904545d6dSBarry Smith       }
1209b54502bSHong Zhang     }
121*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatFactorGetError(((PC_Factor*)dir)->fact,&err));
12200e125f8SBarry Smith     if (err) { /* FactorSymbolic() fails */
12300e125f8SBarry Smith       pc->failedreason = (PCFailedReason)err;
1246baea169SHong Zhang       PetscFunctionReturn(0);
1256baea169SHong Zhang     }
1266baea169SHong Zhang 
127*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCholeskyFactorNumeric(((PC_Factor*)dir)->fact,pc->pmat,&((PC_Factor*)dir)->info));
128*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatFactorGetError(((PC_Factor*)dir)->fact,&err));
12900e125f8SBarry Smith     if (err) { /* FactorNumeric() fails */
13000e125f8SBarry Smith       pc->failedreason = (PCFailedReason)err;
1316baea169SHong Zhang     }
1329b54502bSHong Zhang   }
13300c67f3bSHong Zhang 
134*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PCFactorGetMatSolverType(pc,&stype));
13500c67f3bSHong Zhang   if (!stype) {
136ea799195SBarry Smith     MatSolverType solverpackage;
137*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatFactorGetSolverType(((PC_Factor*)dir)->fact,&solverpackage));
138*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PCFactorSetMatSolverType(pc,solverpackage));
13900c67f3bSHong Zhang   }
1409b54502bSHong Zhang   PetscFunctionReturn(0);
1419b54502bSHong Zhang }
1429b54502bSHong Zhang 
14369d2c0f9SBarry Smith static PetscErrorCode PCReset_Cholesky(PC pc)
1449b54502bSHong Zhang {
1459b54502bSHong Zhang   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;
1469b54502bSHong Zhang 
1479b54502bSHong Zhang   PetscFunctionBegin;
148*5f80ce2aSJacob Faibussowitsch   if (!dir->hdr.inplace && ((PC_Factor*)dir)->fact) CHKERRQ(MatDestroy(&((PC_Factor*)dir)->fact));
149*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&dir->row));
150*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&dir->col));
15169d2c0f9SBarry Smith   PetscFunctionReturn(0);
15269d2c0f9SBarry Smith }
15369d2c0f9SBarry Smith 
15469d2c0f9SBarry Smith static PetscErrorCode PCDestroy_Cholesky(PC pc)
15569d2c0f9SBarry Smith {
15669d2c0f9SBarry Smith   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;
15769d2c0f9SBarry Smith 
15869d2c0f9SBarry Smith   PetscFunctionBegin;
159*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PCReset_Cholesky(pc));
160*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(((PC_Factor*)dir)->ordering));
161*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(((PC_Factor*)dir)->solvertype));
162*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(pc->data));
1639b54502bSHong Zhang   PetscFunctionReturn(0);
1649b54502bSHong Zhang }
1659b54502bSHong Zhang 
1669b54502bSHong Zhang static PetscErrorCode PCApply_Cholesky(PC pc,Vec x,Vec y)
1679b54502bSHong Zhang {
1689b54502bSHong Zhang   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;
1699b54502bSHong Zhang 
1709b54502bSHong Zhang   PetscFunctionBegin;
1713d1c1ea0SBarry Smith   if (dir->hdr.inplace) {
172*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSolve(pc->pmat,x,y));
1732fa5cd67SKarl Rupp   } else {
174*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSolve(((PC_Factor*)dir)->fact,x,y));
1752fa5cd67SKarl Rupp   }
1769b54502bSHong Zhang   PetscFunctionReturn(0);
1779b54502bSHong Zhang }
1789b54502bSHong Zhang 
1797b6e2003SPierre Jolivet static PetscErrorCode PCMatApply_Cholesky(PC pc,Mat X,Mat Y)
1807b6e2003SPierre Jolivet {
1817b6e2003SPierre Jolivet   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;
1827b6e2003SPierre Jolivet 
1837b6e2003SPierre Jolivet   PetscFunctionBegin;
1847b6e2003SPierre Jolivet   if (dir->hdr.inplace) {
185*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMatSolve(pc->pmat,X,Y));
1867b6e2003SPierre Jolivet   } else {
187*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMatSolve(((PC_Factor*)dir)->fact,X,Y));
1887b6e2003SPierre Jolivet   }
1897b6e2003SPierre Jolivet   PetscFunctionReturn(0);
1907b6e2003SPierre Jolivet }
1917b6e2003SPierre Jolivet 
1923d72fe1eSJed Brown static PetscErrorCode PCApplySymmetricLeft_Cholesky(PC pc,Vec x,Vec y)
1933d72fe1eSJed Brown {
1943d72fe1eSJed Brown   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;
1953d72fe1eSJed Brown 
1963d72fe1eSJed Brown   PetscFunctionBegin;
1973d72fe1eSJed Brown   if (dir->hdr.inplace) {
198*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatForwardSolve(pc->pmat,x,y));
1993d72fe1eSJed Brown   } else {
200*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatForwardSolve(((PC_Factor*)dir)->fact,x,y));
2013d72fe1eSJed Brown   }
2023d72fe1eSJed Brown   PetscFunctionReturn(0);
2033d72fe1eSJed Brown }
2043d72fe1eSJed Brown 
2053d72fe1eSJed Brown static PetscErrorCode PCApplySymmetricRight_Cholesky(PC pc,Vec x,Vec y)
2063d72fe1eSJed Brown {
2073d72fe1eSJed Brown   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;
2083d72fe1eSJed Brown 
2093d72fe1eSJed Brown   PetscFunctionBegin;
2103d72fe1eSJed Brown   if (dir->hdr.inplace) {
211*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatBackwardSolve(pc->pmat,x,y));
2123d72fe1eSJed Brown   } else {
213*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatBackwardSolve(((PC_Factor*)dir)->fact,x,y));
2143d72fe1eSJed Brown   }
2153d72fe1eSJed Brown   PetscFunctionReturn(0);
2163d72fe1eSJed Brown }
2173d72fe1eSJed Brown 
2189b54502bSHong Zhang static PetscErrorCode PCApplyTranspose_Cholesky(PC pc,Vec x,Vec y)
2199b54502bSHong Zhang {
2209b54502bSHong Zhang   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;
2219b54502bSHong Zhang 
2229b54502bSHong Zhang   PetscFunctionBegin;
2233d1c1ea0SBarry Smith   if (dir->hdr.inplace) {
224*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSolveTranspose(pc->pmat,x,y));
2252fa5cd67SKarl Rupp   } else {
226*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSolveTranspose(((PC_Factor*)dir)->fact,x,y));
2272fa5cd67SKarl Rupp   }
2289b54502bSHong Zhang   PetscFunctionReturn(0);
2299b54502bSHong Zhang }
2309b54502bSHong Zhang 
2319b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/
2329b54502bSHong Zhang 
2339b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/
2349b54502bSHong Zhang 
2359b54502bSHong Zhang /*@
2362401956bSBarry Smith    PCFactorSetReuseOrdering - When similar matrices are factored, this
2379b54502bSHong Zhang    causes the ordering computed in the first factor to be used for all
2389b54502bSHong Zhang    following factors.
2399b54502bSHong Zhang 
240ad4df100SBarry Smith    Logically Collective on PC
2419b54502bSHong Zhang 
2429b54502bSHong Zhang    Input Parameters:
2439b54502bSHong Zhang +  pc - the preconditioner context
2449b54502bSHong Zhang -  flag - PETSC_TRUE to reuse else PETSC_FALSE
2459b54502bSHong Zhang 
2469b54502bSHong Zhang    Options Database Key:
2472401956bSBarry Smith .  -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering()
2489b54502bSHong Zhang 
2499b54502bSHong Zhang    Level: intermediate
2509b54502bSHong Zhang 
2512401956bSBarry Smith .seealso: PCFactorSetReuseFill()
2529b54502bSHong Zhang @*/
2537087cfbeSBarry Smith PetscErrorCode  PCFactorSetReuseOrdering(PC pc,PetscBool flag)
2549b54502bSHong Zhang {
2559b54502bSHong Zhang   PetscFunctionBegin;
2560700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
257acfcf0e5SJed Brown   PetscValidLogicalCollectiveBool(pc,flag,2);
258*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscTryMethod(pc,"PCFactorSetReuseOrdering_C",(PC,PetscBool),(pc,flag)));
2599b54502bSHong Zhang   PetscFunctionReturn(0);
2609b54502bSHong Zhang }
2619b54502bSHong Zhang 
2629b54502bSHong Zhang /*MC
26396fc60bcSBarry Smith    PCCHOLESKY - Uses a direct solver, based on Cholesky factorization, as a preconditioner
2649b54502bSHong Zhang 
2659b54502bSHong Zhang    Options Database Keys:
2662401956bSBarry Smith +  -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering()
2673ca39a21SBarry Smith .  -pc_factor_mat_solver_type - Actives PCFactorSetMatSolverType() to choose the direct solver, like superlu
2682401956bSBarry Smith .  -pc_factor_reuse_fill - Activates PCFactorSetReuseFill()
26955ba2a51SBarry Smith .  -pc_factor_fill <fill> - Sets fill amount
2702401956bSBarry Smith .  -pc_factor_in_place - Activates in-place factorization
271145b9266SHong Zhang -  -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine
2729b54502bSHong Zhang 
27395452b02SPatrick Sanan    Notes:
27495452b02SPatrick Sanan     Not all options work for all matrix formats
2759b54502bSHong Zhang 
2769b54502bSHong Zhang    Level: beginner
2779b54502bSHong Zhang 
27895452b02SPatrick Sanan    Notes:
27995452b02SPatrick Sanan     Usually this will compute an "exact" solution in one iteration and does
2809b54502bSHong Zhang           not need a Krylov method (i.e. you can use -ksp_type preonly, or
2819b54502bSHong Zhang           KSPSetType(ksp,KSPPREONLY) for the Krylov method
2829b54502bSHong Zhang 
2839b54502bSHong Zhang .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
284a4fd02acSBarry Smith            PCILU, PCLU, PCICC, PCFactorSetReuseOrdering(), PCFactorSetReuseFill(), PCFactorGetMatrix(),
285145b9266SHong Zhang            PCFactorSetFill(), PCFactorSetShiftNonzero(), PCFactorSetShiftType(), PCFactorSetShiftAmount()
2868e37d05fSBarry Smith            PCFactorSetUseInPlace(), PCFactorGetUseInPlace(), PCFactorSetMatOrderingType()
2879b54502bSHong Zhang 
2889b54502bSHong Zhang M*/
2899b54502bSHong Zhang 
2908cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_Cholesky(PC pc)
2919b54502bSHong Zhang {
2929b54502bSHong Zhang   PC_Cholesky    *dir;
2939b54502bSHong Zhang 
2949b54502bSHong Zhang   PetscFunctionBegin;
295*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscNewLog(pc,&dir));
2963d1c1ea0SBarry Smith   pc->data = (void*)dir;
297*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PCFactorInitialize(pc,MAT_FACTOR_CHOLESKY));
2982fa5cd67SKarl Rupp 
299075768bcSBarry Smith   ((PC_Factor*)dir)->info.fill  = 5.0;
3002fa5cd67SKarl Rupp 
3019b54502bSHong Zhang   pc->ops->destroy             = PCDestroy_Cholesky;
30269d2c0f9SBarry Smith   pc->ops->reset               = PCReset_Cholesky;
3039b54502bSHong Zhang   pc->ops->apply               = PCApply_Cholesky;
3047b6e2003SPierre Jolivet   pc->ops->matapply            = PCMatApply_Cholesky;
3053d72fe1eSJed Brown   pc->ops->applysymmetricleft  = PCApplySymmetricLeft_Cholesky;
3063d72fe1eSJed Brown   pc->ops->applysymmetricright = PCApplySymmetricRight_Cholesky;
3079b54502bSHong Zhang   pc->ops->applytranspose      = PCApplyTranspose_Cholesky;
3089b54502bSHong Zhang   pc->ops->setup               = PCSetUp_Cholesky;
3099b54502bSHong Zhang   pc->ops->setfromoptions      = PCSetFromOptions_Cholesky;
31092e08861SBarry Smith   pc->ops->view                = PCView_Factor;
3110a545947SLisandro Dalcin   pc->ops->applyrichardson     = NULL;
3129b54502bSHong Zhang   PetscFunctionReturn(0);
3139b54502bSHong Zhang }
314