xref: /petsc/src/ksp/pc/impls/factor/factor.c (revision 03e5aca47aa9bd9d5a48c628e4311d43ab2119ff)
1 
2 #include <../src/ksp/pc/impls/factor/factor.h> /*I "petscpc.h" I*/
3 #include <petsc/private/matimpl.h>
4 
5 /*
6     If an ordering is not yet set and the matrix is available determine a default ordering
7 */
8 PetscErrorCode PCFactorSetDefaultOrdering_Factor(PC pc)
9 {
10   PetscBool   foundmtype, flg;
11   const char *prefix;
12 
13   PetscFunctionBegin;
14   if (pc->pmat) {
15     PetscCall(PCGetOptionsPrefix(pc, &prefix));
16     PetscCall(MatSetOptionsPrefixFactor(pc->pmat, prefix));
17     PC_Factor *fact = (PC_Factor *)pc->data;
18     PetscCall(MatSolverTypeGet(fact->solvertype, ((PetscObject)pc->pmat)->type_name, fact->factortype, NULL, &foundmtype, NULL));
19     if (foundmtype) {
20       if (!fact->fact) { PetscCall(MatGetFactor(pc->pmat, fact->solvertype, fact->factortype, &fact->fact)); }
21       if (!fact->fact) PetscFunctionReturn(PETSC_SUCCESS);
22       if (!fact->fact->assembled) {
23         PetscCall(PetscStrcmp(fact->solvertype, fact->fact->solvertype, &flg));
24         if (!flg) {
25           Mat B;
26           PetscCall(MatGetFactor(pc->pmat, fact->solvertype, fact->factortype, &B));
27           PetscCall(MatHeaderReplace(fact->fact, &B));
28         }
29       }
30       if (!fact->ordering) {
31         PetscBool       canuseordering;
32         MatOrderingType otype;
33 
34         PetscCall(MatFactorGetCanUseOrdering(fact->fact, &canuseordering));
35         if (canuseordering) {
36           PetscCall(MatFactorGetPreferredOrdering(fact->fact, fact->factortype, &otype));
37         } else otype = MATORDERINGEXTERNAL;
38         PetscCall(PetscStrallocpy(otype, (char **)&fact->ordering));
39       }
40     }
41   }
42   PetscFunctionReturn(PETSC_SUCCESS);
43 }
44 
45 static PetscErrorCode PCFactorSetReuseOrdering_Factor(PC pc, PetscBool flag)
46 {
47   PC_Factor *lu = (PC_Factor *)pc->data;
48 
49   PetscFunctionBegin;
50   lu->reuseordering = flag;
51   PetscFunctionReturn(PETSC_SUCCESS);
52 }
53 
54 static PetscErrorCode PCFactorSetReuseFill_Factor(PC pc, PetscBool flag)
55 {
56   PC_Factor *lu = (PC_Factor *)pc->data;
57 
58   PetscFunctionBegin;
59   lu->reusefill = flag;
60   PetscFunctionReturn(PETSC_SUCCESS);
61 }
62 
63 static PetscErrorCode PCFactorSetUseInPlace_Factor(PC pc, PetscBool flg)
64 {
65   PC_Factor *dir = (PC_Factor *)pc->data;
66 
67   PetscFunctionBegin;
68   dir->inplace = flg;
69   PetscFunctionReturn(PETSC_SUCCESS);
70 }
71 
72 static PetscErrorCode PCFactorGetUseInPlace_Factor(PC pc, PetscBool *flg)
73 {
74   PC_Factor *dir = (PC_Factor *)pc->data;
75 
76   PetscFunctionBegin;
77   *flg = dir->inplace;
78   PetscFunctionReturn(PETSC_SUCCESS);
79 }
80 
81 /*@
82     PCFactorSetUpMatSolverType - Can be called after `KSPSetOperators()` or `PCSetOperators()`, causes `MatGetFactor()` to be called so then one may
83        set the options for that particular factorization object.
84 
85   Input Parameter:
86 .  pc  - the preconditioner context
87 
88   Note:
89   After you have called this function (which has to be after the `KSPSetOperators()` or `PCSetOperators()`) you can call `PCFactorGetMatrix()` and then set factor options on that matrix.
90   This function raises an error if the requested combination of solver package and matrix type is not supported.
91 
92   Level: intermediate
93 
94 .seealso: `PCCHOLESKY`, `PCLU`, `PCFactorSetMatSolverType()`, `PCFactorGetMatrix()`
95 @*/
96 PetscErrorCode PCFactorSetUpMatSolverType(PC pc)
97 {
98   PetscFunctionBegin;
99   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
100   PetscTryMethod(pc, "PCFactorSetUpMatSolverType_C", (PC), (pc));
101   PetscFunctionReturn(PETSC_SUCCESS);
102 }
103 
104 /*@
105    PCFactorSetZeroPivot - Sets the size at which smaller pivots are declared to be zero
106 
107    Logically Collective
108 
109    Input Parameters:
110 +  pc - the preconditioner context
111 -  zero - all pivots smaller than this will be considered zero
112 
113    Options Database Key:
114 .  -pc_factor_zeropivot <zero> - Sets tolerance for what is considered a zero pivot
115 
116    Level: intermediate
117 
118 .seealso: `PCCHOLESKY`, `PCLU`, `PCFactorSetShiftType()`, `PCFactorSetShiftAmount()`
119 @*/
120 PetscErrorCode PCFactorSetZeroPivot(PC pc, PetscReal zero)
121 {
122   PetscFunctionBegin;
123   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
124   PetscValidLogicalCollectiveReal(pc, zero, 2);
125   PetscTryMethod(pc, "PCFactorSetZeroPivot_C", (PC, PetscReal), (pc, zero));
126   PetscFunctionReturn(PETSC_SUCCESS);
127 }
128 
129 /*@
130    PCFactorSetShiftType - adds a particular type of quantity to the diagonal of the matrix during
131      numerical factorization, thus the matrix has nonzero pivots
132 
133    Logically Collective
134 
135    Input Parameters:
136 +  pc - the preconditioner context
137 -  shifttype - type of shift; one of `MAT_SHIFT_NONE`, `MAT_SHIFT_NONZERO`, `MAT_SHIFT_POSITIVE_DEFINITE`, `MAT_SHIFT_INBLOCKS`
138 
139    Options Database Key:
140 .  -pc_factor_shift_type <shifttype> - Sets shift type; use '-help' for a list of available types
141 
142    Level: intermediate
143 
144 .seealso: `PCCHOLESKY`, `PCLU`, `PCFactorSetZeroPivot()`, `PCFactorSetShiftAmount()`
145 @*/
146 PetscErrorCode PCFactorSetShiftType(PC pc, MatFactorShiftType shifttype)
147 {
148   PetscFunctionBegin;
149   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
150   PetscValidLogicalCollectiveEnum(pc, shifttype, 2);
151   PetscTryMethod(pc, "PCFactorSetShiftType_C", (PC, MatFactorShiftType), (pc, shifttype));
152   PetscFunctionReturn(PETSC_SUCCESS);
153 }
154 
155 /*@
156    PCFactorSetShiftAmount - adds a quantity to the diagonal of the matrix during
157      numerical factorization, thus the matrix has nonzero pivots
158 
159    Logically Collective
160 
161    Input Parameters:
162 +  pc - the preconditioner context
163 -  shiftamount - amount of shift or `PETSC_DECIDE` for the default
164 
165    Options Database Key:
166 .  -pc_factor_shift_amount <shiftamount> - Sets shift amount or -1 for the default
167 
168    Level: intermediate
169 
170 .seealso: `PCCHOLESKY`, `PCLU`, ``PCFactorSetZeroPivot()`, `PCFactorSetShiftType()`
171 @*/
172 PetscErrorCode PCFactorSetShiftAmount(PC pc, PetscReal shiftamount)
173 {
174   PetscFunctionBegin;
175   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
176   PetscValidLogicalCollectiveReal(pc, shiftamount, 2);
177   PetscTryMethod(pc, "PCFactorSetShiftAmount_C", (PC, PetscReal), (pc, shiftamount));
178   PetscFunctionReturn(PETSC_SUCCESS);
179 }
180 
181 /*@
182    PCFactorSetDropTolerance - The preconditioner will use an `PCILU`
183    based on a drop tolerance.
184 
185    Logically Collective
186 
187    Input Parameters:
188 +  pc - the preconditioner context
189 .  dt - the drop tolerance, try from 1.e-10 to .1
190 .  dtcol - tolerance for column pivot, good values [0.1 to 0.01]
191 -  maxrowcount - the max number of nonzeros allowed in a row, best value
192                  depends on the number of nonzeros in row of original matrix
193 
194    Options Database Key:
195 .  -pc_factor_drop_tolerance <dt,dtcol,maxrowcount> - Sets drop tolerance
196 
197    Level: intermediate
198 
199    Note:
200    There are NO default values for the 3 parameters, you must set them with reasonable values for your
201    matrix. We don't know how to compute reasonable values.
202 
203 .seealso: `PCILU`
204 @*/
205 PetscErrorCode PCFactorSetDropTolerance(PC pc, PetscReal dt, PetscReal dtcol, PetscInt maxrowcount)
206 {
207   PetscFunctionBegin;
208   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
209   PetscValidLogicalCollectiveReal(pc, dtcol, 3);
210   PetscValidLogicalCollectiveInt(pc, maxrowcount, 4);
211   PetscTryMethod(pc, "PCFactorSetDropTolerance_C", (PC, PetscReal, PetscReal, PetscInt), (pc, dt, dtcol, maxrowcount));
212   PetscFunctionReturn(PETSC_SUCCESS);
213 }
214 
215 /*@
216    PCFactorGetZeroPivot - Gets the tolerance used to define a zero privot
217 
218    Not Collective
219 
220    Input Parameter:
221 .  pc - the preconditioner context
222 
223    Output Parameter:
224 .  pivot - the tolerance
225 
226    Level: intermediate
227 
228 .seealso: `PCLU`, `PCCHOLESKY`, `PCFactorSetZeroPivot()`
229 @*/
230 PetscErrorCode PCFactorGetZeroPivot(PC pc, PetscReal *pivot)
231 {
232   PetscFunctionBegin;
233   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
234   PetscUseMethod(pc, "PCFactorGetZeroPivot_C", (PC, PetscReal *), (pc, pivot));
235   PetscFunctionReturn(PETSC_SUCCESS);
236 }
237 
238 /*@
239    PCFactorGetShiftAmount - Gets the tolerance used to define a zero privot
240 
241    Not Collective
242 
243    Input Parameter:
244 .  pc - the preconditioner context
245 
246    Output Parameter:
247 .  shift - how much to shift the diagonal entry
248 
249    Level: intermediate
250 
251 .seealso: `PCLU`, `PCCHOLESKY`, `PCFactorSetShiftAmount()`, `PCFactorSetShiftType()`, `PCFactorGetShiftType()`
252 @*/
253 PetscErrorCode PCFactorGetShiftAmount(PC pc, PetscReal *shift)
254 {
255   PetscFunctionBegin;
256   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
257   PetscUseMethod(pc, "PCFactorGetShiftAmount_C", (PC, PetscReal *), (pc, shift));
258   PetscFunctionReturn(PETSC_SUCCESS);
259 }
260 
261 /*@
262    PCFactorGetShiftType - Gets the type of shift, if any, done when a zero pivot is detected
263 
264    Not Collective
265 
266    Input Parameter:
267 .  pc - the preconditioner context
268 
269    Output Parameter:
270 .  type - one of `MAT_SHIFT_NONE`, `MAT_SHIFT_NONZERO`, `MAT_SHIFT_POSITIVE_DEFINITE`, or `MAT_SHIFT_INBLOCKS`
271 
272    Level: intermediate
273 
274 .seealso: `PCLU`, `PCCHOLESKY`, `PCFactorSetShiftType()`, `MatFactorShiftType`, `PCFactorSetShiftAmount()`, `PCFactorGetShiftAmount()`
275 @*/
276 PetscErrorCode PCFactorGetShiftType(PC pc, MatFactorShiftType *type)
277 {
278   PetscFunctionBegin;
279   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
280   PetscUseMethod(pc, "PCFactorGetShiftType_C", (PC, MatFactorShiftType *), (pc, type));
281   PetscFunctionReturn(PETSC_SUCCESS);
282 }
283 
284 /*@
285    PCFactorGetLevels - Gets the number of levels of fill to use.
286 
287    Logically Collective
288 
289    Input Parameter:
290 .  pc - the preconditioner context
291 
292    Output Parameter:
293 .  levels - number of levels of fill
294 
295    Level: intermediate
296 
297 .seealso: `PCILU`, `PCICC`, `PCFactorSetLevels()`
298 @*/
299 PetscErrorCode PCFactorGetLevels(PC pc, PetscInt *levels)
300 {
301   PetscFunctionBegin;
302   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
303   PetscUseMethod(pc, "PCFactorGetLevels_C", (PC, PetscInt *), (pc, levels));
304   PetscFunctionReturn(PETSC_SUCCESS);
305 }
306 
307 /*@
308    PCFactorSetLevels - Sets the number of levels of fill to use.
309 
310    Logically Collective
311 
312    Input Parameters:
313 +  pc - the preconditioner context
314 -  levels - number of levels of fill
315 
316    Options Database Key:
317 .  -pc_factor_levels <levels> - Sets fill level
318 
319    Level: intermediate
320 
321 .seealso: `PCILU`, `PCICC`, `PCFactorGetLevels()`
322 @*/
323 PetscErrorCode PCFactorSetLevels(PC pc, PetscInt levels)
324 {
325   PetscFunctionBegin;
326   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
327   PetscCheck(levels >= 0, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "negative levels");
328   PetscValidLogicalCollectiveInt(pc, levels, 2);
329   PetscTryMethod(pc, "PCFactorSetLevels_C", (PC, PetscInt), (pc, levels));
330   PetscFunctionReturn(PETSC_SUCCESS);
331 }
332 
333 /*@
334    PCFactorSetAllowDiagonalFill - Causes all diagonal matrix entries to be
335    treated as level 0 fill even if there is no non-zero location.
336 
337    Logically Collective
338 
339    Input Parameters:
340 +  pc - the preconditioner context
341 -  flg - `PETSC_TRUE` to turn on, `PETSC_FALSE` to turn off
342 
343    Options Database Key:
344 .  -pc_factor_diagonal_fill <bool> - allow the diagonal fill
345 
346    Note:
347    Does not apply with 0 fill.
348 
349    Level: intermediate
350 
351 .seealso: `PCILU`, `PCICC`, `PCFactorGetAllowDiagonalFill()`
352 @*/
353 PetscErrorCode PCFactorSetAllowDiagonalFill(PC pc, PetscBool flg)
354 {
355   PetscFunctionBegin;
356   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
357   PetscTryMethod(pc, "PCFactorSetAllowDiagonalFill_C", (PC, PetscBool), (pc, flg));
358   PetscFunctionReturn(PETSC_SUCCESS);
359 }
360 
361 /*@
362    PCFactorGetAllowDiagonalFill - Determines if all diagonal matrix entries are
363        treated as level 0 fill even if there is no non-zero location.
364 
365    Logically Collective
366 
367    Input Parameter:
368 .  pc - the preconditioner context
369 
370    Output Parameter:
371 .   flg - `PETSC_TRUE` to turn on, `PETSC_FALSE` to turn off
372 
373    Note:
374    Does not apply with 0 fill.
375 
376    Level: intermediate
377 
378 .seealso:  `PCILU`, `PCICC`, `PCFactorSetAllowDiagonalFill()`
379 @*/
380 PetscErrorCode PCFactorGetAllowDiagonalFill(PC pc, PetscBool *flg)
381 {
382   PetscFunctionBegin;
383   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
384   PetscUseMethod(pc, "PCFactorGetAllowDiagonalFill_C", (PC, PetscBool *), (pc, flg));
385   PetscFunctionReturn(PETSC_SUCCESS);
386 }
387 
388 /*@
389    PCFactorReorderForNonzeroDiagonal - reorders rows/columns of matrix to remove zeros from diagonal
390 
391    Logically Collective
392 
393    Input Parameters:
394 +  pc - the preconditioner context
395 -  tol - diagonal entries smaller than this in absolute value are considered zero
396 
397    Options Database Key:
398 .  -pc_factor_nonzeros_along_diagonal <tol> - perform the reordering with the given tolerance
399 
400    Level: intermediate
401 
402 .seealso:  `PCILU`, `PCICC`, `PCFactorSetFill()`, `PCFactorSetShiftNonzero()`, `PCFactorSetZeroPivot()`, `MatReorderForNonzeroDiagonal()`
403 @*/
404 PetscErrorCode PCFactorReorderForNonzeroDiagonal(PC pc, PetscReal rtol)
405 {
406   PetscFunctionBegin;
407   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
408   PetscValidLogicalCollectiveReal(pc, rtol, 2);
409   PetscTryMethod(pc, "PCFactorReorderForNonzeroDiagonal_C", (PC, PetscReal), (pc, rtol));
410   PetscFunctionReturn(PETSC_SUCCESS);
411 }
412 
413 /*@C
414    PCFactorSetMatSolverType - sets the solver package that is used to perform the factorization
415 
416    Logically Collective
417 
418    Input Parameters:
419 +  pc - the preconditioner context
420 -  stype - for example, `MATSOLVERSUPERLU`, `MATSOLVERSUPERLU_DIST`, `MATSOLVERMUMPS`
421 
422    Options Database Key:
423 .  -pc_factor_mat_solver_type <stype> - petsc, superlu, superlu_dist, mumps, cusparse
424 
425    Level: intermediate
426 
427    Note:
428    By default this will use the PETSc factorization if it exists
429 
430 .seealso: `PCLU`, `PCCHOLESKY`, `MatGetFactor()`, `MatSolverType`, `PCFactorGetMatSolverType()`, `MatSolverType`,
431           `MATSOLVERSUPERLU`, `MATSOLVERSUPERLU_DIST`, `MATSOLVERMUMPS`
432 @*/
433 PetscErrorCode PCFactorSetMatSolverType(PC pc, MatSolverType stype)
434 {
435   PetscFunctionBegin;
436   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
437   PetscTryMethod(pc, "PCFactorSetMatSolverType_C", (PC, MatSolverType), (pc, stype));
438   PetscFunctionReturn(PETSC_SUCCESS);
439 }
440 
441 /*@C
442    PCFactorGetMatSolverType - gets the solver package that is used to perform the factorization
443 
444    Not Collective
445 
446    Input Parameter:
447 .  pc - the preconditioner context
448 
449    Output Parameter:
450 .   stype - for example, `MATSOLVERSUPERLU`, `MATSOLVERSUPERLU_DIST`, `MATSOLVERMUMPS`
451 
452    Level: intermediate
453 
454 .seealso: `PCLU`, `PCCHOLESKY`, `MatGetFactor()`, `MatSolverType`, `PCFactorGetMatSolverType()`, `MatSolverType`,
455           `MATSOLVERSUPERLU`, `MATSOLVERSUPERLU_DIST`, `MATSOLVERMUMPS`
456 @*/
457 PetscErrorCode PCFactorGetMatSolverType(PC pc, MatSolverType *stype)
458 {
459   PetscErrorCode (*f)(PC, MatSolverType *);
460 
461   PetscFunctionBegin;
462   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
463   PetscValidPointer(stype, 2);
464   PetscCall(PetscObjectQueryFunction((PetscObject)pc, "PCFactorGetMatSolverType_C", &f));
465   if (f) PetscCall((*f)(pc, stype));
466   else *stype = NULL;
467   PetscFunctionReturn(PETSC_SUCCESS);
468 }
469 
470 /*@
471    PCFactorSetFill - Indicate the amount of fill you expect in the factored matrix,
472    fill = number nonzeros in factor/number nonzeros in original matrix.
473 
474    Not Collective, each process can expect a different amount of fill
475 
476    Input Parameters:
477 +  pc - the preconditioner context
478 -  fill - amount of expected fill
479 
480    Options Database Key:
481 .  -pc_factor_fill <fill> - Sets fill amount
482 
483    Level: intermediate
484 
485    Notes:
486    For sparse matrix factorizations it is difficult to predict how much
487    fill to expect. By running with the option -info PETSc will print the
488    actual amount of fill used; allowing you to set the value accurately for
489    future runs. Default PETSc uses a value of 5.0
490 
491    This is ignored for most solver packages
492 
493    This parameter has NOTHING to do with the levels-of-fill of ILU(). That is set with `PCFactorSetLevels()` or -pc_factor_levels.
494 
495 .seealso: `PCLU`, `PCCHOLESKY`, `PCILU`, `PCICC`, `PCFactorSetReuseFill()`
496 @*/
497 PetscErrorCode PCFactorSetFill(PC pc, PetscReal fill)
498 {
499   PetscFunctionBegin;
500   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
501   PetscCheck(fill >= 1.0, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Fill factor cannot be less then 1.0");
502   PetscTryMethod(pc, "PCFactorSetFill_C", (PC, PetscReal), (pc, fill));
503   PetscFunctionReturn(PETSC_SUCCESS);
504 }
505 
506 /*@
507    PCFactorSetUseInPlace - Tells the system to do an in-place factorization.
508    For dense matrices, this enables the solution of much larger problems.
509    For sparse matrices the factorization cannot be done truly in-place
510    so this does not save memory during the factorization, but after the matrix
511    is factored, the original unfactored matrix is freed, thus recovering that
512    space. For ICC(0) and ILU(0) with the default natural ordering the factorization is done efficiently in-place.
513 
514    Logically Collective
515 
516    Input Parameters:
517 +  pc - the preconditioner context
518 -  flg - `PETSC_TRUE` to enable, `PETSC_FALSE` to disable
519 
520    Options Database Key:
521 .  -pc_factor_in_place <true,false>- Activate/deactivate in-place factorization
522 
523    Note:
524    `PCFactorSetUseInplace()` can only be used with the `KSP` method `KSPPREONLY` or when
525    a different matrix is provided for the multiply and the preconditioner in
526    a call to `KSPSetOperators()`.
527    This is because the Krylov space methods require an application of the
528    matrix multiplication, which is not possible here because the matrix has
529    been factored in-place, replacing the original matrix.
530 
531    Level: intermediate
532 
533 .seealso: `PCLU`, `PCCHOLESKY`, `PCILU`, `PCICC`, `PCFactorGetUseInPlace()`
534 @*/
535 PetscErrorCode PCFactorSetUseInPlace(PC pc, PetscBool flg)
536 {
537   PetscFunctionBegin;
538   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
539   PetscTryMethod(pc, "PCFactorSetUseInPlace_C", (PC, PetscBool), (pc, flg));
540   PetscFunctionReturn(PETSC_SUCCESS);
541 }
542 
543 /*@
544    PCFactorGetUseInPlace - Determines if an in-place factorization is being used.
545 
546    Logically Collective
547 
548    Input Parameter:
549 .  pc - the preconditioner context
550 
551    Output Parameter:
552 .  flg - `PETSC_TRUE` to enable, `PETSC_FALSE` to disable
553 
554    Level: intermediate
555 
556 .seealso: `PCLU`, `PCCHOLESKY`, `PCILU`, `PCICC`, `PCFactorSetUseInPlace()`
557 @*/
558 PetscErrorCode PCFactorGetUseInPlace(PC pc, PetscBool *flg)
559 {
560   PetscFunctionBegin;
561   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
562   PetscUseMethod(pc, "PCFactorGetUseInPlace_C", (PC, PetscBool *), (pc, flg));
563   PetscFunctionReturn(PETSC_SUCCESS);
564 }
565 
566 /*@C
567     PCFactorSetMatOrderingType - Sets the ordering routine (to reduce fill) to
568     be used in the `PCLU`, `PCCHOLESKY`, `PCILU`,  or `PCICC` preconditioners
569 
570     Logically Collective
571 
572     Input Parameters:
573 +   pc - the preconditioner context
574 -   ordering - the matrix ordering name, for example, `MATORDERINGND` or `MATORDERINGRCM`
575 
576     Options Database Key:
577 .   -pc_factor_mat_ordering_type <nd,rcm,...,external> - Sets ordering routine
578 
579     Level: intermediate
580 
581     Notes:
582       Nested dissection is used by default for some of PETSc's sparse matrix formats
583 
584      For `PCCHOLESKY` and `PCICC` and the `MATSBAIJ` format the only reordering available is natural since only the upper half of the matrix is stored
585      and reordering this matrix is very expensive.
586 
587     You can use a `MATSEQAIJ` matrix with Cholesky and ICC and use any ordering.
588 
589     `MATORDERINGEXTERNAL` means PETSc will not compute an ordering and the package will use its own ordering, usable with `MATSOLVERCHOLMOD`, `MATSOLVERUMFPACK`, and others.
590 
591 .seealso: `PCLU`, `PCCHOLESKY`, `PCILU`, `PCICC`, `MatOrderingType`, `MATORDERINGEXTERNAL`, `MATORDERINGND`, `MATORDERINGRCM`
592 @*/
593 PetscErrorCode PCFactorSetMatOrderingType(PC pc, MatOrderingType ordering)
594 {
595   PetscFunctionBegin;
596   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
597   PetscTryMethod(pc, "PCFactorSetMatOrderingType_C", (PC, MatOrderingType), (pc, ordering));
598   PetscFunctionReturn(PETSC_SUCCESS);
599 }
600 
601 /*@
602     PCFactorSetColumnPivot - Determines when column pivoting is done during matrix factorization.
603       For PETSc dense matrices column pivoting is always done, for PETSc sparse matrices
604       it is never done. For the MATLAB and `MATSOLVERSUPERLU` factorization this is used.
605 
606     Logically Collective
607 
608     Input Parameters:
609 +   pc - the preconditioner context
610 -   dtcol - 0.0 implies no pivoting, 1.0 complete pivoting (slower, requires more memory but more stable)
611 
612     Options Database Key:
613 .   -pc_factor_pivoting <dtcol> - perform the pivoting with the given tolerance
614 
615     Level: intermediate
616 
617 .seealso: `PCLU`, `PCCHOLESKY`, `PCILU`, `PCICC`, `PCILUSetMatOrdering()`, `PCFactorSetPivotInBlocks()`
618 @*/
619 PetscErrorCode PCFactorSetColumnPivot(PC pc, PetscReal dtcol)
620 {
621   PetscFunctionBegin;
622   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
623   PetscValidLogicalCollectiveReal(pc, dtcol, 2);
624   PetscTryMethod(pc, "PCFactorSetColumnPivot_C", (PC, PetscReal), (pc, dtcol));
625   PetscFunctionReturn(PETSC_SUCCESS);
626 }
627 
628 /*@
629     PCFactorSetPivotInBlocks - Determines if pivoting is done while factoring each block
630       with `MATBAIJ` or `MATSBAIJ` matrices
631 
632     Logically Collective
633 
634     Input Parameters:
635 +   pc - the preconditioner context
636 -   pivot - `PETSC_TRUE` or `PETSC_FALSE`
637 
638     Options Database Key:
639 .   -pc_factor_pivot_in_blocks <true,false> - Pivot inside matrix dense blocks for `MATBAIJ` and `MATSBAIJ`
640 
641     Level: intermediate
642 
643 .seealso: `PCLU`, `PCCHOLESKY`, `PCILU`, `PCICC`, `PCILUSetMatOrdering()`, `PCFactorSetColumnPivot()`
644 @*/
645 PetscErrorCode PCFactorSetPivotInBlocks(PC pc, PetscBool pivot)
646 {
647   PetscFunctionBegin;
648   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
649   PetscValidLogicalCollectiveBool(pc, pivot, 2);
650   PetscTryMethod(pc, "PCFactorSetPivotInBlocks_C", (PC, PetscBool), (pc, pivot));
651   PetscFunctionReturn(PETSC_SUCCESS);
652 }
653 
654 /*@
655    PCFactorSetReuseFill - When matrices with different nonzero structure are factored,
656    this causes later ones to use the fill ratio computed in the initial factorization.
657 
658    Logically Collective
659 
660    Input Parameters:
661 +  pc - the preconditioner context
662 -  flag - `PETSC_TRUE` to reuse else `PETSC_FALSE`
663 
664    Options Database Key:
665 .  -pc_factor_reuse_fill - Activates `PCFactorSetReuseFill()`
666 
667    Level: intermediate
668 
669 .seealso: `PCLU`, `PCCHOLESKY`, `PCILU`, `PCICC`, `PCFactorSetReuseOrdering()`, `PCFactorSetFill()`
670 @*/
671 PetscErrorCode PCFactorSetReuseFill(PC pc, PetscBool flag)
672 {
673   PetscFunctionBegin;
674   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
675   PetscValidLogicalCollectiveBool(pc, flag, 2);
676   PetscTryMethod(pc, "PCFactorSetReuseFill_C", (PC, PetscBool), (pc, flag));
677   PetscFunctionReturn(PETSC_SUCCESS);
678 }
679 
680 PetscErrorCode PCFactorInitialize(PC pc, MatFactorType ftype)
681 {
682   PC_Factor *fact = (PC_Factor *)pc->data;
683 
684   PetscFunctionBegin;
685   PetscCall(MatFactorInfoInitialize(&fact->info));
686   fact->factortype           = ftype;
687   fact->info.shifttype       = (PetscReal)MAT_SHIFT_NONE;
688   fact->info.shiftamount     = 100.0 * PETSC_MACHINE_EPSILON;
689   fact->info.zeropivot       = 100.0 * PETSC_MACHINE_EPSILON;
690   fact->info.pivotinblocks   = 1.0;
691   pc->ops->getfactoredmatrix = PCFactorGetMatrix_Factor;
692 
693   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetZeroPivot_C", PCFactorSetZeroPivot_Factor));
694   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetZeroPivot_C", PCFactorGetZeroPivot_Factor));
695   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetShiftType_C", PCFactorSetShiftType_Factor));
696   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetShiftType_C", PCFactorGetShiftType_Factor));
697   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetShiftAmount_C", PCFactorSetShiftAmount_Factor));
698   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetShiftAmount_C", PCFactorGetShiftAmount_Factor));
699   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetMatSolverType_C", PCFactorGetMatSolverType_Factor));
700   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetMatSolverType_C", PCFactorSetMatSolverType_Factor));
701   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetUpMatSolverType_C", PCFactorSetUpMatSolverType_Factor));
702   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetFill_C", PCFactorSetFill_Factor));
703   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetMatOrderingType_C", PCFactorSetMatOrderingType_Factor));
704   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetLevels_C", PCFactorSetLevels_Factor));
705   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetLevels_C", PCFactorGetLevels_Factor));
706   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetAllowDiagonalFill_C", PCFactorSetAllowDiagonalFill_Factor));
707   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetAllowDiagonalFill_C", PCFactorGetAllowDiagonalFill_Factor));
708   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetPivotInBlocks_C", PCFactorSetPivotInBlocks_Factor));
709   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetUseInPlace_C", PCFactorSetUseInPlace_Factor));
710   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetUseInPlace_C", PCFactorGetUseInPlace_Factor));
711   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetReuseOrdering_C", PCFactorSetReuseOrdering_Factor));
712   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetReuseFill_C", PCFactorSetReuseFill_Factor));
713   PetscFunctionReturn(PETSC_SUCCESS);
714 }
715 
716 PetscErrorCode PCFactorClearComposedFunctions(PC pc)
717 {
718   PetscFunctionBegin;
719   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetZeroPivot_C", NULL));
720   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetZeroPivot_C", NULL));
721   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetShiftType_C", NULL));
722   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetShiftType_C", NULL));
723   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetShiftAmount_C", NULL));
724   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetShiftAmount_C", NULL));
725   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetMatSolverType_C", NULL));
726   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetMatSolverType_C", NULL));
727   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetUpMatSolverType_C", NULL));
728   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetFill_C", NULL));
729   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetMatOrderingType_C", NULL));
730   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetLevels_C", NULL));
731   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetLevels_C", NULL));
732   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetAllowDiagonalFill_C", NULL));
733   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetAllowDiagonalFill_C", NULL));
734   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetPivotInBlocks_C", NULL));
735   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetUseInPlace_C", NULL));
736   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetUseInPlace_C", NULL));
737   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetReuseOrdering_C", NULL));
738   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetReuseFill_C", NULL));
739   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorReorderForNonzeroDiagonal_C", NULL));
740   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetDropTolerance_C", NULL));
741   PetscFunctionReturn(PETSC_SUCCESS);
742 }
743