xref: /petsc/src/ksp/pc/impls/factor/cholesky/cholesky.c (revision 145b9266b3dcbabdc57db56abc48d2ce1118ea0a)
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 {
64*145b9266SHong Zhang   PC_Cholesky    *chol = (PC_Cholesky*)pc->data;
659b54502bSHong Zhang   PetscErrorCode ierr;
669b54502bSHong Zhang   PetscTruth     iascii,isstring;
679b54502bSHong Zhang 
689b54502bSHong Zhang   PetscFunctionBegin;
699b54502bSHong Zhang   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
709b54502bSHong Zhang   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr);
719b54502bSHong Zhang   if (iascii) {
72*145b9266SHong Zhang     if (((PC_Factor*)chol)->info.shifttype==MAT_SHIFT_POSITIVE_DEFINITE) {
73*145b9266SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  Cholesky: using Manteuffel shift\n");CHKERRQ(ierr);
74*145b9266SHong Zhang     }
75*145b9266SHong Zhang     if (((PC_Factor*)chol)->info.shifttype==MAT_SHIFT_NONZERO) {
76*145b9266SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  Cholesky: using diagonal shift to prevent zero pivot\n");CHKERRQ(ierr);
77*145b9266SHong Zhang     }
78*145b9266SHong Zhang     if (((PC_Factor*)chol)->info.shifttype==MAT_SHIFT_INBLOCKS) {
79*145b9266SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  Cholesky: using diagonal shift on blocks to prevent zero pivot\n");CHKERRQ(ierr);
80*145b9266SHong Zhang     }
819b54502bSHong Zhang 
82*145b9266SHong Zhang     if (chol->inplace) {ierr = PetscViewerASCIIPrintf(viewer,"  Cholesky: in-place factorization\n");CHKERRQ(ierr);}
839b54502bSHong Zhang     else             {ierr = PetscViewerASCIIPrintf(viewer,"  Cholesky: out-of-place factorization\n");CHKERRQ(ierr);}
84*145b9266SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"    matrix ordering: %s\n",((PC_Factor*)chol)->ordering);CHKERRQ(ierr);
85*145b9266SHong Zhang     if (((PC_Factor*)chol)->fact) {
86*145b9266SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  Cholesky: factor fill ratio needed %G\n",chol->actualfill);CHKERRQ(ierr);
87f3a39becSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"       Factored matrix follows\n");CHKERRQ(ierr);
88f3a39becSBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
89f3a39becSBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
90f3a39becSBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
91f3a39becSBarry Smith       ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
92*145b9266SHong Zhang       ierr = MatView(((PC_Factor*)chol)->fact,viewer);CHKERRQ(ierr);
939b54502bSHong Zhang       ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
94f3a39becSBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
95f3a39becSBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
96f3a39becSBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
979b54502bSHong Zhang     }
98*145b9266SHong Zhang     if (chol->reusefill)    {ierr = PetscViewerASCIIPrintf(viewer,"       Reusing fill from past factorization\n");CHKERRQ(ierr);}
99*145b9266SHong Zhang     if (chol->reuseordering) {ierr = PetscViewerASCIIPrintf(viewer,"       Reusing reordering from past factorization\n");CHKERRQ(ierr);}
1009b54502bSHong Zhang   } else if (isstring) {
101*145b9266SHong Zhang     ierr = PetscViewerStringSPrintf(viewer," order=%s",((PC_Factor*)chol)->ordering);CHKERRQ(ierr);CHKERRQ(ierr);
1029b54502bSHong Zhang   } else {
1032401956bSBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCCHOLESKY",((PetscObject)viewer)->type_name);
1049b54502bSHong Zhang   }
1059b54502bSHong Zhang   PetscFunctionReturn(0);
1069b54502bSHong Zhang }
1079b54502bSHong Zhang 
1089b54502bSHong Zhang 
1099b54502bSHong Zhang #undef __FUNCT__
1109b54502bSHong Zhang #define __FUNCT__ "PCSetUp_Cholesky"
1119b54502bSHong Zhang static PetscErrorCode PCSetUp_Cholesky(PC pc)
1129b54502bSHong Zhang {
1139b54502bSHong Zhang   PetscErrorCode ierr;
1149b54502bSHong Zhang   PetscTruth     flg;
1159b54502bSHong Zhang   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;
1169b54502bSHong Zhang 
1179b54502bSHong Zhang   PetscFunctionBegin;
118075768bcSBarry Smith   if (dir->reusefill && pc->setupcalled) ((PC_Factor*)dir)->info.fill = dir->actualfill;
1199b54502bSHong Zhang 
1209b54502bSHong Zhang   if (dir->inplace) {
1219b54502bSHong Zhang     if (dir->row && dir->col && (dir->row != dir->col)) {
1229b54502bSHong Zhang       ierr = ISDestroy(dir->row);CHKERRQ(ierr);
1239b54502bSHong Zhang       dir->row = 0;
1249b54502bSHong Zhang     }
1259b54502bSHong Zhang     if (dir->col) {
1269b54502bSHong Zhang       ierr = ISDestroy(dir->col);CHKERRQ(ierr);
1279b54502bSHong Zhang       dir->col = 0;
1289b54502bSHong Zhang     }
129075768bcSBarry Smith     ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr);
1309b54502bSHong Zhang     if (dir->col && (dir->row != dir->col)) {  /* only use row ordering for SBAIJ */
1319b54502bSHong Zhang       ierr = ISDestroy(dir->col);CHKERRQ(ierr);
1329b54502bSHong Zhang       dir->col=0;
1339b54502bSHong Zhang     }
134efee365bSSatish Balay     if (dir->row) {ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr);}
135075768bcSBarry Smith     ierr = MatCholeskyFactor(pc->pmat,dir->row,&((PC_Factor*)dir)->info);CHKERRQ(ierr);
136075768bcSBarry Smith     ((PC_Factor*)dir)->fact = pc->pmat;
1379b54502bSHong Zhang   } else {
1389b54502bSHong Zhang     MatInfo info;
1399b54502bSHong Zhang     if (!pc->setupcalled) {
140075768bcSBarry Smith       ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr);
1419bfd6278SHong Zhang       /* check if dir->row == dir->col */
1429bfd6278SHong Zhang       ierr = ISEqual(dir->row,dir->col,&flg);CHKERRQ(ierr);
1439bfd6278SHong Zhang       if (!flg) SETERRQ(PETSC_ERR_ARG_INCOMP,"row and column permutations must equal");
1449bfd6278SHong Zhang       ierr = ISDestroy(dir->col);CHKERRQ(ierr); /* only pass one ordering into CholeskyFactor */
1459b54502bSHong Zhang       dir->col=0;
1469bfd6278SHong Zhang 
14790d69ab7SBarry Smith       flg  = PETSC_FALSE;
14890d69ab7SBarry Smith       ierr = PetscOptionsGetTruth(((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&flg,PETSC_NULL);CHKERRQ(ierr);
1499b54502bSHong Zhang       if (flg) {
1509b54502bSHong Zhang         PetscReal tol = 1.e-10;
1517adad957SLisandro Dalcin         ierr = PetscOptionsGetReal(((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&tol,PETSC_NULL);CHKERRQ(ierr);
1529b54502bSHong Zhang         ierr = MatReorderForNonzeroDiagonal(pc->pmat,tol,dir->row,dir->row);CHKERRQ(ierr);
1539b54502bSHong Zhang       }
154efee365bSSatish Balay       if (dir->row) {ierr = PetscLogObjectParent(pc,dir->row);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     } else if (pc->flag != SAME_NONZERO_PATTERN) {
1619b54502bSHong Zhang       if (!dir->reuseordering) {
1629b54502bSHong Zhang         if (dir->row && dir->col && (dir->row != dir->col)) {
1639b54502bSHong Zhang           ierr = ISDestroy(dir->row);CHKERRQ(ierr);
1649b54502bSHong Zhang           dir->row = 0;
1659b54502bSHong Zhang         }
1669b54502bSHong Zhang         if (dir->col) {
1679b54502bSHong Zhang           ierr = ISDestroy(dir->col);CHKERRQ(ierr);
1689b54502bSHong Zhang           dir->col =0;
1699b54502bSHong Zhang         }
170075768bcSBarry Smith         ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr);
1719b54502bSHong Zhang         if (dir->col && (dir->row != dir->col)) {  /* only use row ordering for SBAIJ */
1729b54502bSHong Zhang           ierr = ISDestroy(dir->col);CHKERRQ(ierr);
1739b54502bSHong Zhang           dir->col=0;
1749b54502bSHong Zhang         }
17590d69ab7SBarry Smith         flg  = PETSC_FALSE;
17690d69ab7SBarry Smith         ierr = PetscOptionsGetTruth(((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&flg,PETSC_NULL);CHKERRQ(ierr);
1779b54502bSHong Zhang         if (flg) {
1789b54502bSHong Zhang           PetscReal tol = 1.e-10;
1797adad957SLisandro Dalcin           ierr = PetscOptionsGetReal(((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&tol,PETSC_NULL);CHKERRQ(ierr);
1809b54502bSHong Zhang           ierr = MatReorderForNonzeroDiagonal(pc->pmat,tol,dir->row,dir->row);CHKERRQ(ierr);
1819b54502bSHong Zhang         }
182efee365bSSatish Balay         if (dir->row) {ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr);}
1839b54502bSHong Zhang       }
184075768bcSBarry Smith       ierr = MatDestroy(((PC_Factor*)dir)->fact);CHKERRQ(ierr);
185075768bcSBarry Smith       ierr = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_CHOLESKY,&((PC_Factor*)dir)->fact);CHKERRQ(ierr);
186075768bcSBarry Smith       ierr = MatCholeskyFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,&((PC_Factor*)dir)->info);CHKERRQ(ierr);
187075768bcSBarry Smith       ierr = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
1889b54502bSHong Zhang       dir->actualfill = info.fill_ratio_needed;
189075768bcSBarry Smith       ierr = PetscLogObjectParent(pc,((PC_Factor*)dir)->fact);CHKERRQ(ierr);
1909b54502bSHong Zhang     }
191075768bcSBarry Smith     ierr = MatCholeskyFactorNumeric(((PC_Factor*)dir)->fact,pc->pmat,&((PC_Factor*)dir)->info);CHKERRQ(ierr);
1929b54502bSHong Zhang   }
1939b54502bSHong Zhang   PetscFunctionReturn(0);
1949b54502bSHong Zhang }
1959b54502bSHong Zhang 
1969b54502bSHong Zhang #undef __FUNCT__
1979b54502bSHong Zhang #define __FUNCT__ "PCDestroy_Cholesky"
1989b54502bSHong Zhang static PetscErrorCode PCDestroy_Cholesky(PC pc)
1999b54502bSHong Zhang {
2009b54502bSHong Zhang   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;
2019b54502bSHong Zhang   PetscErrorCode ierr;
2029b54502bSHong Zhang 
2039b54502bSHong Zhang   PetscFunctionBegin;
204075768bcSBarry Smith   if (!dir->inplace && ((PC_Factor*)dir)->fact) {ierr = MatDestroy(((PC_Factor*)dir)->fact);CHKERRQ(ierr);}
2059b54502bSHong Zhang   if (dir->row) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);}
2069b54502bSHong Zhang   if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);}
207075768bcSBarry Smith   ierr = PetscStrfree(((PC_Factor*)dir)->ordering);CHKERRQ(ierr);
208075768bcSBarry Smith   ierr = PetscStrfree(((PC_Factor*)dir)->solvertype);CHKERRQ(ierr);
2099b54502bSHong Zhang   ierr = PetscFree(dir);CHKERRQ(ierr);
2109b54502bSHong Zhang   PetscFunctionReturn(0);
2119b54502bSHong Zhang }
2129b54502bSHong Zhang 
2139b54502bSHong Zhang #undef __FUNCT__
2149b54502bSHong Zhang #define __FUNCT__ "PCApply_Cholesky"
2159b54502bSHong Zhang static PetscErrorCode PCApply_Cholesky(PC pc,Vec x,Vec y)
2169b54502bSHong Zhang {
2179b54502bSHong Zhang   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;
2189b54502bSHong Zhang   PetscErrorCode ierr;
2199b54502bSHong Zhang 
2209b54502bSHong Zhang   PetscFunctionBegin;
2219b54502bSHong Zhang   if (dir->inplace) {ierr = MatSolve(pc->pmat,x,y);CHKERRQ(ierr);}
222075768bcSBarry Smith   else              {ierr = MatSolve(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr);}
2239b54502bSHong Zhang   PetscFunctionReturn(0);
2249b54502bSHong Zhang }
2259b54502bSHong Zhang 
2269b54502bSHong Zhang #undef __FUNCT__
2279b54502bSHong Zhang #define __FUNCT__ "PCApplyTranspose_Cholesky"
2289b54502bSHong Zhang static PetscErrorCode PCApplyTranspose_Cholesky(PC pc,Vec x,Vec y)
2299b54502bSHong Zhang {
2309b54502bSHong Zhang   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;
2319b54502bSHong Zhang   PetscErrorCode ierr;
2329b54502bSHong Zhang 
2339b54502bSHong Zhang   PetscFunctionBegin;
2349b54502bSHong Zhang   if (dir->inplace) {ierr = MatSolveTranspose(pc->pmat,x,y);CHKERRQ(ierr);}
235075768bcSBarry Smith   else              {ierr = MatSolveTranspose(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr);}
2369b54502bSHong Zhang   PetscFunctionReturn(0);
2379b54502bSHong Zhang }
2389b54502bSHong Zhang 
2399b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/
2409b54502bSHong Zhang 
2419b54502bSHong Zhang EXTERN_C_BEGIN
2429b54502bSHong Zhang #undef __FUNCT__
2432401956bSBarry Smith #define __FUNCT__ "PCFactorSetUseInPlace_Cholesky"
2442401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseInPlace_Cholesky(PC pc)
2459b54502bSHong Zhang {
2469b54502bSHong Zhang   PC_Cholesky *dir;
2479b54502bSHong Zhang 
2489b54502bSHong Zhang   PetscFunctionBegin;
2499b54502bSHong Zhang   dir = (PC_Cholesky*)pc->data;
2509b54502bSHong Zhang   dir->inplace = PETSC_TRUE;
2519b54502bSHong Zhang   PetscFunctionReturn(0);
2529b54502bSHong Zhang }
2539b54502bSHong Zhang EXTERN_C_END
2549b54502bSHong Zhang 
2559b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/
2569b54502bSHong Zhang 
2579b54502bSHong Zhang #undef __FUNCT__
2582401956bSBarry Smith #define __FUNCT__ "PCFactorSetReuseOrdering"
2599b54502bSHong Zhang /*@
2602401956bSBarry Smith    PCFactorSetReuseOrdering - When similar matrices are factored, this
2619b54502bSHong Zhang    causes the ordering computed in the first factor to be used for all
2629b54502bSHong Zhang    following factors.
2639b54502bSHong Zhang 
2649b54502bSHong Zhang    Collective on PC
2659b54502bSHong Zhang 
2669b54502bSHong Zhang    Input Parameters:
2679b54502bSHong Zhang +  pc - the preconditioner context
2689b54502bSHong Zhang -  flag - PETSC_TRUE to reuse else PETSC_FALSE
2699b54502bSHong Zhang 
2709b54502bSHong Zhang    Options Database Key:
2712401956bSBarry Smith .  -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering()
2729b54502bSHong Zhang 
2739b54502bSHong Zhang    Level: intermediate
2749b54502bSHong Zhang 
2759b54502bSHong Zhang .keywords: PC, levels, reordering, factorization, incomplete, LU
2769b54502bSHong Zhang 
2772401956bSBarry Smith .seealso: PCFactorSetReuseFill()
2789b54502bSHong Zhang @*/
2792401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseOrdering(PC pc,PetscTruth flag)
2809b54502bSHong Zhang {
2819b54502bSHong Zhang   PetscErrorCode ierr,(*f)(PC,PetscTruth);
2829b54502bSHong Zhang 
2839b54502bSHong Zhang   PetscFunctionBegin;
2849b54502bSHong Zhang   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
2852401956bSBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetReuseOrdering_C",(void (**)(void))&f);CHKERRQ(ierr);
2869b54502bSHong Zhang   if (f) {
2879b54502bSHong Zhang     ierr = (*f)(pc,flag);CHKERRQ(ierr);
2889b54502bSHong Zhang   }
2899b54502bSHong Zhang   PetscFunctionReturn(0);
2909b54502bSHong Zhang }
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()
298c7393fdbSBarry Smith .  -pc_factor_mat_solver_package - Actives PCFactorSetMatSolverPackage() to choose the direct solver, like spooles
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
302*145b9266SHong 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(),
316*145b9266SHong Zhang            PCFactorSetFill(), PCFactorSetShiftNonzero(), PCFactorSetShiftType(), PCFactorSetShiftAmount()
3178ff23777SHong Zhang 	   PCFactorSetUseInPlace(), PCFactorSetMatOrderingType()
3189b54502bSHong Zhang 
3199b54502bSHong Zhang M*/
3209b54502bSHong Zhang 
3219b54502bSHong Zhang EXTERN_C_BEGIN
3229b54502bSHong Zhang #undef __FUNCT__
3239b54502bSHong Zhang #define __FUNCT__ "PCCreate_Cholesky"
324dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_Cholesky(PC pc)
3259b54502bSHong Zhang {
3269b54502bSHong Zhang   PetscErrorCode ierr;
3279b54502bSHong Zhang   PC_Cholesky    *dir;
3289b54502bSHong Zhang 
3299b54502bSHong Zhang   PetscFunctionBegin;
33038f2d2fdSLisandro Dalcin   ierr = PetscNewLog(pc,PC_Cholesky,&dir);CHKERRQ(ierr);
3319b54502bSHong Zhang 
332075768bcSBarry Smith   ((PC_Factor*)dir)->fact                   = 0;
3339b54502bSHong Zhang   dir->inplace                = PETSC_FALSE;
334075768bcSBarry Smith   ierr = MatFactorInfoInitialize(&((PC_Factor*)dir)->info);CHKERRQ(ierr);
335075768bcSBarry Smith   ((PC_Factor*)dir)->info.fill          = 5.0;
336915743fcSHong Zhang   ((PC_Factor*)dir)->info.shifttype     = MAT_SHIFT_NONE;
337915743fcSHong Zhang   ((PC_Factor*)dir)->info.shiftamount   = 0.0;
338075768bcSBarry Smith   ((PC_Factor*)dir)->info.pivotinblocks = 1.0;
3399b54502bSHong Zhang   dir->col                    = 0;
3409b54502bSHong Zhang   dir->row                    = 0;
341075768bcSBarry Smith   ierr = PetscStrallocpy(MATORDERING_NATURAL,&((PC_Factor*)dir)->ordering);CHKERRQ(ierr);
342075768bcSBarry Smith   ierr = PetscStrallocpy(MAT_SOLVER_PETSC,&((PC_Factor*)dir)->solvertype);CHKERRQ(ierr);
3439b54502bSHong Zhang   dir->reusefill        = PETSC_FALSE;
3449b54502bSHong Zhang   dir->reuseordering    = PETSC_FALSE;
3459b54502bSHong Zhang   pc->data              = (void*)dir;
3469b54502bSHong Zhang 
3479b54502bSHong Zhang   pc->ops->destroy           = PCDestroy_Cholesky;
3489b54502bSHong Zhang   pc->ops->apply             = PCApply_Cholesky;
3499b54502bSHong Zhang   pc->ops->applytranspose    = PCApplyTranspose_Cholesky;
3509b54502bSHong Zhang   pc->ops->setup             = PCSetUp_Cholesky;
3519b54502bSHong Zhang   pc->ops->setfromoptions    = PCSetFromOptions_Cholesky;
3529b54502bSHong Zhang   pc->ops->view              = PCView_Cholesky;
3539b54502bSHong Zhang   pc->ops->applyrichardson   = 0;
35485317021SBarry Smith   pc->ops->getfactoredmatrix = PCFactorGetMatrix_Factor;
3559b54502bSHong Zhang 
35685317021SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatSolverPackage_C","PCFactorSetMatSolverPackage_Factor",
35785317021SBarry Smith                     PCFactorSetMatSolverPackage_Factor);CHKERRQ(ierr);
3587112b564SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorGetMatSolverPackage_C","PCFactorGetMatSolverPackage_Factor",
3597112b564SBarry Smith                     PCFactorGetMatSolverPackage_Factor);CHKERRQ(ierr);
36085317021SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_Factor",
36185317021SBarry Smith                     PCFactorSetZeroPivot_Factor);CHKERRQ(ierr);
362d90ac83dSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftType_C","PCFactorSetShiftType_Factor",
363d90ac83dSHong Zhang                     PCFactorSetShiftType_Factor);CHKERRQ(ierr);
364d90ac83dSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftAmount_C","PCFactorSetShiftAmount_Factor",
365d90ac83dSHong Zhang                     PCFactorSetShiftAmount_Factor);CHKERRQ(ierr);
36685317021SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetFill_C","PCFactorSetFill_Factor",
36785317021SBarry Smith                     PCFactorSetFill_Factor);CHKERRQ(ierr);
3682401956bSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUseInPlace_C","PCFactorSetUseInPlace_Cholesky",
3692401956bSBarry Smith                     PCFactorSetUseInPlace_Cholesky);CHKERRQ(ierr);
37085317021SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatOrderingType_C","PCFactorSetMatOrderingType_Factor",
37185317021SBarry Smith                     PCFactorSetMatOrderingType_Factor);CHKERRQ(ierr);
3722401956bSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseOrdering_C","PCFactorSetReuseOrdering_Cholesky",
3732401956bSBarry Smith                     PCFactorSetReuseOrdering_Cholesky);CHKERRQ(ierr);
3742401956bSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseFill_C","PCFactorSetReuseFill_Cholesky",
3752401956bSBarry Smith                     PCFactorSetReuseFill_Cholesky);CHKERRQ(ierr);
3769b54502bSHong Zhang   PetscFunctionReturn(0);
3779b54502bSHong Zhang }
3789b54502bSHong Zhang EXTERN_C_END
379