xref: /petsc/src/ksp/pc/impls/factor/factimpl.c (revision 7c4f633dc6bb6149cca88d301ead35a99e103cbb)
185317021SBarry Smith #define PETSCKSP_DLL
285317021SBarry Smith 
3*7c4f633dSBarry 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 
5585317021SBarry Smith EXTERN_C_BEGIN
5685317021SBarry Smith #undef __FUNCT__
5785317021SBarry Smith #define __FUNCT__ "PCFactorSetUseDropTolerance_Factor"
5885317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseDropTolerance_Factor(PC pc,PetscReal dt,PetscReal dtcol,PetscInt dtcount)
5985317021SBarry Smith {
6085317021SBarry Smith   PC_Factor         *ilu = (PC_Factor*)pc->data;
6185317021SBarry Smith 
6285317021SBarry Smith   PetscFunctionBegin;
6385317021SBarry 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)) {
6485317021SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Cannot change tolerance after use");
6585317021SBarry Smith   }
6685317021SBarry Smith   ilu->info.usedt   = PETSC_TRUE;
6785317021SBarry Smith   ilu->info.dt      = dt;
6885317021SBarry Smith   ilu->info.dtcol   = dtcol;
6985317021SBarry Smith   ilu->info.dtcount = dtcount;
7085317021SBarry Smith   ilu->info.fill    = PETSC_DEFAULT;
7185317021SBarry Smith   PetscFunctionReturn(0);
7285317021SBarry Smith }
7385317021SBarry Smith EXTERN_C_END
7485317021SBarry Smith 
7585317021SBarry Smith EXTERN_C_BEGIN
7685317021SBarry Smith #undef __FUNCT__
7785317021SBarry Smith #define __FUNCT__ "PCFactorSetFill_Factor"
7885317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetFill_Factor(PC pc,PetscReal fill)
7985317021SBarry Smith {
8085317021SBarry Smith   PC_Factor *dir = (PC_Factor*)pc->data;
8185317021SBarry Smith 
8285317021SBarry Smith   PetscFunctionBegin;
8385317021SBarry Smith   dir->info.fill = fill;
8485317021SBarry Smith   PetscFunctionReturn(0);
8585317021SBarry Smith }
8685317021SBarry Smith EXTERN_C_END
8785317021SBarry Smith 
8885317021SBarry Smith EXTERN_C_BEGIN
8985317021SBarry Smith #undef __FUNCT__
9085317021SBarry Smith #define __FUNCT__ "PCFactorSetMatOrderingType_Factor"
9185317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatOrderingType_Factor(PC pc,const MatOrderingType ordering)
9285317021SBarry Smith {
9385317021SBarry Smith   PC_Factor      *dir = (PC_Factor*)pc->data;
9485317021SBarry Smith   PetscErrorCode ierr;
9585317021SBarry Smith   PetscTruth     flg;
9685317021SBarry Smith 
9785317021SBarry Smith   PetscFunctionBegin;
9885317021SBarry Smith   if (!pc->setupcalled) {
9985317021SBarry Smith      ierr = PetscStrfree(dir->ordering);CHKERRQ(ierr);
10085317021SBarry Smith      ierr = PetscStrallocpy(ordering,&dir->ordering);CHKERRQ(ierr);
10185317021SBarry Smith   } else {
10285317021SBarry Smith     ierr = PetscStrcmp(dir->ordering,ordering,&flg);CHKERRQ(ierr);
10385317021SBarry Smith     if (!flg) {
10485317021SBarry Smith       SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Cannot change ordering after use");
10585317021SBarry Smith     }
10685317021SBarry Smith   }
10785317021SBarry Smith   PetscFunctionReturn(0);
10885317021SBarry Smith }
10985317021SBarry Smith EXTERN_C_END
11085317021SBarry Smith 
11185317021SBarry Smith EXTERN_C_BEGIN
11285317021SBarry Smith #undef __FUNCT__
11385317021SBarry Smith #define __FUNCT__ "PCFactorSetLevels_Factor"
11485317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetLevels_Factor(PC pc,PetscInt levels)
11585317021SBarry Smith {
11685317021SBarry Smith   PC_Factor      *ilu = (PC_Factor*)pc->data;
11785317021SBarry Smith 
11885317021SBarry Smith   PetscFunctionBegin;
11985317021SBarry Smith   if (!pc->setupcalled) {
12085317021SBarry Smith     ilu->info.levels = levels;
12185317021SBarry Smith   } else if (ilu->info.usedt || ilu->info.levels != levels) {
12285317021SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Cannot change levels after use");
12385317021SBarry Smith   }
12485317021SBarry Smith   PetscFunctionReturn(0);
12585317021SBarry Smith }
12685317021SBarry Smith EXTERN_C_END
12785317021SBarry Smith 
12885317021SBarry Smith EXTERN_C_BEGIN
12985317021SBarry Smith #undef __FUNCT__
13085317021SBarry Smith #define __FUNCT__ "PCFactorSetAllowDiagonalFill_Factor"
13185317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetAllowDiagonalFill_Factor(PC pc)
13285317021SBarry Smith {
13385317021SBarry Smith   PC_Factor *dir = (PC_Factor*)pc->data;
13485317021SBarry Smith 
13585317021SBarry Smith   PetscFunctionBegin;
13685317021SBarry Smith   dir->info.diagonal_fill = 1;
13785317021SBarry Smith   PetscFunctionReturn(0);
13885317021SBarry Smith }
13985317021SBarry Smith EXTERN_C_END
14085317021SBarry Smith 
14185317021SBarry Smith 
14285317021SBarry Smith /* ------------------------------------------------------------------------------------------*/
14385317021SBarry Smith 
14485317021SBarry Smith EXTERN_C_BEGIN
14585317021SBarry Smith #undef __FUNCT__
14685317021SBarry Smith #define __FUNCT__ "PCFactorSetPivotInBlocks_Factor"
14785317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivotInBlocks_Factor(PC pc,PetscTruth pivot)
14885317021SBarry Smith {
14985317021SBarry Smith   PC_Factor *dir = (PC_Factor*)pc->data;
15085317021SBarry Smith 
15185317021SBarry Smith   PetscFunctionBegin;
15285317021SBarry Smith   dir->info.pivotinblocks = pivot ? 1.0 : 0.0;
15385317021SBarry Smith   PetscFunctionReturn(0);
15485317021SBarry Smith }
15585317021SBarry Smith EXTERN_C_END
15685317021SBarry Smith 
15785317021SBarry Smith EXTERN_C_BEGIN
15885317021SBarry Smith #undef __FUNCT__
15985317021SBarry Smith #define __FUNCT__ "PCFactorSetShiftInBlocks_Factor"
16085317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftInBlocks_Factor(PC pc,PetscReal shift)
16185317021SBarry Smith {
16285317021SBarry Smith   PC_Factor *dir = (PC_Factor*)pc->data;
16385317021SBarry Smith 
16485317021SBarry Smith   PetscFunctionBegin;
16585317021SBarry Smith   if (shift == PETSC_DEFAULT) {
16685317021SBarry Smith     dir->info.shiftinblocks = 1.e-12;
16785317021SBarry Smith   } else {
16885317021SBarry Smith     dir->info.shiftinblocks = shift;
16985317021SBarry Smith   }
17085317021SBarry Smith   PetscFunctionReturn(0);
17185317021SBarry Smith }
17285317021SBarry Smith EXTERN_C_END
17385317021SBarry Smith 
17485317021SBarry Smith 
17585317021SBarry Smith EXTERN_C_BEGIN
17685317021SBarry Smith #undef __FUNCT__
17785317021SBarry Smith #define __FUNCT__ "PCFactorGetMatrix_Factor"
17885317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorGetMatrix_Factor(PC pc,Mat *mat)
17985317021SBarry Smith {
18085317021SBarry Smith   PC_Factor *ilu = (PC_Factor*)pc->data;
18185317021SBarry Smith 
18285317021SBarry Smith   PetscFunctionBegin;
18385317021SBarry Smith   if (!ilu->fact) SETERRQ(PETSC_ERR_ORDER,"Matrix not yet factored; call after KSPSetUp() or PCSetUp()");
18485317021SBarry Smith   *mat = ilu->fact;
18585317021SBarry Smith   PetscFunctionReturn(0);
18685317021SBarry Smith }
18785317021SBarry Smith EXTERN_C_END
18885317021SBarry Smith 
18985317021SBarry Smith EXTERN_C_BEGIN
19085317021SBarry Smith #undef __FUNCT__
19185317021SBarry Smith #define __FUNCT__ "PCFactorSetMatSolverPackage_Factor"
19285317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatSolverPackage_Factor(PC pc,const MatSolverPackage stype)
19385317021SBarry Smith {
19485317021SBarry Smith   PetscErrorCode ierr;
19585317021SBarry Smith   PC_Factor      *lu = (PC_Factor*)pc->data;
19685317021SBarry Smith 
19785317021SBarry Smith   PetscFunctionBegin;
1987112b564SBarry Smith   if (lu->fact) {
19985317021SBarry Smith     const MatSolverPackage ltype;
20085317021SBarry Smith     PetscTruth             flg;
20185317021SBarry Smith     ierr = MatFactorGetSolverPackage(lu->fact,&ltype);CHKERRQ(ierr);
20285317021SBarry Smith     ierr = PetscStrcmp(stype,ltype,&flg);CHKERRQ(ierr);
20385317021SBarry Smith     if (!flg) {
20485317021SBarry Smith       SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Cannot change solver matrix package after PC has been setup or used");
20585317021SBarry Smith     }
20685317021SBarry Smith   }
20785317021SBarry Smith   ierr = PetscStrfree(lu->solvertype);CHKERRQ(ierr);
20885317021SBarry Smith   ierr = PetscStrallocpy(stype,&lu->solvertype);CHKERRQ(ierr);
20985317021SBarry Smith   PetscFunctionReturn(0);
21085317021SBarry Smith }
21185317021SBarry Smith EXTERN_C_END
21285317021SBarry Smith 
21385317021SBarry Smith EXTERN_C_BEGIN
21485317021SBarry Smith #undef __FUNCT__
2157112b564SBarry Smith #define __FUNCT__ "PCFactorGetMatSolverPackage_Factor"
2167112b564SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorGetMatSolverPackage_Factor(PC pc,const MatSolverPackage *stype)
2177112b564SBarry Smith {
2187112b564SBarry Smith   PC_Factor      *lu = (PC_Factor*)pc->data;
2197112b564SBarry Smith 
2207112b564SBarry Smith   PetscFunctionBegin;
2217112b564SBarry Smith   *stype = lu->solvertype;
2227112b564SBarry Smith   PetscFunctionReturn(0);
2237112b564SBarry Smith }
2247112b564SBarry Smith EXTERN_C_END
2257112b564SBarry Smith 
2267112b564SBarry Smith EXTERN_C_BEGIN
2277112b564SBarry Smith #undef __FUNCT__
22885317021SBarry Smith #define __FUNCT__ "PCFactorSetPivoting_Factor"
22985317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivoting_Factor(PC pc,PetscReal dtcol)
23085317021SBarry Smith {
23185317021SBarry Smith   PC_Factor *dir = (PC_Factor*)pc->data;
23285317021SBarry Smith 
23385317021SBarry Smith   PetscFunctionBegin;
23485317021SBarry 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);
23585317021SBarry Smith   dir->info.dtcol = dtcol;
23685317021SBarry Smith   PetscFunctionReturn(0);
23785317021SBarry Smith }
23885317021SBarry Smith EXTERN_C_END
239