xref: /petsc/src/ksp/pc/impls/factor/cholesky/cholesky.c (revision 3d72fe1e4a0c20b4668a5ec3aff87a83aa5ba620)
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   PetscErrorCode ierr;
179b54502bSHong Zhang 
189b54502bSHong Zhang   PetscFunctionBegin;
19e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Cholesky options");CHKERRQ(ierr);
20e55864a3SBarry Smith   ierr = PCSetFromOptions_Factor(PetscOptionsObject,pc);CHKERRQ(ierr);
219b54502bSHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
229b54502bSHong Zhang   PetscFunctionReturn(0);
239b54502bSHong Zhang }
249b54502bSHong Zhang 
259b54502bSHong Zhang static PetscErrorCode PCView_Cholesky(PC pc,PetscViewer viewer)
269b54502bSHong Zhang {
279b54502bSHong Zhang   PetscErrorCode ierr;
28ace3abfcSBarry Smith   PetscBool      iascii;
299b54502bSHong Zhang 
309b54502bSHong Zhang   PetscFunctionBegin;
31251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
32914a5d51SHong Zhang   ierr = PCView_Factor(pc,viewer);CHKERRQ(ierr);
339b54502bSHong Zhang   PetscFunctionReturn(0);
349b54502bSHong Zhang }
359b54502bSHong Zhang 
369b54502bSHong Zhang static PetscErrorCode PCSetUp_Cholesky(PC pc)
379b54502bSHong Zhang {
389b54502bSHong Zhang   PetscErrorCode         ierr;
39ace3abfcSBarry Smith   PetscBool              flg;
409b54502bSHong Zhang   PC_Cholesky            *dir = (PC_Cholesky*)pc->data;
41ea799195SBarry Smith   MatSolverType          stype;
4200e125f8SBarry Smith   MatFactorError         err;
439b54502bSHong Zhang 
449b54502bSHong Zhang   PetscFunctionBegin;
45c6e4fdc6SHong Zhang   pc->failedreason = PC_NOERROR;
463d1c1ea0SBarry Smith   if (dir->hdr.reusefill && pc->setupcalled) ((PC_Factor*)dir)->info.fill = dir->hdr.actualfill;
479b54502bSHong Zhang 
4884d44b13SHong Zhang   ierr = MatSetErrorIfFailure(pc->pmat,pc->erroriffailure);CHKERRQ(ierr);
493d1c1ea0SBarry Smith   if (dir->hdr.inplace) {
509b54502bSHong Zhang     if (dir->row && dir->col && (dir->row != dir->col)) {
51fcfd50ebSBarry Smith       ierr = ISDestroy(&dir->row);CHKERRQ(ierr);
529b54502bSHong Zhang     }
53fcfd50ebSBarry Smith     ierr = ISDestroy(&dir->col);CHKERRQ(ierr);
54075768bcSBarry Smith     ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr);
559b54502bSHong Zhang     if (dir->col && (dir->row != dir->col)) {  /* only use row ordering for SBAIJ */
56fcfd50ebSBarry Smith       ierr = ISDestroy(&dir->col);CHKERRQ(ierr);
579b54502bSHong Zhang     }
583bb1ff40SBarry Smith     if (dir->row) {ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);CHKERRQ(ierr);}
59075768bcSBarry Smith     ierr = MatCholeskyFactor(pc->pmat,dir->row,&((PC_Factor*)dir)->info);CHKERRQ(ierr);
6000e125f8SBarry Smith     ierr = MatFactorGetError(pc->pmat,&err);CHKERRQ(ierr);
6100e125f8SBarry Smith     if (err) { /* Factor() fails */
6200e125f8SBarry Smith       pc->failedreason = (PCFailedReason)err;
636baea169SHong Zhang       PetscFunctionReturn(0);
646baea169SHong Zhang     }
652fa5cd67SKarl Rupp 
66075768bcSBarry Smith     ((PC_Factor*)dir)->fact = pc->pmat;
679b54502bSHong Zhang   } else {
689b54502bSHong Zhang     MatInfo info;
6900e125f8SBarry Smith 
709b54502bSHong Zhang     if (!pc->setupcalled) {
71075768bcSBarry Smith       ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr);
729bfd6278SHong Zhang       /* check if dir->row == dir->col */
739bfd6278SHong Zhang       ierr = ISEqual(dir->row,dir->col,&flg);CHKERRQ(ierr);
74e32f2f54SBarry Smith       if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"row and column permutations must equal");
75fcfd50ebSBarry Smith       ierr = ISDestroy(&dir->col);CHKERRQ(ierr); /* only pass one ordering into CholeskyFactor */
769bfd6278SHong Zhang 
7790d69ab7SBarry Smith       flg  = PETSC_FALSE;
78c5929fdfSBarry Smith       ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&flg,NULL);CHKERRQ(ierr);
799b54502bSHong Zhang       if (flg) {
809b54502bSHong Zhang         PetscReal tol = 1.e-10;
81c5929fdfSBarry Smith         ierr = PetscOptionsGetReal(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&tol,NULL);CHKERRQ(ierr);
829b54502bSHong Zhang         ierr = MatReorderForNonzeroDiagonal(pc->pmat,tol,dir->row,dir->row);CHKERRQ(ierr);
839b54502bSHong Zhang       }
843bb1ff40SBarry Smith       if (dir->row) {ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);CHKERRQ(ierr);}
85d09a07f4SBarry Smith       if (!((PC_Factor*)dir)->fact) {
86075768bcSBarry Smith         ierr = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_CHOLESKY,&((PC_Factor*)dir)->fact);CHKERRQ(ierr);
87a1f19f5aSHong Zhang       }
88075768bcSBarry Smith       ierr                = MatCholeskyFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,&((PC_Factor*)dir)->info);CHKERRQ(ierr);
89075768bcSBarry Smith       ierr                = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
903d1c1ea0SBarry Smith       dir->hdr.actualfill = info.fill_ratio_needed;
913bb1ff40SBarry Smith       ierr                = PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact);CHKERRQ(ierr);
929b54502bSHong Zhang     } else if (pc->flag != SAME_NONZERO_PATTERN) {
933d1c1ea0SBarry Smith       if (!dir->hdr.reuseordering) {
94fcfd50ebSBarry Smith         ierr = ISDestroy(&dir->row);CHKERRQ(ierr);
95075768bcSBarry Smith         ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr);
969d61c4f5SHong Zhang         ierr = ISDestroy(&dir->col);CHKERRQ(ierr); /* only use dir->row ordering in CholeskyFactor */
979d61c4f5SHong Zhang 
9890d69ab7SBarry Smith         flg  = PETSC_FALSE;
99c5929fdfSBarry Smith         ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&flg,NULL);CHKERRQ(ierr);
1009b54502bSHong Zhang         if (flg) {
1019b54502bSHong Zhang           PetscReal tol = 1.e-10;
102c5929fdfSBarry Smith           ierr = PetscOptionsGetReal(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&tol,NULL);CHKERRQ(ierr);
1039b54502bSHong Zhang           ierr = MatReorderForNonzeroDiagonal(pc->pmat,tol,dir->row,dir->row);CHKERRQ(ierr);
1049b54502bSHong Zhang         }
1053bb1ff40SBarry Smith         if (dir->row) {ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);CHKERRQ(ierr);}
1069b54502bSHong Zhang       }
1076bf464f9SBarry Smith       ierr                = MatDestroy(&((PC_Factor*)dir)->fact);CHKERRQ(ierr);
108075768bcSBarry Smith       ierr                = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_CHOLESKY,&((PC_Factor*)dir)->fact);CHKERRQ(ierr);
109075768bcSBarry Smith       ierr                = MatCholeskyFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,&((PC_Factor*)dir)->info);CHKERRQ(ierr);
110075768bcSBarry Smith       ierr                = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
1113d1c1ea0SBarry Smith       dir->hdr.actualfill = info.fill_ratio_needed;
1123bb1ff40SBarry Smith       ierr                = PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact);CHKERRQ(ierr);
11304545d6dSBarry Smith     } else {
114b8b68cfdSBarry Smith       ierr = MatFactorGetError(((PC_Factor*)dir)->fact,&err);CHKERRQ(ierr);
115160a8794SBarry Smith       if (err == MAT_FACTOR_NUMERIC_ZEROPIVOT) {
116b8b68cfdSBarry Smith         ierr = MatFactorClearError(((PC_Factor*)dir)->fact);CHKERRQ(ierr);
117b8b68cfdSBarry Smith         pc->failedreason = PC_NOERROR;
11804545d6dSBarry Smith       }
1199b54502bSHong Zhang     }
12000e125f8SBarry Smith     ierr = MatFactorGetError(((PC_Factor*)dir)->fact,&err);CHKERRQ(ierr);
12100e125f8SBarry Smith     if (err) { /* FactorSymbolic() fails */
12200e125f8SBarry Smith       pc->failedreason = (PCFailedReason)err;
1236baea169SHong Zhang       PetscFunctionReturn(0);
1246baea169SHong Zhang     }
1256baea169SHong Zhang 
126075768bcSBarry Smith     ierr = MatCholeskyFactorNumeric(((PC_Factor*)dir)->fact,pc->pmat,&((PC_Factor*)dir)->info);CHKERRQ(ierr);
12700e125f8SBarry Smith     ierr = MatFactorGetError(((PC_Factor*)dir)->fact,&err);CHKERRQ(ierr);
12800e125f8SBarry Smith     if (err) { /* FactorNumeric() fails */
12900e125f8SBarry Smith       pc->failedreason = (PCFailedReason)err;
1306baea169SHong Zhang     }
1319b54502bSHong Zhang   }
13200c67f3bSHong Zhang 
1333ca39a21SBarry Smith   ierr = PCFactorGetMatSolverType(pc,&stype);CHKERRQ(ierr);
13400c67f3bSHong Zhang   if (!stype) {
135ea799195SBarry Smith     MatSolverType solverpackage;
1363ca39a21SBarry Smith     ierr = MatFactorGetSolverType(((PC_Factor*)dir)->fact,&solverpackage);CHKERRQ(ierr);
1373ca39a21SBarry Smith     ierr = PCFactorSetMatSolverType(pc,solverpackage);CHKERRQ(ierr);
13800c67f3bSHong Zhang   }
1399b54502bSHong Zhang   PetscFunctionReturn(0);
1409b54502bSHong Zhang }
1419b54502bSHong Zhang 
14269d2c0f9SBarry Smith static PetscErrorCode PCReset_Cholesky(PC pc)
1439b54502bSHong Zhang {
1449b54502bSHong Zhang   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;
1459b54502bSHong Zhang   PetscErrorCode ierr;
1469b54502bSHong Zhang 
1479b54502bSHong Zhang   PetscFunctionBegin;
1483d1c1ea0SBarry Smith   if (!dir->hdr.inplace && ((PC_Factor*)dir)->fact) {ierr = MatDestroy(&((PC_Factor*)dir)->fact);CHKERRQ(ierr);}
149fcfd50ebSBarry Smith   ierr = ISDestroy(&dir->row);CHKERRQ(ierr);
150fcfd50ebSBarry Smith   ierr = ISDestroy(&dir->col);CHKERRQ(ierr);
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   PetscErrorCode ierr;
15869d2c0f9SBarry Smith 
15969d2c0f9SBarry Smith   PetscFunctionBegin;
16069d2c0f9SBarry Smith   ierr = PCReset_Cholesky(pc);CHKERRQ(ierr);
161503cfb0cSBarry Smith   ierr = PetscFree(((PC_Factor*)dir)->ordering);CHKERRQ(ierr);
162503cfb0cSBarry Smith   ierr = PetscFree(((PC_Factor*)dir)->solvertype);CHKERRQ(ierr);
163c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
1649b54502bSHong Zhang   PetscFunctionReturn(0);
1659b54502bSHong Zhang }
1669b54502bSHong Zhang 
1679b54502bSHong Zhang static PetscErrorCode PCApply_Cholesky(PC pc,Vec x,Vec y)
1689b54502bSHong Zhang {
1699b54502bSHong Zhang   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;
1709b54502bSHong Zhang   PetscErrorCode ierr;
1719b54502bSHong Zhang 
1729b54502bSHong Zhang   PetscFunctionBegin;
1733d1c1ea0SBarry Smith   if (dir->hdr.inplace) {
1742fa5cd67SKarl Rupp     ierr = MatSolve(pc->pmat,x,y);CHKERRQ(ierr);
1752fa5cd67SKarl Rupp   } else {
1762fa5cd67SKarl Rupp     ierr = MatSolve(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr);
1772fa5cd67SKarl Rupp   }
1789b54502bSHong Zhang   PetscFunctionReturn(0);
1799b54502bSHong Zhang }
1809b54502bSHong Zhang 
181*3d72fe1eSJed Brown static PetscErrorCode PCApplySymmetricLeft_Cholesky(PC pc,Vec x,Vec y)
182*3d72fe1eSJed Brown {
183*3d72fe1eSJed Brown   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;
184*3d72fe1eSJed Brown   PetscErrorCode ierr;
185*3d72fe1eSJed Brown 
186*3d72fe1eSJed Brown   PetscFunctionBegin;
187*3d72fe1eSJed Brown   if (dir->hdr.inplace) {
188*3d72fe1eSJed Brown     ierr = MatForwardSolve(pc->pmat,x,y);CHKERRQ(ierr);
189*3d72fe1eSJed Brown   } else {
190*3d72fe1eSJed Brown     ierr = MatForwardSolve(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr);
191*3d72fe1eSJed Brown   }
192*3d72fe1eSJed Brown   PetscFunctionReturn(0);
193*3d72fe1eSJed Brown }
194*3d72fe1eSJed Brown 
195*3d72fe1eSJed Brown static PetscErrorCode PCApplySymmetricRight_Cholesky(PC pc,Vec x,Vec y)
196*3d72fe1eSJed Brown {
197*3d72fe1eSJed Brown   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;
198*3d72fe1eSJed Brown   PetscErrorCode ierr;
199*3d72fe1eSJed Brown 
200*3d72fe1eSJed Brown   PetscFunctionBegin;
201*3d72fe1eSJed Brown   if (dir->hdr.inplace) {
202*3d72fe1eSJed Brown     ierr = MatBackwardSolve(pc->pmat,x,y);CHKERRQ(ierr);
203*3d72fe1eSJed Brown   } else {
204*3d72fe1eSJed Brown     ierr = MatBackwardSolve(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr);
205*3d72fe1eSJed Brown   }
206*3d72fe1eSJed Brown   PetscFunctionReturn(0);
207*3d72fe1eSJed Brown }
208*3d72fe1eSJed Brown 
2099b54502bSHong Zhang static PetscErrorCode PCApplyTranspose_Cholesky(PC pc,Vec x,Vec y)
2109b54502bSHong Zhang {
2119b54502bSHong Zhang   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;
2129b54502bSHong Zhang   PetscErrorCode ierr;
2139b54502bSHong Zhang 
2149b54502bSHong Zhang   PetscFunctionBegin;
2153d1c1ea0SBarry Smith   if (dir->hdr.inplace) {
2162fa5cd67SKarl Rupp     ierr = MatSolveTranspose(pc->pmat,x,y);CHKERRQ(ierr);
2172fa5cd67SKarl Rupp   } else {
2182fa5cd67SKarl Rupp     ierr = MatSolveTranspose(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr);
2192fa5cd67SKarl Rupp   }
2209b54502bSHong Zhang   PetscFunctionReturn(0);
2219b54502bSHong Zhang }
2229b54502bSHong Zhang 
2239b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/
2249b54502bSHong Zhang 
2259b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/
2269b54502bSHong Zhang 
2279b54502bSHong Zhang /*@
2282401956bSBarry Smith    PCFactorSetReuseOrdering - When similar matrices are factored, this
2299b54502bSHong Zhang    causes the ordering computed in the first factor to be used for all
2309b54502bSHong Zhang    following factors.
2319b54502bSHong Zhang 
232ad4df100SBarry Smith    Logically Collective on PC
2339b54502bSHong Zhang 
2349b54502bSHong Zhang    Input Parameters:
2359b54502bSHong Zhang +  pc - the preconditioner context
2369b54502bSHong Zhang -  flag - PETSC_TRUE to reuse else PETSC_FALSE
2379b54502bSHong Zhang 
2389b54502bSHong Zhang    Options Database Key:
2392401956bSBarry Smith .  -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering()
2409b54502bSHong Zhang 
2419b54502bSHong Zhang    Level: intermediate
2429b54502bSHong Zhang 
2439b54502bSHong Zhang .keywords: PC, levels, reordering, factorization, incomplete, LU
2449b54502bSHong Zhang 
2452401956bSBarry Smith .seealso: PCFactorSetReuseFill()
2469b54502bSHong Zhang @*/
2477087cfbeSBarry Smith PetscErrorCode  PCFactorSetReuseOrdering(PC pc,PetscBool flag)
2489b54502bSHong Zhang {
2494ac538c5SBarry Smith   PetscErrorCode ierr;
2509b54502bSHong Zhang 
2519b54502bSHong Zhang   PetscFunctionBegin;
2520700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
253acfcf0e5SJed Brown   PetscValidLogicalCollectiveBool(pc,flag,2);
2544ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFactorSetReuseOrdering_C",(PC,PetscBool),(pc,flag));CHKERRQ(ierr);
2559b54502bSHong Zhang   PetscFunctionReturn(0);
2569b54502bSHong Zhang }
2579b54502bSHong Zhang 
2589b54502bSHong Zhang /*MC
25996fc60bcSBarry Smith    PCCHOLESKY - Uses a direct solver, based on Cholesky factorization, as a preconditioner
2609b54502bSHong Zhang 
2619b54502bSHong Zhang    Options Database Keys:
2622401956bSBarry Smith +  -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering()
2633ca39a21SBarry Smith .  -pc_factor_mat_solver_type - Actives PCFactorSetMatSolverType() to choose the direct solver, like superlu
2642401956bSBarry Smith .  -pc_factor_reuse_fill - Activates PCFactorSetReuseFill()
26555ba2a51SBarry Smith .  -pc_factor_fill <fill> - Sets fill amount
2662401956bSBarry Smith .  -pc_factor_in_place - Activates in-place factorization
267145b9266SHong Zhang -  -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine
2689b54502bSHong Zhang 
26995452b02SPatrick Sanan    Notes:
27095452b02SPatrick Sanan     Not all options work for all matrix formats
2719b54502bSHong Zhang 
2729b54502bSHong Zhang    Level: beginner
2739b54502bSHong Zhang 
2749b54502bSHong Zhang    Concepts: Cholesky factorization, direct solver
2759b54502bSHong Zhang 
27695452b02SPatrick Sanan    Notes:
27795452b02SPatrick Sanan     Usually this will compute an "exact" solution in one iteration and does
2789b54502bSHong Zhang           not need a Krylov method (i.e. you can use -ksp_type preonly, or
2799b54502bSHong Zhang           KSPSetType(ksp,KSPPREONLY) for the Krylov method
2809b54502bSHong Zhang 
2819b54502bSHong Zhang .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
282a4fd02acSBarry Smith            PCILU, PCLU, PCICC, PCFactorSetReuseOrdering(), PCFactorSetReuseFill(), PCFactorGetMatrix(),
283145b9266SHong Zhang            PCFactorSetFill(), PCFactorSetShiftNonzero(), PCFactorSetShiftType(), PCFactorSetShiftAmount()
2848e37d05fSBarry Smith            PCFactorSetUseInPlace(), PCFactorGetUseInPlace(), PCFactorSetMatOrderingType()
2859b54502bSHong Zhang 
2869b54502bSHong Zhang M*/
2879b54502bSHong Zhang 
2888cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_Cholesky(PC pc)
2899b54502bSHong Zhang {
2909b54502bSHong Zhang   PetscErrorCode ierr;
2919b54502bSHong Zhang   PC_Cholesky    *dir;
2929b54502bSHong Zhang 
2939b54502bSHong Zhang   PetscFunctionBegin;
294b00a9115SJed Brown   ierr     = PetscNewLog(pc,&dir);CHKERRQ(ierr);
2953d1c1ea0SBarry Smith   pc->data = (void*)dir;
2963d1c1ea0SBarry Smith   ierr     = PCFactorInitialize(pc);CHKERRQ(ierr);
2972fa5cd67SKarl Rupp 
298879e8a4dSBarry Smith   ((PC_Factor*)dir)->factortype         = MAT_FACTOR_CHOLESKY;
299075768bcSBarry Smith   ((PC_Factor*)dir)->info.fill          = 5.0;
3002fa5cd67SKarl Rupp 
3019b54502bSHong Zhang   dir->col = 0;
3029b54502bSHong Zhang   dir->row = 0;
3032fa5cd67SKarl Rupp 
30419fd82e9SBarry Smith   ierr = PetscStrallocpy(MATORDERINGNATURAL,(char**)&((PC_Factor*)dir)->ordering);CHKERRQ(ierr);
3059b54502bSHong Zhang 
3069b54502bSHong Zhang   pc->ops->destroy             = PCDestroy_Cholesky;
30769d2c0f9SBarry Smith   pc->ops->reset               = PCReset_Cholesky;
3089b54502bSHong Zhang   pc->ops->apply               = PCApply_Cholesky;
309*3d72fe1eSJed Brown   pc->ops->applysymmetricleft  = PCApplySymmetricLeft_Cholesky;
310*3d72fe1eSJed Brown   pc->ops->applysymmetricright = PCApplySymmetricRight_Cholesky;
3119b54502bSHong Zhang   pc->ops->applytranspose      = PCApplyTranspose_Cholesky;
3129b54502bSHong Zhang   pc->ops->setup               = PCSetUp_Cholesky;
3139b54502bSHong Zhang   pc->ops->setfromoptions      = PCSetFromOptions_Cholesky;
3149b54502bSHong Zhang   pc->ops->view                = PCView_Cholesky;
3159b54502bSHong Zhang   pc->ops->applyrichardson     = 0;
3169b54502bSHong Zhang   PetscFunctionReturn(0);
3179b54502bSHong Zhang }
318