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; 85*00e125f8SBarry 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); 102*00e125f8SBarry Smith ierr = MatFactorGetError(pc->pmat,&err);CHKERRQ(ierr); 103*00e125f8SBarry Smith if (err) { /* Factor() fails */ 104*00e125f8SBarry 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; 111*00e125f8SBarry 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); 1559b54502bSHong Zhang } 156*00e125f8SBarry Smith ierr = MatFactorGetError(((PC_Factor*)dir)->fact,&err);CHKERRQ(ierr); 157*00e125f8SBarry Smith if (err) { /* FactorSymbolic() fails */ 158*00e125f8SBarry Smith pc->failedreason = (PCFailedReason)err; 1596baea169SHong Zhang PetscFunctionReturn(0); 1606baea169SHong Zhang } 1616baea169SHong Zhang 162075768bcSBarry Smith ierr = MatCholeskyFactorNumeric(((PC_Factor*)dir)->fact,pc->pmat,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 163*00e125f8SBarry Smith ierr = MatFactorGetError(((PC_Factor*)dir)->fact,&err);CHKERRQ(ierr); 164*00e125f8SBarry Smith if (err) { /* FactorNumeric() fails */ 165*00e125f8SBarry Smith pc->failedreason = (PCFailedReason)err; 1666baea169SHong Zhang } 1679b54502bSHong Zhang } 16800c67f3bSHong Zhang 16900c67f3bSHong Zhang ierr = PCFactorGetMatSolverPackage(pc,&stype);CHKERRQ(ierr); 17000c67f3bSHong Zhang if (!stype) { 171*00e125f8SBarry Smith const MatSolverPackage solverpackage; 172*00e125f8SBarry Smith ierr = MatFactorGetSolverPackage(((PC_Factor*)dir)->fact,&solverpackage);CHKERRQ(ierr); 173*00e125f8SBarry Smith ierr = PCFactorSetMatSolverPackage(pc,solverpackage);CHKERRQ(ierr); 17400c67f3bSHong Zhang } 1759b54502bSHong Zhang PetscFunctionReturn(0); 1769b54502bSHong Zhang } 1779b54502bSHong Zhang 1789b54502bSHong Zhang #undef __FUNCT__ 17969d2c0f9SBarry Smith #define __FUNCT__ "PCReset_Cholesky" 18069d2c0f9SBarry Smith static PetscErrorCode PCReset_Cholesky(PC pc) 1819b54502bSHong Zhang { 1829b54502bSHong Zhang PC_Cholesky *dir = (PC_Cholesky*)pc->data; 1839b54502bSHong Zhang PetscErrorCode ierr; 1849b54502bSHong Zhang 1859b54502bSHong Zhang PetscFunctionBegin; 1866bf464f9SBarry Smith if (!dir->inplace && ((PC_Factor*)dir)->fact) {ierr = MatDestroy(&((PC_Factor*)dir)->fact);CHKERRQ(ierr);} 187fcfd50ebSBarry Smith ierr = ISDestroy(&dir->row);CHKERRQ(ierr); 188fcfd50ebSBarry Smith ierr = ISDestroy(&dir->col);CHKERRQ(ierr); 18969d2c0f9SBarry Smith PetscFunctionReturn(0); 19069d2c0f9SBarry Smith } 19169d2c0f9SBarry Smith 19269d2c0f9SBarry Smith #undef __FUNCT__ 19369d2c0f9SBarry Smith #define __FUNCT__ "PCDestroy_Cholesky" 19469d2c0f9SBarry Smith static PetscErrorCode PCDestroy_Cholesky(PC pc) 19569d2c0f9SBarry Smith { 19669d2c0f9SBarry Smith PC_Cholesky *dir = (PC_Cholesky*)pc->data; 19769d2c0f9SBarry Smith PetscErrorCode ierr; 19869d2c0f9SBarry Smith 19969d2c0f9SBarry Smith PetscFunctionBegin; 20069d2c0f9SBarry Smith ierr = PCReset_Cholesky(pc);CHKERRQ(ierr); 201503cfb0cSBarry Smith ierr = PetscFree(((PC_Factor*)dir)->ordering);CHKERRQ(ierr); 202503cfb0cSBarry Smith ierr = PetscFree(((PC_Factor*)dir)->solvertype);CHKERRQ(ierr); 203c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 2049b54502bSHong Zhang PetscFunctionReturn(0); 2059b54502bSHong Zhang } 2069b54502bSHong Zhang 2079b54502bSHong Zhang #undef __FUNCT__ 2089b54502bSHong Zhang #define __FUNCT__ "PCApply_Cholesky" 2099b54502bSHong Zhang static PetscErrorCode PCApply_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; 2152fa5cd67SKarl Rupp if (dir->inplace) { 2162fa5cd67SKarl Rupp ierr = MatSolve(pc->pmat,x,y);CHKERRQ(ierr); 2172fa5cd67SKarl Rupp } else { 2182fa5cd67SKarl Rupp ierr = MatSolve(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr); 2192fa5cd67SKarl Rupp } 2209b54502bSHong Zhang PetscFunctionReturn(0); 2219b54502bSHong Zhang } 2229b54502bSHong Zhang 2239b54502bSHong Zhang #undef __FUNCT__ 2249b54502bSHong Zhang #define __FUNCT__ "PCApplyTranspose_Cholesky" 2259b54502bSHong Zhang static PetscErrorCode PCApplyTranspose_Cholesky(PC pc,Vec x,Vec y) 2269b54502bSHong Zhang { 2279b54502bSHong Zhang PC_Cholesky *dir = (PC_Cholesky*)pc->data; 2289b54502bSHong Zhang PetscErrorCode ierr; 2299b54502bSHong Zhang 2309b54502bSHong Zhang PetscFunctionBegin; 2312fa5cd67SKarl Rupp if (dir->inplace) { 2322fa5cd67SKarl Rupp ierr = MatSolveTranspose(pc->pmat,x,y);CHKERRQ(ierr); 2332fa5cd67SKarl Rupp } else { 2342fa5cd67SKarl Rupp ierr = MatSolveTranspose(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr); 2352fa5cd67SKarl Rupp } 2369b54502bSHong Zhang PetscFunctionReturn(0); 2379b54502bSHong Zhang } 2389b54502bSHong Zhang 2399b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/ 2409b54502bSHong Zhang 2419b54502bSHong Zhang #undef __FUNCT__ 2422401956bSBarry Smith #define __FUNCT__ "PCFactorSetUseInPlace_Cholesky" 2438e37d05fSBarry Smith static PetscErrorCode PCFactorSetUseInPlace_Cholesky(PC pc,PetscBool flg) 2449b54502bSHong Zhang { 2458e37d05fSBarry Smith PC_Cholesky *dir = (PC_Cholesky*)pc->data; 2469b54502bSHong Zhang 2479b54502bSHong Zhang PetscFunctionBegin; 2488e37d05fSBarry Smith dir->inplace = flg; 2498e37d05fSBarry Smith PetscFunctionReturn(0); 2508e37d05fSBarry Smith } 2518e37d05fSBarry Smith 2528e37d05fSBarry Smith #undef __FUNCT__ 2538e37d05fSBarry Smith #define __FUNCT__ "PCFactorGetUseInPlace_Cholesky" 2548e37d05fSBarry Smith static PetscErrorCode PCFactorGetUseInPlace_Cholesky(PC pc,PetscBool *flg) 2558e37d05fSBarry Smith { 2568e37d05fSBarry Smith PC_Cholesky *dir = (PC_Cholesky*)pc->data; 2578e37d05fSBarry Smith 2588e37d05fSBarry Smith PetscFunctionBegin; 2598e37d05fSBarry Smith *flg = dir->inplace; 2609b54502bSHong Zhang PetscFunctionReturn(0); 2619b54502bSHong Zhang } 2629b54502bSHong Zhang 2639b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/ 2649b54502bSHong Zhang 2659b54502bSHong Zhang #undef __FUNCT__ 2662401956bSBarry Smith #define __FUNCT__ "PCFactorSetReuseOrdering" 2679b54502bSHong Zhang /*@ 2682401956bSBarry Smith PCFactorSetReuseOrdering - When similar matrices are factored, this 2699b54502bSHong Zhang causes the ordering computed in the first factor to be used for all 2709b54502bSHong Zhang following factors. 2719b54502bSHong Zhang 272ad4df100SBarry Smith Logically Collective on PC 2739b54502bSHong Zhang 2749b54502bSHong Zhang Input Parameters: 2759b54502bSHong Zhang + pc - the preconditioner context 2769b54502bSHong Zhang - flag - PETSC_TRUE to reuse else PETSC_FALSE 2779b54502bSHong Zhang 2789b54502bSHong Zhang Options Database Key: 2792401956bSBarry Smith . -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering() 2809b54502bSHong Zhang 2819b54502bSHong Zhang Level: intermediate 2829b54502bSHong Zhang 2839b54502bSHong Zhang .keywords: PC, levels, reordering, factorization, incomplete, LU 2849b54502bSHong Zhang 2852401956bSBarry Smith .seealso: PCFactorSetReuseFill() 2869b54502bSHong Zhang @*/ 2877087cfbeSBarry Smith PetscErrorCode PCFactorSetReuseOrdering(PC pc,PetscBool flag) 2889b54502bSHong Zhang { 2894ac538c5SBarry Smith PetscErrorCode ierr; 2909b54502bSHong Zhang 2919b54502bSHong Zhang PetscFunctionBegin; 2920700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 293acfcf0e5SJed Brown PetscValidLogicalCollectiveBool(pc,flag,2); 2944ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFactorSetReuseOrdering_C",(PC,PetscBool),(pc,flag));CHKERRQ(ierr); 2959b54502bSHong Zhang PetscFunctionReturn(0); 2969b54502bSHong Zhang } 2979b54502bSHong Zhang 2989b54502bSHong Zhang /*MC 29996fc60bcSBarry Smith PCCHOLESKY - Uses a direct solver, based on Cholesky factorization, as a preconditioner 3009b54502bSHong Zhang 3019b54502bSHong Zhang Options Database Keys: 3022401956bSBarry Smith + -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering() 303f60c3dc2SHong Zhang . -pc_factor_mat_solver_package - Actives PCFactorSetMatSolverPackage() to choose the direct solver, like superlu 3042401956bSBarry Smith . -pc_factor_reuse_fill - Activates PCFactorSetReuseFill() 30555ba2a51SBarry Smith . -pc_factor_fill <fill> - Sets fill amount 3062401956bSBarry Smith . -pc_factor_in_place - Activates in-place factorization 307145b9266SHong Zhang - -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine 3089b54502bSHong Zhang 3099b54502bSHong Zhang Notes: Not all options work for all matrix formats 3109b54502bSHong Zhang 3119b54502bSHong Zhang Level: beginner 3129b54502bSHong Zhang 3139b54502bSHong Zhang Concepts: Cholesky factorization, direct solver 3149b54502bSHong Zhang 3159b54502bSHong Zhang Notes: Usually this will compute an "exact" solution in one iteration and does 3169b54502bSHong Zhang not need a Krylov method (i.e. you can use -ksp_type preonly, or 3179b54502bSHong Zhang KSPSetType(ksp,KSPPREONLY) for the Krylov method 3189b54502bSHong Zhang 3199b54502bSHong Zhang .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 320a4fd02acSBarry Smith PCILU, PCLU, PCICC, PCFactorSetReuseOrdering(), PCFactorSetReuseFill(), PCFactorGetMatrix(), 321145b9266SHong Zhang PCFactorSetFill(), PCFactorSetShiftNonzero(), PCFactorSetShiftType(), PCFactorSetShiftAmount() 3228e37d05fSBarry Smith PCFactorSetUseInPlace(), PCFactorGetUseInPlace(), PCFactorSetMatOrderingType() 3239b54502bSHong Zhang 3249b54502bSHong Zhang M*/ 3259b54502bSHong Zhang 3269b54502bSHong Zhang #undef __FUNCT__ 3279b54502bSHong Zhang #define __FUNCT__ "PCCreate_Cholesky" 3288cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_Cholesky(PC pc) 3299b54502bSHong Zhang { 3309b54502bSHong Zhang PetscErrorCode ierr; 3319b54502bSHong Zhang PC_Cholesky *dir; 3329b54502bSHong Zhang 3339b54502bSHong Zhang PetscFunctionBegin; 334b00a9115SJed Brown ierr = PetscNewLog(pc,&dir);CHKERRQ(ierr); 3359b54502bSHong Zhang 336075768bcSBarry Smith ((PC_Factor*)dir)->fact = 0; 3379b54502bSHong Zhang dir->inplace = PETSC_FALSE; 3382fa5cd67SKarl Rupp 339075768bcSBarry Smith ierr = MatFactorInfoInitialize(&((PC_Factor*)dir)->info);CHKERRQ(ierr); 3402fa5cd67SKarl Rupp 341879e8a4dSBarry Smith ((PC_Factor*)dir)->factortype = MAT_FACTOR_CHOLESKY; 342075768bcSBarry Smith ((PC_Factor*)dir)->info.fill = 5.0; 343f4db908eSBarry Smith ((PC_Factor*)dir)->info.shifttype = (PetscReal) MAT_SHIFT_NONE; 344915743fcSHong Zhang ((PC_Factor*)dir)->info.shiftamount = 0.0; 3450ed735ceSHong Zhang ((PC_Factor*)dir)->info.zeropivot = 100.0*PETSC_MACHINE_EPSILON; 346075768bcSBarry Smith ((PC_Factor*)dir)->info.pivotinblocks = 1.0; 3472fa5cd67SKarl Rupp 3489b54502bSHong Zhang dir->col = 0; 3499b54502bSHong Zhang dir->row = 0; 3502fa5cd67SKarl Rupp 35119fd82e9SBarry Smith ierr = PetscStrallocpy(MATORDERINGNATURAL,(char**)&((PC_Factor*)dir)->ordering);CHKERRQ(ierr); 3529b54502bSHong Zhang dir->reusefill = PETSC_FALSE; 3539b54502bSHong Zhang dir->reuseordering = PETSC_FALSE; 3549b54502bSHong Zhang pc->data = (void*)dir; 3559b54502bSHong Zhang 3569b54502bSHong Zhang pc->ops->destroy = PCDestroy_Cholesky; 35769d2c0f9SBarry Smith pc->ops->reset = PCReset_Cholesky; 3589b54502bSHong Zhang pc->ops->apply = PCApply_Cholesky; 3599b54502bSHong Zhang pc->ops->applytranspose = PCApplyTranspose_Cholesky; 3609b54502bSHong Zhang pc->ops->setup = PCSetUp_Cholesky; 3619b54502bSHong Zhang pc->ops->setfromoptions = PCSetFromOptions_Cholesky; 3629b54502bSHong Zhang pc->ops->view = PCView_Cholesky; 3639b54502bSHong Zhang pc->ops->applyrichardson = 0; 36485317021SBarry Smith pc->ops->getfactoredmatrix = PCFactorGetMatrix_Factor; 3659b54502bSHong Zhang 366bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetUpMatSolverPackage_C",PCFactorSetUpMatSolverPackage_Factor);CHKERRQ(ierr); 367bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetMatSolverPackage_C",PCFactorSetMatSolverPackage_Factor);CHKERRQ(ierr); 368bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetMatSolverPackage_C",PCFactorGetMatSolverPackage_Factor);CHKERRQ(ierr); 369bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetZeroPivot_C",PCFactorSetZeroPivot_Factor);CHKERRQ(ierr); 370bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetShiftType_C",PCFactorSetShiftType_Factor);CHKERRQ(ierr); 371bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetShiftAmount_C",PCFactorSetShiftAmount_Factor);CHKERRQ(ierr); 372bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetFill_C",PCFactorSetFill_Factor);CHKERRQ(ierr); 373bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetUseInPlace_C",PCFactorSetUseInPlace_Cholesky);CHKERRQ(ierr); 3748e37d05fSBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetUseInPlace_C",PCFactorGetUseInPlace_Cholesky);CHKERRQ(ierr); 375bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetMatOrderingType_C",PCFactorSetMatOrderingType_Factor);CHKERRQ(ierr); 376bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetReuseOrdering_C",PCFactorSetReuseOrdering_Cholesky);CHKERRQ(ierr); 377bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetReuseFill_C",PCFactorSetReuseFill_Cholesky);CHKERRQ(ierr); 3789b54502bSHong Zhang PetscFunctionReturn(0); 3799b54502bSHong Zhang } 380