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