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 */ 444ac6704cSBarry Smith ierr = PCFactorSetDefaultOrdering_Factor(pc);CHKERRQ(ierr); 45075768bcSBarry Smith ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 469b54502bSHong Zhang if (dir->col && (dir->row != dir->col)) { /* only use row ordering for SBAIJ */ 47fcfd50ebSBarry Smith ierr = ISDestroy(&dir->col);CHKERRQ(ierr); 489b54502bSHong Zhang } 493bb1ff40SBarry Smith if (dir->row) {ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);CHKERRQ(ierr);} 50075768bcSBarry Smith ierr = MatCholeskyFactor(pc->pmat,dir->row,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 5100e125f8SBarry Smith ierr = MatFactorGetError(pc->pmat,&err);CHKERRQ(ierr); 5200e125f8SBarry Smith if (err) { /* Factor() fails */ 5300e125f8SBarry Smith pc->failedreason = (PCFailedReason)err; 546baea169SHong Zhang PetscFunctionReturn(0); 556baea169SHong Zhang } 562fa5cd67SKarl Rupp 57075768bcSBarry Smith ((PC_Factor*)dir)->fact = pc->pmat; 589b54502bSHong Zhang } else { 599b54502bSHong Zhang MatInfo info; 6000e125f8SBarry Smith 619b54502bSHong Zhang if (!pc->setupcalled) { 62f73b0415SBarry Smith PetscBool canuseordering; 632c7c0729SBarry Smith if (!((PC_Factor*)dir)->fact) { 642c7c0729SBarry Smith ierr = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_CHOLESKY,&((PC_Factor*)dir)->fact);CHKERRQ(ierr); 652c7c0729SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact);CHKERRQ(ierr); 662c7c0729SBarry Smith } 67f73b0415SBarry Smith ierr = MatFactorGetCanUseOrdering(((PC_Factor*)dir)->fact,&canuseordering);CHKERRQ(ierr); 68f73b0415SBarry Smith if (canuseordering) { 694ac6704cSBarry Smith ierr = PCFactorSetDefaultOrdering_Factor(pc);CHKERRQ(ierr); 70075768bcSBarry Smith ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 719bfd6278SHong Zhang /* check if dir->row == dir->col */ 724ac6704cSBarry Smith if (dir->row) { 739bfd6278SHong Zhang ierr = ISEqual(dir->row,dir->col,&flg);CHKERRQ(ierr); 74*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(!flg,PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"row and column permutations must be equal"); 754ac6704cSBarry Smith } 76fcfd50ebSBarry Smith ierr = ISDestroy(&dir->col);CHKERRQ(ierr); /* only pass one ordering into CholeskyFactor */ 779bfd6278SHong Zhang 7890d69ab7SBarry Smith flg = PETSC_FALSE; 79c5929fdfSBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&flg,NULL);CHKERRQ(ierr); 809b54502bSHong Zhang if (flg) { 819b54502bSHong Zhang PetscReal tol = 1.e-10; 82c5929fdfSBarry Smith ierr = PetscOptionsGetReal(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&tol,NULL);CHKERRQ(ierr); 839b54502bSHong Zhang ierr = MatReorderForNonzeroDiagonal(pc->pmat,tol,dir->row,dir->row);CHKERRQ(ierr); 849b54502bSHong Zhang } 852c7c0729SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);CHKERRQ(ierr); 86a1f19f5aSHong Zhang } 87075768bcSBarry Smith ierr = MatCholeskyFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 88075768bcSBarry Smith ierr = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 893d1c1ea0SBarry Smith dir->hdr.actualfill = info.fill_ratio_needed; 909b54502bSHong Zhang } else if (pc->flag != SAME_NONZERO_PATTERN) { 913d1c1ea0SBarry Smith if (!dir->hdr.reuseordering) { 92f73b0415SBarry Smith PetscBool canuseordering; 932c7c0729SBarry Smith ierr = MatDestroy(&((PC_Factor*)dir)->fact);CHKERRQ(ierr); 942c7c0729SBarry Smith ierr = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_CHOLESKY,&((PC_Factor*)dir)->fact);CHKERRQ(ierr); 952c7c0729SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact);CHKERRQ(ierr); 96f73b0415SBarry Smith ierr = MatFactorGetCanUseOrdering(((PC_Factor*)dir)->fact,&canuseordering);CHKERRQ(ierr); 97f73b0415SBarry Smith if (canuseordering) { 98fcfd50ebSBarry Smith ierr = ISDestroy(&dir->row);CHKERRQ(ierr); 994ac6704cSBarry Smith ierr = PCFactorSetDefaultOrdering_Factor(pc);CHKERRQ(ierr); 100075768bcSBarry Smith ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 1019d61c4f5SHong Zhang ierr = ISDestroy(&dir->col);CHKERRQ(ierr); /* only use dir->row ordering in CholeskyFactor */ 1029d61c4f5SHong Zhang 10390d69ab7SBarry Smith flg = PETSC_FALSE; 104c5929fdfSBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&flg,NULL);CHKERRQ(ierr); 1059b54502bSHong Zhang if (flg) { 1069b54502bSHong Zhang PetscReal tol = 1.e-10; 107c5929fdfSBarry Smith ierr = PetscOptionsGetReal(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&tol,NULL);CHKERRQ(ierr); 1089b54502bSHong Zhang ierr = MatReorderForNonzeroDiagonal(pc->pmat,tol,dir->row,dir->row);CHKERRQ(ierr); 1099b54502bSHong Zhang } 1102c7c0729SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);CHKERRQ(ierr); 1119b54502bSHong Zhang } 1122c7c0729SBarry Smith } 113075768bcSBarry Smith ierr = MatCholeskyFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 114075768bcSBarry Smith ierr = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 1153d1c1ea0SBarry Smith dir->hdr.actualfill = info.fill_ratio_needed; 1163bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact);CHKERRQ(ierr); 11704545d6dSBarry Smith } else { 118b8b68cfdSBarry Smith ierr = MatFactorGetError(((PC_Factor*)dir)->fact,&err);CHKERRQ(ierr); 119160a8794SBarry Smith if (err == MAT_FACTOR_NUMERIC_ZEROPIVOT) { 120b8b68cfdSBarry Smith ierr = MatFactorClearError(((PC_Factor*)dir)->fact);CHKERRQ(ierr); 121b8b68cfdSBarry Smith pc->failedreason = PC_NOERROR; 12204545d6dSBarry Smith } 1239b54502bSHong Zhang } 12400e125f8SBarry Smith ierr = MatFactorGetError(((PC_Factor*)dir)->fact,&err);CHKERRQ(ierr); 12500e125f8SBarry Smith if (err) { /* FactorSymbolic() fails */ 12600e125f8SBarry Smith pc->failedreason = (PCFailedReason)err; 1276baea169SHong Zhang PetscFunctionReturn(0); 1286baea169SHong Zhang } 1296baea169SHong Zhang 130075768bcSBarry Smith ierr = MatCholeskyFactorNumeric(((PC_Factor*)dir)->fact,pc->pmat,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 13100e125f8SBarry Smith ierr = MatFactorGetError(((PC_Factor*)dir)->fact,&err);CHKERRQ(ierr); 13200e125f8SBarry Smith if (err) { /* FactorNumeric() fails */ 13300e125f8SBarry Smith pc->failedreason = (PCFailedReason)err; 1346baea169SHong Zhang } 1359b54502bSHong Zhang } 13600c67f3bSHong Zhang 1373ca39a21SBarry Smith ierr = PCFactorGetMatSolverType(pc,&stype);CHKERRQ(ierr); 13800c67f3bSHong Zhang if (!stype) { 139ea799195SBarry Smith MatSolverType solverpackage; 1403ca39a21SBarry Smith ierr = MatFactorGetSolverType(((PC_Factor*)dir)->fact,&solverpackage);CHKERRQ(ierr); 1413ca39a21SBarry Smith ierr = PCFactorSetMatSolverType(pc,solverpackage);CHKERRQ(ierr); 14200c67f3bSHong Zhang } 1439b54502bSHong Zhang PetscFunctionReturn(0); 1449b54502bSHong Zhang } 1459b54502bSHong Zhang 14669d2c0f9SBarry Smith static PetscErrorCode PCReset_Cholesky(PC pc) 1479b54502bSHong Zhang { 1489b54502bSHong Zhang PC_Cholesky *dir = (PC_Cholesky*)pc->data; 1499b54502bSHong Zhang PetscErrorCode ierr; 1509b54502bSHong Zhang 1519b54502bSHong Zhang PetscFunctionBegin; 1523d1c1ea0SBarry Smith if (!dir->hdr.inplace && ((PC_Factor*)dir)->fact) {ierr = MatDestroy(&((PC_Factor*)dir)->fact);CHKERRQ(ierr);} 153fcfd50ebSBarry Smith ierr = ISDestroy(&dir->row);CHKERRQ(ierr); 154fcfd50ebSBarry Smith ierr = ISDestroy(&dir->col);CHKERRQ(ierr); 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 PetscErrorCode ierr; 16269d2c0f9SBarry Smith 16369d2c0f9SBarry Smith PetscFunctionBegin; 16469d2c0f9SBarry Smith ierr = PCReset_Cholesky(pc);CHKERRQ(ierr); 165503cfb0cSBarry Smith ierr = PetscFree(((PC_Factor*)dir)->ordering);CHKERRQ(ierr); 166503cfb0cSBarry Smith ierr = PetscFree(((PC_Factor*)dir)->solvertype);CHKERRQ(ierr); 167c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 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 PetscErrorCode ierr; 1759b54502bSHong Zhang 1769b54502bSHong Zhang PetscFunctionBegin; 1773d1c1ea0SBarry Smith if (dir->hdr.inplace) { 1782fa5cd67SKarl Rupp ierr = MatSolve(pc->pmat,x,y);CHKERRQ(ierr); 1792fa5cd67SKarl Rupp } else { 1802fa5cd67SKarl Rupp ierr = MatSolve(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr); 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 PetscErrorCode ierr; 1897b6e2003SPierre Jolivet 1907b6e2003SPierre Jolivet PetscFunctionBegin; 1917b6e2003SPierre Jolivet if (dir->hdr.inplace) { 1927b6e2003SPierre Jolivet ierr = MatMatSolve(pc->pmat,X,Y);CHKERRQ(ierr); 1937b6e2003SPierre Jolivet } else { 1947b6e2003SPierre Jolivet ierr = MatMatSolve(((PC_Factor*)dir)->fact,X,Y);CHKERRQ(ierr); 1957b6e2003SPierre Jolivet } 1967b6e2003SPierre Jolivet PetscFunctionReturn(0); 1977b6e2003SPierre Jolivet } 1987b6e2003SPierre Jolivet 1993d72fe1eSJed Brown static PetscErrorCode PCApplySymmetricLeft_Cholesky(PC pc,Vec x,Vec y) 2003d72fe1eSJed Brown { 2013d72fe1eSJed Brown PC_Cholesky *dir = (PC_Cholesky*)pc->data; 2023d72fe1eSJed Brown PetscErrorCode ierr; 2033d72fe1eSJed Brown 2043d72fe1eSJed Brown PetscFunctionBegin; 2053d72fe1eSJed Brown if (dir->hdr.inplace) { 2063d72fe1eSJed Brown ierr = MatForwardSolve(pc->pmat,x,y);CHKERRQ(ierr); 2073d72fe1eSJed Brown } else { 2083d72fe1eSJed Brown ierr = MatForwardSolve(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr); 2093d72fe1eSJed Brown } 2103d72fe1eSJed Brown PetscFunctionReturn(0); 2113d72fe1eSJed Brown } 2123d72fe1eSJed Brown 2133d72fe1eSJed Brown static PetscErrorCode PCApplySymmetricRight_Cholesky(PC pc,Vec x,Vec y) 2143d72fe1eSJed Brown { 2153d72fe1eSJed Brown PC_Cholesky *dir = (PC_Cholesky*)pc->data; 2163d72fe1eSJed Brown PetscErrorCode ierr; 2173d72fe1eSJed Brown 2183d72fe1eSJed Brown PetscFunctionBegin; 2193d72fe1eSJed Brown if (dir->hdr.inplace) { 2203d72fe1eSJed Brown ierr = MatBackwardSolve(pc->pmat,x,y);CHKERRQ(ierr); 2213d72fe1eSJed Brown } else { 2223d72fe1eSJed Brown ierr = MatBackwardSolve(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr); 2233d72fe1eSJed Brown } 2243d72fe1eSJed Brown PetscFunctionReturn(0); 2253d72fe1eSJed Brown } 2263d72fe1eSJed Brown 2279b54502bSHong Zhang static PetscErrorCode PCApplyTranspose_Cholesky(PC pc,Vec x,Vec y) 2289b54502bSHong Zhang { 2299b54502bSHong Zhang PC_Cholesky *dir = (PC_Cholesky*)pc->data; 2309b54502bSHong Zhang PetscErrorCode ierr; 2319b54502bSHong Zhang 2329b54502bSHong Zhang PetscFunctionBegin; 2333d1c1ea0SBarry Smith if (dir->hdr.inplace) { 2342fa5cd67SKarl Rupp ierr = MatSolveTranspose(pc->pmat,x,y);CHKERRQ(ierr); 2352fa5cd67SKarl Rupp } else { 2362fa5cd67SKarl Rupp ierr = MatSolveTranspose(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr); 2372fa5cd67SKarl Rupp } 2389b54502bSHong Zhang PetscFunctionReturn(0); 2399b54502bSHong Zhang } 2409b54502bSHong Zhang 2419b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/ 2429b54502bSHong Zhang 2439b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/ 2449b54502bSHong Zhang 2459b54502bSHong Zhang /*@ 2462401956bSBarry Smith PCFactorSetReuseOrdering - When similar matrices are factored, this 2479b54502bSHong Zhang causes the ordering computed in the first factor to be used for all 2489b54502bSHong Zhang following factors. 2499b54502bSHong Zhang 250ad4df100SBarry Smith Logically Collective on PC 2519b54502bSHong Zhang 2529b54502bSHong Zhang Input Parameters: 2539b54502bSHong Zhang + pc - the preconditioner context 2549b54502bSHong Zhang - flag - PETSC_TRUE to reuse else PETSC_FALSE 2559b54502bSHong Zhang 2569b54502bSHong Zhang Options Database Key: 2572401956bSBarry Smith . -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering() 2589b54502bSHong Zhang 2599b54502bSHong Zhang Level: intermediate 2609b54502bSHong Zhang 2612401956bSBarry Smith .seealso: PCFactorSetReuseFill() 2629b54502bSHong Zhang @*/ 2637087cfbeSBarry Smith PetscErrorCode PCFactorSetReuseOrdering(PC pc,PetscBool flag) 2649b54502bSHong Zhang { 2654ac538c5SBarry Smith PetscErrorCode ierr; 2669b54502bSHong Zhang 2679b54502bSHong Zhang PetscFunctionBegin; 2680700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 269acfcf0e5SJed Brown PetscValidLogicalCollectiveBool(pc,flag,2); 2704ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFactorSetReuseOrdering_C",(PC,PetscBool),(pc,flag));CHKERRQ(ierr); 2719b54502bSHong Zhang PetscFunctionReturn(0); 2729b54502bSHong Zhang } 2739b54502bSHong Zhang 2749b54502bSHong Zhang /*MC 27596fc60bcSBarry Smith PCCHOLESKY - Uses a direct solver, based on Cholesky factorization, as a preconditioner 2769b54502bSHong Zhang 2779b54502bSHong Zhang Options Database Keys: 2782401956bSBarry Smith + -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering() 2793ca39a21SBarry Smith . -pc_factor_mat_solver_type - Actives PCFactorSetMatSolverType() to choose the direct solver, like superlu 2802401956bSBarry Smith . -pc_factor_reuse_fill - Activates PCFactorSetReuseFill() 28155ba2a51SBarry Smith . -pc_factor_fill <fill> - Sets fill amount 2822401956bSBarry Smith . -pc_factor_in_place - Activates in-place factorization 283145b9266SHong Zhang - -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine 2849b54502bSHong Zhang 28595452b02SPatrick Sanan Notes: 28695452b02SPatrick Sanan Not all options work for all matrix formats 2879b54502bSHong Zhang 2889b54502bSHong Zhang Level: beginner 2899b54502bSHong Zhang 29095452b02SPatrick Sanan Notes: 29195452b02SPatrick Sanan Usually this will compute an "exact" solution in one iteration and does 2929b54502bSHong Zhang not need a Krylov method (i.e. you can use -ksp_type preonly, or 2939b54502bSHong Zhang KSPSetType(ksp,KSPPREONLY) for the Krylov method 2949b54502bSHong Zhang 2959b54502bSHong Zhang .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 296a4fd02acSBarry Smith PCILU, PCLU, PCICC, PCFactorSetReuseOrdering(), PCFactorSetReuseFill(), PCFactorGetMatrix(), 297145b9266SHong Zhang PCFactorSetFill(), PCFactorSetShiftNonzero(), PCFactorSetShiftType(), PCFactorSetShiftAmount() 2988e37d05fSBarry Smith PCFactorSetUseInPlace(), PCFactorGetUseInPlace(), PCFactorSetMatOrderingType() 2999b54502bSHong Zhang 3009b54502bSHong Zhang M*/ 3019b54502bSHong Zhang 3028cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_Cholesky(PC pc) 3039b54502bSHong Zhang { 3049b54502bSHong Zhang PetscErrorCode ierr; 3059b54502bSHong Zhang PC_Cholesky *dir; 3069b54502bSHong Zhang 3079b54502bSHong Zhang PetscFunctionBegin; 308b00a9115SJed Brown ierr = PetscNewLog(pc,&dir);CHKERRQ(ierr); 3093d1c1ea0SBarry Smith pc->data = (void*)dir; 3104ac6704cSBarry Smith ierr = PCFactorInitialize(pc,MAT_FACTOR_CHOLESKY);CHKERRQ(ierr); 3112fa5cd67SKarl Rupp 312075768bcSBarry Smith ((PC_Factor*)dir)->info.fill = 5.0; 3132fa5cd67SKarl Rupp 3149b54502bSHong Zhang pc->ops->destroy = PCDestroy_Cholesky; 31569d2c0f9SBarry Smith pc->ops->reset = PCReset_Cholesky; 3169b54502bSHong Zhang pc->ops->apply = PCApply_Cholesky; 3177b6e2003SPierre Jolivet pc->ops->matapply = PCMatApply_Cholesky; 3183d72fe1eSJed Brown pc->ops->applysymmetricleft = PCApplySymmetricLeft_Cholesky; 3193d72fe1eSJed Brown pc->ops->applysymmetricright = PCApplySymmetricRight_Cholesky; 3209b54502bSHong Zhang pc->ops->applytranspose = PCApplyTranspose_Cholesky; 3219b54502bSHong Zhang pc->ops->setup = PCSetUp_Cholesky; 3229b54502bSHong Zhang pc->ops->setfromoptions = PCSetFromOptions_Cholesky; 32392e08861SBarry Smith pc->ops->view = PCView_Factor; 3240a545947SLisandro Dalcin pc->ops->applyrichardson = NULL; 3259b54502bSHong Zhang PetscFunctionReturn(0); 3269b54502bSHong Zhang } 327