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 { 229b54502bSHong Zhang PC_Cholesky *lu; 239b54502bSHong Zhang 249b54502bSHong Zhang PetscFunctionBegin; 259b54502bSHong Zhang lu = (PC_Cholesky*)pc->data; 269b54502bSHong Zhang lu->reuseordering = flag; 279b54502bSHong Zhang PetscFunctionReturn(0); 289b54502bSHong Zhang } 299b54502bSHong Zhang 309b54502bSHong Zhang #undef __FUNCT__ 312401956bSBarry Smith #define __FUNCT__ "PCFactorSetReuseFill_Cholesky" 32b2573a8aSBarry Smith static PetscErrorCode PCFactorSetReuseFill_Cholesky(PC pc,PetscBool flag) 339b54502bSHong Zhang { 349b54502bSHong Zhang PC_Cholesky *lu; 359b54502bSHong Zhang 369b54502bSHong Zhang PetscFunctionBegin; 379b54502bSHong Zhang lu = (PC_Cholesky*)pc->data; 389b54502bSHong Zhang lu->reusefill = flag; 399b54502bSHong Zhang PetscFunctionReturn(0); 409b54502bSHong Zhang } 419b54502bSHong Zhang 429b54502bSHong Zhang #undef __FUNCT__ 439b54502bSHong Zhang #define __FUNCT__ "PCSetFromOptions_Cholesky" 449b54502bSHong Zhang static PetscErrorCode PCSetFromOptions_Cholesky(PC pc) 459b54502bSHong Zhang { 469b54502bSHong Zhang PetscErrorCode ierr; 479b54502bSHong Zhang 489b54502bSHong Zhang PetscFunctionBegin; 499b54502bSHong Zhang ierr = PetscOptionsHead("Cholesky options");CHKERRQ(ierr); 508ff23777SHong Zhang ierr = PCSetFromOptions_Factor(pc);CHKERRQ(ierr); 519b54502bSHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 529b54502bSHong Zhang PetscFunctionReturn(0); 539b54502bSHong Zhang } 549b54502bSHong Zhang 559b54502bSHong Zhang #undef __FUNCT__ 569b54502bSHong Zhang #define __FUNCT__ "PCView_Cholesky" 579b54502bSHong Zhang static PetscErrorCode PCView_Cholesky(PC pc,PetscViewer viewer) 589b54502bSHong Zhang { 59145b9266SHong Zhang PC_Cholesky *chol = (PC_Cholesky*)pc->data; 609b54502bSHong Zhang PetscErrorCode ierr; 61ace3abfcSBarry Smith PetscBool iascii; 629b54502bSHong Zhang 639b54502bSHong Zhang PetscFunctionBegin; 64251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 659b54502bSHong Zhang if (iascii) { 66914a5d51SHong Zhang if (chol->inplace) { 67914a5d51SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Cholesky: in-place factorization\n");CHKERRQ(ierr); 68914a5d51SHong Zhang } else { 69914a5d51SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Cholesky: out-of-place factorization\n");CHKERRQ(ierr); 70145b9266SHong Zhang } 719b54502bSHong Zhang 72145b9266SHong Zhang if (chol->reusefill) {ierr = PetscViewerASCIIPrintf(viewer," Reusing fill from past factorization\n");CHKERRQ(ierr);} 73145b9266SHong Zhang if (chol->reuseordering) {ierr = PetscViewerASCIIPrintf(viewer," Reusing reordering from past factorization\n");CHKERRQ(ierr);} 749b54502bSHong Zhang } 75914a5d51SHong Zhang ierr = PCView_Factor(pc,viewer);CHKERRQ(ierr); 769b54502bSHong Zhang PetscFunctionReturn(0); 779b54502bSHong Zhang } 789b54502bSHong Zhang 799b54502bSHong Zhang 809b54502bSHong Zhang #undef __FUNCT__ 819b54502bSHong Zhang #define __FUNCT__ "PCSetUp_Cholesky" 829b54502bSHong Zhang static PetscErrorCode PCSetUp_Cholesky(PC pc) 839b54502bSHong Zhang { 849b54502bSHong Zhang PetscErrorCode ierr; 85ace3abfcSBarry Smith PetscBool flg; 869b54502bSHong Zhang PC_Cholesky *dir = (PC_Cholesky*)pc->data; 879b54502bSHong Zhang 889b54502bSHong Zhang PetscFunctionBegin; 89075768bcSBarry Smith if (dir->reusefill && pc->setupcalled) ((PC_Factor*)dir)->info.fill = dir->actualfill; 909b54502bSHong Zhang 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 } 100*3bb1ff40SBarry 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); 1022fa5cd67SKarl Rupp 103075768bcSBarry Smith ((PC_Factor*)dir)->fact = pc->pmat; 1049b54502bSHong Zhang } else { 1059b54502bSHong Zhang MatInfo info; 1069b54502bSHong Zhang if (!pc->setupcalled) { 107075768bcSBarry Smith ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 1089bfd6278SHong Zhang /* check if dir->row == dir->col */ 1099bfd6278SHong Zhang ierr = ISEqual(dir->row,dir->col,&flg);CHKERRQ(ierr); 110e32f2f54SBarry Smith if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"row and column permutations must equal"); 111fcfd50ebSBarry Smith ierr = ISDestroy(&dir->col);CHKERRQ(ierr); /* only pass one ordering into CholeskyFactor */ 1129bfd6278SHong Zhang 11390d69ab7SBarry Smith flg = PETSC_FALSE; 1140298fd71SBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&flg,NULL);CHKERRQ(ierr); 1159b54502bSHong Zhang if (flg) { 1169b54502bSHong Zhang PetscReal tol = 1.e-10; 1170298fd71SBarry Smith ierr = PetscOptionsGetReal(((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&tol,NULL);CHKERRQ(ierr); 1189b54502bSHong Zhang ierr = MatReorderForNonzeroDiagonal(pc->pmat,tol,dir->row,dir->row);CHKERRQ(ierr); 1199b54502bSHong Zhang } 120*3bb1ff40SBarry Smith if (dir->row) {ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);CHKERRQ(ierr);} 121d09a07f4SBarry Smith if (!((PC_Factor*)dir)->fact) { 122075768bcSBarry Smith ierr = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_CHOLESKY,&((PC_Factor*)dir)->fact);CHKERRQ(ierr); 123a1f19f5aSHong Zhang } 124075768bcSBarry Smith ierr = MatCholeskyFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 125075768bcSBarry Smith ierr = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 1269b54502bSHong Zhang dir->actualfill = info.fill_ratio_needed; 127*3bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact);CHKERRQ(ierr); 1289b54502bSHong Zhang } else if (pc->flag != SAME_NONZERO_PATTERN) { 1299b54502bSHong Zhang if (!dir->reuseordering) { 130fcfd50ebSBarry Smith ierr = ISDestroy(&dir->row);CHKERRQ(ierr); 131075768bcSBarry Smith ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 1329d61c4f5SHong Zhang ierr = ISDestroy(&dir->col);CHKERRQ(ierr); /* only use dir->row ordering in CholeskyFactor */ 1339d61c4f5SHong Zhang 13490d69ab7SBarry Smith flg = PETSC_FALSE; 1350298fd71SBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&flg,NULL);CHKERRQ(ierr); 1369b54502bSHong Zhang if (flg) { 1379b54502bSHong Zhang PetscReal tol = 1.e-10; 1380298fd71SBarry Smith ierr = PetscOptionsGetReal(((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&tol,NULL);CHKERRQ(ierr); 1399b54502bSHong Zhang ierr = MatReorderForNonzeroDiagonal(pc->pmat,tol,dir->row,dir->row);CHKERRQ(ierr); 1409b54502bSHong Zhang } 141*3bb1ff40SBarry Smith if (dir->row) {ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);CHKERRQ(ierr);} 1429b54502bSHong Zhang } 1436bf464f9SBarry Smith ierr = MatDestroy(&((PC_Factor*)dir)->fact);CHKERRQ(ierr); 144075768bcSBarry Smith ierr = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_CHOLESKY,&((PC_Factor*)dir)->fact);CHKERRQ(ierr); 145075768bcSBarry Smith ierr = MatCholeskyFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 146075768bcSBarry Smith ierr = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 1479b54502bSHong Zhang dir->actualfill = info.fill_ratio_needed; 148*3bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact);CHKERRQ(ierr); 1499b54502bSHong Zhang } 150075768bcSBarry Smith ierr = MatCholeskyFactorNumeric(((PC_Factor*)dir)->fact,pc->pmat,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 1519b54502bSHong Zhang } 1529b54502bSHong Zhang PetscFunctionReturn(0); 1539b54502bSHong Zhang } 1549b54502bSHong Zhang 1559b54502bSHong Zhang #undef __FUNCT__ 15669d2c0f9SBarry Smith #define __FUNCT__ "PCReset_Cholesky" 15769d2c0f9SBarry Smith static PetscErrorCode PCReset_Cholesky(PC pc) 1589b54502bSHong Zhang { 1599b54502bSHong Zhang PC_Cholesky *dir = (PC_Cholesky*)pc->data; 1609b54502bSHong Zhang PetscErrorCode ierr; 1619b54502bSHong Zhang 1629b54502bSHong Zhang PetscFunctionBegin; 1636bf464f9SBarry Smith if (!dir->inplace && ((PC_Factor*)dir)->fact) {ierr = MatDestroy(&((PC_Factor*)dir)->fact);CHKERRQ(ierr);} 164fcfd50ebSBarry Smith ierr = ISDestroy(&dir->row);CHKERRQ(ierr); 165fcfd50ebSBarry Smith ierr = ISDestroy(&dir->col);CHKERRQ(ierr); 16669d2c0f9SBarry Smith PetscFunctionReturn(0); 16769d2c0f9SBarry Smith } 16869d2c0f9SBarry Smith 16969d2c0f9SBarry Smith #undef __FUNCT__ 17069d2c0f9SBarry Smith #define __FUNCT__ "PCDestroy_Cholesky" 17169d2c0f9SBarry Smith static PetscErrorCode PCDestroy_Cholesky(PC pc) 17269d2c0f9SBarry Smith { 17369d2c0f9SBarry Smith PC_Cholesky *dir = (PC_Cholesky*)pc->data; 17469d2c0f9SBarry Smith PetscErrorCode ierr; 17569d2c0f9SBarry Smith 17669d2c0f9SBarry Smith PetscFunctionBegin; 17769d2c0f9SBarry Smith ierr = PCReset_Cholesky(pc);CHKERRQ(ierr); 178503cfb0cSBarry Smith ierr = PetscFree(((PC_Factor*)dir)->ordering);CHKERRQ(ierr); 179503cfb0cSBarry Smith ierr = PetscFree(((PC_Factor*)dir)->solvertype);CHKERRQ(ierr); 180c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 1819b54502bSHong Zhang PetscFunctionReturn(0); 1829b54502bSHong Zhang } 1839b54502bSHong Zhang 1849b54502bSHong Zhang #undef __FUNCT__ 1859b54502bSHong Zhang #define __FUNCT__ "PCApply_Cholesky" 1869b54502bSHong Zhang static PetscErrorCode PCApply_Cholesky(PC pc,Vec x,Vec y) 1879b54502bSHong Zhang { 1889b54502bSHong Zhang PC_Cholesky *dir = (PC_Cholesky*)pc->data; 1899b54502bSHong Zhang PetscErrorCode ierr; 1909b54502bSHong Zhang 1919b54502bSHong Zhang PetscFunctionBegin; 1922fa5cd67SKarl Rupp if (dir->inplace) { 1932fa5cd67SKarl Rupp ierr = MatSolve(pc->pmat,x,y);CHKERRQ(ierr); 1942fa5cd67SKarl Rupp } else { 1952fa5cd67SKarl Rupp ierr = MatSolve(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr); 1962fa5cd67SKarl Rupp } 1979b54502bSHong Zhang PetscFunctionReturn(0); 1989b54502bSHong Zhang } 1999b54502bSHong Zhang 2009b54502bSHong Zhang #undef __FUNCT__ 2019b54502bSHong Zhang #define __FUNCT__ "PCApplyTranspose_Cholesky" 2029b54502bSHong Zhang static PetscErrorCode PCApplyTranspose_Cholesky(PC pc,Vec x,Vec y) 2039b54502bSHong Zhang { 2049b54502bSHong Zhang PC_Cholesky *dir = (PC_Cholesky*)pc->data; 2059b54502bSHong Zhang PetscErrorCode ierr; 2069b54502bSHong Zhang 2079b54502bSHong Zhang PetscFunctionBegin; 2082fa5cd67SKarl Rupp if (dir->inplace) { 2092fa5cd67SKarl Rupp ierr = MatSolveTranspose(pc->pmat,x,y);CHKERRQ(ierr); 2102fa5cd67SKarl Rupp } else { 2112fa5cd67SKarl Rupp ierr = MatSolveTranspose(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr); 2122fa5cd67SKarl Rupp } 2139b54502bSHong Zhang PetscFunctionReturn(0); 2149b54502bSHong Zhang } 2159b54502bSHong Zhang 2169b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/ 2179b54502bSHong Zhang 2189b54502bSHong Zhang #undef __FUNCT__ 2192401956bSBarry Smith #define __FUNCT__ "PCFactorSetUseInPlace_Cholesky" 220b2573a8aSBarry Smith static PetscErrorCode 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 2309b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/ 2319b54502bSHong Zhang 2329b54502bSHong Zhang #undef __FUNCT__ 2332401956bSBarry Smith #define __FUNCT__ "PCFactorSetReuseOrdering" 2349b54502bSHong Zhang /*@ 2352401956bSBarry Smith PCFactorSetReuseOrdering - When similar matrices are factored, this 2369b54502bSHong Zhang causes the ordering computed in the first factor to be used for all 2379b54502bSHong Zhang following factors. 2389b54502bSHong Zhang 239ad4df100SBarry Smith Logically Collective on PC 2409b54502bSHong Zhang 2419b54502bSHong Zhang Input Parameters: 2429b54502bSHong Zhang + pc - the preconditioner context 2439b54502bSHong Zhang - flag - PETSC_TRUE to reuse else PETSC_FALSE 2449b54502bSHong Zhang 2459b54502bSHong Zhang Options Database Key: 2462401956bSBarry Smith . -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering() 2479b54502bSHong Zhang 2489b54502bSHong Zhang Level: intermediate 2499b54502bSHong Zhang 2509b54502bSHong Zhang .keywords: PC, levels, reordering, factorization, incomplete, LU 2519b54502bSHong Zhang 2522401956bSBarry Smith .seealso: PCFactorSetReuseFill() 2539b54502bSHong Zhang @*/ 2547087cfbeSBarry Smith PetscErrorCode PCFactorSetReuseOrdering(PC pc,PetscBool flag) 2559b54502bSHong Zhang { 2564ac538c5SBarry Smith PetscErrorCode ierr; 2579b54502bSHong Zhang 2589b54502bSHong Zhang PetscFunctionBegin; 2590700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 260acfcf0e5SJed Brown PetscValidLogicalCollectiveBool(pc,flag,2); 2614ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFactorSetReuseOrdering_C",(PC,PetscBool),(pc,flag));CHKERRQ(ierr); 2629b54502bSHong Zhang PetscFunctionReturn(0); 2639b54502bSHong Zhang } 2649b54502bSHong Zhang 2659b54502bSHong Zhang 2669b54502bSHong Zhang /*MC 26796fc60bcSBarry Smith PCCHOLESKY - Uses a direct solver, based on Cholesky factorization, as a preconditioner 2689b54502bSHong Zhang 2699b54502bSHong Zhang Options Database Keys: 2702401956bSBarry Smith + -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering() 271f60c3dc2SHong Zhang . -pc_factor_mat_solver_package - Actives PCFactorSetMatSolverPackage() to choose the direct solver, like superlu 2722401956bSBarry Smith . -pc_factor_reuse_fill - Activates PCFactorSetReuseFill() 27355ba2a51SBarry Smith . -pc_factor_fill <fill> - Sets fill amount 2742401956bSBarry Smith . -pc_factor_in_place - Activates in-place factorization 275145b9266SHong Zhang - -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine 2769b54502bSHong Zhang 2779b54502bSHong Zhang Notes: Not all options work for all matrix formats 2789b54502bSHong Zhang 2799b54502bSHong Zhang Level: beginner 2809b54502bSHong Zhang 2819b54502bSHong Zhang Concepts: Cholesky factorization, direct solver 2829b54502bSHong Zhang 2839b54502bSHong Zhang Notes: Usually this will compute an "exact" solution in one iteration and does 2849b54502bSHong Zhang not need a Krylov method (i.e. you can use -ksp_type preonly, or 2859b54502bSHong Zhang KSPSetType(ksp,KSPPREONLY) for the Krylov method 2869b54502bSHong Zhang 2879b54502bSHong Zhang .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 288a4fd02acSBarry Smith PCILU, PCLU, PCICC, PCFactorSetReuseOrdering(), PCFactorSetReuseFill(), PCFactorGetMatrix(), 289145b9266SHong Zhang PCFactorSetFill(), PCFactorSetShiftNonzero(), PCFactorSetShiftType(), PCFactorSetShiftAmount() 2908ff23777SHong Zhang PCFactorSetUseInPlace(), PCFactorSetMatOrderingType() 2919b54502bSHong Zhang 2929b54502bSHong Zhang M*/ 2939b54502bSHong Zhang 2949b54502bSHong Zhang #undef __FUNCT__ 2959b54502bSHong Zhang #define __FUNCT__ "PCCreate_Cholesky" 2968cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_Cholesky(PC pc) 2979b54502bSHong Zhang { 2989b54502bSHong Zhang PetscErrorCode ierr; 2999b54502bSHong Zhang PC_Cholesky *dir; 3009b54502bSHong Zhang 3019b54502bSHong Zhang PetscFunctionBegin; 30238f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_Cholesky,&dir);CHKERRQ(ierr); 3039b54502bSHong Zhang 304075768bcSBarry Smith ((PC_Factor*)dir)->fact = 0; 3059b54502bSHong Zhang dir->inplace = PETSC_FALSE; 3062fa5cd67SKarl Rupp 307075768bcSBarry Smith ierr = MatFactorInfoInitialize(&((PC_Factor*)dir)->info);CHKERRQ(ierr); 3082fa5cd67SKarl Rupp 309879e8a4dSBarry Smith ((PC_Factor*)dir)->factortype = MAT_FACTOR_CHOLESKY; 310075768bcSBarry Smith ((PC_Factor*)dir)->info.fill = 5.0; 311f4db908eSBarry Smith ((PC_Factor*)dir)->info.shifttype = (PetscReal) MAT_SHIFT_NONE; 312915743fcSHong Zhang ((PC_Factor*)dir)->info.shiftamount = 0.0; 3130ed735ceSHong Zhang ((PC_Factor*)dir)->info.zeropivot = 100.0*PETSC_MACHINE_EPSILON; 314075768bcSBarry Smith ((PC_Factor*)dir)->info.pivotinblocks = 1.0; 3152fa5cd67SKarl Rupp 3169b54502bSHong Zhang dir->col = 0; 3179b54502bSHong Zhang dir->row = 0; 3182fa5cd67SKarl Rupp 31919fd82e9SBarry Smith ierr = PetscStrallocpy(MATORDERINGNATURAL,(char**)&((PC_Factor*)dir)->ordering);CHKERRQ(ierr); 3202692d6eeSBarry Smith ierr = PetscStrallocpy(MATSOLVERPETSC,&((PC_Factor*)dir)->solvertype);CHKERRQ(ierr); 3212fa5cd67SKarl Rupp 3229b54502bSHong Zhang dir->reusefill = PETSC_FALSE; 3239b54502bSHong Zhang dir->reuseordering = PETSC_FALSE; 3249b54502bSHong Zhang pc->data = (void*)dir; 3259b54502bSHong Zhang 3269b54502bSHong Zhang pc->ops->destroy = PCDestroy_Cholesky; 32769d2c0f9SBarry Smith pc->ops->reset = PCReset_Cholesky; 3289b54502bSHong Zhang pc->ops->apply = PCApply_Cholesky; 3299b54502bSHong Zhang pc->ops->applytranspose = PCApplyTranspose_Cholesky; 3309b54502bSHong Zhang pc->ops->setup = PCSetUp_Cholesky; 3319b54502bSHong Zhang pc->ops->setfromoptions = PCSetFromOptions_Cholesky; 3329b54502bSHong Zhang pc->ops->view = PCView_Cholesky; 3339b54502bSHong Zhang pc->ops->applyrichardson = 0; 33485317021SBarry Smith pc->ops->getfactoredmatrix = PCFactorGetMatrix_Factor; 3359b54502bSHong Zhang 336bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetUpMatSolverPackage_C",PCFactorSetUpMatSolverPackage_Factor);CHKERRQ(ierr); 337bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetMatSolverPackage_C",PCFactorSetMatSolverPackage_Factor);CHKERRQ(ierr); 338bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetMatSolverPackage_C",PCFactorGetMatSolverPackage_Factor);CHKERRQ(ierr); 339bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetZeroPivot_C",PCFactorSetZeroPivot_Factor);CHKERRQ(ierr); 340bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetShiftType_C",PCFactorSetShiftType_Factor);CHKERRQ(ierr); 341bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetShiftAmount_C",PCFactorSetShiftAmount_Factor);CHKERRQ(ierr); 342bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetFill_C",PCFactorSetFill_Factor);CHKERRQ(ierr); 343bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetUseInPlace_C",PCFactorSetUseInPlace_Cholesky);CHKERRQ(ierr); 344bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetMatOrderingType_C",PCFactorSetMatOrderingType_Factor);CHKERRQ(ierr); 345bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetReuseOrdering_C",PCFactorSetReuseOrdering_Cholesky);CHKERRQ(ierr); 346bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetReuseFill_C",PCFactorSetReuseFill_Cholesky);CHKERRQ(ierr); 3479b54502bSHong Zhang PetscFunctionReturn(0); 3489b54502bSHong Zhang } 349