xref: /petsc/src/ksp/pc/impls/factor/factimpl.c (revision b7c853c42524e99073a4955bd48eeb29623ae6ef)
185317021SBarry Smith #define PETSCKSP_DLL
285317021SBarry Smith 
37c4f633dSBarry Smith #include "../src/ksp/pc/impls/factor/factor.h"     /*I "petscpc.h"  I*/
485317021SBarry Smith 
585317021SBarry Smith /* ------------------------------------------------------------------------------------------*/
685317021SBarry Smith 
785317021SBarry Smith EXTERN_C_BEGIN
885317021SBarry Smith #undef __FUNCT__
985317021SBarry Smith #define __FUNCT__ "PCFactorSetZeroPivot_Factor"
1085317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetZeroPivot_Factor(PC pc,PetscReal z)
1185317021SBarry Smith {
1285317021SBarry Smith   PC_Factor *ilu = (PC_Factor*)pc->data;
1385317021SBarry Smith 
1485317021SBarry Smith   PetscFunctionBegin;
1585317021SBarry Smith   ilu->info.zeropivot = z;
1685317021SBarry Smith   PetscFunctionReturn(0);
1785317021SBarry Smith }
1885317021SBarry Smith EXTERN_C_END
1985317021SBarry Smith 
2085317021SBarry Smith EXTERN_C_BEGIN
2185317021SBarry Smith #undef __FUNCT__
2285317021SBarry Smith #define __FUNCT__ "PCFactorSetShiftNonzero_Factor"
2385317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftNonzero_Factor(PC pc,PetscReal shift)
2485317021SBarry Smith {
2585317021SBarry Smith   PC_Factor *dir = (PC_Factor*)pc->data;
2685317021SBarry Smith 
2785317021SBarry Smith   PetscFunctionBegin;
2885317021SBarry Smith   if (shift == (PetscReal) PETSC_DECIDE) {
2985317021SBarry Smith     dir->info.shiftnz = 1.e-12;
3085317021SBarry Smith   } else {
3185317021SBarry Smith     dir->info.shiftnz = shift;
3285317021SBarry Smith   }
3385317021SBarry Smith   PetscFunctionReturn(0);
3485317021SBarry Smith }
3585317021SBarry Smith EXTERN_C_END
3685317021SBarry Smith 
3785317021SBarry Smith EXTERN_C_BEGIN
3885317021SBarry Smith #undef __FUNCT__
3985317021SBarry Smith #define __FUNCT__ "PCFactorSetShiftPd_Factor"
4085317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftPd_Factor(PC pc,PetscTruth shift)
4185317021SBarry Smith {
4285317021SBarry Smith   PC_Factor *dir = (PC_Factor*)pc->data;
4385317021SBarry Smith 
4485317021SBarry Smith   PetscFunctionBegin;
4585317021SBarry Smith   if (shift) {
4685317021SBarry Smith     dir->info.shiftpd = 1.0;
4785317021SBarry Smith   } else {
4885317021SBarry Smith     dir->info.shiftpd = 0.0;
4985317021SBarry Smith   }
5085317021SBarry Smith   PetscFunctionReturn(0);
5185317021SBarry Smith }
5285317021SBarry Smith EXTERN_C_END
5385317021SBarry Smith 
5485317021SBarry Smith EXTERN_C_BEGIN
5585317021SBarry Smith #undef __FUNCT__
56*b7c853c4SBarry Smith #define __FUNCT__ "PCFactorSetDropTolerance_Factor"
57*b7c853c4SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetDropTolerance_Factor(PC pc,PetscReal dt,PetscReal dtcol,PetscInt dtcount)
5885317021SBarry Smith {
5985317021SBarry Smith   PC_Factor         *ilu = (PC_Factor*)pc->data;
6085317021SBarry Smith 
6185317021SBarry Smith   PetscFunctionBegin;
6285317021SBarry Smith   if (pc->setupcalled && (!ilu->info.usedt || ((PC_Factor*)ilu)->info.dt != dt || ((PC_Factor*)ilu)->info.dtcol != dtcol || ((PC_Factor*)ilu)->info.dtcount != dtcount)) {
6385317021SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Cannot change tolerance after use");
6485317021SBarry Smith   }
6585317021SBarry Smith   ilu->info.usedt   = PETSC_TRUE;
6685317021SBarry Smith   ilu->info.dt      = dt;
6785317021SBarry Smith   ilu->info.dtcol   = dtcol;
6885317021SBarry Smith   ilu->info.dtcount = dtcount;
6985317021SBarry Smith   ilu->info.fill    = PETSC_DEFAULT;
7085317021SBarry Smith   PetscFunctionReturn(0);
7185317021SBarry Smith }
7285317021SBarry Smith EXTERN_C_END
7385317021SBarry Smith 
7485317021SBarry Smith EXTERN_C_BEGIN
7585317021SBarry Smith #undef __FUNCT__
7685317021SBarry Smith #define __FUNCT__ "PCFactorSetFill_Factor"
7785317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetFill_Factor(PC pc,PetscReal fill)
7885317021SBarry Smith {
7985317021SBarry Smith   PC_Factor *dir = (PC_Factor*)pc->data;
8085317021SBarry Smith 
8185317021SBarry Smith   PetscFunctionBegin;
8285317021SBarry Smith   dir->info.fill = fill;
8385317021SBarry Smith   PetscFunctionReturn(0);
8485317021SBarry Smith }
8585317021SBarry Smith EXTERN_C_END
8685317021SBarry Smith 
8785317021SBarry Smith EXTERN_C_BEGIN
8885317021SBarry Smith #undef __FUNCT__
8985317021SBarry Smith #define __FUNCT__ "PCFactorSetMatOrderingType_Factor"
9085317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatOrderingType_Factor(PC pc,const MatOrderingType ordering)
9185317021SBarry Smith {
9285317021SBarry Smith   PC_Factor      *dir = (PC_Factor*)pc->data;
9385317021SBarry Smith   PetscErrorCode ierr;
9485317021SBarry Smith   PetscTruth     flg;
9585317021SBarry Smith 
9685317021SBarry Smith   PetscFunctionBegin;
9785317021SBarry Smith   if (!pc->setupcalled) {
9885317021SBarry Smith      ierr = PetscStrfree(dir->ordering);CHKERRQ(ierr);
9985317021SBarry Smith      ierr = PetscStrallocpy(ordering,&dir->ordering);CHKERRQ(ierr);
10085317021SBarry Smith   } else {
10185317021SBarry Smith     ierr = PetscStrcmp(dir->ordering,ordering,&flg);CHKERRQ(ierr);
10285317021SBarry Smith     if (!flg) {
10385317021SBarry Smith       SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Cannot change ordering after use");
10485317021SBarry Smith     }
10585317021SBarry Smith   }
10685317021SBarry Smith   PetscFunctionReturn(0);
10785317021SBarry Smith }
10885317021SBarry Smith EXTERN_C_END
10985317021SBarry Smith 
11085317021SBarry Smith EXTERN_C_BEGIN
11185317021SBarry Smith #undef __FUNCT__
11285317021SBarry Smith #define __FUNCT__ "PCFactorSetLevels_Factor"
11385317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetLevels_Factor(PC pc,PetscInt levels)
11485317021SBarry Smith {
11585317021SBarry Smith   PC_Factor      *ilu = (PC_Factor*)pc->data;
11685317021SBarry Smith 
11785317021SBarry Smith   PetscFunctionBegin;
11885317021SBarry Smith   if (!pc->setupcalled) {
11985317021SBarry Smith     ilu->info.levels = levels;
12085317021SBarry Smith   } else if (ilu->info.usedt || ilu->info.levels != levels) {
12185317021SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Cannot change levels after use");
12285317021SBarry Smith   }
12385317021SBarry Smith   PetscFunctionReturn(0);
12485317021SBarry Smith }
12585317021SBarry Smith EXTERN_C_END
12685317021SBarry Smith 
12785317021SBarry Smith EXTERN_C_BEGIN
12885317021SBarry Smith #undef __FUNCT__
12985317021SBarry Smith #define __FUNCT__ "PCFactorSetAllowDiagonalFill_Factor"
13085317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetAllowDiagonalFill_Factor(PC pc)
13185317021SBarry Smith {
13285317021SBarry Smith   PC_Factor *dir = (PC_Factor*)pc->data;
13385317021SBarry Smith 
13485317021SBarry Smith   PetscFunctionBegin;
13585317021SBarry Smith   dir->info.diagonal_fill = 1;
13685317021SBarry Smith   PetscFunctionReturn(0);
13785317021SBarry Smith }
13885317021SBarry Smith EXTERN_C_END
13985317021SBarry Smith 
14085317021SBarry Smith 
14185317021SBarry Smith /* ------------------------------------------------------------------------------------------*/
14285317021SBarry Smith 
14385317021SBarry Smith EXTERN_C_BEGIN
14485317021SBarry Smith #undef __FUNCT__
14585317021SBarry Smith #define __FUNCT__ "PCFactorSetPivotInBlocks_Factor"
14685317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivotInBlocks_Factor(PC pc,PetscTruth pivot)
14785317021SBarry Smith {
14885317021SBarry Smith   PC_Factor *dir = (PC_Factor*)pc->data;
14985317021SBarry Smith 
15085317021SBarry Smith   PetscFunctionBegin;
15185317021SBarry Smith   dir->info.pivotinblocks = pivot ? 1.0 : 0.0;
15285317021SBarry Smith   PetscFunctionReturn(0);
15385317021SBarry Smith }
15485317021SBarry Smith EXTERN_C_END
15585317021SBarry Smith 
15685317021SBarry Smith EXTERN_C_BEGIN
15785317021SBarry Smith #undef __FUNCT__
15885317021SBarry Smith #define __FUNCT__ "PCFactorSetShiftInBlocks_Factor"
15985317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftInBlocks_Factor(PC pc,PetscReal shift)
16085317021SBarry Smith {
16185317021SBarry Smith   PC_Factor *dir = (PC_Factor*)pc->data;
16285317021SBarry Smith 
16385317021SBarry Smith   PetscFunctionBegin;
16485317021SBarry Smith   if (shift == PETSC_DEFAULT) {
16585317021SBarry Smith     dir->info.shiftinblocks = 1.e-12;
16685317021SBarry Smith   } else {
16785317021SBarry Smith     dir->info.shiftinblocks = shift;
16885317021SBarry Smith   }
16985317021SBarry Smith   PetscFunctionReturn(0);
17085317021SBarry Smith }
17185317021SBarry Smith EXTERN_C_END
17285317021SBarry Smith 
17385317021SBarry Smith 
17485317021SBarry Smith EXTERN_C_BEGIN
17585317021SBarry Smith #undef __FUNCT__
17685317021SBarry Smith #define __FUNCT__ "PCFactorGetMatrix_Factor"
17785317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorGetMatrix_Factor(PC pc,Mat *mat)
17885317021SBarry Smith {
17985317021SBarry Smith   PC_Factor *ilu = (PC_Factor*)pc->data;
18085317021SBarry Smith 
18185317021SBarry Smith   PetscFunctionBegin;
18285317021SBarry Smith   if (!ilu->fact) SETERRQ(PETSC_ERR_ORDER,"Matrix not yet factored; call after KSPSetUp() or PCSetUp()");
18385317021SBarry Smith   *mat = ilu->fact;
18485317021SBarry Smith   PetscFunctionReturn(0);
18585317021SBarry Smith }
18685317021SBarry Smith EXTERN_C_END
18785317021SBarry Smith 
18885317021SBarry Smith EXTERN_C_BEGIN
18985317021SBarry Smith #undef __FUNCT__
19085317021SBarry Smith #define __FUNCT__ "PCFactorSetMatSolverPackage_Factor"
19185317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatSolverPackage_Factor(PC pc,const MatSolverPackage stype)
19285317021SBarry Smith {
19385317021SBarry Smith   PetscErrorCode ierr;
19485317021SBarry Smith   PC_Factor      *lu = (PC_Factor*)pc->data;
19585317021SBarry Smith 
19685317021SBarry Smith   PetscFunctionBegin;
1977112b564SBarry Smith   if (lu->fact) {
19885317021SBarry Smith     const MatSolverPackage ltype;
19985317021SBarry Smith     PetscTruth             flg;
20085317021SBarry Smith     ierr = MatFactorGetSolverPackage(lu->fact,&ltype);CHKERRQ(ierr);
20185317021SBarry Smith     ierr = PetscStrcmp(stype,ltype,&flg);CHKERRQ(ierr);
20285317021SBarry Smith     if (!flg) {
20385317021SBarry Smith       SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Cannot change solver matrix package after PC has been setup or used");
20485317021SBarry Smith     }
20585317021SBarry Smith   }
20685317021SBarry Smith   ierr = PetscStrfree(lu->solvertype);CHKERRQ(ierr);
20785317021SBarry Smith   ierr = PetscStrallocpy(stype,&lu->solvertype);CHKERRQ(ierr);
20885317021SBarry Smith   PetscFunctionReturn(0);
20985317021SBarry Smith }
21085317021SBarry Smith EXTERN_C_END
21185317021SBarry Smith 
21285317021SBarry Smith EXTERN_C_BEGIN
21385317021SBarry Smith #undef __FUNCT__
2147112b564SBarry Smith #define __FUNCT__ "PCFactorGetMatSolverPackage_Factor"
2157112b564SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorGetMatSolverPackage_Factor(PC pc,const MatSolverPackage *stype)
2167112b564SBarry Smith {
2177112b564SBarry Smith   PC_Factor      *lu = (PC_Factor*)pc->data;
2187112b564SBarry Smith 
2197112b564SBarry Smith   PetscFunctionBegin;
2207112b564SBarry Smith   *stype = lu->solvertype;
2217112b564SBarry Smith   PetscFunctionReturn(0);
2227112b564SBarry Smith }
2237112b564SBarry Smith EXTERN_C_END
2247112b564SBarry Smith 
2257112b564SBarry Smith EXTERN_C_BEGIN
2267112b564SBarry Smith #undef __FUNCT__
2278ff23777SHong Zhang #define __FUNCT__ "PCFactorSetColumnPivot_Factor"
2288ff23777SHong Zhang PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetColumnPivot_Factor(PC pc,PetscReal dtcol)
22985317021SBarry Smith {
23085317021SBarry Smith   PC_Factor       *dir = (PC_Factor*)pc->data;
23185317021SBarry Smith 
23285317021SBarry Smith   PetscFunctionBegin;
23385317021SBarry Smith   if (dtcol < 0.0 || dtcol > 1.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Column pivot tolerance is %G must be between 0 and 1",dtcol);
23485317021SBarry Smith   dir->info.dtcol = dtcol;
23585317021SBarry Smith   PetscFunctionReturn(0);
23685317021SBarry Smith }
23785317021SBarry Smith EXTERN_C_END
2388ff23777SHong Zhang 
2398ff23777SHong Zhang EXTERN_C_BEGIN
2408ff23777SHong Zhang #undef __FUNCT__
2418ff23777SHong Zhang #define __FUNCT__ "PCSetFromOptions_Factor"
2428ff23777SHong Zhang PetscErrorCode PETSCKSP_DLLEXPORT PCSetFromOptions_Factor(PC pc)
2438ff23777SHong Zhang {
2448ff23777SHong Zhang   PC_Factor       *factor = (PC_Factor*)pc->data;
2458ff23777SHong Zhang   PetscErrorCode  ierr;
2468ff23777SHong Zhang   PetscTruth      flg = PETSC_FALSE,set;
2478ff23777SHong Zhang   char            tname[256], solvertype[64];
2488ff23777SHong Zhang   PetscFList      ordlist;
2498ff23777SHong Zhang 
2508ff23777SHong Zhang   PetscFunctionBegin;
2518ff23777SHong Zhang   if (!MatOrderingRegisterAllCalled) {ierr = MatOrderingRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
2528ff23777SHong Zhang   ierr = PetscOptionsTruth("-pc_factor_in_place","Form factored matrix in the same memory as the matrix","PCFactorSetUseInPlace",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
2538ff23777SHong Zhang     if (flg) {
2548ff23777SHong Zhang       ierr = PCFactorSetUseInPlace(pc);CHKERRQ(ierr);
2558ff23777SHong Zhang     }
2568ff23777SHong Zhang     ierr = PetscOptionsReal("-pc_factor_fill","Expected non-zeros in factored matrix","PCFactorSetFill",((PC_Factor*)factor)->info.fill,&((PC_Factor*)factor)->info.fill,0);CHKERRQ(ierr);
2578ff23777SHong Zhang 
2588ff23777SHong Zhang     flg  = PETSC_FALSE;
2598ff23777SHong Zhang     ierr = PetscOptionsName("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",&flg);CHKERRQ(ierr);
2608ff23777SHong Zhang     if (flg) {
2618ff23777SHong Zhang       ierr = PCFactorSetShiftNonzero(pc,(PetscReal) PETSC_DECIDE);CHKERRQ(ierr);
2628ff23777SHong Zhang     }
2638ff23777SHong Zhang     ierr = PetscOptionsReal("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",((PC_Factor*)factor)->info.shiftnz,&((PC_Factor*)factor)->info.shiftnz,0);CHKERRQ(ierr);
2648ff23777SHong Zhang     flg  = PETSC_FALSE;
2658ff23777SHong Zhang     ierr = PetscOptionsTruth("-pc_factor_shift_positive_definite","Manteuffel shift applied to diagonal","PCFactorSetShiftPd",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
2668ff23777SHong Zhang     if (flg) {
2678ff23777SHong Zhang       ierr = PCFactorSetShiftPd(pc,PETSC_TRUE);CHKERRQ(ierr);
2688ff23777SHong Zhang     }
2698ff23777SHong Zhang     flg  = PETSC_FALSE;
2708ff23777SHong Zhang     ierr = PetscOptionsName("-pc_factor_shift_in_blocks","Shift added to diagonal dense blocks","PCFactorSetShiftInBlocks",&flg);CHKERRQ(ierr);
2718ff23777SHong Zhang     if (flg) {
2728ff23777SHong Zhang       ierr = PCFactorSetShiftInBlocks(pc,(PetscReal) PETSC_DEFAULT);CHKERRQ(ierr);
2738ff23777SHong Zhang     }
2748ff23777SHong Zhang     ierr = PetscOptionsReal("-pc_factor_shift_in_blocks","Shift added to diagonal dense blocks","PCFactorSetShiftInBlocks",((PC_Factor*)factor)->info.shiftinblocks,&((PC_Factor*)factor)->info.shiftinblocks,0);CHKERRQ(ierr);
2758ff23777SHong Zhang 
2768ff23777SHong Zhang     ierr = PetscOptionsReal("-pc_factor_zeropivot","Pivot is considered zero if less than","PCFactorSetZeroPivot",((PC_Factor*)factor)->info.zeropivot,&((PC_Factor*)factor)->info.zeropivot,0);CHKERRQ(ierr);
2778ff23777SHong Zhang     ierr = PetscOptionsReal("-pc_factor_column_pivot","Column pivot tolerance (used only for some factorization)","PCFactorSetColumnPivot",((PC_Factor*)factor)->info.dtcol,&((PC_Factor*)factor)->info.dtcol,&flg);CHKERRQ(ierr);
2788ff23777SHong Zhang 
2798ff23777SHong Zhang     flg = ((PC_Factor*)factor)->info.pivotinblocks ? PETSC_TRUE : PETSC_FALSE;
2808ff23777SHong Zhang     ierr = PetscOptionsTruth("-pc_factor_pivot_in_blocks","Pivot inside matrix dense blocks for BAIJ and SBAIJ","PCFactorSetPivotInBlocks",flg,&flg,&set);CHKERRQ(ierr);
2818ff23777SHong Zhang     if (set) {
2828ff23777SHong Zhang       ierr = PCFactorSetPivotInBlocks(pc,flg);CHKERRQ(ierr);
2838ff23777SHong Zhang     }
2848ff23777SHong Zhang 
2858ff23777SHong Zhang     flg  = PETSC_FALSE;
2868ff23777SHong Zhang     ierr = PetscOptionsTruth("-pc_factor_reuse_fill","Use fill from previous factorization","PCFactorSetReuseFill",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
2878ff23777SHong Zhang     if (flg) {
2888ff23777SHong Zhang       ierr = PCFactorSetReuseFill(pc,PETSC_TRUE);CHKERRQ(ierr);
2898ff23777SHong Zhang     }
2908ff23777SHong Zhang     flg  = PETSC_FALSE;
2918ff23777SHong Zhang     ierr = PetscOptionsTruth("-pc_factor_reuse_ordering","Reuse ordering from previous factorization","PCFactorSetReuseOrdering",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
2928ff23777SHong Zhang     if (flg) {
2938ff23777SHong Zhang       ierr = PCFactorSetReuseOrdering(pc,PETSC_TRUE);CHKERRQ(ierr);
2948ff23777SHong Zhang     }
2958ff23777SHong Zhang 
2968ff23777SHong Zhang     ierr = MatGetOrderingList(&ordlist);CHKERRQ(ierr);
2978ff23777SHong Zhang     ierr = PetscOptionsList("-pc_factor_mat_ordering_type","Reordering to reduce nonzeros in factored matrix","PCFactorSetMatOrderingType",ordlist,((PC_Factor*)factor)->ordering,tname,256,&flg);CHKERRQ(ierr);
2988ff23777SHong Zhang     if (flg) {
2998ff23777SHong Zhang       ierr = PCFactorSetMatOrderingType(pc,tname);CHKERRQ(ierr);
3008ff23777SHong Zhang     }
3018ff23777SHong Zhang 
3028ff23777SHong Zhang     /* maybe should have MatGetSolverTypes(Mat,&list) like the ordering list */
3038ff23777SHong Zhang     ierr = PetscOptionsString("-pc_factor_mat_solver_package","Specific direct solver to use","MatGetFactor",((PC_Factor*)factor)->solvertype,solvertype,64,&flg);CHKERRQ(ierr);
3048ff23777SHong Zhang     if (flg) {
3058ff23777SHong Zhang       ierr = PCFactorSetMatSolverPackage(pc,solvertype);CHKERRQ(ierr);
3068ff23777SHong Zhang     }
3078ff23777SHong Zhang 
3088ff23777SHong Zhang     /*
3098ff23777SHong Zhang     ierr = PetscOptionsReal("-pc_factor_column_pivot","Column pivot tolerance (used only for some factorization)","PCFactorSetColumnPivot",((PC_Factor*)factor)->info.dtcol,&((PC_Factor*)factor)->info.dtcol,&flg);CHKERRQ(ierr);
3108ff23777SHong Zhang 
3118ff23777SHong Zhang     flg = ((PC_Factor*)factor)->info.pivotinblocks ? PETSC_TRUE : PETSC_FALSE;
3128ff23777SHong Zhang     ierr = PetscOptionsTruth("-pc_factor_pivot_in_blocks","Pivot inside matrix dense blocks for BAIJ and SBAIJ","PCFactorSetPivotInBlocks",flg,&flg,&set);CHKERRQ(ierr);
3138ff23777SHong Zhang     if (set) {
3148ff23777SHong Zhang       ierr = PCFactorSetPivotInBlocks(pc,flg);CHKERRQ(ierr);
3158ff23777SHong Zhang     }
3168ff23777SHong Zhang     */
3178ff23777SHong Zhang   PetscFunctionReturn(0);
3188ff23777SHong Zhang }
3198ff23777SHong Zhang EXTERN_C_END
320