xref: /petsc/src/ksp/pc/impls/factor/factimpl.c (revision 853170213946f745f70941b14c22455e5c54846f)
1*85317021SBarry Smith #define PETSCKSP_DLL
2*85317021SBarry Smith 
3*85317021SBarry Smith #include "src/ksp/pc/impls/factor/factor.h"     /*I "petscpc.h"  I*/
4*85317021SBarry Smith 
5*85317021SBarry Smith /* ------------------------------------------------------------------------------------------*/
6*85317021SBarry Smith 
7*85317021SBarry Smith EXTERN_C_BEGIN
8*85317021SBarry Smith #undef __FUNCT__
9*85317021SBarry Smith #define __FUNCT__ "PCFactorSetZeroPivot_Factor"
10*85317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetZeroPivot_Factor(PC pc,PetscReal z)
11*85317021SBarry Smith {
12*85317021SBarry Smith   PC_Factor *ilu = (PC_Factor*)pc->data;
13*85317021SBarry Smith 
14*85317021SBarry Smith   PetscFunctionBegin;
15*85317021SBarry Smith   ilu->info.zeropivot = z;
16*85317021SBarry Smith   PetscFunctionReturn(0);
17*85317021SBarry Smith }
18*85317021SBarry Smith EXTERN_C_END
19*85317021SBarry Smith 
20*85317021SBarry Smith EXTERN_C_BEGIN
21*85317021SBarry Smith #undef __FUNCT__
22*85317021SBarry Smith #define __FUNCT__ "PCFactorSetShiftNonzero_Factor"
23*85317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftNonzero_Factor(PC pc,PetscReal shift)
24*85317021SBarry Smith {
25*85317021SBarry Smith   PC_Factor *dir = (PC_Factor*)pc->data;
26*85317021SBarry Smith 
27*85317021SBarry Smith   PetscFunctionBegin;
28*85317021SBarry Smith   if (shift == (PetscReal) PETSC_DECIDE) {
29*85317021SBarry Smith     dir->info.shiftnz = 1.e-12;
30*85317021SBarry Smith   } else {
31*85317021SBarry Smith     dir->info.shiftnz = shift;
32*85317021SBarry Smith   }
33*85317021SBarry Smith   PetscFunctionReturn(0);
34*85317021SBarry Smith }
35*85317021SBarry Smith EXTERN_C_END
36*85317021SBarry Smith 
37*85317021SBarry Smith EXTERN_C_BEGIN
38*85317021SBarry Smith #undef __FUNCT__
39*85317021SBarry Smith #define __FUNCT__ "PCFactorSetShiftPd_Factor"
40*85317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftPd_Factor(PC pc,PetscTruth shift)
41*85317021SBarry Smith {
42*85317021SBarry Smith   PC_Factor *dir = (PC_Factor*)pc->data;
43*85317021SBarry Smith 
44*85317021SBarry Smith   PetscFunctionBegin;
45*85317021SBarry Smith   if (shift) {
46*85317021SBarry Smith     dir->info.shiftpd = 1.0;
47*85317021SBarry Smith   } else {
48*85317021SBarry Smith     dir->info.shiftpd = 0.0;
49*85317021SBarry Smith   }
50*85317021SBarry Smith   PetscFunctionReturn(0);
51*85317021SBarry Smith }
52*85317021SBarry Smith EXTERN_C_END
53*85317021SBarry Smith 
54*85317021SBarry Smith 
55*85317021SBarry Smith EXTERN_C_BEGIN
56*85317021SBarry Smith #undef __FUNCT__
57*85317021SBarry Smith #define __FUNCT__ "PCFactorSetUseDropTolerance_Factor"
58*85317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseDropTolerance_Factor(PC pc,PetscReal dt,PetscReal dtcol,PetscInt dtcount)
59*85317021SBarry Smith {
60*85317021SBarry Smith   PC_Factor         *ilu = (PC_Factor*)pc->data;
61*85317021SBarry Smith 
62*85317021SBarry Smith   PetscFunctionBegin;
63*85317021SBarry 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)) {
64*85317021SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Cannot change tolerance after use");
65*85317021SBarry Smith   }
66*85317021SBarry Smith   ilu->info.usedt   = PETSC_TRUE;
67*85317021SBarry Smith   ilu->info.dt      = dt;
68*85317021SBarry Smith   ilu->info.dtcol   = dtcol;
69*85317021SBarry Smith   ilu->info.dtcount = dtcount;
70*85317021SBarry Smith   ilu->info.fill    = PETSC_DEFAULT;
71*85317021SBarry Smith   PetscFunctionReturn(0);
72*85317021SBarry Smith }
73*85317021SBarry Smith EXTERN_C_END
74*85317021SBarry Smith 
75*85317021SBarry Smith EXTERN_C_BEGIN
76*85317021SBarry Smith #undef __FUNCT__
77*85317021SBarry Smith #define __FUNCT__ "PCFactorSetFill_Factor"
78*85317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetFill_Factor(PC pc,PetscReal fill)
79*85317021SBarry Smith {
80*85317021SBarry Smith   PC_Factor *dir = (PC_Factor*)pc->data;
81*85317021SBarry Smith 
82*85317021SBarry Smith   PetscFunctionBegin;
83*85317021SBarry Smith   dir->info.fill = fill;
84*85317021SBarry Smith   PetscFunctionReturn(0);
85*85317021SBarry Smith }
86*85317021SBarry Smith EXTERN_C_END
87*85317021SBarry Smith 
88*85317021SBarry Smith EXTERN_C_BEGIN
89*85317021SBarry Smith #undef __FUNCT__
90*85317021SBarry Smith #define __FUNCT__ "PCFactorSetMatOrderingType_Factor"
91*85317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatOrderingType_Factor(PC pc,const MatOrderingType ordering)
92*85317021SBarry Smith {
93*85317021SBarry Smith   PC_Factor      *dir = (PC_Factor*)pc->data;
94*85317021SBarry Smith   PetscErrorCode ierr;
95*85317021SBarry Smith   PetscTruth     flg;
96*85317021SBarry Smith 
97*85317021SBarry Smith   PetscFunctionBegin;
98*85317021SBarry Smith   if (!pc->setupcalled) {
99*85317021SBarry Smith      ierr = PetscStrfree(dir->ordering);CHKERRQ(ierr);
100*85317021SBarry Smith      ierr = PetscStrallocpy(ordering,&dir->ordering);CHKERRQ(ierr);
101*85317021SBarry Smith   } else {
102*85317021SBarry Smith     ierr = PetscStrcmp(dir->ordering,ordering,&flg);CHKERRQ(ierr);
103*85317021SBarry Smith     if (!flg) {
104*85317021SBarry Smith       SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Cannot change ordering after use");
105*85317021SBarry Smith     }
106*85317021SBarry Smith   }
107*85317021SBarry Smith   PetscFunctionReturn(0);
108*85317021SBarry Smith }
109*85317021SBarry Smith EXTERN_C_END
110*85317021SBarry Smith 
111*85317021SBarry Smith EXTERN_C_BEGIN
112*85317021SBarry Smith #undef __FUNCT__
113*85317021SBarry Smith #define __FUNCT__ "PCFactorSetLevels_Factor"
114*85317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetLevels_Factor(PC pc,PetscInt levels)
115*85317021SBarry Smith {
116*85317021SBarry Smith   PC_Factor      *ilu = (PC_Factor*)pc->data;
117*85317021SBarry Smith 
118*85317021SBarry Smith   PetscFunctionBegin;
119*85317021SBarry Smith   if (!pc->setupcalled) {
120*85317021SBarry Smith     ilu->info.levels = levels;
121*85317021SBarry Smith   } else if (ilu->info.usedt || ilu->info.levels != levels) {
122*85317021SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Cannot change levels after use");
123*85317021SBarry Smith   }
124*85317021SBarry Smith   PetscFunctionReturn(0);
125*85317021SBarry Smith }
126*85317021SBarry Smith EXTERN_C_END
127*85317021SBarry Smith 
128*85317021SBarry Smith EXTERN_C_BEGIN
129*85317021SBarry Smith #undef __FUNCT__
130*85317021SBarry Smith #define __FUNCT__ "PCFactorSetAllowDiagonalFill_Factor"
131*85317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetAllowDiagonalFill_Factor(PC pc)
132*85317021SBarry Smith {
133*85317021SBarry Smith   PC_Factor *dir = (PC_Factor*)pc->data;
134*85317021SBarry Smith 
135*85317021SBarry Smith   PetscFunctionBegin;
136*85317021SBarry Smith   dir->info.diagonal_fill = 1;
137*85317021SBarry Smith   PetscFunctionReturn(0);
138*85317021SBarry Smith }
139*85317021SBarry Smith EXTERN_C_END
140*85317021SBarry Smith 
141*85317021SBarry Smith 
142*85317021SBarry Smith /* ------------------------------------------------------------------------------------------*/
143*85317021SBarry Smith 
144*85317021SBarry Smith EXTERN_C_BEGIN
145*85317021SBarry Smith #undef __FUNCT__
146*85317021SBarry Smith #define __FUNCT__ "PCFactorSetPivotInBlocks_Factor"
147*85317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivotInBlocks_Factor(PC pc,PetscTruth pivot)
148*85317021SBarry Smith {
149*85317021SBarry Smith   PC_Factor *dir = (PC_Factor*)pc->data;
150*85317021SBarry Smith 
151*85317021SBarry Smith   PetscFunctionBegin;
152*85317021SBarry Smith   dir->info.pivotinblocks = pivot ? 1.0 : 0.0;
153*85317021SBarry Smith   PetscFunctionReturn(0);
154*85317021SBarry Smith }
155*85317021SBarry Smith EXTERN_C_END
156*85317021SBarry Smith 
157*85317021SBarry Smith EXTERN_C_BEGIN
158*85317021SBarry Smith #undef __FUNCT__
159*85317021SBarry Smith #define __FUNCT__ "PCFactorSetShiftInBlocks_Factor"
160*85317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftInBlocks_Factor(PC pc,PetscReal shift)
161*85317021SBarry Smith {
162*85317021SBarry Smith   PC_Factor *dir = (PC_Factor*)pc->data;
163*85317021SBarry Smith 
164*85317021SBarry Smith   PetscFunctionBegin;
165*85317021SBarry Smith   if (shift == PETSC_DEFAULT) {
166*85317021SBarry Smith     dir->info.shiftinblocks = 1.e-12;
167*85317021SBarry Smith   } else {
168*85317021SBarry Smith     dir->info.shiftinblocks = shift;
169*85317021SBarry Smith   }
170*85317021SBarry Smith   PetscFunctionReturn(0);
171*85317021SBarry Smith }
172*85317021SBarry Smith EXTERN_C_END
173*85317021SBarry Smith 
174*85317021SBarry Smith 
175*85317021SBarry Smith EXTERN_C_BEGIN
176*85317021SBarry Smith #undef __FUNCT__
177*85317021SBarry Smith #define __FUNCT__ "PCFactorGetMatrix_Factor"
178*85317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorGetMatrix_Factor(PC pc,Mat *mat)
179*85317021SBarry Smith {
180*85317021SBarry Smith   PC_Factor *ilu = (PC_Factor*)pc->data;
181*85317021SBarry Smith 
182*85317021SBarry Smith   PetscFunctionBegin;
183*85317021SBarry Smith   if (!ilu->fact) SETERRQ(PETSC_ERR_ORDER,"Matrix not yet factored; call after KSPSetUp() or PCSetUp()");
184*85317021SBarry Smith   *mat = ilu->fact;
185*85317021SBarry Smith   PetscFunctionReturn(0);
186*85317021SBarry Smith }
187*85317021SBarry Smith EXTERN_C_END
188*85317021SBarry Smith 
189*85317021SBarry Smith EXTERN_C_BEGIN
190*85317021SBarry Smith #undef __FUNCT__
191*85317021SBarry Smith #define __FUNCT__ "PCFactorSetMatSolverPackage_Factor"
192*85317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatSolverPackage_Factor(PC pc,const MatSolverPackage stype)
193*85317021SBarry Smith {
194*85317021SBarry Smith   PetscErrorCode ierr;
195*85317021SBarry Smith   PC_Factor      *lu = (PC_Factor*)pc->data;
196*85317021SBarry Smith 
197*85317021SBarry Smith   PetscFunctionBegin;
198*85317021SBarry Smith   if (((PC_Factor*)lu)->fact) {
199*85317021SBarry Smith     const MatSolverPackage ltype;
200*85317021SBarry Smith     PetscTruth             flg;
201*85317021SBarry Smith     ierr = MatFactorGetSolverPackage(lu->fact,&ltype);CHKERRQ(ierr);
202*85317021SBarry Smith     ierr = PetscStrcmp(stype,ltype,&flg);CHKERRQ(ierr);
203*85317021SBarry Smith     if (!flg) {
204*85317021SBarry Smith       SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Cannot change solver matrix package after PC has been setup or used");
205*85317021SBarry Smith     }
206*85317021SBarry Smith   }
207*85317021SBarry Smith   ierr = PetscStrfree(lu->solvertype);CHKERRQ(ierr);
208*85317021SBarry Smith   ierr = PetscStrallocpy(stype,&lu->solvertype);CHKERRQ(ierr);
209*85317021SBarry Smith   PetscFunctionReturn(0);
210*85317021SBarry Smith }
211*85317021SBarry Smith EXTERN_C_END
212*85317021SBarry Smith 
213*85317021SBarry Smith EXTERN_C_BEGIN
214*85317021SBarry Smith #undef __FUNCT__
215*85317021SBarry Smith #define __FUNCT__ "PCFactorSetPivoting_Factor"
216*85317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivoting_Factor(PC pc,PetscReal dtcol)
217*85317021SBarry Smith {
218*85317021SBarry Smith   PC_Factor *dir = (PC_Factor*)pc->data;
219*85317021SBarry Smith 
220*85317021SBarry Smith   PetscFunctionBegin;
221*85317021SBarry 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);
222*85317021SBarry Smith   dir->info.dtcol = dtcol;
223*85317021SBarry Smith   PetscFunctionReturn(0);
224*85317021SBarry Smith }
225*85317021SBarry Smith EXTERN_C_END
226