1dba47a55SKris Buschelman #define PETSCKSP_DLL 29b54502bSHong Zhang 39b54502bSHong Zhang /* 49b54502bSHong Zhang Defines a direct factorization preconditioner for any Mat implementation 59b54502bSHong Zhang Note: this need not be consided a preconditioner since it supplies 69b54502bSHong Zhang a direct solver. 79b54502bSHong Zhang */ 87c4f633dSBarry Smith #include "../src/ksp/pc/impls/factor/factor.h" /*I "petscpc.h" I*/ 99b54502bSHong Zhang 109b54502bSHong Zhang typedef struct { 11075768bcSBarry Smith PC_Factor hdr; 129b54502bSHong Zhang PetscReal actualfill; /* actual fill in factor */ 139b54502bSHong Zhang PetscTruth inplace; /* flag indicating in-place factorization */ 149b54502bSHong Zhang IS row,col; /* index sets used for reordering */ 159b54502bSHong Zhang PetscTruth reuseordering; /* reuses previous reordering computed */ 169b54502bSHong Zhang PetscTruth reusefill; /* reuse fill from previous Cholesky */ 179b54502bSHong Zhang } PC_Cholesky; 189b54502bSHong Zhang 199b54502bSHong Zhang EXTERN_C_BEGIN 209b54502bSHong Zhang #undef __FUNCT__ 212401956bSBarry Smith #define __FUNCT__ "PCFactorSetReuseOrdering_Cholesky" 222401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseOrdering_Cholesky(PC pc,PetscTruth flag) 239b54502bSHong Zhang { 249b54502bSHong Zhang PC_Cholesky *lu; 259b54502bSHong Zhang 269b54502bSHong Zhang PetscFunctionBegin; 279b54502bSHong Zhang lu = (PC_Cholesky*)pc->data; 289b54502bSHong Zhang lu->reuseordering = flag; 299b54502bSHong Zhang PetscFunctionReturn(0); 309b54502bSHong Zhang } 319b54502bSHong Zhang EXTERN_C_END 329b54502bSHong Zhang 339b54502bSHong Zhang EXTERN_C_BEGIN 349b54502bSHong Zhang #undef __FUNCT__ 352401956bSBarry Smith #define __FUNCT__ "PCFactorSetReuseFill_Cholesky" 362401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseFill_Cholesky(PC pc,PetscTruth flag) 379b54502bSHong Zhang { 389b54502bSHong Zhang PC_Cholesky *lu; 399b54502bSHong Zhang 409b54502bSHong Zhang PetscFunctionBegin; 419b54502bSHong Zhang lu = (PC_Cholesky*)pc->data; 429b54502bSHong Zhang lu->reusefill = flag; 439b54502bSHong Zhang PetscFunctionReturn(0); 449b54502bSHong Zhang } 459b54502bSHong Zhang EXTERN_C_END 469b54502bSHong Zhang 479b54502bSHong Zhang #undef __FUNCT__ 489b54502bSHong Zhang #define __FUNCT__ "PCSetFromOptions_Cholesky" 499b54502bSHong Zhang static PetscErrorCode PCSetFromOptions_Cholesky(PC pc) 509b54502bSHong Zhang { 519b54502bSHong Zhang PetscErrorCode ierr; 529b54502bSHong Zhang 539b54502bSHong Zhang PetscFunctionBegin; 549b54502bSHong Zhang ierr = PetscOptionsHead("Cholesky options");CHKERRQ(ierr); 558ff23777SHong Zhang ierr = PCSetFromOptions_Factor(pc);CHKERRQ(ierr); 569b54502bSHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 579b54502bSHong Zhang PetscFunctionReturn(0); 589b54502bSHong Zhang } 599b54502bSHong Zhang 609b54502bSHong Zhang #undef __FUNCT__ 619b54502bSHong Zhang #define __FUNCT__ "PCView_Cholesky" 629b54502bSHong Zhang static PetscErrorCode PCView_Cholesky(PC pc,PetscViewer viewer) 639b54502bSHong Zhang { 64145b9266SHong Zhang PC_Cholesky *chol = (PC_Cholesky*)pc->data; 659b54502bSHong Zhang PetscErrorCode ierr; 66914a5d51SHong Zhang PetscTruth iascii; 679b54502bSHong Zhang 689b54502bSHong Zhang PetscFunctionBegin; 699b54502bSHong Zhang ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 709b54502bSHong Zhang if (iascii) { 71914a5d51SHong Zhang if (chol->inplace) { 72914a5d51SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Cholesky: in-place factorization\n");CHKERRQ(ierr); 73914a5d51SHong Zhang } else { 74914a5d51SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Cholesky: out-of-place factorization\n");CHKERRQ(ierr); 75145b9266SHong Zhang } 769b54502bSHong Zhang 77145b9266SHong Zhang if (chol->reusefill) {ierr = PetscViewerASCIIPrintf(viewer," Reusing fill from past factorization\n");CHKERRQ(ierr);} 78145b9266SHong Zhang if (chol->reuseordering) {ierr = PetscViewerASCIIPrintf(viewer," Reusing reordering from past factorization\n");CHKERRQ(ierr);} 799b54502bSHong Zhang } 80914a5d51SHong Zhang ierr = PCView_Factor(pc,viewer);CHKERRQ(ierr); 819b54502bSHong Zhang PetscFunctionReturn(0); 829b54502bSHong Zhang } 839b54502bSHong Zhang 849b54502bSHong Zhang 859b54502bSHong Zhang #undef __FUNCT__ 869b54502bSHong Zhang #define __FUNCT__ "PCSetUp_Cholesky" 879b54502bSHong Zhang static PetscErrorCode PCSetUp_Cholesky(PC pc) 889b54502bSHong Zhang { 899b54502bSHong Zhang PetscErrorCode ierr; 909b54502bSHong Zhang PetscTruth flg; 919b54502bSHong Zhang PC_Cholesky *dir = (PC_Cholesky*)pc->data; 929b54502bSHong Zhang 939b54502bSHong Zhang PetscFunctionBegin; 94075768bcSBarry Smith if (dir->reusefill && pc->setupcalled) ((PC_Factor*)dir)->info.fill = dir->actualfill; 959b54502bSHong Zhang 969b54502bSHong Zhang if (dir->inplace) { 979b54502bSHong Zhang if (dir->row && dir->col && (dir->row != dir->col)) { 989b54502bSHong Zhang ierr = ISDestroy(dir->row);CHKERRQ(ierr); 999b54502bSHong Zhang dir->row = 0; 1009b54502bSHong Zhang } 1019b54502bSHong Zhang if (dir->col) { 1029b54502bSHong Zhang ierr = ISDestroy(dir->col);CHKERRQ(ierr); 1039b54502bSHong Zhang dir->col = 0; 1049b54502bSHong Zhang } 105075768bcSBarry Smith ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 1069b54502bSHong Zhang if (dir->col && (dir->row != dir->col)) { /* only use row ordering for SBAIJ */ 1079b54502bSHong Zhang ierr = ISDestroy(dir->col);CHKERRQ(ierr); 1089b54502bSHong Zhang dir->col=0; 1099b54502bSHong Zhang } 110efee365bSSatish Balay if (dir->row) {ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr);} 111075768bcSBarry Smith ierr = MatCholeskyFactor(pc->pmat,dir->row,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 112075768bcSBarry Smith ((PC_Factor*)dir)->fact = pc->pmat; 1139b54502bSHong Zhang } else { 1149b54502bSHong Zhang MatInfo info; 1159b54502bSHong Zhang if (!pc->setupcalled) { 116075768bcSBarry Smith ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 1179bfd6278SHong Zhang /* check if dir->row == dir->col */ 1189bfd6278SHong Zhang ierr = ISEqual(dir->row,dir->col,&flg);CHKERRQ(ierr); 1199bfd6278SHong Zhang if (!flg) SETERRQ(PETSC_ERR_ARG_INCOMP,"row and column permutations must equal"); 1209bfd6278SHong Zhang ierr = ISDestroy(dir->col);CHKERRQ(ierr); /* only pass one ordering into CholeskyFactor */ 1219b54502bSHong Zhang dir->col=0; 1229bfd6278SHong Zhang 12390d69ab7SBarry Smith flg = PETSC_FALSE; 12490d69ab7SBarry Smith ierr = PetscOptionsGetTruth(((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&flg,PETSC_NULL);CHKERRQ(ierr); 1259b54502bSHong Zhang if (flg) { 1269b54502bSHong Zhang PetscReal tol = 1.e-10; 1277adad957SLisandro Dalcin ierr = PetscOptionsGetReal(((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&tol,PETSC_NULL);CHKERRQ(ierr); 1289b54502bSHong Zhang ierr = MatReorderForNonzeroDiagonal(pc->pmat,tol,dir->row,dir->row);CHKERRQ(ierr); 1299b54502bSHong Zhang } 130efee365bSSatish Balay if (dir->row) {ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr);} 131075768bcSBarry Smith ierr = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_CHOLESKY,&((PC_Factor*)dir)->fact);CHKERRQ(ierr); 132075768bcSBarry Smith ierr = MatCholeskyFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 133075768bcSBarry Smith ierr = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 1349b54502bSHong Zhang dir->actualfill = info.fill_ratio_needed; 135075768bcSBarry Smith ierr = PetscLogObjectParent(pc,((PC_Factor*)dir)->fact);CHKERRQ(ierr); 1369b54502bSHong Zhang } else if (pc->flag != SAME_NONZERO_PATTERN) { 1379b54502bSHong Zhang if (!dir->reuseordering) { 1389b54502bSHong Zhang if (dir->row && dir->col && (dir->row != dir->col)) { 1399b54502bSHong Zhang ierr = ISDestroy(dir->row);CHKERRQ(ierr); 1409b54502bSHong Zhang dir->row = 0; 1419b54502bSHong Zhang } 1429b54502bSHong Zhang if (dir->col) { 1439b54502bSHong Zhang ierr = ISDestroy(dir->col);CHKERRQ(ierr); 1449b54502bSHong Zhang dir->col =0; 1459b54502bSHong Zhang } 146075768bcSBarry Smith ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 1479b54502bSHong Zhang if (dir->col && (dir->row != dir->col)) { /* only use row ordering for SBAIJ */ 1489b54502bSHong Zhang ierr = ISDestroy(dir->col);CHKERRQ(ierr); 1499b54502bSHong Zhang dir->col=0; 1509b54502bSHong Zhang } 15190d69ab7SBarry Smith flg = PETSC_FALSE; 15290d69ab7SBarry Smith ierr = PetscOptionsGetTruth(((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&flg,PETSC_NULL);CHKERRQ(ierr); 1539b54502bSHong Zhang if (flg) { 1549b54502bSHong Zhang PetscReal tol = 1.e-10; 1557adad957SLisandro Dalcin ierr = PetscOptionsGetReal(((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&tol,PETSC_NULL);CHKERRQ(ierr); 1569b54502bSHong Zhang ierr = MatReorderForNonzeroDiagonal(pc->pmat,tol,dir->row,dir->row);CHKERRQ(ierr); 1579b54502bSHong Zhang } 158efee365bSSatish Balay if (dir->row) {ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr);} 1599b54502bSHong Zhang } 160075768bcSBarry Smith ierr = MatDestroy(((PC_Factor*)dir)->fact);CHKERRQ(ierr); 161075768bcSBarry Smith ierr = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_CHOLESKY,&((PC_Factor*)dir)->fact);CHKERRQ(ierr); 162075768bcSBarry Smith ierr = MatCholeskyFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 163075768bcSBarry Smith ierr = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 1649b54502bSHong Zhang dir->actualfill = info.fill_ratio_needed; 165075768bcSBarry Smith ierr = PetscLogObjectParent(pc,((PC_Factor*)dir)->fact);CHKERRQ(ierr); 1669b54502bSHong Zhang } 167075768bcSBarry Smith ierr = MatCholeskyFactorNumeric(((PC_Factor*)dir)->fact,pc->pmat,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 1689b54502bSHong Zhang } 1699b54502bSHong Zhang PetscFunctionReturn(0); 1709b54502bSHong Zhang } 1719b54502bSHong Zhang 1729b54502bSHong Zhang #undef __FUNCT__ 1739b54502bSHong Zhang #define __FUNCT__ "PCDestroy_Cholesky" 1749b54502bSHong Zhang static PetscErrorCode PCDestroy_Cholesky(PC pc) 1759b54502bSHong Zhang { 1769b54502bSHong Zhang PC_Cholesky *dir = (PC_Cholesky*)pc->data; 1779b54502bSHong Zhang PetscErrorCode ierr; 1789b54502bSHong Zhang 1799b54502bSHong Zhang PetscFunctionBegin; 180075768bcSBarry Smith if (!dir->inplace && ((PC_Factor*)dir)->fact) {ierr = MatDestroy(((PC_Factor*)dir)->fact);CHKERRQ(ierr);} 1819b54502bSHong Zhang if (dir->row) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);} 1829b54502bSHong Zhang if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);} 183075768bcSBarry Smith ierr = PetscStrfree(((PC_Factor*)dir)->ordering);CHKERRQ(ierr); 184075768bcSBarry Smith ierr = PetscStrfree(((PC_Factor*)dir)->solvertype);CHKERRQ(ierr); 1859b54502bSHong Zhang ierr = PetscFree(dir);CHKERRQ(ierr); 1869b54502bSHong Zhang PetscFunctionReturn(0); 1879b54502bSHong Zhang } 1889b54502bSHong Zhang 1899b54502bSHong Zhang #undef __FUNCT__ 1909b54502bSHong Zhang #define __FUNCT__ "PCApply_Cholesky" 1919b54502bSHong Zhang static PetscErrorCode PCApply_Cholesky(PC pc,Vec x,Vec y) 1929b54502bSHong Zhang { 1939b54502bSHong Zhang PC_Cholesky *dir = (PC_Cholesky*)pc->data; 1949b54502bSHong Zhang PetscErrorCode ierr; 1959b54502bSHong Zhang 1969b54502bSHong Zhang PetscFunctionBegin; 1979b54502bSHong Zhang if (dir->inplace) {ierr = MatSolve(pc->pmat,x,y);CHKERRQ(ierr);} 198075768bcSBarry Smith else {ierr = MatSolve(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr);} 1999b54502bSHong Zhang PetscFunctionReturn(0); 2009b54502bSHong Zhang } 2019b54502bSHong Zhang 2029b54502bSHong Zhang #undef __FUNCT__ 2039b54502bSHong Zhang #define __FUNCT__ "PCApplyTranspose_Cholesky" 2049b54502bSHong Zhang static PetscErrorCode PCApplyTranspose_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; 2109b54502bSHong Zhang if (dir->inplace) {ierr = MatSolveTranspose(pc->pmat,x,y);CHKERRQ(ierr);} 211075768bcSBarry Smith else {ierr = MatSolveTranspose(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr);} 2129b54502bSHong Zhang PetscFunctionReturn(0); 2139b54502bSHong Zhang } 2149b54502bSHong Zhang 2159b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/ 2169b54502bSHong Zhang 2179b54502bSHong Zhang EXTERN_C_BEGIN 2189b54502bSHong Zhang #undef __FUNCT__ 2192401956bSBarry Smith #define __FUNCT__ "PCFactorSetUseInPlace_Cholesky" 2202401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseInPlace_Cholesky(PC pc) 2219b54502bSHong Zhang { 2229b54502bSHong Zhang PC_Cholesky *dir; 2239b54502bSHong Zhang 2249b54502bSHong Zhang PetscFunctionBegin; 2259b54502bSHong Zhang dir = (PC_Cholesky*)pc->data; 2269b54502bSHong Zhang dir->inplace = PETSC_TRUE; 2279b54502bSHong Zhang PetscFunctionReturn(0); 2289b54502bSHong Zhang } 2299b54502bSHong Zhang EXTERN_C_END 2309b54502bSHong Zhang 2319b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/ 2329b54502bSHong Zhang 2339b54502bSHong Zhang #undef __FUNCT__ 2342401956bSBarry Smith #define __FUNCT__ "PCFactorSetReuseOrdering" 2359b54502bSHong Zhang /*@ 2362401956bSBarry Smith PCFactorSetReuseOrdering - When similar matrices are factored, this 2379b54502bSHong Zhang causes the ordering computed in the first factor to be used for all 2389b54502bSHong Zhang following factors. 2399b54502bSHong Zhang 2409b54502bSHong Zhang Collective on PC 2419b54502bSHong Zhang 2429b54502bSHong Zhang Input Parameters: 2439b54502bSHong Zhang + pc - the preconditioner context 2449b54502bSHong Zhang - flag - PETSC_TRUE to reuse else PETSC_FALSE 2459b54502bSHong Zhang 2469b54502bSHong Zhang Options Database Key: 2472401956bSBarry Smith . -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering() 2489b54502bSHong Zhang 2499b54502bSHong Zhang Level: intermediate 2509b54502bSHong Zhang 2519b54502bSHong Zhang .keywords: PC, levels, reordering, factorization, incomplete, LU 2529b54502bSHong Zhang 2532401956bSBarry Smith .seealso: PCFactorSetReuseFill() 2549b54502bSHong Zhang @*/ 2552401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseOrdering(PC pc,PetscTruth flag) 2569b54502bSHong Zhang { 2579b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscTruth); 2589b54502bSHong Zhang 2599b54502bSHong Zhang PetscFunctionBegin; 260*0700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 2612401956bSBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetReuseOrdering_C",(void (**)(void))&f);CHKERRQ(ierr); 2629b54502bSHong Zhang if (f) { 2639b54502bSHong Zhang ierr = (*f)(pc,flag);CHKERRQ(ierr); 2649b54502bSHong Zhang } 2659b54502bSHong Zhang PetscFunctionReturn(0); 2669b54502bSHong Zhang } 2679b54502bSHong Zhang 2689b54502bSHong Zhang 2699b54502bSHong Zhang /*MC 27096fc60bcSBarry Smith PCCHOLESKY - Uses a direct solver, based on Cholesky factorization, as a preconditioner 2719b54502bSHong Zhang 2729b54502bSHong Zhang Options Database Keys: 2732401956bSBarry Smith + -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering() 274c7393fdbSBarry Smith . -pc_factor_mat_solver_package - Actives PCFactorSetMatSolverPackage() to choose the direct solver, like spooles 2752401956bSBarry Smith . -pc_factor_reuse_fill - Activates PCFactorSetReuseFill() 27655ba2a51SBarry Smith . -pc_factor_fill <fill> - Sets fill amount 2772401956bSBarry Smith . -pc_factor_in_place - Activates in-place factorization 278145b9266SHong Zhang - -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine 2799b54502bSHong Zhang 2809b54502bSHong Zhang Notes: Not all options work for all matrix formats 2819b54502bSHong Zhang 2829b54502bSHong Zhang Level: beginner 2839b54502bSHong Zhang 2849b54502bSHong Zhang Concepts: Cholesky factorization, direct solver 2859b54502bSHong Zhang 2869b54502bSHong Zhang Notes: Usually this will compute an "exact" solution in one iteration and does 2879b54502bSHong Zhang not need a Krylov method (i.e. you can use -ksp_type preonly, or 2889b54502bSHong Zhang KSPSetType(ksp,KSPPREONLY) for the Krylov method 2899b54502bSHong Zhang 2909b54502bSHong Zhang .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 291a4fd02acSBarry Smith PCILU, PCLU, PCICC, PCFactorSetReuseOrdering(), PCFactorSetReuseFill(), PCFactorGetMatrix(), 292145b9266SHong Zhang PCFactorSetFill(), PCFactorSetShiftNonzero(), PCFactorSetShiftType(), PCFactorSetShiftAmount() 2938ff23777SHong Zhang PCFactorSetUseInPlace(), PCFactorSetMatOrderingType() 2949b54502bSHong Zhang 2959b54502bSHong Zhang M*/ 2969b54502bSHong Zhang 2979b54502bSHong Zhang EXTERN_C_BEGIN 2989b54502bSHong Zhang #undef __FUNCT__ 2999b54502bSHong Zhang #define __FUNCT__ "PCCreate_Cholesky" 300dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_Cholesky(PC pc) 3019b54502bSHong Zhang { 3029b54502bSHong Zhang PetscErrorCode ierr; 3039b54502bSHong Zhang PC_Cholesky *dir; 3049b54502bSHong Zhang 3059b54502bSHong Zhang PetscFunctionBegin; 30638f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_Cholesky,&dir);CHKERRQ(ierr); 3079b54502bSHong Zhang 308075768bcSBarry Smith ((PC_Factor*)dir)->fact = 0; 3099b54502bSHong Zhang dir->inplace = PETSC_FALSE; 310075768bcSBarry Smith ierr = MatFactorInfoInitialize(&((PC_Factor*)dir)->info);CHKERRQ(ierr); 311879e8a4dSBarry Smith ((PC_Factor*)dir)->factortype = MAT_FACTOR_CHOLESKY; 312075768bcSBarry Smith ((PC_Factor*)dir)->info.fill = 5.0; 313f4db908eSBarry Smith ((PC_Factor*)dir)->info.shifttype = (PetscReal) MAT_SHIFT_NONE; 314915743fcSHong Zhang ((PC_Factor*)dir)->info.shiftamount = 0.0; 315075768bcSBarry Smith ((PC_Factor*)dir)->info.pivotinblocks = 1.0; 3169b54502bSHong Zhang dir->col = 0; 3179b54502bSHong Zhang dir->row = 0; 318075768bcSBarry Smith ierr = PetscStrallocpy(MATORDERING_NATURAL,&((PC_Factor*)dir)->ordering);CHKERRQ(ierr); 319075768bcSBarry Smith ierr = PetscStrallocpy(MAT_SOLVER_PETSC,&((PC_Factor*)dir)->solvertype);CHKERRQ(ierr); 3209b54502bSHong Zhang dir->reusefill = PETSC_FALSE; 3219b54502bSHong Zhang dir->reuseordering = PETSC_FALSE; 3229b54502bSHong Zhang pc->data = (void*)dir; 3239b54502bSHong Zhang 3249b54502bSHong Zhang pc->ops->destroy = PCDestroy_Cholesky; 3259b54502bSHong Zhang pc->ops->apply = PCApply_Cholesky; 3269b54502bSHong Zhang pc->ops->applytranspose = PCApplyTranspose_Cholesky; 3279b54502bSHong Zhang pc->ops->setup = PCSetUp_Cholesky; 3289b54502bSHong Zhang pc->ops->setfromoptions = PCSetFromOptions_Cholesky; 3299b54502bSHong Zhang pc->ops->view = PCView_Cholesky; 3309b54502bSHong Zhang pc->ops->applyrichardson = 0; 33185317021SBarry Smith pc->ops->getfactoredmatrix = PCFactorGetMatrix_Factor; 3329b54502bSHong Zhang 33385317021SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatSolverPackage_C","PCFactorSetMatSolverPackage_Factor", 33485317021SBarry Smith PCFactorSetMatSolverPackage_Factor);CHKERRQ(ierr); 3357112b564SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorGetMatSolverPackage_C","PCFactorGetMatSolverPackage_Factor", 3367112b564SBarry Smith PCFactorGetMatSolverPackage_Factor);CHKERRQ(ierr); 33785317021SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_Factor", 33885317021SBarry Smith PCFactorSetZeroPivot_Factor);CHKERRQ(ierr); 339d90ac83dSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftType_C","PCFactorSetShiftType_Factor", 340d90ac83dSHong Zhang PCFactorSetShiftType_Factor);CHKERRQ(ierr); 341d90ac83dSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftAmount_C","PCFactorSetShiftAmount_Factor", 342d90ac83dSHong Zhang PCFactorSetShiftAmount_Factor);CHKERRQ(ierr); 34385317021SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetFill_C","PCFactorSetFill_Factor", 34485317021SBarry Smith PCFactorSetFill_Factor);CHKERRQ(ierr); 3452401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUseInPlace_C","PCFactorSetUseInPlace_Cholesky", 3462401956bSBarry Smith PCFactorSetUseInPlace_Cholesky);CHKERRQ(ierr); 34785317021SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatOrderingType_C","PCFactorSetMatOrderingType_Factor", 34885317021SBarry Smith PCFactorSetMatOrderingType_Factor);CHKERRQ(ierr); 3492401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseOrdering_C","PCFactorSetReuseOrdering_Cholesky", 3502401956bSBarry Smith PCFactorSetReuseOrdering_Cholesky);CHKERRQ(ierr); 3512401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseFill_C","PCFactorSetReuseFill_Cholesky", 3522401956bSBarry Smith PCFactorSetReuseFill_Cholesky);CHKERRQ(ierr); 3539b54502bSHong Zhang PetscFunctionReturn(0); 3549b54502bSHong Zhang } 3559b54502bSHong Zhang EXTERN_C_END 356