xref: /petsc/src/ksp/pc/impls/factor/cholesky/cholesky.c (revision f73b041577f311cd4e522d727240702cd4268ffe)
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 PCSetUp_Cholesky(PC pc)
269b54502bSHong Zhang {
279b54502bSHong Zhang   PetscErrorCode         ierr;
28ace3abfcSBarry Smith   PetscBool              flg;
299b54502bSHong Zhang   PC_Cholesky            *dir = (PC_Cholesky*)pc->data;
30ea799195SBarry Smith   MatSolverType          stype;
3100e125f8SBarry Smith   MatFactorError         err;
329b54502bSHong Zhang 
339b54502bSHong Zhang   PetscFunctionBegin;
34c6e4fdc6SHong Zhang   pc->failedreason = PC_NOERROR;
353d1c1ea0SBarry Smith   if (dir->hdr.reusefill && pc->setupcalled) ((PC_Factor*)dir)->info.fill = dir->hdr.actualfill;
369b54502bSHong Zhang 
3784d44b13SHong Zhang   ierr = MatSetErrorIfFailure(pc->pmat,pc->erroriffailure);CHKERRQ(ierr);
383d1c1ea0SBarry Smith   if (dir->hdr.inplace) {
399b54502bSHong Zhang     if (dir->row && dir->col && (dir->row != dir->col)) {
40fcfd50ebSBarry Smith       ierr = ISDestroy(&dir->row);CHKERRQ(ierr);
419b54502bSHong Zhang     }
42fcfd50ebSBarry Smith     ierr = ISDestroy(&dir->col);CHKERRQ(ierr);
432c7c0729SBarry Smith     /* should only get reordering if the factor matrix uses it but cannot determine because MatGetFactor() not called */
44075768bcSBarry Smith     ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr);
459b54502bSHong Zhang     if (dir->col && (dir->row != dir->col)) {  /* only use row ordering for SBAIJ */
46fcfd50ebSBarry Smith       ierr = ISDestroy(&dir->col);CHKERRQ(ierr);
479b54502bSHong Zhang     }
483bb1ff40SBarry Smith     if (dir->row) {ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);CHKERRQ(ierr);}
49075768bcSBarry Smith     ierr = MatCholeskyFactor(pc->pmat,dir->row,&((PC_Factor*)dir)->info);CHKERRQ(ierr);
5000e125f8SBarry Smith     ierr = MatFactorGetError(pc->pmat,&err);CHKERRQ(ierr);
5100e125f8SBarry Smith     if (err) { /* Factor() fails */
5200e125f8SBarry Smith       pc->failedreason = (PCFailedReason)err;
536baea169SHong Zhang       PetscFunctionReturn(0);
546baea169SHong Zhang     }
552fa5cd67SKarl Rupp 
56075768bcSBarry Smith     ((PC_Factor*)dir)->fact = pc->pmat;
579b54502bSHong Zhang   } else {
589b54502bSHong Zhang     MatInfo info;
5900e125f8SBarry Smith 
609b54502bSHong Zhang     if (!pc->setupcalled) {
61*f73b0415SBarry Smith       PetscBool canuseordering;
622c7c0729SBarry Smith       if (!((PC_Factor*)dir)->fact) {
632c7c0729SBarry Smith         ierr = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_CHOLESKY,&((PC_Factor*)dir)->fact);CHKERRQ(ierr);
642c7c0729SBarry Smith         ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact);CHKERRQ(ierr);
652c7c0729SBarry Smith       }
66*f73b0415SBarry Smith       ierr = MatFactorGetCanUseOrdering(((PC_Factor*)dir)->fact,&canuseordering);CHKERRQ(ierr);
67*f73b0415SBarry Smith       if (canuseordering) {
68075768bcSBarry Smith         ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr);
699bfd6278SHong Zhang         /* check if dir->row == dir->col */
709bfd6278SHong Zhang         ierr = ISEqual(dir->row,dir->col,&flg);CHKERRQ(ierr);
71e32f2f54SBarry Smith         if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"row and column permutations must equal");
72fcfd50ebSBarry Smith         ierr = ISDestroy(&dir->col);CHKERRQ(ierr); /* only pass one ordering into CholeskyFactor */
739bfd6278SHong Zhang 
7490d69ab7SBarry Smith         flg  = PETSC_FALSE;
75c5929fdfSBarry Smith         ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&flg,NULL);CHKERRQ(ierr);
769b54502bSHong Zhang         if (flg) {
779b54502bSHong Zhang           PetscReal tol = 1.e-10;
78c5929fdfSBarry Smith           ierr = PetscOptionsGetReal(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&tol,NULL);CHKERRQ(ierr);
799b54502bSHong Zhang           ierr = MatReorderForNonzeroDiagonal(pc->pmat,tol,dir->row,dir->row);CHKERRQ(ierr);
809b54502bSHong Zhang         }
812c7c0729SBarry Smith         ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);CHKERRQ(ierr);
82a1f19f5aSHong Zhang       }
83075768bcSBarry Smith       ierr                = MatCholeskyFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,&((PC_Factor*)dir)->info);CHKERRQ(ierr);
84075768bcSBarry Smith       ierr                = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
853d1c1ea0SBarry Smith       dir->hdr.actualfill = info.fill_ratio_needed;
869b54502bSHong Zhang     } else if (pc->flag != SAME_NONZERO_PATTERN) {
873d1c1ea0SBarry Smith       if (!dir->hdr.reuseordering) {
88*f73b0415SBarry Smith         PetscBool canuseordering;
892c7c0729SBarry Smith         ierr = MatDestroy(&((PC_Factor*)dir)->fact);CHKERRQ(ierr);
902c7c0729SBarry Smith         ierr = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_CHOLESKY,&((PC_Factor*)dir)->fact);CHKERRQ(ierr);
912c7c0729SBarry Smith         ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact);CHKERRQ(ierr);
92*f73b0415SBarry Smith         ierr = MatFactorGetCanUseOrdering(((PC_Factor*)dir)->fact,&canuseordering);CHKERRQ(ierr);
93*f73b0415SBarry Smith         if (canuseordering) {
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           }
1052c7c0729SBarry Smith           ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);CHKERRQ(ierr);
1069b54502bSHong Zhang         }
1072c7c0729SBarry Smith       }
108075768bcSBarry Smith       ierr                = MatCholeskyFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,&((PC_Factor*)dir)->info);CHKERRQ(ierr);
109075768bcSBarry Smith       ierr                = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
1103d1c1ea0SBarry Smith       dir->hdr.actualfill = info.fill_ratio_needed;
1113bb1ff40SBarry Smith       ierr                = PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact);CHKERRQ(ierr);
11204545d6dSBarry Smith     } else {
113b8b68cfdSBarry Smith       ierr = MatFactorGetError(((PC_Factor*)dir)->fact,&err);CHKERRQ(ierr);
114160a8794SBarry Smith       if (err == MAT_FACTOR_NUMERIC_ZEROPIVOT) {
115b8b68cfdSBarry Smith         ierr = MatFactorClearError(((PC_Factor*)dir)->fact);CHKERRQ(ierr);
116b8b68cfdSBarry Smith         pc->failedreason = PC_NOERROR;
11704545d6dSBarry Smith       }
1189b54502bSHong Zhang     }
11900e125f8SBarry Smith     ierr = MatFactorGetError(((PC_Factor*)dir)->fact,&err);CHKERRQ(ierr);
12000e125f8SBarry Smith     if (err) { /* FactorSymbolic() fails */
12100e125f8SBarry Smith       pc->failedreason = (PCFailedReason)err;
1226baea169SHong Zhang       PetscFunctionReturn(0);
1236baea169SHong Zhang     }
1246baea169SHong Zhang 
125075768bcSBarry Smith     ierr = MatCholeskyFactorNumeric(((PC_Factor*)dir)->fact,pc->pmat,&((PC_Factor*)dir)->info);CHKERRQ(ierr);
12600e125f8SBarry Smith     ierr = MatFactorGetError(((PC_Factor*)dir)->fact,&err);CHKERRQ(ierr);
12700e125f8SBarry Smith     if (err) { /* FactorNumeric() fails */
12800e125f8SBarry Smith       pc->failedreason = (PCFailedReason)err;
1296baea169SHong Zhang     }
1309b54502bSHong Zhang   }
13100c67f3bSHong Zhang 
1323ca39a21SBarry Smith   ierr = PCFactorGetMatSolverType(pc,&stype);CHKERRQ(ierr);
13300c67f3bSHong Zhang   if (!stype) {
134ea799195SBarry Smith     MatSolverType solverpackage;
1353ca39a21SBarry Smith     ierr = MatFactorGetSolverType(((PC_Factor*)dir)->fact,&solverpackage);CHKERRQ(ierr);
1363ca39a21SBarry Smith     ierr = PCFactorSetMatSolverType(pc,solverpackage);CHKERRQ(ierr);
13700c67f3bSHong Zhang   }
1389b54502bSHong Zhang   PetscFunctionReturn(0);
1399b54502bSHong Zhang }
1409b54502bSHong Zhang 
14169d2c0f9SBarry Smith static PetscErrorCode PCReset_Cholesky(PC pc)
1429b54502bSHong Zhang {
1439b54502bSHong Zhang   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;
1449b54502bSHong Zhang   PetscErrorCode ierr;
1459b54502bSHong Zhang 
1469b54502bSHong Zhang   PetscFunctionBegin;
1473d1c1ea0SBarry Smith   if (!dir->hdr.inplace && ((PC_Factor*)dir)->fact) {ierr = MatDestroy(&((PC_Factor*)dir)->fact);CHKERRQ(ierr);}
148fcfd50ebSBarry Smith   ierr = ISDestroy(&dir->row);CHKERRQ(ierr);
149fcfd50ebSBarry Smith   ierr = ISDestroy(&dir->col);CHKERRQ(ierr);
15069d2c0f9SBarry Smith   PetscFunctionReturn(0);
15169d2c0f9SBarry Smith }
15269d2c0f9SBarry Smith 
15369d2c0f9SBarry Smith static PetscErrorCode PCDestroy_Cholesky(PC pc)
15469d2c0f9SBarry Smith {
15569d2c0f9SBarry Smith   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;
15669d2c0f9SBarry Smith   PetscErrorCode ierr;
15769d2c0f9SBarry Smith 
15869d2c0f9SBarry Smith   PetscFunctionBegin;
15969d2c0f9SBarry Smith   ierr = PCReset_Cholesky(pc);CHKERRQ(ierr);
160503cfb0cSBarry Smith   ierr = PetscFree(((PC_Factor*)dir)->ordering);CHKERRQ(ierr);
161503cfb0cSBarry Smith   ierr = PetscFree(((PC_Factor*)dir)->solvertype);CHKERRQ(ierr);
162c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
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   PetscErrorCode ierr;
1709b54502bSHong Zhang 
1719b54502bSHong Zhang   PetscFunctionBegin;
1723d1c1ea0SBarry Smith   if (dir->hdr.inplace) {
1732fa5cd67SKarl Rupp     ierr = MatSolve(pc->pmat,x,y);CHKERRQ(ierr);
1742fa5cd67SKarl Rupp   } else {
1752fa5cd67SKarl Rupp     ierr = MatSolve(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr);
1762fa5cd67SKarl Rupp   }
1779b54502bSHong Zhang   PetscFunctionReturn(0);
1789b54502bSHong Zhang }
1799b54502bSHong Zhang 
1807b6e2003SPierre Jolivet static PetscErrorCode PCMatApply_Cholesky(PC pc,Mat X,Mat Y)
1817b6e2003SPierre Jolivet {
1827b6e2003SPierre Jolivet   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;
1837b6e2003SPierre Jolivet   PetscErrorCode ierr;
1847b6e2003SPierre Jolivet 
1857b6e2003SPierre Jolivet   PetscFunctionBegin;
1867b6e2003SPierre Jolivet   if (dir->hdr.inplace) {
1877b6e2003SPierre Jolivet     ierr = MatMatSolve(pc->pmat,X,Y);CHKERRQ(ierr);
1887b6e2003SPierre Jolivet   } else {
1897b6e2003SPierre Jolivet     ierr = MatMatSolve(((PC_Factor*)dir)->fact,X,Y);CHKERRQ(ierr);
1907b6e2003SPierre Jolivet   }
1917b6e2003SPierre Jolivet   PetscFunctionReturn(0);
1927b6e2003SPierre Jolivet }
1937b6e2003SPierre Jolivet 
1943d72fe1eSJed Brown static PetscErrorCode PCApplySymmetricLeft_Cholesky(PC pc,Vec x,Vec y)
1953d72fe1eSJed Brown {
1963d72fe1eSJed Brown   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;
1973d72fe1eSJed Brown   PetscErrorCode ierr;
1983d72fe1eSJed Brown 
1993d72fe1eSJed Brown   PetscFunctionBegin;
2003d72fe1eSJed Brown   if (dir->hdr.inplace) {
2013d72fe1eSJed Brown     ierr = MatForwardSolve(pc->pmat,x,y);CHKERRQ(ierr);
2023d72fe1eSJed Brown   } else {
2033d72fe1eSJed Brown     ierr = MatForwardSolve(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr);
2043d72fe1eSJed Brown   }
2053d72fe1eSJed Brown   PetscFunctionReturn(0);
2063d72fe1eSJed Brown }
2073d72fe1eSJed Brown 
2083d72fe1eSJed Brown static PetscErrorCode PCApplySymmetricRight_Cholesky(PC pc,Vec x,Vec y)
2093d72fe1eSJed Brown {
2103d72fe1eSJed Brown   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;
2113d72fe1eSJed Brown   PetscErrorCode ierr;
2123d72fe1eSJed Brown 
2133d72fe1eSJed Brown   PetscFunctionBegin;
2143d72fe1eSJed Brown   if (dir->hdr.inplace) {
2153d72fe1eSJed Brown     ierr = MatBackwardSolve(pc->pmat,x,y);CHKERRQ(ierr);
2163d72fe1eSJed Brown   } else {
2173d72fe1eSJed Brown     ierr = MatBackwardSolve(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr);
2183d72fe1eSJed Brown   }
2193d72fe1eSJed Brown   PetscFunctionReturn(0);
2203d72fe1eSJed Brown }
2213d72fe1eSJed Brown 
2229b54502bSHong Zhang static PetscErrorCode PCApplyTranspose_Cholesky(PC pc,Vec x,Vec y)
2239b54502bSHong Zhang {
2249b54502bSHong Zhang   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;
2259b54502bSHong Zhang   PetscErrorCode ierr;
2269b54502bSHong Zhang 
2279b54502bSHong Zhang   PetscFunctionBegin;
2283d1c1ea0SBarry Smith   if (dir->hdr.inplace) {
2292fa5cd67SKarl Rupp     ierr = MatSolveTranspose(pc->pmat,x,y);CHKERRQ(ierr);
2302fa5cd67SKarl Rupp   } else {
2312fa5cd67SKarl Rupp     ierr = MatSolveTranspose(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr);
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 
2562401956bSBarry Smith .seealso: PCFactorSetReuseFill()
2579b54502bSHong Zhang @*/
2587087cfbeSBarry Smith PetscErrorCode  PCFactorSetReuseOrdering(PC pc,PetscBool flag)
2599b54502bSHong Zhang {
2604ac538c5SBarry Smith   PetscErrorCode ierr;
2619b54502bSHong Zhang 
2629b54502bSHong Zhang   PetscFunctionBegin;
2630700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
264acfcf0e5SJed Brown   PetscValidLogicalCollectiveBool(pc,flag,2);
2654ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFactorSetReuseOrdering_C",(PC,PetscBool),(pc,flag));CHKERRQ(ierr);
2669b54502bSHong Zhang   PetscFunctionReturn(0);
2679b54502bSHong Zhang }
2689b54502bSHong Zhang 
2699b54502bSHong Zhang /*MC
27096fc60bcSBarry Smith    PCCHOLESKY - Uses a direct solver, based on Cholesky factorization, as a preconditioner
2719b54502bSHong Zhang 
2729b54502bSHong Zhang    Options Database Keys:
2732401956bSBarry Smith +  -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering()
2743ca39a21SBarry Smith .  -pc_factor_mat_solver_type - Actives PCFactorSetMatSolverType() to choose the direct solver, like superlu
2752401956bSBarry Smith .  -pc_factor_reuse_fill - Activates PCFactorSetReuseFill()
27655ba2a51SBarry Smith .  -pc_factor_fill <fill> - Sets fill amount
2772401956bSBarry Smith .  -pc_factor_in_place - Activates in-place factorization
278145b9266SHong Zhang -  -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine
2799b54502bSHong Zhang 
28095452b02SPatrick Sanan    Notes:
28195452b02SPatrick Sanan     Not all options work for all matrix formats
2829b54502bSHong Zhang 
2839b54502bSHong Zhang    Level: beginner
2849b54502bSHong Zhang 
28595452b02SPatrick Sanan    Notes:
28695452b02SPatrick Sanan     Usually this will compute an "exact" solution in one iteration and does
2879b54502bSHong Zhang           not need a Krylov method (i.e. you can use -ksp_type preonly, or
2889b54502bSHong Zhang           KSPSetType(ksp,KSPPREONLY) for the Krylov method
2899b54502bSHong Zhang 
2909b54502bSHong Zhang .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
291a4fd02acSBarry Smith            PCILU, PCLU, PCICC, PCFactorSetReuseOrdering(), PCFactorSetReuseFill(), PCFactorGetMatrix(),
292145b9266SHong Zhang            PCFactorSetFill(), PCFactorSetShiftNonzero(), PCFactorSetShiftType(), PCFactorSetShiftAmount()
2938e37d05fSBarry Smith            PCFactorSetUseInPlace(), PCFactorGetUseInPlace(), PCFactorSetMatOrderingType()
2949b54502bSHong Zhang 
2959b54502bSHong Zhang M*/
2969b54502bSHong Zhang 
2978cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_Cholesky(PC pc)
2989b54502bSHong Zhang {
2999b54502bSHong Zhang   PetscErrorCode ierr;
3009b54502bSHong Zhang   PC_Cholesky    *dir;
3019b54502bSHong Zhang 
3029b54502bSHong Zhang   PetscFunctionBegin;
303b00a9115SJed Brown   ierr     = PetscNewLog(pc,&dir);CHKERRQ(ierr);
3043d1c1ea0SBarry Smith   pc->data = (void*)dir;
3053d1c1ea0SBarry Smith   ierr     = PCFactorInitialize(pc);CHKERRQ(ierr);
3062fa5cd67SKarl Rupp 
307879e8a4dSBarry Smith   ((PC_Factor*)dir)->factortype         = MAT_FACTOR_CHOLESKY;
308075768bcSBarry Smith   ((PC_Factor*)dir)->info.fill          = 5.0;
3092fa5cd67SKarl Rupp 
3100a545947SLisandro Dalcin   dir->col = NULL;
3110a545947SLisandro Dalcin   dir->row = NULL;
3122fa5cd67SKarl Rupp 
3139bd791bbSBarry Smith   /* MATORDERINGNATURAL_OR_ND allows selecting type based on matrix type sbaij or aij */
3149bd791bbSBarry Smith   ierr = PetscStrallocpy(MATORDERINGNATURAL_OR_ND,(char**)&((PC_Factor*)dir)->ordering);CHKERRQ(ierr);
3159b54502bSHong Zhang 
3169b54502bSHong Zhang   pc->ops->destroy             = PCDestroy_Cholesky;
31769d2c0f9SBarry Smith   pc->ops->reset               = PCReset_Cholesky;
3189b54502bSHong Zhang   pc->ops->apply               = PCApply_Cholesky;
3197b6e2003SPierre Jolivet   pc->ops->matapply            = PCMatApply_Cholesky;
3203d72fe1eSJed Brown   pc->ops->applysymmetricleft  = PCApplySymmetricLeft_Cholesky;
3213d72fe1eSJed Brown   pc->ops->applysymmetricright = PCApplySymmetricRight_Cholesky;
3229b54502bSHong Zhang   pc->ops->applytranspose      = PCApplyTranspose_Cholesky;
3239b54502bSHong Zhang   pc->ops->setup               = PCSetUp_Cholesky;
3249b54502bSHong Zhang   pc->ops->setfromoptions      = PCSetFromOptions_Cholesky;
32592e08861SBarry Smith   pc->ops->view                = PCView_Factor;
3260a545947SLisandro Dalcin   pc->ops->applyrichardson     = NULL;
3279b54502bSHong Zhang   PetscFunctionReturn(0);
3289b54502bSHong Zhang }
329