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 PetscReal actualfill; /* actual fill in factor */ 12ace3abfcSBarry Smith PetscBool inplace; /* flag indicating in-place factorization */ 139b54502bSHong Zhang IS row,col; /* index sets used for reordering */ 14ace3abfcSBarry Smith PetscBool reuseordering; /* reuses previous reordering computed */ 15ace3abfcSBarry Smith PetscBool reusefill; /* reuse fill from previous Cholesky */ 169b54502bSHong Zhang } PC_Cholesky; 179b54502bSHong Zhang 189b54502bSHong Zhang #undef __FUNCT__ 192401956bSBarry Smith #define __FUNCT__ "PCFactorSetReuseOrdering_Cholesky" 20b2573a8aSBarry Smith static PetscErrorCode PCFactorSetReuseOrdering_Cholesky(PC pc,PetscBool flag) 219b54502bSHong Zhang { 2285a85e9cSBarry Smith PC_Cholesky *lu = (PC_Cholesky*)pc->data; 239b54502bSHong Zhang 249b54502bSHong Zhang PetscFunctionBegin; 259b54502bSHong Zhang lu->reuseordering = flag; 269b54502bSHong Zhang PetscFunctionReturn(0); 279b54502bSHong Zhang } 289b54502bSHong Zhang 299b54502bSHong Zhang #undef __FUNCT__ 302401956bSBarry Smith #define __FUNCT__ "PCFactorSetReuseFill_Cholesky" 31b2573a8aSBarry Smith static PetscErrorCode PCFactorSetReuseFill_Cholesky(PC pc,PetscBool flag) 329b54502bSHong Zhang { 3385a85e9cSBarry Smith PC_Cholesky *lu = (PC_Cholesky*)pc->data; 349b54502bSHong Zhang 359b54502bSHong Zhang PetscFunctionBegin; 369b54502bSHong Zhang lu->reusefill = flag; 379b54502bSHong Zhang PetscFunctionReturn(0); 389b54502bSHong Zhang } 399b54502bSHong Zhang 409b54502bSHong Zhang #undef __FUNCT__ 419b54502bSHong Zhang #define __FUNCT__ "PCSetFromOptions_Cholesky" 424416b707SBarry Smith static PetscErrorCode PCSetFromOptions_Cholesky(PetscOptionItems *PetscOptionsObject,PC pc) 439b54502bSHong Zhang { 449b54502bSHong Zhang PetscErrorCode ierr; 459b54502bSHong Zhang 469b54502bSHong Zhang PetscFunctionBegin; 47e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"Cholesky options");CHKERRQ(ierr); 48e55864a3SBarry Smith ierr = PCSetFromOptions_Factor(PetscOptionsObject,pc);CHKERRQ(ierr); 499b54502bSHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 509b54502bSHong Zhang PetscFunctionReturn(0); 519b54502bSHong Zhang } 529b54502bSHong Zhang 539b54502bSHong Zhang #undef __FUNCT__ 549b54502bSHong Zhang #define __FUNCT__ "PCView_Cholesky" 559b54502bSHong Zhang static PetscErrorCode PCView_Cholesky(PC pc,PetscViewer viewer) 569b54502bSHong Zhang { 57145b9266SHong Zhang PC_Cholesky *chol = (PC_Cholesky*)pc->data; 589b54502bSHong Zhang PetscErrorCode ierr; 59ace3abfcSBarry Smith PetscBool iascii; 609b54502bSHong Zhang 619b54502bSHong Zhang PetscFunctionBegin; 62251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 639b54502bSHong Zhang if (iascii) { 64914a5d51SHong Zhang if (chol->inplace) { 65914a5d51SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Cholesky: in-place factorization\n");CHKERRQ(ierr); 66914a5d51SHong Zhang } else { 67914a5d51SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Cholesky: out-of-place factorization\n");CHKERRQ(ierr); 68145b9266SHong Zhang } 699b54502bSHong Zhang 70145b9266SHong Zhang if (chol->reusefill) {ierr = PetscViewerASCIIPrintf(viewer," Reusing fill from past factorization\n");CHKERRQ(ierr);} 71145b9266SHong Zhang if (chol->reuseordering) {ierr = PetscViewerASCIIPrintf(viewer," Reusing reordering from past factorization\n");CHKERRQ(ierr);} 729b54502bSHong Zhang } 73914a5d51SHong Zhang ierr = PCView_Factor(pc,viewer);CHKERRQ(ierr); 749b54502bSHong Zhang PetscFunctionReturn(0); 759b54502bSHong Zhang } 769b54502bSHong Zhang 779b54502bSHong Zhang 789b54502bSHong Zhang #undef __FUNCT__ 799b54502bSHong Zhang #define __FUNCT__ "PCSetUp_Cholesky" 809b54502bSHong Zhang static PetscErrorCode PCSetUp_Cholesky(PC pc) 819b54502bSHong Zhang { 829b54502bSHong Zhang PetscErrorCode ierr; 83ace3abfcSBarry Smith PetscBool flg; 849b54502bSHong Zhang PC_Cholesky *dir = (PC_Cholesky*)pc->data; 859b54502bSHong Zhang 869b54502bSHong Zhang PetscFunctionBegin; 87075768bcSBarry Smith if (dir->reusefill && pc->setupcalled) ((PC_Factor*)dir)->info.fill = dir->actualfill; 889b54502bSHong Zhang 899b54502bSHong Zhang if (dir->inplace) { 909b54502bSHong Zhang if (dir->row && dir->col && (dir->row != dir->col)) { 91fcfd50ebSBarry Smith ierr = ISDestroy(&dir->row);CHKERRQ(ierr); 929b54502bSHong Zhang } 93fcfd50ebSBarry Smith ierr = ISDestroy(&dir->col);CHKERRQ(ierr); 94075768bcSBarry Smith ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 959b54502bSHong Zhang if (dir->col && (dir->row != dir->col)) { /* only use row ordering for SBAIJ */ 96fcfd50ebSBarry Smith ierr = ISDestroy(&dir->col);CHKERRQ(ierr); 979b54502bSHong Zhang } 983bb1ff40SBarry Smith if (dir->row) {ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);CHKERRQ(ierr);} 99075768bcSBarry Smith ierr = MatCholeskyFactor(pc->pmat,dir->row,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 100*6baea169SHong Zhang if (((PC_Factor*)dir)->info.errortype) { /* Factor() fails */ 101*6baea169SHong Zhang MatFactorInfo factinfo=((PC_Factor*)dir)->info; 102*6baea169SHong Zhang pc->failedreason = (PCFailedReason)factinfo.errortype; 103*6baea169SHong Zhang PetscFunctionReturn(0); 104*6baea169SHong Zhang } 1052fa5cd67SKarl Rupp 106075768bcSBarry Smith ((PC_Factor*)dir)->fact = pc->pmat; 1079b54502bSHong Zhang } else { 1089b54502bSHong Zhang MatInfo info; 1099b54502bSHong Zhang if (!pc->setupcalled) { 110075768bcSBarry Smith ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 1119bfd6278SHong Zhang /* check if dir->row == dir->col */ 1129bfd6278SHong Zhang ierr = ISEqual(dir->row,dir->col,&flg);CHKERRQ(ierr); 113e32f2f54SBarry Smith if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"row and column permutations must equal"); 114fcfd50ebSBarry Smith ierr = ISDestroy(&dir->col);CHKERRQ(ierr); /* only pass one ordering into CholeskyFactor */ 1159bfd6278SHong Zhang 11690d69ab7SBarry Smith flg = PETSC_FALSE; 117c5929fdfSBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&flg,NULL);CHKERRQ(ierr); 1189b54502bSHong Zhang if (flg) { 1199b54502bSHong Zhang PetscReal tol = 1.e-10; 120c5929fdfSBarry Smith ierr = PetscOptionsGetReal(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&tol,NULL);CHKERRQ(ierr); 1219b54502bSHong Zhang ierr = MatReorderForNonzeroDiagonal(pc->pmat,tol,dir->row,dir->row);CHKERRQ(ierr); 1229b54502bSHong Zhang } 1233bb1ff40SBarry Smith if (dir->row) {ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);CHKERRQ(ierr);} 124d09a07f4SBarry Smith if (!((PC_Factor*)dir)->fact) { 125075768bcSBarry Smith ierr = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_CHOLESKY,&((PC_Factor*)dir)->fact);CHKERRQ(ierr); 126a1f19f5aSHong Zhang } 127075768bcSBarry Smith ierr = MatCholeskyFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 128075768bcSBarry Smith ierr = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 1299b54502bSHong Zhang dir->actualfill = info.fill_ratio_needed; 1303bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact);CHKERRQ(ierr); 1319b54502bSHong Zhang } else if (pc->flag != SAME_NONZERO_PATTERN) { 1329b54502bSHong Zhang if (!dir->reuseordering) { 133fcfd50ebSBarry Smith ierr = ISDestroy(&dir->row);CHKERRQ(ierr); 134075768bcSBarry Smith ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 1359d61c4f5SHong Zhang ierr = ISDestroy(&dir->col);CHKERRQ(ierr); /* only use dir->row ordering in CholeskyFactor */ 1369d61c4f5SHong Zhang 13790d69ab7SBarry Smith flg = PETSC_FALSE; 138c5929fdfSBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&flg,NULL);CHKERRQ(ierr); 1399b54502bSHong Zhang if (flg) { 1409b54502bSHong Zhang PetscReal tol = 1.e-10; 141c5929fdfSBarry Smith ierr = PetscOptionsGetReal(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&tol,NULL);CHKERRQ(ierr); 1429b54502bSHong Zhang ierr = MatReorderForNonzeroDiagonal(pc->pmat,tol,dir->row,dir->row);CHKERRQ(ierr); 1439b54502bSHong Zhang } 1443bb1ff40SBarry Smith if (dir->row) {ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);CHKERRQ(ierr);} 1459b54502bSHong Zhang } 1466bf464f9SBarry Smith ierr = MatDestroy(&((PC_Factor*)dir)->fact);CHKERRQ(ierr); 147075768bcSBarry Smith ierr = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_CHOLESKY,&((PC_Factor*)dir)->fact);CHKERRQ(ierr); 148075768bcSBarry Smith ierr = MatCholeskyFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 149075768bcSBarry Smith ierr = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 1509b54502bSHong Zhang dir->actualfill = info.fill_ratio_needed; 1513bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact);CHKERRQ(ierr); 1529b54502bSHong Zhang } 153*6baea169SHong Zhang if (((PC_Factor*)dir)->info.errortype) { /* FactorSymbolic() fails */ 154*6baea169SHong Zhang MatFactorInfo factinfo=((PC_Factor*)dir)->info; 155*6baea169SHong Zhang pc->failedreason = (PCFailedReason)factinfo.errortype; 156*6baea169SHong Zhang PetscFunctionReturn(0); 157*6baea169SHong Zhang } 158*6baea169SHong Zhang 159075768bcSBarry Smith ierr = MatCholeskyFactorNumeric(((PC_Factor*)dir)->fact,pc->pmat,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 160*6baea169SHong Zhang if (((PC_Factor*)dir)->info.errortype) { /* FactorNumeric() fails */ 161*6baea169SHong Zhang MatFactorInfo factinfo=((PC_Factor*)dir)->info; 162*6baea169SHong Zhang pc->failedreason = (PCFailedReason)factinfo.errortype; 163*6baea169SHong Zhang } 1649b54502bSHong Zhang } 1659b54502bSHong Zhang PetscFunctionReturn(0); 1669b54502bSHong Zhang } 1679b54502bSHong Zhang 1689b54502bSHong Zhang #undef __FUNCT__ 16969d2c0f9SBarry Smith #define __FUNCT__ "PCReset_Cholesky" 17069d2c0f9SBarry Smith static PetscErrorCode PCReset_Cholesky(PC pc) 1719b54502bSHong Zhang { 1729b54502bSHong Zhang PC_Cholesky *dir = (PC_Cholesky*)pc->data; 1739b54502bSHong Zhang PetscErrorCode ierr; 1749b54502bSHong Zhang 1759b54502bSHong Zhang PetscFunctionBegin; 1766bf464f9SBarry Smith if (!dir->inplace && ((PC_Factor*)dir)->fact) {ierr = MatDestroy(&((PC_Factor*)dir)->fact);CHKERRQ(ierr);} 177fcfd50ebSBarry Smith ierr = ISDestroy(&dir->row);CHKERRQ(ierr); 178fcfd50ebSBarry Smith ierr = ISDestroy(&dir->col);CHKERRQ(ierr); 17969d2c0f9SBarry Smith PetscFunctionReturn(0); 18069d2c0f9SBarry Smith } 18169d2c0f9SBarry Smith 18269d2c0f9SBarry Smith #undef __FUNCT__ 18369d2c0f9SBarry Smith #define __FUNCT__ "PCDestroy_Cholesky" 18469d2c0f9SBarry Smith static PetscErrorCode PCDestroy_Cholesky(PC pc) 18569d2c0f9SBarry Smith { 18669d2c0f9SBarry Smith PC_Cholesky *dir = (PC_Cholesky*)pc->data; 18769d2c0f9SBarry Smith PetscErrorCode ierr; 18869d2c0f9SBarry Smith 18969d2c0f9SBarry Smith PetscFunctionBegin; 19069d2c0f9SBarry Smith ierr = PCReset_Cholesky(pc);CHKERRQ(ierr); 191503cfb0cSBarry Smith ierr = PetscFree(((PC_Factor*)dir)->ordering);CHKERRQ(ierr); 192503cfb0cSBarry Smith ierr = PetscFree(((PC_Factor*)dir)->solvertype);CHKERRQ(ierr); 193c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 1949b54502bSHong Zhang PetscFunctionReturn(0); 1959b54502bSHong Zhang } 1969b54502bSHong Zhang 1979b54502bSHong Zhang #undef __FUNCT__ 1989b54502bSHong Zhang #define __FUNCT__ "PCApply_Cholesky" 1999b54502bSHong Zhang static PetscErrorCode PCApply_Cholesky(PC pc,Vec x,Vec y) 2009b54502bSHong Zhang { 2019b54502bSHong Zhang PC_Cholesky *dir = (PC_Cholesky*)pc->data; 2029b54502bSHong Zhang PetscErrorCode ierr; 2039b54502bSHong Zhang 2049b54502bSHong Zhang PetscFunctionBegin; 2052fa5cd67SKarl Rupp if (dir->inplace) { 2062fa5cd67SKarl Rupp ierr = MatSolve(pc->pmat,x,y);CHKERRQ(ierr); 2072fa5cd67SKarl Rupp } else { 2082fa5cd67SKarl Rupp ierr = MatSolve(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr); 2092fa5cd67SKarl Rupp } 2109b54502bSHong Zhang PetscFunctionReturn(0); 2119b54502bSHong Zhang } 2129b54502bSHong Zhang 2139b54502bSHong Zhang #undef __FUNCT__ 2149b54502bSHong Zhang #define __FUNCT__ "PCApplyTranspose_Cholesky" 2159b54502bSHong Zhang static PetscErrorCode PCApplyTranspose_Cholesky(PC pc,Vec x,Vec y) 2169b54502bSHong Zhang { 2179b54502bSHong Zhang PC_Cholesky *dir = (PC_Cholesky*)pc->data; 2189b54502bSHong Zhang PetscErrorCode ierr; 2199b54502bSHong Zhang 2209b54502bSHong Zhang PetscFunctionBegin; 2212fa5cd67SKarl Rupp if (dir->inplace) { 2222fa5cd67SKarl Rupp ierr = MatSolveTranspose(pc->pmat,x,y);CHKERRQ(ierr); 2232fa5cd67SKarl Rupp } else { 2242fa5cd67SKarl Rupp ierr = MatSolveTranspose(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr); 2252fa5cd67SKarl Rupp } 2269b54502bSHong Zhang PetscFunctionReturn(0); 2279b54502bSHong Zhang } 2289b54502bSHong Zhang 2299b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/ 2309b54502bSHong Zhang 2319b54502bSHong Zhang #undef __FUNCT__ 2322401956bSBarry Smith #define __FUNCT__ "PCFactorSetUseInPlace_Cholesky" 2338e37d05fSBarry Smith static PetscErrorCode PCFactorSetUseInPlace_Cholesky(PC pc,PetscBool flg) 2349b54502bSHong Zhang { 2358e37d05fSBarry Smith PC_Cholesky *dir = (PC_Cholesky*)pc->data; 2369b54502bSHong Zhang 2379b54502bSHong Zhang PetscFunctionBegin; 2388e37d05fSBarry Smith dir->inplace = flg; 2398e37d05fSBarry Smith PetscFunctionReturn(0); 2408e37d05fSBarry Smith } 2418e37d05fSBarry Smith 2428e37d05fSBarry Smith #undef __FUNCT__ 2438e37d05fSBarry Smith #define __FUNCT__ "PCFactorGetUseInPlace_Cholesky" 2448e37d05fSBarry Smith static PetscErrorCode PCFactorGetUseInPlace_Cholesky(PC pc,PetscBool *flg) 2458e37d05fSBarry Smith { 2468e37d05fSBarry Smith PC_Cholesky *dir = (PC_Cholesky*)pc->data; 2478e37d05fSBarry Smith 2488e37d05fSBarry Smith PetscFunctionBegin; 2498e37d05fSBarry Smith *flg = dir->inplace; 2509b54502bSHong Zhang PetscFunctionReturn(0); 2519b54502bSHong Zhang } 2529b54502bSHong Zhang 2539b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/ 2549b54502bSHong Zhang 2559b54502bSHong Zhang #undef __FUNCT__ 2562401956bSBarry Smith #define __FUNCT__ "PCFactorSetReuseOrdering" 2579b54502bSHong Zhang /*@ 2582401956bSBarry Smith PCFactorSetReuseOrdering - When similar matrices are factored, this 2599b54502bSHong Zhang causes the ordering computed in the first factor to be used for all 2609b54502bSHong Zhang following factors. 2619b54502bSHong Zhang 262ad4df100SBarry Smith Logically Collective on PC 2639b54502bSHong Zhang 2649b54502bSHong Zhang Input Parameters: 2659b54502bSHong Zhang + pc - the preconditioner context 2669b54502bSHong Zhang - flag - PETSC_TRUE to reuse else PETSC_FALSE 2679b54502bSHong Zhang 2689b54502bSHong Zhang Options Database Key: 2692401956bSBarry Smith . -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering() 2709b54502bSHong Zhang 2719b54502bSHong Zhang Level: intermediate 2729b54502bSHong Zhang 2739b54502bSHong Zhang .keywords: PC, levels, reordering, factorization, incomplete, LU 2749b54502bSHong Zhang 2752401956bSBarry Smith .seealso: PCFactorSetReuseFill() 2769b54502bSHong Zhang @*/ 2777087cfbeSBarry Smith PetscErrorCode PCFactorSetReuseOrdering(PC pc,PetscBool flag) 2789b54502bSHong Zhang { 2794ac538c5SBarry Smith PetscErrorCode ierr; 2809b54502bSHong Zhang 2819b54502bSHong Zhang PetscFunctionBegin; 2820700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 283acfcf0e5SJed Brown PetscValidLogicalCollectiveBool(pc,flag,2); 2844ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFactorSetReuseOrdering_C",(PC,PetscBool),(pc,flag));CHKERRQ(ierr); 2859b54502bSHong Zhang PetscFunctionReturn(0); 2869b54502bSHong Zhang } 2879b54502bSHong Zhang 2889b54502bSHong Zhang /*MC 28996fc60bcSBarry Smith PCCHOLESKY - Uses a direct solver, based on Cholesky factorization, as a preconditioner 2909b54502bSHong Zhang 2919b54502bSHong Zhang Options Database Keys: 2922401956bSBarry Smith + -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering() 293f60c3dc2SHong Zhang . -pc_factor_mat_solver_package - Actives PCFactorSetMatSolverPackage() to choose the direct solver, like superlu 2942401956bSBarry Smith . -pc_factor_reuse_fill - Activates PCFactorSetReuseFill() 29555ba2a51SBarry Smith . -pc_factor_fill <fill> - Sets fill amount 2962401956bSBarry Smith . -pc_factor_in_place - Activates in-place factorization 297145b9266SHong Zhang - -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine 2989b54502bSHong Zhang 2999b54502bSHong Zhang Notes: Not all options work for all matrix formats 3009b54502bSHong Zhang 3019b54502bSHong Zhang Level: beginner 3029b54502bSHong Zhang 3039b54502bSHong Zhang Concepts: Cholesky factorization, direct solver 3049b54502bSHong Zhang 3059b54502bSHong Zhang Notes: Usually this will compute an "exact" solution in one iteration and does 3069b54502bSHong Zhang not need a Krylov method (i.e. you can use -ksp_type preonly, or 3079b54502bSHong Zhang KSPSetType(ksp,KSPPREONLY) for the Krylov method 3089b54502bSHong Zhang 3099b54502bSHong Zhang .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 310a4fd02acSBarry Smith PCILU, PCLU, PCICC, PCFactorSetReuseOrdering(), PCFactorSetReuseFill(), PCFactorGetMatrix(), 311145b9266SHong Zhang PCFactorSetFill(), PCFactorSetShiftNonzero(), PCFactorSetShiftType(), PCFactorSetShiftAmount() 3128e37d05fSBarry Smith PCFactorSetUseInPlace(), PCFactorGetUseInPlace(), PCFactorSetMatOrderingType() 3139b54502bSHong Zhang 3149b54502bSHong Zhang M*/ 3159b54502bSHong Zhang 3169b54502bSHong Zhang #undef __FUNCT__ 3179b54502bSHong Zhang #define __FUNCT__ "PCCreate_Cholesky" 3188cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_Cholesky(PC pc) 3199b54502bSHong Zhang { 3209b54502bSHong Zhang PetscErrorCode ierr; 3219b54502bSHong Zhang PC_Cholesky *dir; 3229b54502bSHong Zhang 3239b54502bSHong Zhang PetscFunctionBegin; 324b00a9115SJed Brown ierr = PetscNewLog(pc,&dir);CHKERRQ(ierr); 3259b54502bSHong Zhang 326075768bcSBarry Smith ((PC_Factor*)dir)->fact = 0; 3279b54502bSHong Zhang dir->inplace = PETSC_FALSE; 3282fa5cd67SKarl Rupp 329075768bcSBarry Smith ierr = MatFactorInfoInitialize(&((PC_Factor*)dir)->info);CHKERRQ(ierr); 3302fa5cd67SKarl Rupp 331879e8a4dSBarry Smith ((PC_Factor*)dir)->factortype = MAT_FACTOR_CHOLESKY; 332075768bcSBarry Smith ((PC_Factor*)dir)->info.fill = 5.0; 333f4db908eSBarry Smith ((PC_Factor*)dir)->info.shifttype = (PetscReal) MAT_SHIFT_NONE; 334915743fcSHong Zhang ((PC_Factor*)dir)->info.shiftamount = 0.0; 3350ed735ceSHong Zhang ((PC_Factor*)dir)->info.zeropivot = 100.0*PETSC_MACHINE_EPSILON; 336075768bcSBarry Smith ((PC_Factor*)dir)->info.pivotinblocks = 1.0; 3372fa5cd67SKarl Rupp 3389b54502bSHong Zhang dir->col = 0; 3399b54502bSHong Zhang dir->row = 0; 3402fa5cd67SKarl Rupp 34119fd82e9SBarry Smith ierr = PetscStrallocpy(MATORDERINGNATURAL,(char**)&((PC_Factor*)dir)->ordering);CHKERRQ(ierr); 3422692d6eeSBarry Smith ierr = PetscStrallocpy(MATSOLVERPETSC,&((PC_Factor*)dir)->solvertype);CHKERRQ(ierr); 3432fa5cd67SKarl Rupp 3449b54502bSHong Zhang dir->reusefill = PETSC_FALSE; 3459b54502bSHong Zhang dir->reuseordering = PETSC_FALSE; 3469b54502bSHong Zhang pc->data = (void*)dir; 3479b54502bSHong Zhang 3489b54502bSHong Zhang pc->ops->destroy = PCDestroy_Cholesky; 34969d2c0f9SBarry Smith pc->ops->reset = PCReset_Cholesky; 3509b54502bSHong Zhang pc->ops->apply = PCApply_Cholesky; 3519b54502bSHong Zhang pc->ops->applytranspose = PCApplyTranspose_Cholesky; 3529b54502bSHong Zhang pc->ops->setup = PCSetUp_Cholesky; 3539b54502bSHong Zhang pc->ops->setfromoptions = PCSetFromOptions_Cholesky; 3549b54502bSHong Zhang pc->ops->view = PCView_Cholesky; 3559b54502bSHong Zhang pc->ops->applyrichardson = 0; 35685317021SBarry Smith pc->ops->getfactoredmatrix = PCFactorGetMatrix_Factor; 3579b54502bSHong Zhang 358bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetUpMatSolverPackage_C",PCFactorSetUpMatSolverPackage_Factor);CHKERRQ(ierr); 359bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetMatSolverPackage_C",PCFactorSetMatSolverPackage_Factor);CHKERRQ(ierr); 360bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetMatSolverPackage_C",PCFactorGetMatSolverPackage_Factor);CHKERRQ(ierr); 361bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetZeroPivot_C",PCFactorSetZeroPivot_Factor);CHKERRQ(ierr); 362bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetShiftType_C",PCFactorSetShiftType_Factor);CHKERRQ(ierr); 363bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetShiftAmount_C",PCFactorSetShiftAmount_Factor);CHKERRQ(ierr); 364bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetFill_C",PCFactorSetFill_Factor);CHKERRQ(ierr); 365bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetUseInPlace_C",PCFactorSetUseInPlace_Cholesky);CHKERRQ(ierr); 3668e37d05fSBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetUseInPlace_C",PCFactorGetUseInPlace_Cholesky);CHKERRQ(ierr); 367bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetMatOrderingType_C",PCFactorSetMatOrderingType_Factor);CHKERRQ(ierr); 368bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetReuseOrdering_C",PCFactorSetReuseOrdering_Cholesky);CHKERRQ(ierr); 369bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetReuseFill_C",PCFactorSetReuseFill_Cholesky);CHKERRQ(ierr); 3709b54502bSHong Zhang PetscFunctionReturn(0); 3719b54502bSHong Zhang } 372