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; 84*00c67f3bSHong Zhang const MatSolverPackage stype; 859b54502bSHong Zhang 869b54502bSHong Zhang PetscFunctionBegin; 87075768bcSBarry Smith if (dir->reusefill && pc->setupcalled) ((PC_Factor*)dir)->info.fill = dir->actualfill; 889b54502bSHong Zhang 8984d44b13SHong Zhang ierr = MatSetErrorIfFailure(pc->pmat,pc->erroriffailure);CHKERRQ(ierr); 909b54502bSHong Zhang if (dir->inplace) { 919b54502bSHong Zhang if (dir->row && dir->col && (dir->row != dir->col)) { 92fcfd50ebSBarry Smith ierr = ISDestroy(&dir->row);CHKERRQ(ierr); 939b54502bSHong Zhang } 94fcfd50ebSBarry Smith ierr = ISDestroy(&dir->col);CHKERRQ(ierr); 95075768bcSBarry Smith ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 969b54502bSHong Zhang if (dir->col && (dir->row != dir->col)) { /* only use row ordering for SBAIJ */ 97fcfd50ebSBarry Smith ierr = ISDestroy(&dir->col);CHKERRQ(ierr); 989b54502bSHong Zhang } 993bb1ff40SBarry Smith if (dir->row) {ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);CHKERRQ(ierr);} 100075768bcSBarry Smith ierr = MatCholeskyFactor(pc->pmat,dir->row,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 101680c5173SHong Zhang if (pc->pmat->errortype) { /* Factor() fails */ 102680c5173SHong Zhang pc->failedreason = (PCFailedReason)pc->pmat->errortype; 1036baea169SHong Zhang PetscFunctionReturn(0); 1046baea169SHong Zhang } 1052fa5cd67SKarl Rupp 106075768bcSBarry Smith ((PC_Factor*)dir)->fact = pc->pmat; 1079b54502bSHong Zhang } else { 1089b54502bSHong Zhang MatInfo info; 109680c5173SHong Zhang Mat F; 1109b54502bSHong Zhang if (!pc->setupcalled) { 111075768bcSBarry Smith ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 1129bfd6278SHong Zhang /* check if dir->row == dir->col */ 1139bfd6278SHong Zhang ierr = ISEqual(dir->row,dir->col,&flg);CHKERRQ(ierr); 114e32f2f54SBarry Smith if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"row and column permutations must equal"); 115fcfd50ebSBarry Smith ierr = ISDestroy(&dir->col);CHKERRQ(ierr); /* only pass one ordering into CholeskyFactor */ 1169bfd6278SHong Zhang 11790d69ab7SBarry Smith flg = PETSC_FALSE; 118c5929fdfSBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&flg,NULL);CHKERRQ(ierr); 1199b54502bSHong Zhang if (flg) { 1209b54502bSHong Zhang PetscReal tol = 1.e-10; 121c5929fdfSBarry Smith ierr = PetscOptionsGetReal(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&tol,NULL);CHKERRQ(ierr); 1229b54502bSHong Zhang ierr = MatReorderForNonzeroDiagonal(pc->pmat,tol,dir->row,dir->row);CHKERRQ(ierr); 1239b54502bSHong Zhang } 1243bb1ff40SBarry Smith if (dir->row) {ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);CHKERRQ(ierr);} 125d09a07f4SBarry Smith if (!((PC_Factor*)dir)->fact) { 126075768bcSBarry Smith ierr = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_CHOLESKY,&((PC_Factor*)dir)->fact);CHKERRQ(ierr); 127a1f19f5aSHong Zhang } 128075768bcSBarry Smith ierr = MatCholeskyFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 129075768bcSBarry Smith ierr = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 1309b54502bSHong Zhang dir->actualfill = info.fill_ratio_needed; 1313bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact);CHKERRQ(ierr); 1329b54502bSHong Zhang } else if (pc->flag != SAME_NONZERO_PATTERN) { 1339b54502bSHong Zhang if (!dir->reuseordering) { 134fcfd50ebSBarry Smith ierr = ISDestroy(&dir->row);CHKERRQ(ierr); 135075768bcSBarry Smith ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 1369d61c4f5SHong Zhang ierr = ISDestroy(&dir->col);CHKERRQ(ierr); /* only use dir->row ordering in CholeskyFactor */ 1379d61c4f5SHong Zhang 13890d69ab7SBarry Smith flg = PETSC_FALSE; 139c5929fdfSBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&flg,NULL);CHKERRQ(ierr); 1409b54502bSHong Zhang if (flg) { 1419b54502bSHong Zhang PetscReal tol = 1.e-10; 142c5929fdfSBarry Smith ierr = PetscOptionsGetReal(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&tol,NULL);CHKERRQ(ierr); 1439b54502bSHong Zhang ierr = MatReorderForNonzeroDiagonal(pc->pmat,tol,dir->row,dir->row);CHKERRQ(ierr); 1449b54502bSHong Zhang } 1453bb1ff40SBarry Smith if (dir->row) {ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);CHKERRQ(ierr);} 1469b54502bSHong Zhang } 1476bf464f9SBarry Smith ierr = MatDestroy(&((PC_Factor*)dir)->fact);CHKERRQ(ierr); 148075768bcSBarry Smith ierr = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_CHOLESKY,&((PC_Factor*)dir)->fact);CHKERRQ(ierr); 149075768bcSBarry Smith ierr = MatCholeskyFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 150075768bcSBarry Smith ierr = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 1519b54502bSHong Zhang dir->actualfill = info.fill_ratio_needed; 1523bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact);CHKERRQ(ierr); 1539b54502bSHong Zhang } 154680c5173SHong Zhang F = ((PC_Factor*)dir)->fact; 155680c5173SHong Zhang if (F->errortype) { /* FactorSymbolic() fails */ 156680c5173SHong Zhang pc->failedreason = (PCFailedReason)F->errortype; 1576baea169SHong Zhang PetscFunctionReturn(0); 1586baea169SHong Zhang } 1596baea169SHong Zhang 160075768bcSBarry Smith ierr = MatCholeskyFactorNumeric(((PC_Factor*)dir)->fact,pc->pmat,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 161680c5173SHong Zhang if (F->errortype) { /* FactorNumeric() fails */ 162680c5173SHong Zhang pc->failedreason = (PCFailedReason)F->errortype; 1636baea169SHong Zhang } 1649b54502bSHong Zhang } 165*00c67f3bSHong Zhang 166*00c67f3bSHong Zhang ierr = PCFactorGetMatSolverPackage(pc,&stype);CHKERRQ(ierr); 167*00c67f3bSHong Zhang if (!stype) { 168*00c67f3bSHong Zhang ierr = PCFactorSetMatSolverPackage(pc,((PC_Factor*)dir)->fact->solvertype);CHKERRQ(ierr); 169*00c67f3bSHong Zhang } 1709b54502bSHong Zhang PetscFunctionReturn(0); 1719b54502bSHong Zhang } 1729b54502bSHong Zhang 1739b54502bSHong Zhang #undef __FUNCT__ 17469d2c0f9SBarry Smith #define __FUNCT__ "PCReset_Cholesky" 17569d2c0f9SBarry Smith static PetscErrorCode PCReset_Cholesky(PC pc) 1769b54502bSHong Zhang { 1779b54502bSHong Zhang PC_Cholesky *dir = (PC_Cholesky*)pc->data; 1789b54502bSHong Zhang PetscErrorCode ierr; 1799b54502bSHong Zhang 1809b54502bSHong Zhang PetscFunctionBegin; 1816bf464f9SBarry Smith if (!dir->inplace && ((PC_Factor*)dir)->fact) {ierr = MatDestroy(&((PC_Factor*)dir)->fact);CHKERRQ(ierr);} 182fcfd50ebSBarry Smith ierr = ISDestroy(&dir->row);CHKERRQ(ierr); 183fcfd50ebSBarry Smith ierr = ISDestroy(&dir->col);CHKERRQ(ierr); 18469d2c0f9SBarry Smith PetscFunctionReturn(0); 18569d2c0f9SBarry Smith } 18669d2c0f9SBarry Smith 18769d2c0f9SBarry Smith #undef __FUNCT__ 18869d2c0f9SBarry Smith #define __FUNCT__ "PCDestroy_Cholesky" 18969d2c0f9SBarry Smith static PetscErrorCode PCDestroy_Cholesky(PC pc) 19069d2c0f9SBarry Smith { 19169d2c0f9SBarry Smith PC_Cholesky *dir = (PC_Cholesky*)pc->data; 19269d2c0f9SBarry Smith PetscErrorCode ierr; 19369d2c0f9SBarry Smith 19469d2c0f9SBarry Smith PetscFunctionBegin; 19569d2c0f9SBarry Smith ierr = PCReset_Cholesky(pc);CHKERRQ(ierr); 196503cfb0cSBarry Smith ierr = PetscFree(((PC_Factor*)dir)->ordering);CHKERRQ(ierr); 197503cfb0cSBarry Smith ierr = PetscFree(((PC_Factor*)dir)->solvertype);CHKERRQ(ierr); 198c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 1999b54502bSHong Zhang PetscFunctionReturn(0); 2009b54502bSHong Zhang } 2019b54502bSHong Zhang 2029b54502bSHong Zhang #undef __FUNCT__ 2039b54502bSHong Zhang #define __FUNCT__ "PCApply_Cholesky" 2049b54502bSHong Zhang static PetscErrorCode PCApply_Cholesky(PC pc,Vec x,Vec y) 2059b54502bSHong Zhang { 2069b54502bSHong Zhang PC_Cholesky *dir = (PC_Cholesky*)pc->data; 2079b54502bSHong Zhang PetscErrorCode ierr; 2089b54502bSHong Zhang 2099b54502bSHong Zhang PetscFunctionBegin; 2102fa5cd67SKarl Rupp if (dir->inplace) { 2112fa5cd67SKarl Rupp ierr = MatSolve(pc->pmat,x,y);CHKERRQ(ierr); 2122fa5cd67SKarl Rupp } else { 2132fa5cd67SKarl Rupp ierr = MatSolve(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr); 2142fa5cd67SKarl Rupp } 2159b54502bSHong Zhang PetscFunctionReturn(0); 2169b54502bSHong Zhang } 2179b54502bSHong Zhang 2189b54502bSHong Zhang #undef __FUNCT__ 2199b54502bSHong Zhang #define __FUNCT__ "PCApplyTranspose_Cholesky" 2209b54502bSHong Zhang static PetscErrorCode PCApplyTranspose_Cholesky(PC pc,Vec x,Vec y) 2219b54502bSHong Zhang { 2229b54502bSHong Zhang PC_Cholesky *dir = (PC_Cholesky*)pc->data; 2239b54502bSHong Zhang PetscErrorCode ierr; 2249b54502bSHong Zhang 2259b54502bSHong Zhang PetscFunctionBegin; 2262fa5cd67SKarl Rupp if (dir->inplace) { 2272fa5cd67SKarl Rupp ierr = MatSolveTranspose(pc->pmat,x,y);CHKERRQ(ierr); 2282fa5cd67SKarl Rupp } else { 2292fa5cd67SKarl Rupp ierr = MatSolveTranspose(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr); 2302fa5cd67SKarl Rupp } 2319b54502bSHong Zhang PetscFunctionReturn(0); 2329b54502bSHong Zhang } 2339b54502bSHong Zhang 2349b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/ 2359b54502bSHong Zhang 2369b54502bSHong Zhang #undef __FUNCT__ 2372401956bSBarry Smith #define __FUNCT__ "PCFactorSetUseInPlace_Cholesky" 2388e37d05fSBarry Smith static PetscErrorCode PCFactorSetUseInPlace_Cholesky(PC pc,PetscBool flg) 2399b54502bSHong Zhang { 2408e37d05fSBarry Smith PC_Cholesky *dir = (PC_Cholesky*)pc->data; 2419b54502bSHong Zhang 2429b54502bSHong Zhang PetscFunctionBegin; 2438e37d05fSBarry Smith dir->inplace = flg; 2448e37d05fSBarry Smith PetscFunctionReturn(0); 2458e37d05fSBarry Smith } 2468e37d05fSBarry Smith 2478e37d05fSBarry Smith #undef __FUNCT__ 2488e37d05fSBarry Smith #define __FUNCT__ "PCFactorGetUseInPlace_Cholesky" 2498e37d05fSBarry Smith static PetscErrorCode PCFactorGetUseInPlace_Cholesky(PC pc,PetscBool *flg) 2508e37d05fSBarry Smith { 2518e37d05fSBarry Smith PC_Cholesky *dir = (PC_Cholesky*)pc->data; 2528e37d05fSBarry Smith 2538e37d05fSBarry Smith PetscFunctionBegin; 2548e37d05fSBarry Smith *flg = dir->inplace; 2559b54502bSHong Zhang PetscFunctionReturn(0); 2569b54502bSHong Zhang } 2579b54502bSHong Zhang 2589b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/ 2599b54502bSHong Zhang 2609b54502bSHong Zhang #undef __FUNCT__ 2612401956bSBarry Smith #define __FUNCT__ "PCFactorSetReuseOrdering" 2629b54502bSHong Zhang /*@ 2632401956bSBarry Smith PCFactorSetReuseOrdering - When similar matrices are factored, this 2649b54502bSHong Zhang causes the ordering computed in the first factor to be used for all 2659b54502bSHong Zhang following factors. 2669b54502bSHong Zhang 267ad4df100SBarry Smith Logically Collective on PC 2689b54502bSHong Zhang 2699b54502bSHong Zhang Input Parameters: 2709b54502bSHong Zhang + pc - the preconditioner context 2719b54502bSHong Zhang - flag - PETSC_TRUE to reuse else PETSC_FALSE 2729b54502bSHong Zhang 2739b54502bSHong Zhang Options Database Key: 2742401956bSBarry Smith . -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering() 2759b54502bSHong Zhang 2769b54502bSHong Zhang Level: intermediate 2779b54502bSHong Zhang 2789b54502bSHong Zhang .keywords: PC, levels, reordering, factorization, incomplete, LU 2799b54502bSHong Zhang 2802401956bSBarry Smith .seealso: PCFactorSetReuseFill() 2819b54502bSHong Zhang @*/ 2827087cfbeSBarry Smith PetscErrorCode PCFactorSetReuseOrdering(PC pc,PetscBool flag) 2839b54502bSHong Zhang { 2844ac538c5SBarry Smith PetscErrorCode ierr; 2859b54502bSHong Zhang 2869b54502bSHong Zhang PetscFunctionBegin; 2870700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 288acfcf0e5SJed Brown PetscValidLogicalCollectiveBool(pc,flag,2); 2894ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFactorSetReuseOrdering_C",(PC,PetscBool),(pc,flag));CHKERRQ(ierr); 2909b54502bSHong Zhang PetscFunctionReturn(0); 2919b54502bSHong Zhang } 2929b54502bSHong Zhang 2939b54502bSHong Zhang /*MC 29496fc60bcSBarry Smith PCCHOLESKY - Uses a direct solver, based on Cholesky factorization, as a preconditioner 2959b54502bSHong Zhang 2969b54502bSHong Zhang Options Database Keys: 2972401956bSBarry Smith + -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering() 298f60c3dc2SHong Zhang . -pc_factor_mat_solver_package - Actives PCFactorSetMatSolverPackage() to choose the direct solver, like superlu 2992401956bSBarry Smith . -pc_factor_reuse_fill - Activates PCFactorSetReuseFill() 30055ba2a51SBarry Smith . -pc_factor_fill <fill> - Sets fill amount 3012401956bSBarry Smith . -pc_factor_in_place - Activates in-place factorization 302145b9266SHong Zhang - -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine 3039b54502bSHong Zhang 3049b54502bSHong Zhang Notes: Not all options work for all matrix formats 3059b54502bSHong Zhang 3069b54502bSHong Zhang Level: beginner 3079b54502bSHong Zhang 3089b54502bSHong Zhang Concepts: Cholesky factorization, direct solver 3099b54502bSHong Zhang 3109b54502bSHong Zhang Notes: Usually this will compute an "exact" solution in one iteration and does 3119b54502bSHong Zhang not need a Krylov method (i.e. you can use -ksp_type preonly, or 3129b54502bSHong Zhang KSPSetType(ksp,KSPPREONLY) for the Krylov method 3139b54502bSHong Zhang 3149b54502bSHong Zhang .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 315a4fd02acSBarry Smith PCILU, PCLU, PCICC, PCFactorSetReuseOrdering(), PCFactorSetReuseFill(), PCFactorGetMatrix(), 316145b9266SHong Zhang PCFactorSetFill(), PCFactorSetShiftNonzero(), PCFactorSetShiftType(), PCFactorSetShiftAmount() 3178e37d05fSBarry Smith PCFactorSetUseInPlace(), PCFactorGetUseInPlace(), PCFactorSetMatOrderingType() 3189b54502bSHong Zhang 3199b54502bSHong Zhang M*/ 3209b54502bSHong Zhang 3219b54502bSHong Zhang #undef __FUNCT__ 3229b54502bSHong Zhang #define __FUNCT__ "PCCreate_Cholesky" 3238cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_Cholesky(PC pc) 3249b54502bSHong Zhang { 3259b54502bSHong Zhang PetscErrorCode ierr; 3269b54502bSHong Zhang PC_Cholesky *dir; 3279b54502bSHong Zhang 3289b54502bSHong Zhang PetscFunctionBegin; 329b00a9115SJed Brown ierr = PetscNewLog(pc,&dir);CHKERRQ(ierr); 3309b54502bSHong Zhang 331075768bcSBarry Smith ((PC_Factor*)dir)->fact = 0; 3329b54502bSHong Zhang dir->inplace = PETSC_FALSE; 3332fa5cd67SKarl Rupp 334075768bcSBarry Smith ierr = MatFactorInfoInitialize(&((PC_Factor*)dir)->info);CHKERRQ(ierr); 3352fa5cd67SKarl Rupp 336879e8a4dSBarry Smith ((PC_Factor*)dir)->factortype = MAT_FACTOR_CHOLESKY; 337075768bcSBarry Smith ((PC_Factor*)dir)->info.fill = 5.0; 338f4db908eSBarry Smith ((PC_Factor*)dir)->info.shifttype = (PetscReal) MAT_SHIFT_NONE; 339915743fcSHong Zhang ((PC_Factor*)dir)->info.shiftamount = 0.0; 3400ed735ceSHong Zhang ((PC_Factor*)dir)->info.zeropivot = 100.0*PETSC_MACHINE_EPSILON; 341075768bcSBarry Smith ((PC_Factor*)dir)->info.pivotinblocks = 1.0; 3422fa5cd67SKarl Rupp 3439b54502bSHong Zhang dir->col = 0; 3449b54502bSHong Zhang dir->row = 0; 3452fa5cd67SKarl Rupp 34619fd82e9SBarry Smith ierr = PetscStrallocpy(MATORDERINGNATURAL,(char**)&((PC_Factor*)dir)->ordering);CHKERRQ(ierr); 3479b54502bSHong Zhang dir->reusefill = PETSC_FALSE; 3489b54502bSHong Zhang dir->reuseordering = PETSC_FALSE; 3499b54502bSHong Zhang pc->data = (void*)dir; 3509b54502bSHong Zhang 3519b54502bSHong Zhang pc->ops->destroy = PCDestroy_Cholesky; 35269d2c0f9SBarry Smith pc->ops->reset = PCReset_Cholesky; 3539b54502bSHong Zhang pc->ops->apply = PCApply_Cholesky; 3549b54502bSHong Zhang pc->ops->applytranspose = PCApplyTranspose_Cholesky; 3559b54502bSHong Zhang pc->ops->setup = PCSetUp_Cholesky; 3569b54502bSHong Zhang pc->ops->setfromoptions = PCSetFromOptions_Cholesky; 3579b54502bSHong Zhang pc->ops->view = PCView_Cholesky; 3589b54502bSHong Zhang pc->ops->applyrichardson = 0; 35985317021SBarry Smith pc->ops->getfactoredmatrix = PCFactorGetMatrix_Factor; 3609b54502bSHong Zhang 361bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetUpMatSolverPackage_C",PCFactorSetUpMatSolverPackage_Factor);CHKERRQ(ierr); 362bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetMatSolverPackage_C",PCFactorSetMatSolverPackage_Factor);CHKERRQ(ierr); 363bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetMatSolverPackage_C",PCFactorGetMatSolverPackage_Factor);CHKERRQ(ierr); 364bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetZeroPivot_C",PCFactorSetZeroPivot_Factor);CHKERRQ(ierr); 365bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetShiftType_C",PCFactorSetShiftType_Factor);CHKERRQ(ierr); 366bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetShiftAmount_C",PCFactorSetShiftAmount_Factor);CHKERRQ(ierr); 367bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetFill_C",PCFactorSetFill_Factor);CHKERRQ(ierr); 368bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetUseInPlace_C",PCFactorSetUseInPlace_Cholesky);CHKERRQ(ierr); 3698e37d05fSBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetUseInPlace_C",PCFactorGetUseInPlace_Cholesky);CHKERRQ(ierr); 370bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetMatOrderingType_C",PCFactorSetMatOrderingType_Factor);CHKERRQ(ierr); 371bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetReuseOrdering_C",PCFactorSetReuseOrdering_Cholesky);CHKERRQ(ierr); 372bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetReuseFill_C",PCFactorSetReuseFill_Cholesky);CHKERRQ(ierr); 3739b54502bSHong Zhang PetscFunctionReturn(0); 3749b54502bSHong Zhang } 375