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 EXTERN_C_BEGIN 199b54502bSHong Zhang #undef __FUNCT__ 202401956bSBarry Smith #define __FUNCT__ "PCFactorSetReuseOrdering_Cholesky" 217087cfbeSBarry Smith PetscErrorCode PCFactorSetReuseOrdering_Cholesky(PC pc,PetscBool flag) 229b54502bSHong Zhang { 239b54502bSHong Zhang PC_Cholesky *lu; 249b54502bSHong Zhang 259b54502bSHong Zhang PetscFunctionBegin; 269b54502bSHong Zhang lu = (PC_Cholesky*)pc->data; 279b54502bSHong Zhang lu->reuseordering = flag; 289b54502bSHong Zhang PetscFunctionReturn(0); 299b54502bSHong Zhang } 309b54502bSHong Zhang EXTERN_C_END 319b54502bSHong Zhang 329b54502bSHong Zhang EXTERN_C_BEGIN 339b54502bSHong Zhang #undef __FUNCT__ 342401956bSBarry Smith #define __FUNCT__ "PCFactorSetReuseFill_Cholesky" 357087cfbeSBarry Smith PetscErrorCode PCFactorSetReuseFill_Cholesky(PC pc,PetscBool flag) 369b54502bSHong Zhang { 379b54502bSHong Zhang PC_Cholesky *lu; 389b54502bSHong Zhang 399b54502bSHong Zhang PetscFunctionBegin; 409b54502bSHong Zhang lu = (PC_Cholesky*)pc->data; 419b54502bSHong Zhang lu->reusefill = flag; 429b54502bSHong Zhang PetscFunctionReturn(0); 439b54502bSHong Zhang } 449b54502bSHong Zhang EXTERN_C_END 459b54502bSHong Zhang 469b54502bSHong Zhang #undef __FUNCT__ 479b54502bSHong Zhang #define __FUNCT__ "PCSetFromOptions_Cholesky" 489b54502bSHong Zhang static PetscErrorCode PCSetFromOptions_Cholesky(PC pc) 499b54502bSHong Zhang { 509b54502bSHong Zhang PetscErrorCode ierr; 519b54502bSHong Zhang 529b54502bSHong Zhang PetscFunctionBegin; 539b54502bSHong Zhang ierr = PetscOptionsHead("Cholesky options");CHKERRQ(ierr); 548ff23777SHong Zhang ierr = PCSetFromOptions_Factor(pc);CHKERRQ(ierr); 559b54502bSHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 569b54502bSHong Zhang PetscFunctionReturn(0); 579b54502bSHong Zhang } 589b54502bSHong Zhang 599b54502bSHong Zhang #undef __FUNCT__ 609b54502bSHong Zhang #define __FUNCT__ "PCView_Cholesky" 619b54502bSHong Zhang static PetscErrorCode PCView_Cholesky(PC pc,PetscViewer viewer) 629b54502bSHong Zhang { 63145b9266SHong Zhang PC_Cholesky *chol = (PC_Cholesky*)pc->data; 649b54502bSHong Zhang PetscErrorCode ierr; 65ace3abfcSBarry Smith PetscBool iascii; 669b54502bSHong Zhang 679b54502bSHong Zhang PetscFunctionBegin; 682692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 699b54502bSHong Zhang if (iascii) { 70914a5d51SHong Zhang if (chol->inplace) { 71914a5d51SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Cholesky: in-place factorization\n");CHKERRQ(ierr); 72914a5d51SHong Zhang } else { 73914a5d51SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Cholesky: out-of-place factorization\n");CHKERRQ(ierr); 74145b9266SHong Zhang } 759b54502bSHong Zhang 76145b9266SHong Zhang if (chol->reusefill) {ierr = PetscViewerASCIIPrintf(viewer," Reusing fill from past factorization\n");CHKERRQ(ierr);} 77145b9266SHong Zhang if (chol->reuseordering) {ierr = PetscViewerASCIIPrintf(viewer," Reusing reordering from past factorization\n");CHKERRQ(ierr);} 789b54502bSHong Zhang } 79914a5d51SHong Zhang ierr = PCView_Factor(pc,viewer);CHKERRQ(ierr); 809b54502bSHong Zhang PetscFunctionReturn(0); 819b54502bSHong Zhang } 829b54502bSHong Zhang 839b54502bSHong Zhang 849b54502bSHong Zhang #undef __FUNCT__ 859b54502bSHong Zhang #define __FUNCT__ "PCSetUp_Cholesky" 869b54502bSHong Zhang static PetscErrorCode PCSetUp_Cholesky(PC pc) 879b54502bSHong Zhang { 889b54502bSHong Zhang PetscErrorCode ierr; 89ace3abfcSBarry Smith PetscBool flg; 909b54502bSHong Zhang PC_Cholesky *dir = (PC_Cholesky*)pc->data; 919b54502bSHong Zhang 929b54502bSHong Zhang PetscFunctionBegin; 93075768bcSBarry Smith if (dir->reusefill && pc->setupcalled) ((PC_Factor*)dir)->info.fill = dir->actualfill; 949b54502bSHong Zhang 959b54502bSHong Zhang if (dir->inplace) { 969b54502bSHong Zhang if (dir->row && dir->col && (dir->row != dir->col)) { 979b54502bSHong Zhang ierr = ISDestroy(dir->row);CHKERRQ(ierr); 989b54502bSHong Zhang } 999b54502bSHong Zhang if (dir->col) { 1009b54502bSHong Zhang ierr = ISDestroy(dir->col);CHKERRQ(ierr); 1019b54502bSHong Zhang } 102075768bcSBarry Smith ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 1039b54502bSHong Zhang if (dir->col && (dir->row != dir->col)) { /* only use row ordering for SBAIJ */ 1049b54502bSHong Zhang ierr = ISDestroy(dir->col);CHKERRQ(ierr); 1059b54502bSHong Zhang } 106efee365bSSatish Balay if (dir->row) {ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr);} 107075768bcSBarry Smith ierr = MatCholeskyFactor(pc->pmat,dir->row,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 108075768bcSBarry Smith ((PC_Factor*)dir)->fact = pc->pmat; 1099b54502bSHong Zhang } else { 1109b54502bSHong Zhang MatInfo info; 1119b54502bSHong Zhang if (!pc->setupcalled) { 112075768bcSBarry Smith ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 1139bfd6278SHong Zhang /* check if dir->row == dir->col */ 1149bfd6278SHong Zhang ierr = ISEqual(dir->row,dir->col,&flg);CHKERRQ(ierr); 115e32f2f54SBarry Smith if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"row and column permutations must equal"); 1169bfd6278SHong Zhang ierr = ISDestroy(dir->col);CHKERRQ(ierr); /* only pass one ordering into CholeskyFactor */ 1179bfd6278SHong Zhang 11890d69ab7SBarry Smith flg = PETSC_FALSE; 119acfcf0e5SJed Brown ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&flg,PETSC_NULL);CHKERRQ(ierr); 1209b54502bSHong Zhang if (flg) { 1219b54502bSHong Zhang PetscReal tol = 1.e-10; 1227adad957SLisandro Dalcin ierr = PetscOptionsGetReal(((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&tol,PETSC_NULL);CHKERRQ(ierr); 1239b54502bSHong Zhang ierr = MatReorderForNonzeroDiagonal(pc->pmat,tol,dir->row,dir->row);CHKERRQ(ierr); 1249b54502bSHong Zhang } 125efee365bSSatish Balay if (dir->row) {ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr);} 126d09a07f4SBarry Smith if (!((PC_Factor*)dir)->fact){ 127075768bcSBarry Smith ierr = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_CHOLESKY,&((PC_Factor*)dir)->fact);CHKERRQ(ierr); 128a1f19f5aSHong Zhang } 129075768bcSBarry Smith ierr = MatCholeskyFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 130075768bcSBarry Smith ierr = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 1319b54502bSHong Zhang dir->actualfill = info.fill_ratio_needed; 132075768bcSBarry Smith ierr = PetscLogObjectParent(pc,((PC_Factor*)dir)->fact);CHKERRQ(ierr); 1339b54502bSHong Zhang } else if (pc->flag != SAME_NONZERO_PATTERN) { 1349b54502bSHong Zhang if (!dir->reuseordering) { 1359b54502bSHong Zhang if (dir->row && dir->col && (dir->row != dir->col)) { 1369b54502bSHong Zhang ierr = ISDestroy(dir->row);CHKERRQ(ierr); 1379b54502bSHong Zhang } 1389b54502bSHong Zhang if (dir->col) { 1399b54502bSHong Zhang ierr = ISDestroy(dir->col);CHKERRQ(ierr); 1409b54502bSHong Zhang } 141075768bcSBarry Smith ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 1429b54502bSHong Zhang if (dir->col && (dir->row != dir->col)) { /* only use row ordering for SBAIJ */ 1439b54502bSHong Zhang ierr = ISDestroy(dir->col);CHKERRQ(ierr); 1449b54502bSHong Zhang } 14590d69ab7SBarry Smith flg = PETSC_FALSE; 146acfcf0e5SJed Brown ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&flg,PETSC_NULL);CHKERRQ(ierr); 1479b54502bSHong Zhang if (flg) { 1489b54502bSHong Zhang PetscReal tol = 1.e-10; 1497adad957SLisandro Dalcin ierr = PetscOptionsGetReal(((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&tol,PETSC_NULL);CHKERRQ(ierr); 1509b54502bSHong Zhang ierr = MatReorderForNonzeroDiagonal(pc->pmat,tol,dir->row,dir->row);CHKERRQ(ierr); 1519b54502bSHong Zhang } 152efee365bSSatish Balay if (dir->row) {ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr);} 1539b54502bSHong Zhang } 154075768bcSBarry Smith ierr = MatDestroy(((PC_Factor*)dir)->fact);CHKERRQ(ierr); 155075768bcSBarry Smith ierr = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_CHOLESKY,&((PC_Factor*)dir)->fact);CHKERRQ(ierr); 156075768bcSBarry Smith ierr = MatCholeskyFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 157075768bcSBarry Smith ierr = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 1589b54502bSHong Zhang dir->actualfill = info.fill_ratio_needed; 159075768bcSBarry Smith ierr = PetscLogObjectParent(pc,((PC_Factor*)dir)->fact);CHKERRQ(ierr); 1609b54502bSHong Zhang } 161075768bcSBarry Smith ierr = MatCholeskyFactorNumeric(((PC_Factor*)dir)->fact,pc->pmat,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 1629b54502bSHong Zhang } 1639b54502bSHong Zhang PetscFunctionReturn(0); 1649b54502bSHong Zhang } 1659b54502bSHong Zhang 1669b54502bSHong Zhang #undef __FUNCT__ 16769d2c0f9SBarry Smith #define __FUNCT__ "PCReset_Cholesky" 16869d2c0f9SBarry Smith static PetscErrorCode PCReset_Cholesky(PC pc) 1699b54502bSHong Zhang { 1709b54502bSHong Zhang PC_Cholesky *dir = (PC_Cholesky*)pc->data; 1719b54502bSHong Zhang PetscErrorCode ierr; 1729b54502bSHong Zhang 1739b54502bSHong Zhang PetscFunctionBegin; 174075768bcSBarry Smith if (!dir->inplace && ((PC_Factor*)dir)->fact) {ierr = MatDestroy(((PC_Factor*)dir)->fact);CHKERRQ(ierr);} 1759b54502bSHong Zhang if (dir->row) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);} 1769b54502bSHong Zhang if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);} 17769d2c0f9SBarry Smith PetscFunctionReturn(0); 17869d2c0f9SBarry Smith } 17969d2c0f9SBarry Smith 18069d2c0f9SBarry Smith #undef __FUNCT__ 18169d2c0f9SBarry Smith #define __FUNCT__ "PCDestroy_Cholesky" 18269d2c0f9SBarry Smith static PetscErrorCode PCDestroy_Cholesky(PC pc) 18369d2c0f9SBarry Smith { 18469d2c0f9SBarry Smith PC_Cholesky *dir = (PC_Cholesky*)pc->data; 18569d2c0f9SBarry Smith PetscErrorCode ierr; 18669d2c0f9SBarry Smith 18769d2c0f9SBarry Smith PetscFunctionBegin; 18869d2c0f9SBarry Smith ierr = PCReset_Cholesky(pc);CHKERRQ(ierr); 189503cfb0cSBarry Smith ierr = PetscFree(((PC_Factor*)dir)->ordering);CHKERRQ(ierr); 190503cfb0cSBarry Smith ierr = PetscFree(((PC_Factor*)dir)->solvertype);CHKERRQ(ierr); 191c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 1929b54502bSHong Zhang PetscFunctionReturn(0); 1939b54502bSHong Zhang } 1949b54502bSHong Zhang 1959b54502bSHong Zhang #undef __FUNCT__ 1969b54502bSHong Zhang #define __FUNCT__ "PCApply_Cholesky" 1979b54502bSHong Zhang static PetscErrorCode PCApply_Cholesky(PC pc,Vec x,Vec y) 1989b54502bSHong Zhang { 1999b54502bSHong Zhang PC_Cholesky *dir = (PC_Cholesky*)pc->data; 2009b54502bSHong Zhang PetscErrorCode ierr; 2019b54502bSHong Zhang 2029b54502bSHong Zhang PetscFunctionBegin; 2039b54502bSHong Zhang if (dir->inplace) {ierr = MatSolve(pc->pmat,x,y);CHKERRQ(ierr);} 204075768bcSBarry Smith else {ierr = MatSolve(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr);} 2059b54502bSHong Zhang PetscFunctionReturn(0); 2069b54502bSHong Zhang } 2079b54502bSHong Zhang 2089b54502bSHong Zhang #undef __FUNCT__ 2099b54502bSHong Zhang #define __FUNCT__ "PCApplyTranspose_Cholesky" 2109b54502bSHong Zhang static PetscErrorCode PCApplyTranspose_Cholesky(PC pc,Vec x,Vec y) 2119b54502bSHong Zhang { 2129b54502bSHong Zhang PC_Cholesky *dir = (PC_Cholesky*)pc->data; 2139b54502bSHong Zhang PetscErrorCode ierr; 2149b54502bSHong Zhang 2159b54502bSHong Zhang PetscFunctionBegin; 2169b54502bSHong Zhang if (dir->inplace) {ierr = MatSolveTranspose(pc->pmat,x,y);CHKERRQ(ierr);} 217075768bcSBarry Smith else {ierr = MatSolveTranspose(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr);} 2189b54502bSHong Zhang PetscFunctionReturn(0); 2199b54502bSHong Zhang } 2209b54502bSHong Zhang 2219b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/ 2229b54502bSHong Zhang 2239b54502bSHong Zhang EXTERN_C_BEGIN 2249b54502bSHong Zhang #undef __FUNCT__ 2252401956bSBarry Smith #define __FUNCT__ "PCFactorSetUseInPlace_Cholesky" 2267087cfbeSBarry Smith PetscErrorCode PCFactorSetUseInPlace_Cholesky(PC pc) 2279b54502bSHong Zhang { 2289b54502bSHong Zhang PC_Cholesky *dir; 2299b54502bSHong Zhang 2309b54502bSHong Zhang PetscFunctionBegin; 2319b54502bSHong Zhang dir = (PC_Cholesky*)pc->data; 2329b54502bSHong Zhang dir->inplace = PETSC_TRUE; 2339b54502bSHong Zhang PetscFunctionReturn(0); 2349b54502bSHong Zhang } 2359b54502bSHong Zhang EXTERN_C_END 2369b54502bSHong Zhang 2379b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/ 2389b54502bSHong Zhang 2399b54502bSHong Zhang #undef __FUNCT__ 2402401956bSBarry Smith #define __FUNCT__ "PCFactorSetReuseOrdering" 2419b54502bSHong Zhang /*@ 2422401956bSBarry Smith PCFactorSetReuseOrdering - When similar matrices are factored, this 2439b54502bSHong Zhang causes the ordering computed in the first factor to be used for all 2449b54502bSHong Zhang following factors. 2459b54502bSHong Zhang 246ad4df100SBarry Smith Logically Collective on PC 2479b54502bSHong Zhang 2489b54502bSHong Zhang Input Parameters: 2499b54502bSHong Zhang + pc - the preconditioner context 2509b54502bSHong Zhang - flag - PETSC_TRUE to reuse else PETSC_FALSE 2519b54502bSHong Zhang 2529b54502bSHong Zhang Options Database Key: 2532401956bSBarry Smith . -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering() 2549b54502bSHong Zhang 2559b54502bSHong Zhang Level: intermediate 2569b54502bSHong Zhang 2579b54502bSHong Zhang .keywords: PC, levels, reordering, factorization, incomplete, LU 2589b54502bSHong Zhang 2592401956bSBarry Smith .seealso: PCFactorSetReuseFill() 2609b54502bSHong Zhang @*/ 2617087cfbeSBarry Smith PetscErrorCode PCFactorSetReuseOrdering(PC pc,PetscBool flag) 2629b54502bSHong Zhang { 2634ac538c5SBarry Smith PetscErrorCode ierr; 2649b54502bSHong Zhang 2659b54502bSHong Zhang PetscFunctionBegin; 2660700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 267acfcf0e5SJed Brown PetscValidLogicalCollectiveBool(pc,flag,2); 2684ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFactorSetReuseOrdering_C",(PC,PetscBool),(pc,flag));CHKERRQ(ierr); 2699b54502bSHong Zhang PetscFunctionReturn(0); 2709b54502bSHong Zhang } 2719b54502bSHong Zhang 2729b54502bSHong Zhang 2739b54502bSHong Zhang /*MC 27496fc60bcSBarry Smith PCCHOLESKY - Uses a direct solver, based on Cholesky factorization, as a preconditioner 2759b54502bSHong Zhang 2769b54502bSHong Zhang Options Database Keys: 2772401956bSBarry Smith + -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering() 278c7393fdbSBarry Smith . -pc_factor_mat_solver_package - Actives PCFactorSetMatSolverPackage() to choose the direct solver, like spooles 2792401956bSBarry Smith . -pc_factor_reuse_fill - Activates PCFactorSetReuseFill() 28055ba2a51SBarry Smith . -pc_factor_fill <fill> - Sets fill amount 2812401956bSBarry Smith . -pc_factor_in_place - Activates in-place factorization 282145b9266SHong Zhang - -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine 2839b54502bSHong Zhang 2849b54502bSHong Zhang Notes: Not all options work for all matrix formats 2859b54502bSHong Zhang 2869b54502bSHong Zhang Level: beginner 2879b54502bSHong Zhang 2889b54502bSHong Zhang Concepts: Cholesky factorization, direct solver 2899b54502bSHong Zhang 2909b54502bSHong Zhang Notes: Usually this will compute an "exact" solution in one iteration and does 2919b54502bSHong Zhang not need a Krylov method (i.e. you can use -ksp_type preonly, or 2929b54502bSHong Zhang KSPSetType(ksp,KSPPREONLY) for the Krylov method 2939b54502bSHong Zhang 2949b54502bSHong Zhang .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 295a4fd02acSBarry Smith PCILU, PCLU, PCICC, PCFactorSetReuseOrdering(), PCFactorSetReuseFill(), PCFactorGetMatrix(), 296145b9266SHong Zhang PCFactorSetFill(), PCFactorSetShiftNonzero(), PCFactorSetShiftType(), PCFactorSetShiftAmount() 2978ff23777SHong Zhang PCFactorSetUseInPlace(), PCFactorSetMatOrderingType() 2989b54502bSHong Zhang 2999b54502bSHong Zhang M*/ 3009b54502bSHong Zhang 3019b54502bSHong Zhang EXTERN_C_BEGIN 3029b54502bSHong Zhang #undef __FUNCT__ 3039b54502bSHong Zhang #define __FUNCT__ "PCCreate_Cholesky" 3047087cfbeSBarry Smith PetscErrorCode PCCreate_Cholesky(PC pc) 3059b54502bSHong Zhang { 3069b54502bSHong Zhang PetscErrorCode ierr; 3079b54502bSHong Zhang PC_Cholesky *dir; 3089b54502bSHong Zhang 3099b54502bSHong Zhang PetscFunctionBegin; 31038f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_Cholesky,&dir);CHKERRQ(ierr); 3119b54502bSHong Zhang 312075768bcSBarry Smith ((PC_Factor*)dir)->fact = 0; 3139b54502bSHong Zhang dir->inplace = PETSC_FALSE; 314075768bcSBarry Smith ierr = MatFactorInfoInitialize(&((PC_Factor*)dir)->info);CHKERRQ(ierr); 315879e8a4dSBarry Smith ((PC_Factor*)dir)->factortype = MAT_FACTOR_CHOLESKY; 316075768bcSBarry Smith ((PC_Factor*)dir)->info.fill = 5.0; 317f4db908eSBarry Smith ((PC_Factor*)dir)->info.shifttype = (PetscReal) MAT_SHIFT_NONE; 318915743fcSHong Zhang ((PC_Factor*)dir)->info.shiftamount = 0.0; 319*0ed735ceSHong Zhang ((PC_Factor*)dir)->info.zeropivot = 100.0*PETSC_MACHINE_EPSILON; 320075768bcSBarry Smith ((PC_Factor*)dir)->info.pivotinblocks = 1.0; 3219b54502bSHong Zhang dir->col = 0; 3229b54502bSHong Zhang dir->row = 0; 3232692d6eeSBarry Smith ierr = PetscStrallocpy(MATORDERINGNATURAL,&((PC_Factor*)dir)->ordering);CHKERRQ(ierr); 3242692d6eeSBarry Smith ierr = PetscStrallocpy(MATSOLVERPETSC,&((PC_Factor*)dir)->solvertype);CHKERRQ(ierr); 3259b54502bSHong Zhang dir->reusefill = PETSC_FALSE; 3269b54502bSHong Zhang dir->reuseordering = PETSC_FALSE; 3279b54502bSHong Zhang pc->data = (void*)dir; 3289b54502bSHong Zhang 3299b54502bSHong Zhang pc->ops->destroy = PCDestroy_Cholesky; 33069d2c0f9SBarry Smith pc->ops->reset = PCReset_Cholesky; 3319b54502bSHong Zhang pc->ops->apply = PCApply_Cholesky; 3329b54502bSHong Zhang pc->ops->applytranspose = PCApplyTranspose_Cholesky; 3339b54502bSHong Zhang pc->ops->setup = PCSetUp_Cholesky; 3349b54502bSHong Zhang pc->ops->setfromoptions = PCSetFromOptions_Cholesky; 3359b54502bSHong Zhang pc->ops->view = PCView_Cholesky; 3369b54502bSHong Zhang pc->ops->applyrichardson = 0; 33785317021SBarry Smith pc->ops->getfactoredmatrix = PCFactorGetMatrix_Factor; 3389b54502bSHong Zhang 339f8260c8fSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUpMatSolverPackage_C","PCFactorSetUpMatSolverPackage_Factor", 340f8260c8fSBarry Smith PCFactorSetUpMatSolverPackage_Factor);CHKERRQ(ierr); 34185317021SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatSolverPackage_C","PCFactorSetMatSolverPackage_Factor", 34285317021SBarry Smith PCFactorSetMatSolverPackage_Factor);CHKERRQ(ierr); 3437112b564SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorGetMatSolverPackage_C","PCFactorGetMatSolverPackage_Factor", 3447112b564SBarry Smith PCFactorGetMatSolverPackage_Factor);CHKERRQ(ierr); 34585317021SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_Factor", 34685317021SBarry Smith PCFactorSetZeroPivot_Factor);CHKERRQ(ierr); 347d90ac83dSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftType_C","PCFactorSetShiftType_Factor", 348d90ac83dSHong Zhang PCFactorSetShiftType_Factor);CHKERRQ(ierr); 349d90ac83dSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftAmount_C","PCFactorSetShiftAmount_Factor", 350d90ac83dSHong Zhang PCFactorSetShiftAmount_Factor);CHKERRQ(ierr); 35185317021SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetFill_C","PCFactorSetFill_Factor", 35285317021SBarry Smith PCFactorSetFill_Factor);CHKERRQ(ierr); 3532401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUseInPlace_C","PCFactorSetUseInPlace_Cholesky", 3542401956bSBarry Smith PCFactorSetUseInPlace_Cholesky);CHKERRQ(ierr); 35585317021SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatOrderingType_C","PCFactorSetMatOrderingType_Factor", 35685317021SBarry Smith PCFactorSetMatOrderingType_Factor);CHKERRQ(ierr); 3572401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseOrdering_C","PCFactorSetReuseOrdering_Cholesky", 3582401956bSBarry Smith PCFactorSetReuseOrdering_Cholesky);CHKERRQ(ierr); 3592401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseFill_C","PCFactorSetReuseFill_Cholesky", 3602401956bSBarry Smith PCFactorSetReuseFill_Cholesky);CHKERRQ(ierr); 3619b54502bSHong Zhang PetscFunctionReturn(0); 3629b54502bSHong Zhang } 3639b54502bSHong Zhang EXTERN_C_END 364