xref: /petsc/src/ksp/pc/impls/factor/factor.c (revision 7c4f633dc6bb6149cca88d301ead35a99e103cbb)
1dba47a55SKris Buschelman #define PETSCKSP_DLL
25e8efad8SHong Zhang 
3*7c4f633dSBarry Smith #include "../src/ksp/pc/impls/factor/factor.h"  /*I "petscpc.h" I*/
45e8efad8SHong Zhang 
5ee45ca4aSHong Zhang #undef __FUNCT__
6ee45ca4aSHong Zhang #define __FUNCT__ "PCFactorSetZeroPivot"
7ee45ca4aSHong Zhang /*@
8ee45ca4aSHong Zhang    PCFactorSetZeroPivot - Sets the size at which smaller pivots are declared to be zero
9ee45ca4aSHong Zhang 
10ee45ca4aSHong Zhang    Collective on PC
11ee45ca4aSHong Zhang 
12ee45ca4aSHong Zhang    Input Parameters:
13afaefe49SHong Zhang +  pc - the preconditioner context
14afaefe49SHong Zhang -  zero - all pivots smaller than this will be considered zero
15ee45ca4aSHong Zhang 
16ee45ca4aSHong Zhang    Options Database Key:
17ee45ca4aSHong Zhang .  -pc_factor_zeropivot <zero> - Sets tolerance for what is considered a zero pivot
18ee45ca4aSHong Zhang 
19ee45ca4aSHong Zhang    Level: intermediate
20ee45ca4aSHong Zhang 
21ee45ca4aSHong Zhang .keywords: PC, set, factorization, direct, fill
22ee45ca4aSHong Zhang 
23ee45ca4aSHong Zhang .seealso: PCFactorSetShiftNonzero(), PCFactorSetShiftPd()
24ee45ca4aSHong Zhang @*/
25dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetZeroPivot(PC pc,PetscReal zero)
26ee45ca4aSHong Zhang {
27afaefe49SHong Zhang   PetscErrorCode ierr,(*f)(PC,PetscReal);
28afaefe49SHong Zhang 
29ee45ca4aSHong Zhang   PetscFunctionBegin;
30afaefe49SHong Zhang   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
31afaefe49SHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetZeroPivot_C",(void (**)(void))&f);CHKERRQ(ierr);
32afaefe49SHong Zhang   if (f) {
33afaefe49SHong Zhang     ierr = (*f)(pc,zero);CHKERRQ(ierr);
34afaefe49SHong Zhang   }
35ee45ca4aSHong Zhang   PetscFunctionReturn(0);
36ee45ca4aSHong Zhang }
37ee45ca4aSHong Zhang 
385e8efad8SHong Zhang #undef __FUNCT__
395e8efad8SHong Zhang #define __FUNCT__ "PCFactorSetShiftNonzero"
405e8efad8SHong Zhang /*@
415e8efad8SHong Zhang    PCFactorSetShiftNonzero - adds this quantity to the diagonal of the matrix during
425e8efad8SHong Zhang      numerical factorization, thus the matrix has nonzero pivots
435e8efad8SHong Zhang 
445e8efad8SHong Zhang    Collective on PC
455e8efad8SHong Zhang 
465e8efad8SHong Zhang    Input Parameters:
47afaefe49SHong Zhang +  pc - the preconditioner context
48afaefe49SHong Zhang -  shift - amount of shift
495e8efad8SHong Zhang 
505e8efad8SHong Zhang    Options Database Key:
519f95998fSHong Zhang .  -pc_factor_shift_nonzero <shift> - Sets shift amount or PETSC_DECIDE for the default
525e8efad8SHong Zhang 
535e8efad8SHong Zhang    Note: If 0.0 is given, then no shift is used. If a diagonal element is classified as a zero
545e8efad8SHong Zhang          pivot, then the shift is doubled until this is alleviated.
555e8efad8SHong Zhang 
565e8efad8SHong Zhang    Level: intermediate
575e8efad8SHong Zhang 
585e8efad8SHong Zhang .keywords: PC, set, factorization, direct, fill
595e8efad8SHong Zhang 
60ee45ca4aSHong Zhang .seealso: PCFactorSetZeroPivot(), PCFactorSetShiftPd()
615e8efad8SHong Zhang @*/
62dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftNonzero(PC pc,PetscReal shift)
635e8efad8SHong Zhang {
64afaefe49SHong Zhang   PetscErrorCode ierr,(*f)(PC,PetscReal);
65afaefe49SHong Zhang 
665e8efad8SHong Zhang   PetscFunctionBegin;
67afaefe49SHong Zhang   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
68afaefe49SHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetShiftNonzero_C",(void (**)(void))&f);CHKERRQ(ierr);
69afaefe49SHong Zhang   if (f) {
70afaefe49SHong Zhang     ierr = (*f)(pc,shift);CHKERRQ(ierr);
715e8efad8SHong Zhang   }
725e8efad8SHong Zhang   PetscFunctionReturn(0);
735e8efad8SHong Zhang }
740a29876aSHong Zhang 
750a29876aSHong Zhang #undef __FUNCT__
760a29876aSHong Zhang #define __FUNCT__ "PCFactorSetShiftPd"
770a29876aSHong Zhang /*@
780a29876aSHong Zhang    PCFactorSetShiftPd - specify whether to use Manteuffel shifting.
790a29876aSHong Zhang    If a matrix factorisation breaks down because of nonpositive pivots,
800a29876aSHong Zhang    adding sufficient identity to the diagonal will remedy this.
810a29876aSHong Zhang    Setting this causes a bisection method to find the minimum shift that
820a29876aSHong Zhang    will lead to a well-defined matrix factor.
830a29876aSHong Zhang 
840a29876aSHong Zhang    Collective on PC
850a29876aSHong Zhang 
860a29876aSHong Zhang    Input parameters:
87afaefe49SHong Zhang +  pc - the preconditioner context
88afaefe49SHong Zhang -  shifting - PETSC_TRUE to set shift else PETSC_FALSE
890a29876aSHong Zhang 
900a29876aSHong Zhang    Options Database Key:
919f95998fSHong Zhang .  -pc_factor_shift_positive_definite [PETSC_TRUE/PETSC_FALSE] - Activate/Deactivate PCFactorSetShiftPd(); the value
929f95998fSHong Zhang    is optional with PETSC_TRUE being the default
930a29876aSHong Zhang 
940a29876aSHong Zhang    Level: intermediate
950a29876aSHong Zhang 
960a29876aSHong Zhang .keywords: PC, indefinite, factorization
970a29876aSHong Zhang 
98ee45ca4aSHong Zhang .seealso: PCFactorSetZeroPivot(), PCFactorSetShiftNonzero()
990a29876aSHong Zhang @*/
100dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftPd(PC pc,PetscTruth shift)
1010a29876aSHong Zhang {
102afaefe49SHong Zhang   PetscErrorCode ierr,(*f)(PC,PetscTruth);
103afaefe49SHong Zhang 
1040a29876aSHong Zhang   PetscFunctionBegin;
105afaefe49SHong Zhang   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
106afaefe49SHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetShiftPd_C",(void (**)(void))&f);CHKERRQ(ierr);
107afaefe49SHong Zhang   if (f) {
108afaefe49SHong Zhang     ierr = (*f)(pc,shift);CHKERRQ(ierr);
109afaefe49SHong Zhang   }
1100a29876aSHong Zhang   PetscFunctionReturn(0);
1110a29876aSHong Zhang }
11285317021SBarry Smith 
11385317021SBarry Smith #undef __FUNCT__
11485317021SBarry Smith #define __FUNCT__ "PCFactorSetUseDropTolerance"
11585317021SBarry Smith /*@
11685317021SBarry Smith    PCFactorSetUseDropTolerance - The preconditioner will use an ILU
11785317021SBarry Smith    based on a drop tolerance.
11885317021SBarry Smith 
11985317021SBarry Smith    Collective on PC
12085317021SBarry Smith 
12185317021SBarry Smith    Input Parameters:
12285317021SBarry Smith +  pc - the preconditioner context
12385317021SBarry Smith .  dt - the drop tolerance, try from 1.e-10 to .1
12485317021SBarry Smith .  dtcol - tolerance for column pivot, good values [0.1 to 0.01]
12585317021SBarry Smith -  maxrowcount - the max number of nonzeros allowed in a row, best value
12685317021SBarry Smith                  depends on the number of nonzeros in row of original matrix
12785317021SBarry Smith 
12885317021SBarry Smith    Options Database Key:
12985317021SBarry Smith .  -pc_factor_use_drop_tolerance <dt,dtcol,maxrowcount> - Sets drop tolerance
13085317021SBarry Smith 
13185317021SBarry Smith    Level: intermediate
13285317021SBarry Smith 
13385317021SBarry Smith     Notes:
13485317021SBarry Smith       This uses the iludt() code of Saad's SPARSKIT package
13585317021SBarry Smith 
13685317021SBarry Smith       There are NO default values for the 3 parameters, you must set them with reasonable values for your
13785317021SBarry Smith       matrix. We don't know how to compute reasonable values.
13885317021SBarry Smith 
13985317021SBarry Smith .keywords: PC, levels, reordering, factorization, incomplete, ILU
14085317021SBarry Smith @*/
14185317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseDropTolerance(PC pc,PetscReal dt,PetscReal dtcol,PetscInt maxrowcount)
14285317021SBarry Smith {
14385317021SBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscReal,PetscReal,PetscInt);
14485317021SBarry Smith 
14585317021SBarry Smith   PetscFunctionBegin;
14685317021SBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
14785317021SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetUseDropTolerance_C",(void (**)(void))&f);CHKERRQ(ierr);
14885317021SBarry Smith   if (f) {
14985317021SBarry Smith     ierr = (*f)(pc,dt,dtcol,maxrowcount);CHKERRQ(ierr);
15085317021SBarry Smith   }
15185317021SBarry Smith   PetscFunctionReturn(0);
15285317021SBarry Smith }
15385317021SBarry Smith 
15485317021SBarry Smith #undef __FUNCT__
15585317021SBarry Smith #define __FUNCT__ "PCFactorSetLevels"
15685317021SBarry Smith /*@
15785317021SBarry Smith    PCFactorSetLevels - Sets the number of levels of fill to use.
15885317021SBarry Smith 
15985317021SBarry Smith    Collective on PC
16085317021SBarry Smith 
16185317021SBarry Smith    Input Parameters:
16285317021SBarry Smith +  pc - the preconditioner context
16385317021SBarry Smith -  levels - number of levels of fill
16485317021SBarry Smith 
16585317021SBarry Smith    Options Database Key:
16685317021SBarry Smith .  -pc_factor_levels <levels> - Sets fill level
16785317021SBarry Smith 
16885317021SBarry Smith    Level: intermediate
16985317021SBarry Smith 
17085317021SBarry Smith .keywords: PC, levels, fill, factorization, incomplete, ILU
17185317021SBarry Smith @*/
17285317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetLevels(PC pc,PetscInt levels)
17385317021SBarry Smith {
17485317021SBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscInt);
17585317021SBarry Smith 
17685317021SBarry Smith   PetscFunctionBegin;
17785317021SBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
17885317021SBarry Smith   if (levels < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"negative levels");
17985317021SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetLevels_C",(void (**)(void))&f);CHKERRQ(ierr);
18085317021SBarry Smith   if (f) {
18185317021SBarry Smith     ierr = (*f)(pc,levels);CHKERRQ(ierr);
18285317021SBarry Smith   }
18385317021SBarry Smith   PetscFunctionReturn(0);
18485317021SBarry Smith }
18585317021SBarry Smith 
18685317021SBarry Smith #undef __FUNCT__
18785317021SBarry Smith #define __FUNCT__ "PCFactorSetAllowDiagonalFill"
18885317021SBarry Smith /*@
18985317021SBarry Smith    PCFactorSetAllowDiagonalFill - Causes all diagonal matrix entries to be
19085317021SBarry Smith    treated as level 0 fill even if there is no non-zero location.
19185317021SBarry Smith 
19285317021SBarry Smith    Collective on PC
19385317021SBarry Smith 
19485317021SBarry Smith    Input Parameters:
19585317021SBarry Smith +  pc - the preconditioner context
19685317021SBarry Smith 
19785317021SBarry Smith    Options Database Key:
19885317021SBarry Smith .  -pc_factor_diagonal_fill
19985317021SBarry Smith 
20085317021SBarry Smith    Notes:
20185317021SBarry Smith    Does not apply with 0 fill.
20285317021SBarry Smith 
20385317021SBarry Smith    Level: intermediate
20485317021SBarry Smith 
20585317021SBarry Smith .keywords: PC, levels, fill, factorization, incomplete, ILU
20685317021SBarry Smith @*/
20785317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetAllowDiagonalFill(PC pc)
20885317021SBarry Smith {
20985317021SBarry Smith   PetscErrorCode ierr,(*f)(PC);
21085317021SBarry Smith 
21185317021SBarry Smith   PetscFunctionBegin;
21285317021SBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
21385317021SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetAllowDiagonalFill_C",(void (**)(void))&f);CHKERRQ(ierr);
21485317021SBarry Smith   if (f) {
21585317021SBarry Smith     ierr = (*f)(pc);CHKERRQ(ierr);
21685317021SBarry Smith   }
21785317021SBarry Smith   PetscFunctionReturn(0);
21885317021SBarry Smith }
21985317021SBarry Smith 
22085317021SBarry Smith #undef __FUNCT__
22185317021SBarry Smith #define __FUNCT__ "PCFactorSetShiftInBlocks"
22285317021SBarry Smith /*@
22385317021SBarry Smith     PCFactorSetShiftInBlocks - Shifts the diagonal of any block if LU on that block detects a zero pivot.
22485317021SBarry Smith 
22585317021SBarry Smith     Collective on PC
22685317021SBarry Smith 
22785317021SBarry Smith     Input Parameters:
22885317021SBarry Smith +   pc - the preconditioner context
22985317021SBarry Smith -   shift - amount to shift or PETSC_DECIDE
23085317021SBarry Smith 
23185317021SBarry Smith     Options Database Key:
23285317021SBarry Smith .   -pc_factor_shift_in_blocks <shift>
23385317021SBarry Smith 
23485317021SBarry Smith     Level: intermediate
23585317021SBarry Smith 
23685317021SBarry Smith .seealso: PCILUSetMatOrdering(), PCFactorSetPivoting(), PCFactorSetShiftInBlocks()
23785317021SBarry Smith @*/
23885317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftInBlocks(PC pc,PetscReal shift)
23985317021SBarry Smith {
24085317021SBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscReal);
24185317021SBarry Smith 
24285317021SBarry Smith   PetscFunctionBegin;
24385317021SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetShiftInBlocks_C",(void (**)(void))&f);CHKERRQ(ierr);
24485317021SBarry Smith   if (f) {
24585317021SBarry Smith     ierr = (*f)(pc,shift);CHKERRQ(ierr);
24685317021SBarry Smith   }
24785317021SBarry Smith   PetscFunctionReturn(0);
24885317021SBarry Smith }
24985317021SBarry Smith 
25085317021SBarry Smith #undef __FUNCT__
25185317021SBarry Smith #define __FUNCT__ "PCFactorReorderForNonzeroDiagonal"
25285317021SBarry Smith /*@
25385317021SBarry Smith    PCFactorReorderForNonzeroDiagonal - reorders rows/columns of matrix to remove zeros from diagonal
25485317021SBarry Smith 
25585317021SBarry Smith    Collective on PC
25685317021SBarry Smith 
25785317021SBarry Smith    Input Parameters:
25885317021SBarry Smith +  pc - the preconditioner context
25985317021SBarry Smith -  tol - diagonal entries smaller than this in absolute value are considered zero
26085317021SBarry Smith 
26185317021SBarry Smith    Options Database Key:
26285317021SBarry Smith .  -pc_factor_nonzeros_along_diagonal
26385317021SBarry Smith 
26485317021SBarry Smith    Level: intermediate
26585317021SBarry Smith 
26685317021SBarry Smith .keywords: PC, set, factorization, direct, fill
26785317021SBarry Smith 
26885317021SBarry Smith .seealso: PCFactorSetFill(), PCFactorSetShiftNonzero(), PCFactorSetZeroPivot(), MatReorderForNonzeroDiagonal()
26985317021SBarry Smith @*/
27085317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorReorderForNonzeroDiagonal(PC pc,PetscReal rtol)
27185317021SBarry Smith {
27285317021SBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscReal);
27385317021SBarry Smith 
27485317021SBarry Smith   PetscFunctionBegin;
27585317021SBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
27685317021SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C",(void (**)(void))&f);CHKERRQ(ierr);
27785317021SBarry Smith   if (f) {
27885317021SBarry Smith     ierr = (*f)(pc,rtol);CHKERRQ(ierr);
27985317021SBarry Smith   }
28085317021SBarry Smith   PetscFunctionReturn(0);
28185317021SBarry Smith }
28285317021SBarry Smith 
28385317021SBarry Smith #undef __FUNCT__
28485317021SBarry Smith #define __FUNCT__ "PCFactorSetMatSolverPackage"
285bf6011e8SBarry Smith /*@C
28685317021SBarry Smith    PCFactorSetMatSolverPackage - sets the software that is used to perform the factorization
28785317021SBarry Smith 
28885317021SBarry Smith    Collective on PC
28985317021SBarry Smith 
29085317021SBarry Smith    Input Parameters:
29185317021SBarry Smith +  pc - the preconditioner context
29285317021SBarry Smith -  stype - for example, spooles, superlu, superlu_dist
29385317021SBarry Smith 
29485317021SBarry Smith    Options Database Key:
29585317021SBarry Smith .  -pc_factor_mat_solver_package <stype> - spooles, petsc, superlu, superlu_dist, mumps
29685317021SBarry Smith 
29785317021SBarry Smith    Level: intermediate
29885317021SBarry Smith 
29985317021SBarry Smith    Note:
30085317021SBarry Smith      By default this will use the PETSc factorization if it exists
30185317021SBarry Smith 
30285317021SBarry Smith 
30385317021SBarry Smith .keywords: PC, set, factorization, direct, fill
30485317021SBarry Smith 
3057112b564SBarry Smith .seealso: MatGetFactor(), MatSolverPackage, PCFactorGetMatSolverPackage()
30685317021SBarry Smith 
30785317021SBarry Smith @*/
30885317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatSolverPackage(PC pc,const MatSolverPackage stype)
30985317021SBarry Smith {
31085317021SBarry Smith   PetscErrorCode ierr,(*f)(PC,const MatSolverPackage);
31185317021SBarry Smith 
31285317021SBarry Smith   PetscFunctionBegin;
31385317021SBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
31485317021SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetMatSolverPackage_C",(void (**)(void))&f);CHKERRQ(ierr);
31585317021SBarry Smith   if (f) {
31685317021SBarry Smith     ierr = (*f)(pc,stype);CHKERRQ(ierr);
31785317021SBarry Smith   }
31885317021SBarry Smith   PetscFunctionReturn(0);
31985317021SBarry Smith }
32085317021SBarry Smith 
32185317021SBarry Smith #undef __FUNCT__
3227112b564SBarry Smith #define __FUNCT__ "PCFactorGetMatSolverPackage"
323bf6011e8SBarry Smith /*@C
3247112b564SBarry Smith    PCFactorGetMatSolverPackage - gets the software that is used to perform the factorization
3257112b564SBarry Smith 
3267112b564SBarry Smith    Collective on PC
3277112b564SBarry Smith 
3287112b564SBarry Smith    Input Parameter:
3297112b564SBarry Smith .  pc - the preconditioner context
3307112b564SBarry Smith 
3317112b564SBarry Smith    Output Parameter:
3327112b564SBarry Smith .   stype - for example, spooles, superlu, superlu_dist
3337112b564SBarry Smith 
3347112b564SBarry Smith    Level: intermediate
3357112b564SBarry Smith 
3367112b564SBarry Smith 
3377112b564SBarry Smith .keywords: PC, set, factorization, direct, fill
3387112b564SBarry Smith 
3397112b564SBarry Smith .seealso: MatGetFactor(), MatSolverPackage, PCFactorGetMatSolverPackage()
3407112b564SBarry Smith 
3417112b564SBarry Smith @*/
3427112b564SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorGetMatSolverPackage(PC pc,const MatSolverPackage *stype)
3437112b564SBarry Smith {
3447112b564SBarry Smith   PetscErrorCode ierr,(*f)(PC,const MatSolverPackage*);
3457112b564SBarry Smith 
3467112b564SBarry Smith   PetscFunctionBegin;
3477112b564SBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
3487112b564SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorGetMatSolverPackage_C",(void (**)(void))&f);CHKERRQ(ierr);
3497112b564SBarry Smith   if (f) {
3507112b564SBarry Smith     ierr = (*f)(pc,stype);CHKERRQ(ierr);
3517112b564SBarry Smith   }
3527112b564SBarry Smith   PetscFunctionReturn(0);
3537112b564SBarry Smith }
3547112b564SBarry Smith 
3557112b564SBarry Smith #undef __FUNCT__
35685317021SBarry Smith #define __FUNCT__ "PCFactorSetFill"
35785317021SBarry Smith /*@
35885317021SBarry Smith    PCFactorSetFill - Indicate the amount of fill you expect in the factored matrix,
35985317021SBarry Smith    fill = number nonzeros in factor/number nonzeros in original matrix.
36085317021SBarry Smith 
36185317021SBarry Smith    Collective on PC
36285317021SBarry Smith 
36385317021SBarry Smith    Input Parameters:
36485317021SBarry Smith +  pc - the preconditioner context
36585317021SBarry Smith -  fill - amount of expected fill
36685317021SBarry Smith 
36785317021SBarry Smith    Options Database Key:
36885317021SBarry Smith .  -pc_factor_fill <fill> - Sets fill amount
36985317021SBarry Smith 
37085317021SBarry Smith    Level: intermediate
37185317021SBarry Smith 
37285317021SBarry Smith    Note:
37385317021SBarry Smith    For sparse matrix factorizations it is difficult to predict how much
37485317021SBarry Smith    fill to expect. By running with the option -info PETSc will print the
37585317021SBarry Smith    actual amount of fill used; allowing you to set the value accurately for
37685317021SBarry Smith    future runs. Default PETSc uses a value of 5.0
37785317021SBarry Smith 
37885317021SBarry Smith .keywords: PC, set, factorization, direct, fill
37985317021SBarry Smith 
38085317021SBarry Smith @*/
38185317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetFill(PC pc,PetscReal fill)
38285317021SBarry Smith {
38385317021SBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscReal);
38485317021SBarry Smith 
38585317021SBarry Smith   PetscFunctionBegin;
38685317021SBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
38785317021SBarry Smith   if (fill < 1.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Fill factor cannot be less then 1.0");
38885317021SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetFill_C",(void (**)(void))&f);CHKERRQ(ierr);
38985317021SBarry Smith   if (f) {
39085317021SBarry Smith     ierr = (*f)(pc,fill);CHKERRQ(ierr);
39185317021SBarry Smith   }
39285317021SBarry Smith   PetscFunctionReturn(0);
39385317021SBarry Smith }
39485317021SBarry Smith 
39585317021SBarry Smith #undef __FUNCT__
39685317021SBarry Smith #define __FUNCT__ "PCFactorSetUseInPlace"
39785317021SBarry Smith /*@
39885317021SBarry Smith    PCFactorSetUseInPlace - Tells the system to do an in-place factorization.
39985317021SBarry Smith    For dense matrices, this enables the solution of much larger problems.
40085317021SBarry Smith    For sparse matrices the factorization cannot be done truly in-place
40185317021SBarry Smith    so this does not save memory during the factorization, but after the matrix
40285317021SBarry Smith    is factored, the original unfactored matrix is freed, thus recovering that
40385317021SBarry Smith    space.
40485317021SBarry Smith 
40585317021SBarry Smith    Collective on PC
40685317021SBarry Smith 
40785317021SBarry Smith    Input Parameters:
40885317021SBarry Smith .  pc - the preconditioner context
40985317021SBarry Smith 
41085317021SBarry Smith    Options Database Key:
41185317021SBarry Smith .  -pc_factor_in_place - Activates in-place factorization
41285317021SBarry Smith 
41385317021SBarry Smith    Notes:
41485317021SBarry Smith    PCFactorSetUseInplace() can only be used with the KSP method KSPPREONLY or when
41585317021SBarry Smith    a different matrix is provided for the multiply and the preconditioner in
41685317021SBarry Smith    a call to KSPSetOperators().
41785317021SBarry Smith    This is because the Krylov space methods require an application of the
41885317021SBarry Smith    matrix multiplication, which is not possible here because the matrix has
41985317021SBarry Smith    been factored in-place, replacing the original matrix.
42085317021SBarry Smith 
42185317021SBarry Smith    Level: intermediate
42285317021SBarry Smith 
42385317021SBarry Smith .keywords: PC, set, factorization, direct, inplace, in-place, LU
42485317021SBarry Smith 
42585317021SBarry Smith .seealso: PCILUSetUseInPlace()
42685317021SBarry Smith @*/
42785317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseInPlace(PC pc)
42885317021SBarry Smith {
42985317021SBarry Smith   PetscErrorCode ierr,(*f)(PC);
43085317021SBarry Smith 
43185317021SBarry Smith   PetscFunctionBegin;
43285317021SBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
43385317021SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetUseInPlace_C",(void (**)(void))&f);CHKERRQ(ierr);
43485317021SBarry Smith   if (f) {
43585317021SBarry Smith     ierr = (*f)(pc);CHKERRQ(ierr);
43685317021SBarry Smith   }
43785317021SBarry Smith   PetscFunctionReturn(0);
43885317021SBarry Smith }
43985317021SBarry Smith 
44085317021SBarry Smith #undef __FUNCT__
44185317021SBarry Smith #define __FUNCT__ "PCFactorSetMatOrderingType"
44285317021SBarry Smith /*@C
44385317021SBarry Smith     PCFactorSetMatOrderingType - Sets the ordering routine (to reduce fill) to
44485317021SBarry Smith     be used in the LU factorization.
44585317021SBarry Smith 
44685317021SBarry Smith     Collective on PC
44785317021SBarry Smith 
44885317021SBarry Smith     Input Parameters:
44985317021SBarry Smith +   pc - the preconditioner context
45085317021SBarry Smith -   ordering - the matrix ordering name, for example, MATORDERING_ND or MATORDERING_RCM
45185317021SBarry Smith 
45285317021SBarry Smith     Options Database Key:
45385317021SBarry Smith .   -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine
45485317021SBarry Smith 
45585317021SBarry Smith     Level: intermediate
45685317021SBarry Smith 
45785317021SBarry Smith     Notes: nested dissection is used by default
45885317021SBarry Smith 
45985317021SBarry Smith     For Cholesky and ICC and the SBAIJ format reorderings are not available,
46085317021SBarry Smith     since only the upper triangular part of the matrix is stored. You can use the
46185317021SBarry Smith     SeqAIJ format in this case to get reorderings.
46285317021SBarry Smith 
46385317021SBarry Smith @*/
46485317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatOrderingType(PC pc,const MatOrderingType ordering)
46585317021SBarry Smith {
46685317021SBarry Smith   PetscErrorCode ierr,(*f)(PC,const MatOrderingType);
46785317021SBarry Smith 
46885317021SBarry Smith   PetscFunctionBegin;
46985317021SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetMatOrderingType_C",(void (**)(void))&f);CHKERRQ(ierr);
47085317021SBarry Smith   if (f) {
47185317021SBarry Smith     ierr = (*f)(pc,ordering);CHKERRQ(ierr);
47285317021SBarry Smith   }
47385317021SBarry Smith   PetscFunctionReturn(0);
47485317021SBarry Smith }
47585317021SBarry Smith 
47685317021SBarry Smith #undef __FUNCT__
47785317021SBarry Smith #define __FUNCT__ "PCFactorSetPivoting"
47885317021SBarry Smith /*@
47985317021SBarry Smith     PCFactorSetPivoting - Determines when pivoting is done during LU.
48085317021SBarry Smith       For PETSc dense matrices column pivoting is always done, for PETSc sparse matrices
48185317021SBarry Smith       it is never done. For the Matlab and SuperLU factorization this is used.
48285317021SBarry Smith 
48385317021SBarry Smith     Collective on PC
48485317021SBarry Smith 
48585317021SBarry Smith     Input Parameters:
48685317021SBarry Smith +   pc - the preconditioner context
48785317021SBarry Smith -   dtcol - 0.0 implies no pivoting, 1.0 complete pivoting (slower, requires more memory but more stable)
48885317021SBarry Smith 
48985317021SBarry Smith     Options Database Key:
49085317021SBarry Smith .   -pc_factor_pivoting <dtcol>
49185317021SBarry Smith 
49285317021SBarry Smith     Level: intermediate
49385317021SBarry Smith 
49485317021SBarry Smith .seealso: PCILUSetMatOrdering(), PCFactorSetPivotInBlocks()
49585317021SBarry Smith @*/
49685317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivoting(PC pc,PetscReal dtcol)
49785317021SBarry Smith {
49885317021SBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscReal);
49985317021SBarry Smith 
50085317021SBarry Smith   PetscFunctionBegin;
50185317021SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetPivoting_C",(void (**)(void))&f);CHKERRQ(ierr);
50285317021SBarry Smith   if (f) {
50385317021SBarry Smith     ierr = (*f)(pc,dtcol);CHKERRQ(ierr);
50485317021SBarry Smith   }
50585317021SBarry Smith   PetscFunctionReturn(0);
50685317021SBarry Smith }
50785317021SBarry Smith 
50885317021SBarry Smith #undef __FUNCT__
50985317021SBarry Smith #define __FUNCT__ "PCFactorSetPivotInBlocks"
51085317021SBarry Smith /*@
51185317021SBarry Smith     PCFactorSetPivotInBlocks - Determines if pivoting is done while factoring each block
51285317021SBarry Smith       with BAIJ or SBAIJ matrices
51385317021SBarry Smith 
51485317021SBarry Smith     Collective on PC
51585317021SBarry Smith 
51685317021SBarry Smith     Input Parameters:
51785317021SBarry Smith +   pc - the preconditioner context
51885317021SBarry Smith -   pivot - PETSC_TRUE or PETSC_FALSE
51985317021SBarry Smith 
52085317021SBarry Smith     Options Database Key:
52185317021SBarry Smith .   -pc_factor_pivot_in_blocks <true,false>
52285317021SBarry Smith 
52385317021SBarry Smith     Level: intermediate
52485317021SBarry Smith 
52585317021SBarry Smith .seealso: PCILUSetMatOrdering(), PCFactorSetPivoting()
52685317021SBarry Smith @*/
52785317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivotInBlocks(PC pc,PetscTruth pivot)
52885317021SBarry Smith {
52985317021SBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscTruth);
53085317021SBarry Smith 
53185317021SBarry Smith   PetscFunctionBegin;
53285317021SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetPivotInBlocks_C",(void (**)(void))&f);CHKERRQ(ierr);
53385317021SBarry Smith   if (f) {
53485317021SBarry Smith     ierr = (*f)(pc,pivot);CHKERRQ(ierr);
53585317021SBarry Smith   }
53685317021SBarry Smith   PetscFunctionReturn(0);
53785317021SBarry Smith }
53885317021SBarry Smith 
53985317021SBarry Smith #undef __FUNCT__
54085317021SBarry Smith #define __FUNCT__ "PCFactorSetReuseFill"
54185317021SBarry Smith /*@
54285317021SBarry Smith    PCFactorSetReuseFill - When matrices with same different nonzero structure are factored,
54385317021SBarry Smith    this causes later ones to use the fill ratio computed in the initial factorization.
54485317021SBarry Smith 
54585317021SBarry Smith    Collective on PC
54685317021SBarry Smith 
54785317021SBarry Smith    Input Parameters:
54885317021SBarry Smith +  pc - the preconditioner context
54985317021SBarry Smith -  flag - PETSC_TRUE to reuse else PETSC_FALSE
55085317021SBarry Smith 
55185317021SBarry Smith    Options Database Key:
55285317021SBarry Smith .  -pc_factor_reuse_fill - Activates PCFactorSetReuseFill()
55385317021SBarry Smith 
55485317021SBarry Smith    Level: intermediate
55585317021SBarry Smith 
55685317021SBarry Smith .keywords: PC, levels, reordering, factorization, incomplete, Cholesky
55785317021SBarry Smith 
55885317021SBarry Smith .seealso: PCFactorSetReuseOrdering()
55985317021SBarry Smith @*/
56085317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseFill(PC pc,PetscTruth flag)
56185317021SBarry Smith {
56285317021SBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscTruth);
56385317021SBarry Smith 
56485317021SBarry Smith   PetscFunctionBegin;
56585317021SBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,2);
56685317021SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetReuseFill_C",(void (**)(void))&f);CHKERRQ(ierr);
56785317021SBarry Smith   if (f) {
56885317021SBarry Smith     ierr = (*f)(pc,flag);CHKERRQ(ierr);
56985317021SBarry Smith   }
57085317021SBarry Smith   PetscFunctionReturn(0);
57185317021SBarry Smith }
572