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 #undef __FUNCT__ 789b54502bSHong Zhang #define __FUNCT__ "PCSetUp_Cholesky" 799b54502bSHong Zhang static PetscErrorCode PCSetUp_Cholesky(PC pc) 809b54502bSHong Zhang { 819b54502bSHong Zhang PetscErrorCode ierr; 82ace3abfcSBarry Smith PetscBool flg; 839b54502bSHong Zhang PC_Cholesky *dir = (PC_Cholesky*)pc->data; 8400c67f3bSHong Zhang const MatSolverPackage stype; 8500e125f8SBarry Smith MatFactorError err; 869b54502bSHong Zhang 879b54502bSHong Zhang PetscFunctionBegin; 88075768bcSBarry Smith if (dir->reusefill && pc->setupcalled) ((PC_Factor*)dir)->info.fill = dir->actualfill; 899b54502bSHong Zhang 9084d44b13SHong Zhang ierr = MatSetErrorIfFailure(pc->pmat,pc->erroriffailure);CHKERRQ(ierr); 919b54502bSHong Zhang if (dir->inplace) { 929b54502bSHong Zhang if (dir->row && dir->col && (dir->row != dir->col)) { 93fcfd50ebSBarry Smith ierr = ISDestroy(&dir->row);CHKERRQ(ierr); 949b54502bSHong Zhang } 95fcfd50ebSBarry Smith ierr = ISDestroy(&dir->col);CHKERRQ(ierr); 96075768bcSBarry Smith ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 979b54502bSHong Zhang if (dir->col && (dir->row != dir->col)) { /* only use row ordering for SBAIJ */ 98fcfd50ebSBarry Smith ierr = ISDestroy(&dir->col);CHKERRQ(ierr); 999b54502bSHong Zhang } 1003bb1ff40SBarry Smith if (dir->row) {ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);CHKERRQ(ierr);} 101075768bcSBarry Smith ierr = MatCholeskyFactor(pc->pmat,dir->row,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 10200e125f8SBarry Smith ierr = MatFactorGetError(pc->pmat,&err);CHKERRQ(ierr); 10300e125f8SBarry Smith if (err) { /* Factor() fails */ 10400e125f8SBarry Smith pc->failedreason = (PCFailedReason)err; 1056baea169SHong Zhang PetscFunctionReturn(0); 1066baea169SHong Zhang } 1072fa5cd67SKarl Rupp 108075768bcSBarry Smith ((PC_Factor*)dir)->fact = pc->pmat; 1099b54502bSHong Zhang } else { 1109b54502bSHong Zhang MatInfo info; 11100e125f8SBarry Smith 1129b54502bSHong Zhang if (!pc->setupcalled) { 113075768bcSBarry Smith ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 1149bfd6278SHong Zhang /* check if dir->row == dir->col */ 1159bfd6278SHong Zhang ierr = ISEqual(dir->row,dir->col,&flg);CHKERRQ(ierr); 116e32f2f54SBarry Smith if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"row and column permutations must equal"); 117fcfd50ebSBarry Smith ierr = ISDestroy(&dir->col);CHKERRQ(ierr); /* only pass one ordering into CholeskyFactor */ 1189bfd6278SHong Zhang 11990d69ab7SBarry Smith flg = PETSC_FALSE; 120c5929fdfSBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&flg,NULL);CHKERRQ(ierr); 1219b54502bSHong Zhang if (flg) { 1229b54502bSHong Zhang PetscReal tol = 1.e-10; 123c5929fdfSBarry Smith ierr = PetscOptionsGetReal(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&tol,NULL);CHKERRQ(ierr); 1249b54502bSHong Zhang ierr = MatReorderForNonzeroDiagonal(pc->pmat,tol,dir->row,dir->row);CHKERRQ(ierr); 1259b54502bSHong Zhang } 1263bb1ff40SBarry Smith if (dir->row) {ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);CHKERRQ(ierr);} 127d09a07f4SBarry Smith if (!((PC_Factor*)dir)->fact) { 128075768bcSBarry Smith ierr = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_CHOLESKY,&((PC_Factor*)dir)->fact);CHKERRQ(ierr); 129a1f19f5aSHong Zhang } 130075768bcSBarry Smith ierr = MatCholeskyFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 131075768bcSBarry Smith ierr = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 1329b54502bSHong Zhang dir->actualfill = info.fill_ratio_needed; 1333bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact);CHKERRQ(ierr); 1349b54502bSHong Zhang } else if (pc->flag != SAME_NONZERO_PATTERN) { 1359b54502bSHong Zhang if (!dir->reuseordering) { 136fcfd50ebSBarry Smith ierr = ISDestroy(&dir->row);CHKERRQ(ierr); 137075768bcSBarry Smith ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 1389d61c4f5SHong Zhang ierr = ISDestroy(&dir->col);CHKERRQ(ierr); /* only use dir->row ordering in CholeskyFactor */ 1399d61c4f5SHong Zhang 14090d69ab7SBarry Smith flg = PETSC_FALSE; 141c5929fdfSBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&flg,NULL);CHKERRQ(ierr); 1429b54502bSHong Zhang if (flg) { 1439b54502bSHong Zhang PetscReal tol = 1.e-10; 144c5929fdfSBarry Smith ierr = PetscOptionsGetReal(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&tol,NULL);CHKERRQ(ierr); 1459b54502bSHong Zhang ierr = MatReorderForNonzeroDiagonal(pc->pmat,tol,dir->row,dir->row);CHKERRQ(ierr); 1469b54502bSHong Zhang } 1473bb1ff40SBarry Smith if (dir->row) {ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);CHKERRQ(ierr);} 1489b54502bSHong Zhang } 1496bf464f9SBarry Smith ierr = MatDestroy(&((PC_Factor*)dir)->fact);CHKERRQ(ierr); 150075768bcSBarry Smith ierr = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_CHOLESKY,&((PC_Factor*)dir)->fact);CHKERRQ(ierr); 151075768bcSBarry Smith ierr = MatCholeskyFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 152075768bcSBarry Smith ierr = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 1539b54502bSHong Zhang dir->actualfill = info.fill_ratio_needed; 1543bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact);CHKERRQ(ierr); 15504545d6dSBarry Smith } else { 156*b8b68cfdSBarry Smith ierr = MatFactorGetError(((PC_Factor*)dir)->fact,&err);CHKERRQ(ierr); 157*b8b68cfdSBarry Smith if (err == PC_FACTOR_NUMERIC_ZEROPIVOT) { 158*b8b68cfdSBarry Smith ierr = MatFactorClearError(((PC_Factor*)dir)->fact);CHKERRQ(ierr); 159*b8b68cfdSBarry Smith pc->failedreason = PC_NOERROR; 16004545d6dSBarry Smith } 1619b54502bSHong Zhang } 16200e125f8SBarry Smith ierr = MatFactorGetError(((PC_Factor*)dir)->fact,&err);CHKERRQ(ierr); 16300e125f8SBarry Smith if (err) { /* FactorSymbolic() fails */ 16400e125f8SBarry Smith pc->failedreason = (PCFailedReason)err; 1656baea169SHong Zhang PetscFunctionReturn(0); 1666baea169SHong Zhang } 1676baea169SHong Zhang 168075768bcSBarry Smith ierr = MatCholeskyFactorNumeric(((PC_Factor*)dir)->fact,pc->pmat,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 16900e125f8SBarry Smith ierr = MatFactorGetError(((PC_Factor*)dir)->fact,&err);CHKERRQ(ierr); 17000e125f8SBarry Smith if (err) { /* FactorNumeric() fails */ 17100e125f8SBarry Smith pc->failedreason = (PCFailedReason)err; 1726baea169SHong Zhang } 1739b54502bSHong Zhang } 17400c67f3bSHong Zhang 17500c67f3bSHong Zhang ierr = PCFactorGetMatSolverPackage(pc,&stype);CHKERRQ(ierr); 17600c67f3bSHong Zhang if (!stype) { 17700e125f8SBarry Smith const MatSolverPackage solverpackage; 17800e125f8SBarry Smith ierr = MatFactorGetSolverPackage(((PC_Factor*)dir)->fact,&solverpackage);CHKERRQ(ierr); 17900e125f8SBarry Smith ierr = PCFactorSetMatSolverPackage(pc,solverpackage);CHKERRQ(ierr); 18000c67f3bSHong Zhang } 1819b54502bSHong Zhang PetscFunctionReturn(0); 1829b54502bSHong Zhang } 1839b54502bSHong Zhang 1849b54502bSHong Zhang #undef __FUNCT__ 18569d2c0f9SBarry Smith #define __FUNCT__ "PCReset_Cholesky" 18669d2c0f9SBarry Smith static PetscErrorCode PCReset_Cholesky(PC pc) 1879b54502bSHong Zhang { 1889b54502bSHong Zhang PC_Cholesky *dir = (PC_Cholesky*)pc->data; 1899b54502bSHong Zhang PetscErrorCode ierr; 1909b54502bSHong Zhang 1919b54502bSHong Zhang PetscFunctionBegin; 1926bf464f9SBarry Smith if (!dir->inplace && ((PC_Factor*)dir)->fact) {ierr = MatDestroy(&((PC_Factor*)dir)->fact);CHKERRQ(ierr);} 193fcfd50ebSBarry Smith ierr = ISDestroy(&dir->row);CHKERRQ(ierr); 194fcfd50ebSBarry Smith ierr = ISDestroy(&dir->col);CHKERRQ(ierr); 19569d2c0f9SBarry Smith PetscFunctionReturn(0); 19669d2c0f9SBarry Smith } 19769d2c0f9SBarry Smith 19869d2c0f9SBarry Smith #undef __FUNCT__ 19969d2c0f9SBarry Smith #define __FUNCT__ "PCDestroy_Cholesky" 20069d2c0f9SBarry Smith static PetscErrorCode PCDestroy_Cholesky(PC pc) 20169d2c0f9SBarry Smith { 20269d2c0f9SBarry Smith PC_Cholesky *dir = (PC_Cholesky*)pc->data; 20369d2c0f9SBarry Smith PetscErrorCode ierr; 20469d2c0f9SBarry Smith 20569d2c0f9SBarry Smith PetscFunctionBegin; 20669d2c0f9SBarry Smith ierr = PCReset_Cholesky(pc);CHKERRQ(ierr); 207503cfb0cSBarry Smith ierr = PetscFree(((PC_Factor*)dir)->ordering);CHKERRQ(ierr); 208503cfb0cSBarry Smith ierr = PetscFree(((PC_Factor*)dir)->solvertype);CHKERRQ(ierr); 209c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 2109b54502bSHong Zhang PetscFunctionReturn(0); 2119b54502bSHong Zhang } 2129b54502bSHong Zhang 2139b54502bSHong Zhang #undef __FUNCT__ 2149b54502bSHong Zhang #define __FUNCT__ "PCApply_Cholesky" 2159b54502bSHong Zhang static PetscErrorCode PCApply_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 = MatSolve(pc->pmat,x,y);CHKERRQ(ierr); 2232fa5cd67SKarl Rupp } else { 2242fa5cd67SKarl Rupp ierr = MatSolve(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr); 2252fa5cd67SKarl Rupp } 2269b54502bSHong Zhang PetscFunctionReturn(0); 2279b54502bSHong Zhang } 2289b54502bSHong Zhang 2299b54502bSHong Zhang #undef __FUNCT__ 2309b54502bSHong Zhang #define __FUNCT__ "PCApplyTranspose_Cholesky" 2319b54502bSHong Zhang static PetscErrorCode PCApplyTranspose_Cholesky(PC pc,Vec x,Vec y) 2329b54502bSHong Zhang { 2339b54502bSHong Zhang PC_Cholesky *dir = (PC_Cholesky*)pc->data; 2349b54502bSHong Zhang PetscErrorCode ierr; 2359b54502bSHong Zhang 2369b54502bSHong Zhang PetscFunctionBegin; 2372fa5cd67SKarl Rupp if (dir->inplace) { 2382fa5cd67SKarl Rupp ierr = MatSolveTranspose(pc->pmat,x,y);CHKERRQ(ierr); 2392fa5cd67SKarl Rupp } else { 2402fa5cd67SKarl Rupp ierr = MatSolveTranspose(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr); 2412fa5cd67SKarl Rupp } 2429b54502bSHong Zhang PetscFunctionReturn(0); 2439b54502bSHong Zhang } 2449b54502bSHong Zhang 2459b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/ 2469b54502bSHong Zhang 2479b54502bSHong Zhang #undef __FUNCT__ 2482401956bSBarry Smith #define __FUNCT__ "PCFactorSetUseInPlace_Cholesky" 2498e37d05fSBarry Smith static PetscErrorCode PCFactorSetUseInPlace_Cholesky(PC pc,PetscBool flg) 2509b54502bSHong Zhang { 2518e37d05fSBarry Smith PC_Cholesky *dir = (PC_Cholesky*)pc->data; 2529b54502bSHong Zhang 2539b54502bSHong Zhang PetscFunctionBegin; 2548e37d05fSBarry Smith dir->inplace = flg; 2558e37d05fSBarry Smith PetscFunctionReturn(0); 2568e37d05fSBarry Smith } 2578e37d05fSBarry Smith 2588e37d05fSBarry Smith #undef __FUNCT__ 2598e37d05fSBarry Smith #define __FUNCT__ "PCFactorGetUseInPlace_Cholesky" 2608e37d05fSBarry Smith static PetscErrorCode PCFactorGetUseInPlace_Cholesky(PC pc,PetscBool *flg) 2618e37d05fSBarry Smith { 2628e37d05fSBarry Smith PC_Cholesky *dir = (PC_Cholesky*)pc->data; 2638e37d05fSBarry Smith 2648e37d05fSBarry Smith PetscFunctionBegin; 2658e37d05fSBarry Smith *flg = dir->inplace; 2669b54502bSHong Zhang PetscFunctionReturn(0); 2679b54502bSHong Zhang } 2689b54502bSHong Zhang 2699b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/ 2709b54502bSHong Zhang 2719b54502bSHong Zhang #undef __FUNCT__ 2722401956bSBarry Smith #define __FUNCT__ "PCFactorSetReuseOrdering" 2739b54502bSHong Zhang /*@ 2742401956bSBarry Smith PCFactorSetReuseOrdering - When similar matrices are factored, this 2759b54502bSHong Zhang causes the ordering computed in the first factor to be used for all 2769b54502bSHong Zhang following factors. 2779b54502bSHong Zhang 278ad4df100SBarry Smith Logically Collective on PC 2799b54502bSHong Zhang 2809b54502bSHong Zhang Input Parameters: 2819b54502bSHong Zhang + pc - the preconditioner context 2829b54502bSHong Zhang - flag - PETSC_TRUE to reuse else PETSC_FALSE 2839b54502bSHong Zhang 2849b54502bSHong Zhang Options Database Key: 2852401956bSBarry Smith . -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering() 2869b54502bSHong Zhang 2879b54502bSHong Zhang Level: intermediate 2889b54502bSHong Zhang 2899b54502bSHong Zhang .keywords: PC, levels, reordering, factorization, incomplete, LU 2909b54502bSHong Zhang 2912401956bSBarry Smith .seealso: PCFactorSetReuseFill() 2929b54502bSHong Zhang @*/ 2937087cfbeSBarry Smith PetscErrorCode PCFactorSetReuseOrdering(PC pc,PetscBool flag) 2949b54502bSHong Zhang { 2954ac538c5SBarry Smith PetscErrorCode ierr; 2969b54502bSHong Zhang 2979b54502bSHong Zhang PetscFunctionBegin; 2980700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 299acfcf0e5SJed Brown PetscValidLogicalCollectiveBool(pc,flag,2); 3004ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFactorSetReuseOrdering_C",(PC,PetscBool),(pc,flag));CHKERRQ(ierr); 3019b54502bSHong Zhang PetscFunctionReturn(0); 3029b54502bSHong Zhang } 3039b54502bSHong Zhang 3049b54502bSHong Zhang /*MC 30596fc60bcSBarry Smith PCCHOLESKY - Uses a direct solver, based on Cholesky factorization, as a preconditioner 3069b54502bSHong Zhang 3079b54502bSHong Zhang Options Database Keys: 3082401956bSBarry Smith + -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering() 309f60c3dc2SHong Zhang . -pc_factor_mat_solver_package - Actives PCFactorSetMatSolverPackage() to choose the direct solver, like superlu 3102401956bSBarry Smith . -pc_factor_reuse_fill - Activates PCFactorSetReuseFill() 31155ba2a51SBarry Smith . -pc_factor_fill <fill> - Sets fill amount 3122401956bSBarry Smith . -pc_factor_in_place - Activates in-place factorization 313145b9266SHong Zhang - -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine 3149b54502bSHong Zhang 3159b54502bSHong Zhang Notes: Not all options work for all matrix formats 3169b54502bSHong Zhang 3179b54502bSHong Zhang Level: beginner 3189b54502bSHong Zhang 3199b54502bSHong Zhang Concepts: Cholesky factorization, direct solver 3209b54502bSHong Zhang 3219b54502bSHong Zhang Notes: Usually this will compute an "exact" solution in one iteration and does 3229b54502bSHong Zhang not need a Krylov method (i.e. you can use -ksp_type preonly, or 3239b54502bSHong Zhang KSPSetType(ksp,KSPPREONLY) for the Krylov method 3249b54502bSHong Zhang 3259b54502bSHong Zhang .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 326a4fd02acSBarry Smith PCILU, PCLU, PCICC, PCFactorSetReuseOrdering(), PCFactorSetReuseFill(), PCFactorGetMatrix(), 327145b9266SHong Zhang PCFactorSetFill(), PCFactorSetShiftNonzero(), PCFactorSetShiftType(), PCFactorSetShiftAmount() 3288e37d05fSBarry Smith PCFactorSetUseInPlace(), PCFactorGetUseInPlace(), PCFactorSetMatOrderingType() 3299b54502bSHong Zhang 3309b54502bSHong Zhang M*/ 3319b54502bSHong Zhang 3329b54502bSHong Zhang #undef __FUNCT__ 3339b54502bSHong Zhang #define __FUNCT__ "PCCreate_Cholesky" 3348cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_Cholesky(PC pc) 3359b54502bSHong Zhang { 3369b54502bSHong Zhang PetscErrorCode ierr; 3379b54502bSHong Zhang PC_Cholesky *dir; 3389b54502bSHong Zhang 3399b54502bSHong Zhang PetscFunctionBegin; 340b00a9115SJed Brown ierr = PetscNewLog(pc,&dir);CHKERRQ(ierr); 3419b54502bSHong Zhang 342075768bcSBarry Smith ((PC_Factor*)dir)->fact = 0; 3439b54502bSHong Zhang dir->inplace = PETSC_FALSE; 3442fa5cd67SKarl Rupp 345075768bcSBarry Smith ierr = MatFactorInfoInitialize(&((PC_Factor*)dir)->info);CHKERRQ(ierr); 3462fa5cd67SKarl Rupp 347879e8a4dSBarry Smith ((PC_Factor*)dir)->factortype = MAT_FACTOR_CHOLESKY; 348075768bcSBarry Smith ((PC_Factor*)dir)->info.fill = 5.0; 349f4db908eSBarry Smith ((PC_Factor*)dir)->info.shifttype = (PetscReal) MAT_SHIFT_NONE; 350915743fcSHong Zhang ((PC_Factor*)dir)->info.shiftamount = 0.0; 3510ed735ceSHong Zhang ((PC_Factor*)dir)->info.zeropivot = 100.0*PETSC_MACHINE_EPSILON; 352075768bcSBarry Smith ((PC_Factor*)dir)->info.pivotinblocks = 1.0; 3532fa5cd67SKarl Rupp 3549b54502bSHong Zhang dir->col = 0; 3559b54502bSHong Zhang dir->row = 0; 3562fa5cd67SKarl Rupp 35719fd82e9SBarry Smith ierr = PetscStrallocpy(MATORDERINGNATURAL,(char**)&((PC_Factor*)dir)->ordering);CHKERRQ(ierr); 3589b54502bSHong Zhang dir->reusefill = PETSC_FALSE; 3599b54502bSHong Zhang dir->reuseordering = PETSC_FALSE; 3609b54502bSHong Zhang pc->data = (void*)dir; 3619b54502bSHong Zhang 3629b54502bSHong Zhang pc->ops->destroy = PCDestroy_Cholesky; 36369d2c0f9SBarry Smith pc->ops->reset = PCReset_Cholesky; 3649b54502bSHong Zhang pc->ops->apply = PCApply_Cholesky; 3659b54502bSHong Zhang pc->ops->applytranspose = PCApplyTranspose_Cholesky; 3669b54502bSHong Zhang pc->ops->setup = PCSetUp_Cholesky; 3679b54502bSHong Zhang pc->ops->setfromoptions = PCSetFromOptions_Cholesky; 3689b54502bSHong Zhang pc->ops->view = PCView_Cholesky; 3699b54502bSHong Zhang pc->ops->applyrichardson = 0; 37085317021SBarry Smith pc->ops->getfactoredmatrix = PCFactorGetMatrix_Factor; 3719b54502bSHong Zhang 372bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetUpMatSolverPackage_C",PCFactorSetUpMatSolverPackage_Factor);CHKERRQ(ierr); 373bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetMatSolverPackage_C",PCFactorSetMatSolverPackage_Factor);CHKERRQ(ierr); 374bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetMatSolverPackage_C",PCFactorGetMatSolverPackage_Factor);CHKERRQ(ierr); 375bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetZeroPivot_C",PCFactorSetZeroPivot_Factor);CHKERRQ(ierr); 376bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetShiftType_C",PCFactorSetShiftType_Factor);CHKERRQ(ierr); 377bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetShiftAmount_C",PCFactorSetShiftAmount_Factor);CHKERRQ(ierr); 378bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetFill_C",PCFactorSetFill_Factor);CHKERRQ(ierr); 379bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetUseInPlace_C",PCFactorSetUseInPlace_Cholesky);CHKERRQ(ierr); 3808e37d05fSBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetUseInPlace_C",PCFactorGetUseInPlace_Cholesky);CHKERRQ(ierr); 381bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetMatOrderingType_C",PCFactorSetMatOrderingType_Factor);CHKERRQ(ierr); 382bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetReuseOrdering_C",PCFactorSetReuseOrdering_Cholesky);CHKERRQ(ierr); 383bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetReuseFill_C",PCFactorSetReuseFill_Cholesky);CHKERRQ(ierr); 3849b54502bSHong Zhang PetscFunctionReturn(0); 3859b54502bSHong Zhang } 386