1 /*
2 The PC (preconditioner) interface routines, callable by users.
3 */
4 #include <petsc/private/pcimpl.h> /*I "petscksp.h" I*/
5 #include <petscdm.h>
6
7 /* Logging support */
8 PetscClassId PC_CLASSID;
9 PetscLogEvent PC_SetUp, PC_SetUpOnBlocks, PC_Apply, PC_MatApply, PC_ApplyCoarse, PC_ApplySymmetricLeft;
10 PetscLogEvent PC_ApplySymmetricRight, PC_ModifySubMatrices, PC_ApplyOnBlocks, PC_ApplyTransposeOnBlocks;
11 PetscInt PetscMGLevelId;
12 PetscLogStage PCMPIStage;
13
PCGetDefaultType_Private(PC pc,const char * type[])14 PETSC_INTERN PetscErrorCode PCGetDefaultType_Private(PC pc, const char *type[])
15 {
16 PetscMPIInt size;
17 PetscBool hasopblock, hasopsolve, flg1, flg2, set, flg3, isnormal;
18
19 PetscFunctionBegin;
20 PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size));
21 if (pc->pmat) {
22 PetscCall(MatHasOperation(pc->pmat, MATOP_GET_DIAGONAL_BLOCK, &hasopblock));
23 PetscCall(MatHasOperation(pc->pmat, MATOP_SOLVE, &hasopsolve));
24 if (size == 1) {
25 PetscCall(MatGetFactorAvailable(pc->pmat, "petsc", MAT_FACTOR_ICC, &flg1));
26 PetscCall(MatGetFactorAvailable(pc->pmat, "petsc", MAT_FACTOR_ILU, &flg2));
27 PetscCall(MatIsSymmetricKnown(pc->pmat, &set, &flg3));
28 PetscCall(PetscObjectTypeCompareAny((PetscObject)pc->pmat, &isnormal, MATNORMAL, MATNORMALHERMITIAN, NULL));
29 if (flg1 && (!flg2 || (set && flg3))) {
30 *type = PCICC;
31 } else if (flg2) {
32 *type = PCILU;
33 } else if (isnormal) {
34 *type = PCNONE;
35 } else if (hasopblock) { /* likely is a parallel matrix run on one processor */
36 *type = PCBJACOBI;
37 } else if (hasopsolve) {
38 *type = PCMAT;
39 } else {
40 *type = PCNONE;
41 }
42 } else {
43 if (hasopblock) {
44 *type = PCBJACOBI;
45 } else if (hasopsolve) {
46 *type = PCMAT;
47 } else {
48 *type = PCNONE;
49 }
50 }
51 } else *type = NULL;
52 PetscFunctionReturn(PETSC_SUCCESS);
53 }
54
55 /* do not log solves, setup and applications of preconditioners while constructing preconditioners; perhaps they should be logged separately from the regular solves */
56 PETSC_EXTERN PetscLogEvent KSP_Solve, KSP_SetUp;
57
PCLogEventsDeactivatePush(void)58 static PetscErrorCode PCLogEventsDeactivatePush(void)
59 {
60 PetscFunctionBegin;
61 PetscCall(KSPInitializePackage());
62 PetscCall(PetscLogEventDeactivatePush(KSP_Solve));
63 PetscCall(PetscLogEventDeactivatePush(KSP_SetUp));
64 PetscCall(PetscLogEventDeactivatePush(PC_Apply));
65 PetscCall(PetscLogEventDeactivatePush(PC_SetUp));
66 PetscCall(PetscLogEventDeactivatePush(PC_SetUpOnBlocks));
67 PetscFunctionReturn(PETSC_SUCCESS);
68 }
69
PCLogEventsDeactivatePop(void)70 static PetscErrorCode PCLogEventsDeactivatePop(void)
71 {
72 PetscFunctionBegin;
73 PetscCall(KSPInitializePackage());
74 PetscCall(PetscLogEventDeactivatePop(KSP_Solve));
75 PetscCall(PetscLogEventDeactivatePop(KSP_SetUp));
76 PetscCall(PetscLogEventDeactivatePop(PC_Apply));
77 PetscCall(PetscLogEventDeactivatePop(PC_SetUp));
78 PetscCall(PetscLogEventDeactivatePop(PC_SetUpOnBlocks));
79 PetscFunctionReturn(PETSC_SUCCESS);
80 }
81
82 /*@
83 PCReset - Resets a `PC` context to the state it was in before `PCSetUp()` was called, and removes any allocated `Vec` and `Mat` from its data structure
84
85 Collective
86
87 Input Parameter:
88 . pc - the `PC` preconditioner context
89
90 Level: developer
91
92 Notes:
93 Any options set, including those set with `KSPSetFromOptions()` remain.
94
95 This allows a `PC` to be reused for a different sized linear system but using the same options that have been previously set in `pc`
96
97 .seealso: [](ch_ksp), `PC`, `PCCreate()`, `PCSetUp()`
98 @*/
PCReset(PC pc)99 PetscErrorCode PCReset(PC pc)
100 {
101 PetscFunctionBegin;
102 PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
103 PetscTryTypeMethod(pc, reset);
104 PetscCall(VecDestroy(&pc->diagonalscaleright));
105 PetscCall(VecDestroy(&pc->diagonalscaleleft));
106 PetscCall(MatDestroy(&pc->pmat));
107 PetscCall(MatDestroy(&pc->mat));
108
109 pc->setupcalled = PETSC_FALSE;
110 PetscFunctionReturn(PETSC_SUCCESS);
111 }
112
113 /*@
114 PCDestroy - Destroys `PC` context that was created with `PCCreate()`.
115
116 Collective
117
118 Input Parameter:
119 . pc - the `PC` preconditioner context
120
121 Level: developer
122
123 .seealso: [](ch_ksp), `PC`, `PCCreate()`, `PCSetUp()`
124 @*/
PCDestroy(PC * pc)125 PetscErrorCode PCDestroy(PC *pc)
126 {
127 PetscFunctionBegin;
128 if (!*pc) PetscFunctionReturn(PETSC_SUCCESS);
129 PetscValidHeaderSpecific(*pc, PC_CLASSID, 1);
130 if (--((PetscObject)*pc)->refct > 0) {
131 *pc = NULL;
132 PetscFunctionReturn(PETSC_SUCCESS);
133 }
134
135 PetscCall(PCReset(*pc));
136
137 /* if memory was published with SAWs then destroy it */
138 PetscCall(PetscObjectSAWsViewOff((PetscObject)*pc));
139 PetscTryTypeMethod(*pc, destroy);
140 PetscCall(DMDestroy(&(*pc)->dm));
141 PetscCall(PetscHeaderDestroy(pc));
142 PetscFunctionReturn(PETSC_SUCCESS);
143 }
144
145 /*@
146 PCGetDiagonalScale - Indicates if the preconditioner applies an additional left and right
147 scaling as needed by certain time-stepping codes.
148
149 Logically Collective
150
151 Input Parameter:
152 . pc - the `PC` preconditioner context
153
154 Output Parameter:
155 . flag - `PETSC_TRUE` if it applies the scaling
156
157 Level: developer
158
159 Note:
160 If this returns `PETSC_TRUE` then the system solved via the Krylov method is, for left and right preconditioning,
161
162 $$
163 \begin{align*}
164 D M A D^{-1} y = D M b \\
165 D A M D^{-1} z = D b.
166 \end{align*}
167 $$
168
169 .seealso: [](ch_ksp), `PC`, `PCCreate()`, `PCSetUp()`, `PCDiagonalScaleLeft()`, `PCDiagonalScaleRight()`, `PCSetDiagonalScale()`
170 @*/
PCGetDiagonalScale(PC pc,PetscBool * flag)171 PetscErrorCode PCGetDiagonalScale(PC pc, PetscBool *flag)
172 {
173 PetscFunctionBegin;
174 PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
175 PetscAssertPointer(flag, 2);
176 *flag = pc->diagonalscale;
177 PetscFunctionReturn(PETSC_SUCCESS);
178 }
179
180 /*@
181 PCSetDiagonalScale - Indicates the left scaling to use to apply an additional left and right
182 scaling as needed by certain time-stepping codes.
183
184 Logically Collective
185
186 Input Parameters:
187 + pc - the `PC` preconditioner context
188 - s - scaling vector
189
190 Level: intermediate
191
192 Notes:
193 The system solved via the Krylov method is, for left and right preconditioning,
194 $$
195 \begin{align*}
196 D M A D^{-1} y = D M b \\
197 D A M D^{-1} z = D b.
198 \end{align*}
199 $$
200
201 `PCDiagonalScaleLeft()` scales a vector by $D$. `PCDiagonalScaleRight()` scales a vector by $D^{-1}$.
202
203 .seealso: [](ch_ksp), `PCCreate()`, `PCSetUp()`, `PCDiagonalScaleLeft()`, `PCDiagonalScaleRight()`, `PCGetDiagonalScale()`
204 @*/
PCSetDiagonalScale(PC pc,Vec s)205 PetscErrorCode PCSetDiagonalScale(PC pc, Vec s)
206 {
207 PetscFunctionBegin;
208 PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
209 PetscValidHeaderSpecific(s, VEC_CLASSID, 2);
210 pc->diagonalscale = PETSC_TRUE;
211
212 PetscCall(PetscObjectReference((PetscObject)s));
213 PetscCall(VecDestroy(&pc->diagonalscaleleft));
214
215 pc->diagonalscaleleft = s;
216
217 PetscCall(VecDuplicate(s, &pc->diagonalscaleright));
218 PetscCall(VecCopy(s, pc->diagonalscaleright));
219 PetscCall(VecReciprocal(pc->diagonalscaleright));
220 PetscFunctionReturn(PETSC_SUCCESS);
221 }
222
223 /*@
224 PCDiagonalScaleLeft - Scales a vector by the left scaling as needed by certain time-stepping codes.
225
226 Logically Collective
227
228 Input Parameters:
229 + pc - the `PC` preconditioner context
230 . in - input vector
231 - out - scaled vector (maybe the same as in)
232
233 Level: intermediate
234
235 Notes:
236 The system solved via the Krylov method is, for left and right preconditioning,
237
238 $$
239 \begin{align*}
240 D M A D^{-1} y = D M b \\
241 D A M D^{-1} z = D b.
242 \end{align*}
243 $$
244
245 `PCDiagonalScaleLeft()` scales a vector by $D$. `PCDiagonalScaleRight()` scales a vector by $D^{-1}$.
246
247 If diagonal scaling is turned off and `in` is not `out` then `in` is copied to `out`
248
249 .seealso: [](ch_ksp), `PCCreate()`, `PCSetUp()`, `PCSetDiagonalScale()`, `PCDiagonalScaleRight()`, `MatDiagonalScale()`
250 @*/
PCDiagonalScaleLeft(PC pc,Vec in,Vec out)251 PetscErrorCode PCDiagonalScaleLeft(PC pc, Vec in, Vec out)
252 {
253 PetscFunctionBegin;
254 PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
255 PetscValidHeaderSpecific(in, VEC_CLASSID, 2);
256 PetscValidHeaderSpecific(out, VEC_CLASSID, 3);
257 if (pc->diagonalscale) {
258 PetscCall(VecPointwiseMult(out, pc->diagonalscaleleft, in));
259 } else if (in != out) {
260 PetscCall(VecCopy(in, out));
261 }
262 PetscFunctionReturn(PETSC_SUCCESS);
263 }
264
265 /*@
266 PCDiagonalScaleRight - Scales a vector by the right scaling as needed by certain time-stepping codes.
267
268 Logically Collective
269
270 Input Parameters:
271 + pc - the `PC` preconditioner context
272 . in - input vector
273 - out - scaled vector (maybe the same as in)
274
275 Level: intermediate
276
277 Notes:
278 The system solved via the Krylov method is, for left and right preconditioning,
279
280 $$
281 \begin{align*}
282 D M A D^{-1} y = D M b \\
283 D A M D^{-1} z = D b.
284 \end{align*}
285 $$
286
287 `PCDiagonalScaleLeft()` scales a vector by $D$. `PCDiagonalScaleRight()` scales a vector by $D^{-1}$.
288
289 If diagonal scaling is turned off and `in` is not `out` then `in` is copied to `out`
290
291 .seealso: [](ch_ksp), `PCCreate()`, `PCSetUp()`, `PCDiagonalScaleLeft()`, `PCSetDiagonalScale()`, `MatDiagonalScale()`
292 @*/
PCDiagonalScaleRight(PC pc,Vec in,Vec out)293 PetscErrorCode PCDiagonalScaleRight(PC pc, Vec in, Vec out)
294 {
295 PetscFunctionBegin;
296 PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
297 PetscValidHeaderSpecific(in, VEC_CLASSID, 2);
298 PetscValidHeaderSpecific(out, VEC_CLASSID, 3);
299 if (pc->diagonalscale) {
300 PetscCall(VecPointwiseMult(out, pc->diagonalscaleright, in));
301 } else if (in != out) {
302 PetscCall(VecCopy(in, out));
303 }
304 PetscFunctionReturn(PETSC_SUCCESS);
305 }
306
307 /*@
308 PCSetUseAmat - Sets a flag to indicate that when the preconditioner needs to apply (part of) the
309 operator during the preconditioning process it applies the Amat provided to `TSSetRHSJacobian()`,
310 `TSSetIJacobian()`, `SNESSetJacobian()`, `KSPSetOperators()` or `PCSetOperators()` not the Pmat.
311
312 Logically Collective
313
314 Input Parameters:
315 + pc - the `PC` preconditioner context
316 - flg - `PETSC_TRUE` to use the Amat, `PETSC_FALSE` to use the Pmat (default is false)
317
318 Options Database Key:
319 . -pc_use_amat <true,false> - use the amat argument to `KSPSetOperators()` or `PCSetOperators()` to apply the operator
320
321 Level: intermediate
322
323 Note:
324 For the common case in which the linear system matrix and the matrix used to construct the
325 preconditioner are identical, this routine has no affect.
326
327 .seealso: [](ch_ksp), `PC`, `PCGetUseAmat()`, `PCBJACOBI`, `PCMG`, `PCFIELDSPLIT`, `PCCOMPOSITE`,
328 `KSPSetOperators()`, `PCSetOperators()`
329 @*/
PCSetUseAmat(PC pc,PetscBool flg)330 PetscErrorCode PCSetUseAmat(PC pc, PetscBool flg)
331 {
332 PetscFunctionBegin;
333 PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
334 pc->useAmat = flg;
335 PetscFunctionReturn(PETSC_SUCCESS);
336 }
337
338 /*@
339 PCSetErrorIfFailure - Causes `PC` to generate an error if a floating point exception, for example a zero pivot, is detected.
340
341 Logically Collective
342
343 Input Parameters:
344 + pc - iterative context obtained from `PCCreate()`
345 - flg - `PETSC_TRUE` indicates you want the error generated
346
347 Level: advanced
348
349 Notes:
350 Normally PETSc continues if a linear solver fails due to a failed setup of a preconditioner, you can call `KSPGetConvergedReason()` after a `KSPSolve()`
351 to determine if it has converged or failed. Or use -ksp_error_if_not_converged to cause the program to terminate as soon as lack of convergence is
352 detected.
353
354 This is propagated into `KSP`s used by this `PC`, which then propagate it into `PC`s used by those `KSP`s
355
356 .seealso: [](ch_ksp), `PC`, `KSPSetErrorIfNotConverged()`, `PCGetInitialGuessNonzero()`, `PCSetInitialGuessKnoll()`, `PCGetInitialGuessKnoll()`
357 @*/
PCSetErrorIfFailure(PC pc,PetscBool flg)358 PetscErrorCode PCSetErrorIfFailure(PC pc, PetscBool flg)
359 {
360 PetscFunctionBegin;
361 PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
362 PetscValidLogicalCollectiveBool(pc, flg, 2);
363 pc->erroriffailure = flg;
364 PetscFunctionReturn(PETSC_SUCCESS);
365 }
366
367 /*@
368 PCGetUseAmat - Gets a flag to indicate that when the preconditioner needs to apply (part of) the
369 operator during the preconditioning process it applies the Amat provided to `TSSetRHSJacobian()`,
370 `TSSetIJacobian()`, `SNESSetJacobian()`, `KSPSetOperators()` or `PCSetOperators()` not the Pmat.
371
372 Logically Collective
373
374 Input Parameter:
375 . pc - the `PC` preconditioner context
376
377 Output Parameter:
378 . flg - `PETSC_TRUE` to use the Amat, `PETSC_FALSE` to use the Pmat (default is false)
379
380 Level: intermediate
381
382 Note:
383 For the common case in which the linear system matrix and the matrix used to construct the
384 preconditioner are identical, this routine is does nothing.
385
386 .seealso: [](ch_ksp), `PC`, `PCSetUseAmat()`, `PCBJACOBI`, `PCMG`, `PCFIELDSPLIT`, `PCCOMPOSITE`
387 @*/
PCGetUseAmat(PC pc,PetscBool * flg)388 PetscErrorCode PCGetUseAmat(PC pc, PetscBool *flg)
389 {
390 PetscFunctionBegin;
391 PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
392 *flg = pc->useAmat;
393 PetscFunctionReturn(PETSC_SUCCESS);
394 }
395
396 /*@
397 PCSetKSPNestLevel - sets the amount of nesting the `KSP` that contains this `PC` has
398
399 Collective
400
401 Input Parameters:
402 + pc - the `PC`
403 - level - the nest level
404
405 Level: developer
406
407 .seealso: [](ch_ksp), `KSPSetUp()`, `KSPSolve()`, `KSPDestroy()`, `KSP`, `KSPGMRES`, `KSPType`, `KSPGetNestLevel()`, `PCGetKSPNestLevel()`, `KSPSetNestLevel()`
408 @*/
PCSetKSPNestLevel(PC pc,PetscInt level)409 PetscErrorCode PCSetKSPNestLevel(PC pc, PetscInt level)
410 {
411 PetscFunctionBegin;
412 PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
413 PetscValidLogicalCollectiveInt(pc, level, 2);
414 pc->kspnestlevel = level;
415 PetscFunctionReturn(PETSC_SUCCESS);
416 }
417
418 /*@
419 PCGetKSPNestLevel - gets the amount of nesting the `KSP` that contains this `PC` has
420
421 Not Collective
422
423 Input Parameter:
424 . pc - the `PC`
425
426 Output Parameter:
427 . level - the nest level
428
429 Level: developer
430
431 .seealso: [](ch_ksp), `KSPSetUp()`, `KSPSolve()`, `KSPDestroy()`, `KSP`, `KSPGMRES`, `KSPType`, `KSPSetNestLevel()`, `PCSetKSPNestLevel()`, `KSPGetNestLevel()`
432 @*/
PCGetKSPNestLevel(PC pc,PetscInt * level)433 PetscErrorCode PCGetKSPNestLevel(PC pc, PetscInt *level)
434 {
435 PetscFunctionBegin;
436 PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
437 PetscAssertPointer(level, 2);
438 *level = pc->kspnestlevel;
439 PetscFunctionReturn(PETSC_SUCCESS);
440 }
441
442 /*@
443 PCCreate - Creates a preconditioner context, `PC`
444
445 Collective
446
447 Input Parameter:
448 . comm - MPI communicator
449
450 Output Parameter:
451 . newpc - location to put the `PC` preconditioner context
452
453 Level: developer
454
455 Notes:
456 This is rarely called directly by users since `KSP` manages the `PC` objects it uses. Use `KSPGetPC()` to access the `PC` used by a `KSP`.
457
458 Use `PCSetType()` or `PCSetFromOptions()` with the option `-pc_type pctype` to set the `PCType` for this `PC`
459
460 The default preconditioner type `PCType` for sparse matrices is `PCILU` or `PCICC` with 0 fill on one process and block Jacobi (`PCBJACOBI`) with `PCILU` or `PCICC`
461 in parallel. For dense matrices it is always `PCNONE`.
462
463 .seealso: [](ch_ksp), `PC`, `PCType`, `PCSetType`, `PCSetUp()`, `PCApply()`, `PCDestroy()`, `KSP`, `KSPGetPC()`
464 @*/
PCCreate(MPI_Comm comm,PC * newpc)465 PetscErrorCode PCCreate(MPI_Comm comm, PC *newpc)
466 {
467 PC pc;
468
469 PetscFunctionBegin;
470 PetscAssertPointer(newpc, 2);
471 PetscCall(PCInitializePackage());
472
473 PetscCall(PetscHeaderCreate(pc, PC_CLASSID, "PC", "Preconditioner", "PC", comm, PCDestroy, PCView));
474 pc->mat = NULL;
475 pc->pmat = NULL;
476 pc->setupcalled = PETSC_FALSE;
477 pc->setfromoptionscalled = 0;
478 pc->data = NULL;
479 pc->diagonalscale = PETSC_FALSE;
480 pc->diagonalscaleleft = NULL;
481 pc->diagonalscaleright = NULL;
482
483 pc->modifysubmatrices = NULL;
484 pc->modifysubmatricesP = NULL;
485
486 *newpc = pc;
487 PetscFunctionReturn(PETSC_SUCCESS);
488 }
489
490 /*@
491 PCApply - Applies the preconditioner to a vector.
492
493 Collective
494
495 Input Parameters:
496 + pc - the `PC` preconditioner context
497 - x - input vector
498
499 Output Parameter:
500 . y - output vector
501
502 Level: developer
503
504 .seealso: [](ch_ksp), `PC`, `PCApplyTranspose()`, `PCApplyBAorAB()`
505 @*/
PCApply(PC pc,Vec x,Vec y)506 PetscErrorCode PCApply(PC pc, Vec x, Vec y)
507 {
508 PetscInt m, n, mv, nv;
509
510 PetscFunctionBegin;
511 PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
512 PetscValidHeaderSpecific(x, VEC_CLASSID, 2);
513 PetscValidHeaderSpecific(y, VEC_CLASSID, 3);
514 PetscCheck(x != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "x and y must be different vectors");
515 if (pc->erroriffailure) PetscCall(VecValidValues_Internal(x, 2, PETSC_TRUE));
516 /* use pmat to check vector sizes since for KSPLSQR the pmat may be of a different size than mat */
517 PetscCall(MatGetLocalSize(pc->pmat, &m, &n));
518 PetscCall(VecGetLocalSize(x, &mv));
519 PetscCall(VecGetLocalSize(y, &nv));
520 /* check pmat * y = x is feasible */
521 PetscCheck(mv == m, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Preconditioner number of local rows %" PetscInt_FMT " does not equal input vector size %" PetscInt_FMT, m, mv);
522 PetscCheck(nv == n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Preconditioner number of local columns %" PetscInt_FMT " does not equal output vector size %" PetscInt_FMT, n, nv);
523 PetscCall(VecSetErrorIfLocked(y, 3));
524
525 PetscCall(PCSetUp(pc));
526 PetscCall(VecLockReadPush(x));
527 PetscCall(PetscLogEventBegin(PC_Apply, pc, x, y, 0));
528 PetscUseTypeMethod(pc, apply, x, y);
529 PetscCall(PetscLogEventEnd(PC_Apply, pc, x, y, 0));
530 if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 3, PETSC_FALSE));
531 PetscCall(VecLockReadPop(x));
532 PetscFunctionReturn(PETSC_SUCCESS);
533 }
534
PCMatApplyTranspose_Private(PC pc,Mat X,Mat Y,PetscBool transpose)535 static PetscErrorCode PCMatApplyTranspose_Private(PC pc, Mat X, Mat Y, PetscBool transpose)
536 {
537 Mat A;
538 Vec cy, cx;
539 PetscInt m1, M1, m2, M2, n1, N1, n2, N2, m3, M3, n3, N3;
540 PetscBool match;
541
542 PetscFunctionBegin;
543 PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
544 PetscValidHeaderSpecific(X, MAT_CLASSID, 2);
545 PetscValidHeaderSpecific(Y, MAT_CLASSID, 3);
546 PetscCheckSameComm(pc, 1, X, 2);
547 PetscCheckSameComm(pc, 1, Y, 3);
548 PetscCheck(Y != X, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "Y and X must be different matrices");
549 PetscCall(PCGetOperators(pc, NULL, &A));
550 PetscCall(MatGetLocalSize(A, &m3, &n3));
551 PetscCall(MatGetLocalSize(X, &m2, &n2));
552 PetscCall(MatGetLocalSize(Y, &m1, &n1));
553 PetscCall(MatGetSize(A, &M3, &N3));
554 PetscCall(MatGetSize(X, &M2, &N2));
555 PetscCall(MatGetSize(Y, &M1, &N1));
556 PetscCheck(n1 == n2 && N1 == N2, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Incompatible number of columns between block of input vectors (n,N) = (%" PetscInt_FMT ",%" PetscInt_FMT ") and block of output vectors (n,N) = (%" PetscInt_FMT ",%" PetscInt_FMT ")", n2, N2, n1, N1);
557 PetscCheck(m2 == m3 && M2 == M3, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Incompatible layout between block of input vectors (m,M) = (%" PetscInt_FMT ",%" PetscInt_FMT ") and Pmat (m,M)x(n,N) = (%" PetscInt_FMT ",%" PetscInt_FMT ")x(%" PetscInt_FMT ",%" PetscInt_FMT ")", m2, M2, m3, M3, n3, N3);
558 PetscCheck(m1 == n3 && M1 == N3, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Incompatible layout between block of output vectors (m,M) = (%" PetscInt_FMT ",%" PetscInt_FMT ") and Pmat (m,M)x(n,N) = (%" PetscInt_FMT ",%" PetscInt_FMT ")x(%" PetscInt_FMT ",%" PetscInt_FMT ")", m1, M1, m3, M3, n3, N3);
559 PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)Y, &match, MATSEQDENSE, MATMPIDENSE, ""));
560 PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Provided block of output vectors not stored in a dense Mat");
561 PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)X, &match, MATSEQDENSE, MATMPIDENSE, ""));
562 PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Provided block of input vectors not stored in a dense Mat");
563 PetscCall(PCSetUp(pc));
564 if (!transpose && pc->ops->matapply) {
565 PetscCall(PetscLogEventBegin(PC_MatApply, pc, X, Y, 0));
566 PetscUseTypeMethod(pc, matapply, X, Y);
567 PetscCall(PetscLogEventEnd(PC_MatApply, pc, X, Y, 0));
568 } else if (transpose && pc->ops->matapplytranspose) {
569 PetscCall(PetscLogEventBegin(PC_MatApply, pc, X, Y, 0));
570 PetscUseTypeMethod(pc, matapplytranspose, X, Y);
571 PetscCall(PetscLogEventEnd(PC_MatApply, pc, X, Y, 0));
572 } else {
573 PetscCall(PetscInfo(pc, "PC type %s applying column by column\n", ((PetscObject)pc)->type_name));
574 for (n1 = 0; n1 < N1; ++n1) {
575 PetscCall(MatDenseGetColumnVecRead(X, n1, &cx));
576 PetscCall(MatDenseGetColumnVecWrite(Y, n1, &cy));
577 if (!transpose) PetscCall(PCApply(pc, cx, cy));
578 else PetscCall(PCApplyTranspose(pc, cx, cy));
579 PetscCall(MatDenseRestoreColumnVecWrite(Y, n1, &cy));
580 PetscCall(MatDenseRestoreColumnVecRead(X, n1, &cx));
581 }
582 }
583 PetscFunctionReturn(PETSC_SUCCESS);
584 }
585
586 /*@
587 PCMatApply - Applies the preconditioner to multiple vectors stored as a `MATDENSE`. Like `PCApply()`, `Y` and `X` must be different matrices.
588
589 Collective
590
591 Input Parameters:
592 + pc - the `PC` preconditioner context
593 - X - block of input vectors
594
595 Output Parameter:
596 . Y - block of output vectors
597
598 Level: developer
599
600 .seealso: [](ch_ksp), `PC`, `PCApply()`, `KSPMatSolve()`
601 @*/
PCMatApply(PC pc,Mat X,Mat Y)602 PetscErrorCode PCMatApply(PC pc, Mat X, Mat Y)
603 {
604 PetscFunctionBegin;
605 PetscCall(PCMatApplyTranspose_Private(pc, X, Y, PETSC_FALSE));
606 PetscFunctionReturn(PETSC_SUCCESS);
607 }
608
609 /*@
610 PCMatApplyTranspose - Applies the transpose of preconditioner to multiple vectors stored as a `MATDENSE`. Like `PCApplyTranspose()`, `Y` and `X` must be different matrices.
611
612 Collective
613
614 Input Parameters:
615 + pc - the `PC` preconditioner context
616 - X - block of input vectors
617
618 Output Parameter:
619 . Y - block of output vectors
620
621 Level: developer
622
623 .seealso: [](ch_ksp), `PC`, `PCApplyTranspose()`, `KSPMatSolveTranspose()`
624 @*/
PCMatApplyTranspose(PC pc,Mat X,Mat Y)625 PetscErrorCode PCMatApplyTranspose(PC pc, Mat X, Mat Y)
626 {
627 PetscFunctionBegin;
628 PetscCall(PCMatApplyTranspose_Private(pc, X, Y, PETSC_TRUE));
629 PetscFunctionReturn(PETSC_SUCCESS);
630 }
631
632 /*@
633 PCApplySymmetricLeft - Applies the left part of a symmetric preconditioner to a vector.
634
635 Collective
636
637 Input Parameters:
638 + pc - the `PC` preconditioner context
639 - x - input vector
640
641 Output Parameter:
642 . y - output vector
643
644 Level: developer
645
646 Note:
647 Currently, this routine is implemented only for `PCICC` and `PCJACOBI` preconditioners.
648
649 .seealso: [](ch_ksp), `PC`, `PCApply()`, `PCApplySymmetricRight()`
650 @*/
PCApplySymmetricLeft(PC pc,Vec x,Vec y)651 PetscErrorCode PCApplySymmetricLeft(PC pc, Vec x, Vec y)
652 {
653 PetscFunctionBegin;
654 PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
655 PetscValidHeaderSpecific(x, VEC_CLASSID, 2);
656 PetscValidHeaderSpecific(y, VEC_CLASSID, 3);
657 PetscCheck(x != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "x and y must be different vectors");
658 if (pc->erroriffailure) PetscCall(VecValidValues_Internal(x, 2, PETSC_TRUE));
659 PetscCall(PCSetUp(pc));
660 PetscCall(VecLockReadPush(x));
661 PetscCall(PetscLogEventBegin(PC_ApplySymmetricLeft, pc, x, y, 0));
662 PetscUseTypeMethod(pc, applysymmetricleft, x, y);
663 PetscCall(PetscLogEventEnd(PC_ApplySymmetricLeft, pc, x, y, 0));
664 PetscCall(VecLockReadPop(x));
665 if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 3, PETSC_FALSE));
666 PetscFunctionReturn(PETSC_SUCCESS);
667 }
668
669 /*@
670 PCApplySymmetricRight - Applies the right part of a symmetric preconditioner to a vector.
671
672 Collective
673
674 Input Parameters:
675 + pc - the `PC` preconditioner context
676 - x - input vector
677
678 Output Parameter:
679 . y - output vector
680
681 Level: developer
682
683 Note:
684 Currently, this routine is implemented only for `PCICC` and `PCJACOBI` preconditioners.
685
686 .seealso: [](ch_ksp), `PC`, `PCApply()`, `PCApplySymmetricLeft()`
687 @*/
PCApplySymmetricRight(PC pc,Vec x,Vec y)688 PetscErrorCode PCApplySymmetricRight(PC pc, Vec x, Vec y)
689 {
690 PetscFunctionBegin;
691 PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
692 PetscValidHeaderSpecific(x, VEC_CLASSID, 2);
693 PetscValidHeaderSpecific(y, VEC_CLASSID, 3);
694 PetscCheck(x != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "x and y must be different vectors");
695 if (pc->erroriffailure) PetscCall(VecValidValues_Internal(x, 2, PETSC_TRUE));
696 PetscCall(PCSetUp(pc));
697 PetscCall(VecLockReadPush(x));
698 PetscCall(PetscLogEventBegin(PC_ApplySymmetricRight, pc, x, y, 0));
699 PetscUseTypeMethod(pc, applysymmetricright, x, y);
700 PetscCall(PetscLogEventEnd(PC_ApplySymmetricRight, pc, x, y, 0));
701 PetscCall(VecLockReadPop(x));
702 if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 3, PETSC_FALSE));
703 PetscFunctionReturn(PETSC_SUCCESS);
704 }
705
706 /*@
707 PCApplyTranspose - Applies the transpose of preconditioner to a vector.
708
709 Collective
710
711 Input Parameters:
712 + pc - the `PC` preconditioner context
713 - x - input vector
714
715 Output Parameter:
716 . y - output vector
717
718 Level: developer
719
720 Note:
721 For complex numbers this applies the non-Hermitian transpose.
722
723 Developer Note:
724 We need to implement a `PCApplyHermitianTranspose()`
725
726 .seealso: [](ch_ksp), `PC`, `PCApply()`, `PCApplyBAorAB()`, `PCApplyBAorABTranspose()`, `PCApplyTransposeExists()`
727 @*/
PCApplyTranspose(PC pc,Vec x,Vec y)728 PetscErrorCode PCApplyTranspose(PC pc, Vec x, Vec y)
729 {
730 PetscFunctionBegin;
731 PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
732 PetscValidHeaderSpecific(x, VEC_CLASSID, 2);
733 PetscValidHeaderSpecific(y, VEC_CLASSID, 3);
734 PetscCheck(x != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "x and y must be different vectors");
735 if (pc->erroriffailure) PetscCall(VecValidValues_Internal(x, 2, PETSC_TRUE));
736 PetscCall(PCSetUp(pc));
737 PetscCall(VecLockReadPush(x));
738 PetscCall(PetscLogEventBegin(PC_Apply, pc, x, y, 0));
739 PetscUseTypeMethod(pc, applytranspose, x, y);
740 PetscCall(PetscLogEventEnd(PC_Apply, pc, x, y, 0));
741 PetscCall(VecLockReadPop(x));
742 if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 3, PETSC_FALSE));
743 PetscFunctionReturn(PETSC_SUCCESS);
744 }
745
746 /*@
747 PCApplyTransposeExists - Test whether the preconditioner has a transpose apply operation
748
749 Collective
750
751 Input Parameter:
752 . pc - the `PC` preconditioner context
753
754 Output Parameter:
755 . flg - `PETSC_TRUE` if a transpose operation is defined
756
757 Level: developer
758
759 .seealso: [](ch_ksp), `PC`, `PCApplyTranspose()`
760 @*/
PCApplyTransposeExists(PC pc,PetscBool * flg)761 PetscErrorCode PCApplyTransposeExists(PC pc, PetscBool *flg)
762 {
763 PetscFunctionBegin;
764 PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
765 PetscAssertPointer(flg, 2);
766 if (pc->ops->applytranspose) *flg = PETSC_TRUE;
767 else *flg = PETSC_FALSE;
768 PetscFunctionReturn(PETSC_SUCCESS);
769 }
770
771 /*@
772 PCApplyBAorAB - Applies the preconditioner and operator to a vector. $y = B*A*x $ or $ y = A*B*x$.
773
774 Collective
775
776 Input Parameters:
777 + pc - the `PC` preconditioner context
778 . side - indicates the preconditioner side, one of `PC_LEFT`, `PC_RIGHT`, or `PC_SYMMETRIC`
779 . x - input vector
780 - work - work vector
781
782 Output Parameter:
783 . y - output vector
784
785 Level: developer
786
787 Note:
788 If the `PC` has had `PCSetDiagonalScale()` set then $ D M A D^{-1} $ for left preconditioning or $ D A M D^{-1} $ is actually applied.
789 The specific `KSPSolve()` method must also be written to handle the post-solve "correction" for the diagonal scaling.
790
791 .seealso: [](ch_ksp), `PC`, `PCApply()`, `PCApplyTranspose()`, `PCApplyBAorABTranspose()`
792 @*/
PCApplyBAorAB(PC pc,PCSide side,Vec x,Vec y,Vec work)793 PetscErrorCode PCApplyBAorAB(PC pc, PCSide side, Vec x, Vec y, Vec work)
794 {
795 PetscFunctionBegin;
796 PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
797 PetscValidLogicalCollectiveEnum(pc, side, 2);
798 PetscValidHeaderSpecific(x, VEC_CLASSID, 3);
799 PetscValidHeaderSpecific(y, VEC_CLASSID, 4);
800 PetscValidHeaderSpecific(work, VEC_CLASSID, 5);
801 PetscCheckSameComm(pc, 1, x, 3);
802 PetscCheckSameComm(pc, 1, y, 4);
803 PetscCheckSameComm(pc, 1, work, 5);
804 PetscCheck(x != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "x and y must be different vectors");
805 PetscCheck(side == PC_LEFT || side == PC_SYMMETRIC || side == PC_RIGHT, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Side must be right, left, or symmetric");
806 PetscCheck(!pc->diagonalscale || side != PC_SYMMETRIC, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Cannot include diagonal scaling with symmetric preconditioner application");
807 if (pc->erroriffailure) PetscCall(VecValidValues_Internal(x, 3, PETSC_TRUE));
808
809 PetscCall(PCSetUp(pc));
810 if (pc->diagonalscale) {
811 if (pc->ops->applyBA) {
812 Vec work2; /* this is expensive, but to fix requires a second work vector argument to PCApplyBAorAB() */
813 PetscCall(VecDuplicate(x, &work2));
814 PetscCall(PCDiagonalScaleRight(pc, x, work2));
815 PetscUseTypeMethod(pc, applyBA, side, work2, y, work);
816 PetscCall(PCDiagonalScaleLeft(pc, y, y));
817 PetscCall(VecDestroy(&work2));
818 } else if (side == PC_RIGHT) {
819 PetscCall(PCDiagonalScaleRight(pc, x, y));
820 PetscCall(PCApply(pc, y, work));
821 PetscCall(MatMult(pc->mat, work, y));
822 PetscCall(PCDiagonalScaleLeft(pc, y, y));
823 } else if (side == PC_LEFT) {
824 PetscCall(PCDiagonalScaleRight(pc, x, y));
825 PetscCall(MatMult(pc->mat, y, work));
826 PetscCall(PCApply(pc, work, y));
827 PetscCall(PCDiagonalScaleLeft(pc, y, y));
828 } else PetscCheck(side != PC_SYMMETRIC, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Cannot provide diagonal scaling with symmetric application of preconditioner");
829 } else {
830 if (pc->ops->applyBA) {
831 PetscUseTypeMethod(pc, applyBA, side, x, y, work);
832 } else if (side == PC_RIGHT) {
833 PetscCall(PCApply(pc, x, work));
834 PetscCall(MatMult(pc->mat, work, y));
835 } else if (side == PC_LEFT) {
836 PetscCall(MatMult(pc->mat, x, work));
837 PetscCall(PCApply(pc, work, y));
838 } else if (side == PC_SYMMETRIC) {
839 /* There's an extra copy here; maybe should provide 2 work vectors instead? */
840 PetscCall(PCApplySymmetricRight(pc, x, work));
841 PetscCall(MatMult(pc->mat, work, y));
842 PetscCall(VecCopy(y, work));
843 PetscCall(PCApplySymmetricLeft(pc, work, y));
844 }
845 }
846 if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 4, PETSC_FALSE));
847 PetscFunctionReturn(PETSC_SUCCESS);
848 }
849
850 /*@
851 PCApplyBAorABTranspose - Applies the transpose of the preconditioner
852 and operator to a vector. That is, applies $B^T * A^T$ with left preconditioning,
853 NOT $(B*A)^T = A^T*B^T$.
854
855 Collective
856
857 Input Parameters:
858 + pc - the `PC` preconditioner context
859 . side - indicates the preconditioner side, one of `PC_LEFT`, `PC_RIGHT`, or `PC_SYMMETRIC`
860 . x - input vector
861 - work - work vector
862
863 Output Parameter:
864 . y - output vector
865
866 Level: developer
867
868 Note:
869 This routine is used internally so that the same Krylov code can be used to solve $A x = b$ and $A^T x = b$, with a preconditioner
870 defined by $B^T$. This is why this has the funny form that it computes $B^T * A^T$
871
872 .seealso: [](ch_ksp), `PC`, `PCApply()`, `PCApplyTranspose()`, `PCApplyBAorAB()`
873 @*/
PCApplyBAorABTranspose(PC pc,PCSide side,Vec x,Vec y,Vec work)874 PetscErrorCode PCApplyBAorABTranspose(PC pc, PCSide side, Vec x, Vec y, Vec work)
875 {
876 PetscFunctionBegin;
877 PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
878 PetscValidHeaderSpecific(x, VEC_CLASSID, 3);
879 PetscValidHeaderSpecific(y, VEC_CLASSID, 4);
880 PetscValidHeaderSpecific(work, VEC_CLASSID, 5);
881 PetscCheck(x != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "x and y must be different vectors");
882 if (pc->erroriffailure) PetscCall(VecValidValues_Internal(x, 3, PETSC_TRUE));
883 if (pc->ops->applyBAtranspose) {
884 PetscUseTypeMethod(pc, applyBAtranspose, side, x, y, work);
885 if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 4, PETSC_FALSE));
886 PetscFunctionReturn(PETSC_SUCCESS);
887 }
888 PetscCheck(side == PC_LEFT || side == PC_RIGHT, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Side must be right or left");
889
890 PetscCall(PCSetUp(pc));
891 if (side == PC_RIGHT) {
892 PetscCall(PCApplyTranspose(pc, x, work));
893 PetscCall(MatMultTranspose(pc->mat, work, y));
894 } else if (side == PC_LEFT) {
895 PetscCall(MatMultTranspose(pc->mat, x, work));
896 PetscCall(PCApplyTranspose(pc, work, y));
897 }
898 /* add support for PC_SYMMETRIC */
899 if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 4, PETSC_FALSE));
900 PetscFunctionReturn(PETSC_SUCCESS);
901 }
902
903 /*@
904 PCApplyRichardsonExists - Determines whether a particular preconditioner has a
905 built-in fast application of Richardson's method.
906
907 Not Collective
908
909 Input Parameter:
910 . pc - the preconditioner
911
912 Output Parameter:
913 . exists - `PETSC_TRUE` or `PETSC_FALSE`
914
915 Level: developer
916
917 .seealso: [](ch_ksp), `PC`, `KSPRICHARDSON`, `PCApplyRichardson()`
918 @*/
PCApplyRichardsonExists(PC pc,PetscBool * exists)919 PetscErrorCode PCApplyRichardsonExists(PC pc, PetscBool *exists)
920 {
921 PetscFunctionBegin;
922 PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
923 PetscAssertPointer(exists, 2);
924 if (pc->ops->applyrichardson) *exists = PETSC_TRUE;
925 else *exists = PETSC_FALSE;
926 PetscFunctionReturn(PETSC_SUCCESS);
927 }
928
929 /*@
930 PCApplyRichardson - Applies several steps of Richardson iteration with
931 the particular preconditioner. This routine is usually used by the
932 Krylov solvers and not the application code directly.
933
934 Collective
935
936 Input Parameters:
937 + pc - the `PC` preconditioner context
938 . b - the right-hand side
939 . w - one work vector
940 . rtol - relative decrease in residual norm convergence criteria
941 . abstol - absolute residual norm convergence criteria
942 . dtol - divergence residual norm increase criteria
943 . its - the number of iterations to apply.
944 - guesszero - if the input x contains nonzero initial guess
945
946 Output Parameters:
947 + outits - number of iterations actually used (for SOR this always equals its)
948 . reason - the reason the apply terminated
949 - y - the solution (also contains initial guess if guesszero is `PETSC_FALSE`
950
951 Level: developer
952
953 Notes:
954 Most preconditioners do not support this function. Use the command
955 `PCApplyRichardsonExists()` to determine if one does.
956
957 Except for the `PCMG` this routine ignores the convergence tolerances
958 and always runs for the number of iterations
959
960 .seealso: [](ch_ksp), `PC`, `PCApplyRichardsonExists()`
961 @*/
PCApplyRichardson(PC pc,Vec b,Vec y,Vec w,PetscReal rtol,PetscReal abstol,PetscReal dtol,PetscInt its,PetscBool guesszero,PetscInt * outits,PCRichardsonConvergedReason * reason)962 PetscErrorCode PCApplyRichardson(PC pc, Vec b, Vec y, Vec w, PetscReal rtol, PetscReal abstol, PetscReal dtol, PetscInt its, PetscBool guesszero, PetscInt *outits, PCRichardsonConvergedReason *reason)
963 {
964 PetscFunctionBegin;
965 PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
966 PetscValidHeaderSpecific(b, VEC_CLASSID, 2);
967 PetscValidHeaderSpecific(y, VEC_CLASSID, 3);
968 PetscValidHeaderSpecific(w, VEC_CLASSID, 4);
969 PetscCheck(b != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "b and y must be different vectors");
970 PetscCall(PCSetUp(pc));
971 PetscUseTypeMethod(pc, applyrichardson, b, y, w, rtol, abstol, dtol, its, guesszero, outits, reason);
972 PetscFunctionReturn(PETSC_SUCCESS);
973 }
974
975 /*@
976 PCSetFailedReason - Sets the reason a `PCSetUp()` failed or `PC_NOERROR` if it did not fail
977
978 Logically Collective
979
980 Input Parameters:
981 + pc - the `PC` preconditioner context
982 - reason - the reason it failed
983
984 Level: advanced
985
986 .seealso: [](ch_ksp), `PC`, `PCCreate()`, `PCApply()`, `PCDestroy()`, `PCFailedReason`
987 @*/
PCSetFailedReason(PC pc,PCFailedReason reason)988 PetscErrorCode PCSetFailedReason(PC pc, PCFailedReason reason)
989 {
990 PetscFunctionBegin;
991 PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
992 pc->failedreason = reason;
993 PetscFunctionReturn(PETSC_SUCCESS);
994 }
995
996 /*@
997 PCGetFailedReason - Gets the reason a `PCSetUp()` failed or `PC_NOERROR` if it did not fail
998
999 Not Collective
1000
1001 Input Parameter:
1002 . pc - the `PC` preconditioner context
1003
1004 Output Parameter:
1005 . reason - the reason it failed
1006
1007 Level: advanced
1008
1009 Note:
1010 After a call to `KSPCheckDot()` or `KSPCheckNorm()` inside a `KSPSolve()` or a call to `PCReduceFailedReason()`
1011 this is the maximum reason over all MPI processes in the `PC` communicator and hence logically collective.
1012 Otherwise it returns the local value.
1013
1014 .seealso: [](ch_ksp), `PC`, `PCCreate()`, `PCApply()`, `PCDestroy()`, `PCSetFailedReason()`, `PCFailedReason`
1015 @*/
PCGetFailedReason(PC pc,PCFailedReason * reason)1016 PetscErrorCode PCGetFailedReason(PC pc, PCFailedReason *reason)
1017 {
1018 PetscFunctionBegin;
1019 PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1020 *reason = pc->failedreason;
1021 PetscFunctionReturn(PETSC_SUCCESS);
1022 }
1023
1024 /*@
1025 PCReduceFailedReason - Reduce the failed reason among the MPI processes that share the `PC`
1026
1027 Collective
1028
1029 Input Parameter:
1030 . pc - the `PC` preconditioner context
1031
1032 Level: advanced
1033
1034 Note:
1035 Different MPI processes may have different reasons or no reason, see `PCGetFailedReason()`. This routine
1036 makes them have a common value (failure if any MPI process had a failure).
1037
1038 .seealso: [](ch_ksp), `PC`, `PCCreate()`, `PCApply()`, `PCDestroy()`, `PCGetFailedReason()`, `PCSetFailedReason()`, `PCFailedReason`
1039 @*/
PCReduceFailedReason(PC pc)1040 PetscErrorCode PCReduceFailedReason(PC pc)
1041 {
1042 PetscInt buf;
1043
1044 PetscFunctionBegin;
1045 PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1046 buf = (PetscInt)pc->failedreason;
1047 PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &buf, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)pc)));
1048 pc->failedreason = (PCFailedReason)buf;
1049 PetscFunctionReturn(PETSC_SUCCESS);
1050 }
1051
1052 /*
1053 a setupcall of 0 indicates never setup,
1054 1 indicates has been previously setup
1055 -1 indicates a PCSetUp() was attempted and failed
1056 */
1057 /*@
1058 PCSetUp - Prepares for the use of a preconditioner. Performs all the one-time operations needed before the preconditioner
1059 can be used with `PCApply()`
1060
1061 Collective
1062
1063 Input Parameter:
1064 . pc - the `PC` preconditioner context
1065
1066 Level: developer
1067
1068 Notes:
1069 For example, for `PCLU` this will compute the factorization.
1070
1071 This is called automatically by `KSPSetUp()` or `PCApply()` so rarely needs to be called directly.
1072
1073 For nested preconditioners, such as `PCFIELDSPLIT` or `PCBJACOBI` this may not finish the construction of the preconditioner
1074 on the inner levels, the routine `PCSetUpOnBlocks()` may compute more of the preconditioner in those situations.
1075
1076 .seealso: [](ch_ksp), `PC`, `PCCreate()`, `PCApply()`, `PCDestroy()`, `KSPSetUp()`, `PCSetUpOnBlocks()`
1077 @*/
PCSetUp(PC pc)1078 PetscErrorCode PCSetUp(PC pc)
1079 {
1080 const char *def;
1081 PetscObjectState matstate, matnonzerostate;
1082
1083 PetscFunctionBegin;
1084 PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1085 PetscCheck(pc->mat, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Matrix must be set first");
1086
1087 if (pc->setupcalled && pc->reusepreconditioner) {
1088 PetscCall(PetscInfo(pc, "Leaving PC with identical preconditioner since reuse preconditioner is set\n"));
1089 PetscFunctionReturn(PETSC_SUCCESS);
1090 }
1091
1092 PetscCall(PetscObjectStateGet((PetscObject)pc->pmat, &matstate));
1093 PetscCall(MatGetNonzeroState(pc->pmat, &matnonzerostate));
1094 if (!pc->setupcalled) {
1095 //PetscCall(PetscInfo(pc, "Setting up PC for first time\n"));
1096 pc->flag = DIFFERENT_NONZERO_PATTERN;
1097 } else if (matstate == pc->matstate) PetscFunctionReturn(PETSC_SUCCESS);
1098 else {
1099 if (matnonzerostate != pc->matnonzerostate) {
1100 PetscCall(PetscInfo(pc, "Setting up PC with different nonzero pattern\n"));
1101 pc->flag = DIFFERENT_NONZERO_PATTERN;
1102 } else {
1103 //PetscCall(PetscInfo(pc, "Setting up PC with same nonzero pattern\n"));
1104 pc->flag = SAME_NONZERO_PATTERN;
1105 }
1106 }
1107 pc->matstate = matstate;
1108 pc->matnonzerostate = matnonzerostate;
1109
1110 if (!((PetscObject)pc)->type_name) {
1111 PetscCall(PCGetDefaultType_Private(pc, &def));
1112 PetscCall(PCSetType(pc, def));
1113 }
1114
1115 PetscCall(MatSetErrorIfFailure(pc->pmat, pc->erroriffailure));
1116 PetscCall(MatSetErrorIfFailure(pc->mat, pc->erroriffailure));
1117 PetscCall(PetscLogEventBegin(PC_SetUp, pc, 0, 0, 0));
1118 if (pc->ops->setup) {
1119 PetscCall(PCLogEventsDeactivatePush());
1120 PetscUseTypeMethod(pc, setup);
1121 PetscCall(PCLogEventsDeactivatePop());
1122 }
1123 PetscCall(PetscLogEventEnd(PC_SetUp, pc, 0, 0, 0));
1124 if (pc->postsetup) PetscCall((*pc->postsetup)(pc));
1125 if (!pc->setupcalled) pc->setupcalled = PETSC_TRUE;
1126 PetscFunctionReturn(PETSC_SUCCESS);
1127 }
1128
1129 /*@
1130 PCSetUpOnBlocks - Sets up the preconditioner for each block in
1131 the block Jacobi, overlapping Schwarz, and fieldsplit methods.
1132
1133 Collective
1134
1135 Input Parameter:
1136 . pc - the `PC` preconditioner context
1137
1138 Level: developer
1139
1140 Notes:
1141 For nested preconditioners such as `PCBJACOBI`, `PCSetUp()` is not called on each sub-`KSP` when `PCSetUp()` is
1142 called on the outer `PC`, this routine ensures it is called.
1143
1144 It calls `PCSetUp()` if not yet called.
1145
1146 .seealso: [](ch_ksp), `PC`, `PCSetUp()`, `PCCreate()`, `PCApply()`, `PCDestroy()`
1147 @*/
PCSetUpOnBlocks(PC pc)1148 PetscErrorCode PCSetUpOnBlocks(PC pc)
1149 {
1150 PetscFunctionBegin;
1151 PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1152 if (!pc->setupcalled) PetscCall(PCSetUp(pc)); /* "if" to prevent -info extra prints */
1153 if (!pc->ops->setuponblocks) PetscFunctionReturn(PETSC_SUCCESS);
1154 PetscCall(MatSetErrorIfFailure(pc->pmat, pc->erroriffailure));
1155 PetscCall(PetscLogEventBegin(PC_SetUpOnBlocks, pc, 0, 0, 0));
1156 PetscCall(PCLogEventsDeactivatePush());
1157 PetscUseTypeMethod(pc, setuponblocks);
1158 PetscCall(PCLogEventsDeactivatePop());
1159 PetscCall(PetscLogEventEnd(PC_SetUpOnBlocks, pc, 0, 0, 0));
1160 PetscFunctionReturn(PETSC_SUCCESS);
1161 }
1162
1163 /*@C
1164 PCSetModifySubMatrices - Sets a user-defined routine for modifying the
1165 submatrices that arise within certain subdomain-based preconditioners such as `PCASM`
1166
1167 Logically Collective
1168
1169 Input Parameters:
1170 + pc - the `PC` preconditioner context
1171 . func - routine for modifying the submatrices, see `PCModifySubMatricesFn`
1172 - ctx - optional user-defined context (may be `NULL`)
1173
1174 Level: advanced
1175
1176 Notes:
1177 The basic submatrices are extracted from the matrix used to construct the preconditioner as
1178 usual; the user can then alter these (for example, to set different boundary
1179 conditions for each submatrix) before they are used for the local solves.
1180
1181 `PCSetModifySubMatrices()` MUST be called before `KSPSetUp()` and
1182 `KSPSolve()`.
1183
1184 A routine set by `PCSetModifySubMatrices()` is currently called within
1185 `PCBJACOBI`, `PCASM`, `PCGASM`, and `PCHPDDM`.
1186 All other preconditioners ignore this routine.
1187
1188 .seealso: [](ch_ksp), `PC`, `PCModifySubMatricesFn`, `PCBJACOBI`, `PCASM`, `PCModifySubMatrices()`
1189 @*/
PCSetModifySubMatrices(PC pc,PCModifySubMatricesFn * func,PetscCtx ctx)1190 PetscErrorCode PCSetModifySubMatrices(PC pc, PCModifySubMatricesFn *func, PetscCtx ctx)
1191 {
1192 PetscFunctionBegin;
1193 PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1194 pc->modifysubmatrices = func;
1195 pc->modifysubmatricesP = ctx;
1196 PetscFunctionReturn(PETSC_SUCCESS);
1197 }
1198
1199 /*@C
1200 PCModifySubMatrices - Calls an optional user-defined routine within
1201 certain preconditioners if one has been set with `PCSetModifySubMatrices()`.
1202
1203 Collective
1204
1205 Input Parameters:
1206 + pc - the `PC` preconditioner context
1207 . nsub - the number of local submatrices
1208 . row - an array of index sets that contain the global row numbers
1209 that comprise each local submatrix
1210 . col - an array of index sets that contain the global column numbers
1211 that comprise each local submatrix
1212 . submat - array of local submatrices
1213 - ctx - optional user-defined context for private data for the
1214 user-defined routine (may be `NULL`)
1215
1216 Output Parameter:
1217 . submat - array of local submatrices (the entries of which may
1218 have been modified)
1219
1220 Level: developer
1221
1222 Note:
1223 The user should NOT generally call this routine, as it will
1224 automatically be called within certain preconditioners.
1225
1226 .seealso: [](ch_ksp), `PC`, `PCModifySubMatricesFn`, `PCSetModifySubMatrices()`
1227 @*/
PCModifySubMatrices(PC pc,PetscInt nsub,const IS row[],const IS col[],Mat submat[],PetscCtx ctx)1228 PetscErrorCode PCModifySubMatrices(PC pc, PetscInt nsub, const IS row[], const IS col[], Mat submat[], PetscCtx ctx)
1229 {
1230 PetscFunctionBegin;
1231 PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1232 if (!pc->modifysubmatrices) PetscFunctionReturn(PETSC_SUCCESS);
1233 PetscCall(PetscLogEventBegin(PC_ModifySubMatrices, pc, 0, 0, 0));
1234 PetscCall((*pc->modifysubmatrices)(pc, nsub, row, col, submat, ctx));
1235 PetscCall(PetscLogEventEnd(PC_ModifySubMatrices, pc, 0, 0, 0));
1236 PetscFunctionReturn(PETSC_SUCCESS);
1237 }
1238
1239 /*@
1240 PCSetOperators - Sets the matrix associated with the linear system and
1241 a (possibly) different one from which the preconditioner will be constructed.
1242
1243 Logically Collective
1244
1245 Input Parameters:
1246 + pc - the `PC` preconditioner context
1247 . Amat - the matrix that defines the linear system
1248 - Pmat - the matrix to be used in constructing the preconditioner, usually the same as Amat.
1249
1250 Level: advanced
1251
1252 Notes:
1253 Using this routine directly is rarely needed, the preferred, and equivalent, usage is `KSPSetOperators()`.
1254
1255 Passing a `NULL` for `Amat` or `Pmat` removes the matrix that is currently used.
1256
1257 If you wish to replace either `Amat` or `Pmat` but leave the other one untouched then
1258 first call `KSPGetOperators()` to get the one you wish to keep, call `PetscObjectReference()`
1259 on it and then pass it back in in your call to `KSPSetOperators()`.
1260
1261 More Notes about Repeated Solution of Linear Systems:
1262 PETSc does NOT reset the matrix entries of either `Amat` or `Pmat`
1263 to zero after a linear solve; the user is completely responsible for
1264 matrix assembly. See the routine `MatZeroEntries()` if desiring to
1265 zero all elements of a matrix.
1266
1267 .seealso: [](ch_ksp), `PC`, `PCGetOperators()`, `MatZeroEntries()`
1268 @*/
PCSetOperators(PC pc,Mat Amat,Mat Pmat)1269 PetscErrorCode PCSetOperators(PC pc, Mat Amat, Mat Pmat)
1270 {
1271 PetscInt m1, n1, m2, n2;
1272
1273 PetscFunctionBegin;
1274 PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1275 if (Amat) PetscValidHeaderSpecific(Amat, MAT_CLASSID, 2);
1276 if (Pmat) PetscValidHeaderSpecific(Pmat, MAT_CLASSID, 3);
1277 if (Amat) PetscCheckSameComm(pc, 1, Amat, 2);
1278 if (Pmat) PetscCheckSameComm(pc, 1, Pmat, 3);
1279 if (pc->setupcalled && pc->mat && pc->pmat && Amat && Pmat) {
1280 PetscCall(MatGetLocalSize(Amat, &m1, &n1));
1281 PetscCall(MatGetLocalSize(pc->mat, &m2, &n2));
1282 PetscCheck(m1 == m2 && n1 == n2, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Cannot change local size of Amat after use old sizes %" PetscInt_FMT " %" PetscInt_FMT " new sizes %" PetscInt_FMT " %" PetscInt_FMT, m2, n2, m1, n1);
1283 PetscCall(MatGetLocalSize(Pmat, &m1, &n1));
1284 PetscCall(MatGetLocalSize(pc->pmat, &m2, &n2));
1285 PetscCheck(m1 == m2 && n1 == n2, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Cannot change local size of Pmat after use old sizes %" PetscInt_FMT " %" PetscInt_FMT " new sizes %" PetscInt_FMT " %" PetscInt_FMT, m2, n2, m1, n1);
1286 }
1287
1288 if (Pmat != pc->pmat) {
1289 /* changing the operator that defines the preconditioner thus reneed to clear current states so new preconditioner is built */
1290 pc->matnonzerostate = -1;
1291 pc->matstate = -1;
1292 }
1293
1294 /* reference first in case the matrices are the same */
1295 if (Amat) PetscCall(PetscObjectReference((PetscObject)Amat));
1296 PetscCall(MatDestroy(&pc->mat));
1297 if (Pmat) PetscCall(PetscObjectReference((PetscObject)Pmat));
1298 PetscCall(MatDestroy(&pc->pmat));
1299 pc->mat = Amat;
1300 pc->pmat = Pmat;
1301 PetscFunctionReturn(PETSC_SUCCESS);
1302 }
1303
1304 /*@
1305 PCSetReusePreconditioner - reuse the current preconditioner even if the operator in the preconditioner `PC` has changed.
1306
1307 Logically Collective
1308
1309 Input Parameters:
1310 + pc - the `PC` preconditioner context
1311 - flag - `PETSC_TRUE` do not compute a new preconditioner, `PETSC_FALSE` do compute a new preconditioner
1312
1313 Level: intermediate
1314
1315 Note:
1316 Normally if a matrix inside a `PC` changes the `PC` automatically updates itself using information from the changed matrix. This option
1317 prevents this.
1318
1319 .seealso: [](ch_ksp), `PC`, `PCGetOperators()`, `MatZeroEntries()`, `PCGetReusePreconditioner()`, `KSPSetReusePreconditioner()`
1320 @*/
PCSetReusePreconditioner(PC pc,PetscBool flag)1321 PetscErrorCode PCSetReusePreconditioner(PC pc, PetscBool flag)
1322 {
1323 PetscFunctionBegin;
1324 PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1325 PetscValidLogicalCollectiveBool(pc, flag, 2);
1326 pc->reusepreconditioner = flag;
1327 PetscTryMethod(pc, "PCSetReusePreconditioner_C", (PC, PetscBool), (pc, flag));
1328 PetscFunctionReturn(PETSC_SUCCESS);
1329 }
1330
1331 /*@
1332 PCGetReusePreconditioner - Determines if the `PC` reuses the current preconditioner even if the operator in the preconditioner has changed.
1333
1334 Not Collective
1335
1336 Input Parameter:
1337 . pc - the `PC` preconditioner context
1338
1339 Output Parameter:
1340 . flag - `PETSC_TRUE` do not compute a new preconditioner, `PETSC_FALSE` do compute a new preconditioner
1341
1342 Level: intermediate
1343
1344 .seealso: [](ch_ksp), `PC`, `PCGetOperators()`, `MatZeroEntries()`, `PCSetReusePreconditioner()`
1345 @*/
PCGetReusePreconditioner(PC pc,PetscBool * flag)1346 PetscErrorCode PCGetReusePreconditioner(PC pc, PetscBool *flag)
1347 {
1348 PetscFunctionBegin;
1349 PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1350 PetscAssertPointer(flag, 2);
1351 *flag = pc->reusepreconditioner;
1352 PetscFunctionReturn(PETSC_SUCCESS);
1353 }
1354
1355 /*@
1356 PCGetOperators - Gets the matrix associated with the linear system and
1357 possibly a different one which is used to construct the preconditioner.
1358
1359 Not Collective, though parallel `Mat`s are returned if `pc` is parallel
1360
1361 Input Parameter:
1362 . pc - the `PC` preconditioner context
1363
1364 Output Parameters:
1365 + Amat - the matrix defining the linear system
1366 - Pmat - the matrix from which the preconditioner is constructed, usually the same as Amat.
1367
1368 Level: intermediate
1369
1370 Note:
1371 Does not increase the reference count of the matrices, so you should not destroy them
1372
1373 Alternative usage: If the operators have NOT been set with `KSPSetOperators()`/`PCSetOperators()` then the operators
1374 are created in `PC` and returned to the user. In this case, if both operators
1375 mat and pmat are requested, two DIFFERENT operators will be returned. If
1376 only one is requested both operators in the PC will be the same (i.e. as
1377 if one had called `KSPSetOperators()`/`PCSetOperators()` with the same argument for both Mats).
1378 The user must set the sizes of the returned matrices and their type etc just
1379 as if the user created them with `MatCreate()`. For example,
1380
1381 .vb
1382 KSP/PCGetOperators(ksp/pc,&Amat,NULL); is equivalent to
1383 set size, type, etc of Amat
1384
1385 MatCreate(comm,&mat);
1386 KSP/PCSetOperators(ksp/pc,Amat,Amat);
1387 PetscObjectDereference((PetscObject)mat);
1388 set size, type, etc of Amat
1389 .ve
1390
1391 and
1392
1393 .vb
1394 KSP/PCGetOperators(ksp/pc,&Amat,&Pmat); is equivalent to
1395 set size, type, etc of Amat and Pmat
1396
1397 MatCreate(comm,&Amat);
1398 MatCreate(comm,&Pmat);
1399 KSP/PCSetOperators(ksp/pc,Amat,Pmat);
1400 PetscObjectDereference((PetscObject)Amat);
1401 PetscObjectDereference((PetscObject)Pmat);
1402 set size, type, etc of Amat and Pmat
1403 .ve
1404
1405 The rationale for this support is so that when creating a `TS`, `SNES`, or `KSP` the hierarchy
1406 of underlying objects (i.e. `SNES`, `KSP`, `PC`, `Mat`) and their lifespans can be completely
1407 managed by the top most level object (i.e. the `TS`, `SNES`, or `KSP`). Another way to look
1408 at this is when you create a `SNES` you do not NEED to create a `KSP` and attach it to
1409 the `SNES` object (the `SNES` object manages it for you). Similarly when you create a KSP
1410 you do not need to attach a `PC` to it (the `KSP` object manages the `PC` object for you).
1411 Thus, why should YOU have to create the `Mat` and attach it to the `SNES`/`KSP`/`PC`, when
1412 it can be created for you?
1413
1414 .seealso: [](ch_ksp), `PC`, `PCSetOperators()`, `KSPGetOperators()`, `KSPSetOperators()`, `PCGetOperatorsSet()`
1415 @*/
PCGetOperators(PC pc,Mat * Amat,Mat * Pmat)1416 PetscErrorCode PCGetOperators(PC pc, Mat *Amat, Mat *Pmat)
1417 {
1418 PetscFunctionBegin;
1419 PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1420 if (Amat) {
1421 if (!pc->mat) {
1422 if (pc->pmat && !Pmat) { /* Pmat has been set, but user did not request it, so use for Amat */
1423 pc->mat = pc->pmat;
1424 PetscCall(PetscObjectReference((PetscObject)pc->mat));
1425 } else { /* both Amat and Pmat are empty */
1426 PetscCall(MatCreate(PetscObjectComm((PetscObject)pc), &pc->mat));
1427 if (!Pmat) { /* user did NOT request Pmat, so make same as Amat */
1428 pc->pmat = pc->mat;
1429 PetscCall(PetscObjectReference((PetscObject)pc->pmat));
1430 }
1431 }
1432 }
1433 *Amat = pc->mat;
1434 }
1435 if (Pmat) {
1436 if (!pc->pmat) {
1437 if (pc->mat && !Amat) { /* Amat has been set but was not requested, so use for pmat */
1438 pc->pmat = pc->mat;
1439 PetscCall(PetscObjectReference((PetscObject)pc->pmat));
1440 } else {
1441 PetscCall(MatCreate(PetscObjectComm((PetscObject)pc), &pc->pmat));
1442 if (!Amat) { /* user did NOT request Amat, so make same as Pmat */
1443 pc->mat = pc->pmat;
1444 PetscCall(PetscObjectReference((PetscObject)pc->mat));
1445 }
1446 }
1447 }
1448 *Pmat = pc->pmat;
1449 }
1450 PetscFunctionReturn(PETSC_SUCCESS);
1451 }
1452
1453 /*@
1454 PCGetOperatorsSet - Determines if the matrix associated with the linear system and
1455 possibly a different one associated with the preconditioner have been set in the `PC`.
1456
1457 Not Collective, though the results on all processes should be the same
1458
1459 Input Parameter:
1460 . pc - the `PC` preconditioner context
1461
1462 Output Parameters:
1463 + mat - the matrix associated with the linear system was set
1464 - pmat - matrix associated with the preconditioner was set, usually the same
1465
1466 Level: intermediate
1467
1468 .seealso: [](ch_ksp), `PC`, `PCSetOperators()`, `KSPGetOperators()`, `KSPSetOperators()`, `PCGetOperators()`
1469 @*/
PCGetOperatorsSet(PC pc,PetscBool * mat,PetscBool * pmat)1470 PetscErrorCode PCGetOperatorsSet(PC pc, PetscBool *mat, PetscBool *pmat)
1471 {
1472 PetscFunctionBegin;
1473 PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1474 if (mat) *mat = (pc->mat) ? PETSC_TRUE : PETSC_FALSE;
1475 if (pmat) *pmat = (pc->pmat) ? PETSC_TRUE : PETSC_FALSE;
1476 PetscFunctionReturn(PETSC_SUCCESS);
1477 }
1478
1479 /*@
1480 PCFactorGetMatrix - Gets the factored matrix from the
1481 preconditioner context. This routine is valid only for the `PCLU`,
1482 `PCILU`, `PCCHOLESKY`, and `PCICC` methods.
1483
1484 Not Collective though `mat` is parallel if `pc` is parallel
1485
1486 Input Parameter:
1487 . pc - the `PC` preconditioner context
1488
1489 Output Parameters:
1490 . mat - the factored matrix
1491
1492 Level: advanced
1493
1494 Note:
1495 Does not increase the reference count for `mat` so DO NOT destroy it
1496
1497 .seealso: [](ch_ksp), `PC`, `PCLU`, `PCILU`, `PCCHOLESKY`, `PCICC`
1498 @*/
PCFactorGetMatrix(PC pc,Mat * mat)1499 PetscErrorCode PCFactorGetMatrix(PC pc, Mat *mat)
1500 {
1501 PetscFunctionBegin;
1502 PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1503 PetscAssertPointer(mat, 2);
1504 PetscCall(PCFactorSetUpMatSolverType(pc));
1505 PetscUseTypeMethod(pc, getfactoredmatrix, mat);
1506 PetscFunctionReturn(PETSC_SUCCESS);
1507 }
1508
1509 /*@
1510 PCSetOptionsPrefix - Sets the prefix used for searching for all
1511 `PC` options in the database.
1512
1513 Logically Collective
1514
1515 Input Parameters:
1516 + pc - the `PC` preconditioner context
1517 - prefix - the prefix string to prepend to all `PC` option requests
1518
1519 Note:
1520 A hyphen (-) must NOT be given at the beginning of the prefix name.
1521 The first character of all runtime options is AUTOMATICALLY the
1522 hyphen.
1523
1524 Level: advanced
1525
1526 .seealso: [](ch_ksp), `PC`, `PCSetFromOptions`, `PCAppendOptionsPrefix()`, `PCGetOptionsPrefix()`
1527 @*/
PCSetOptionsPrefix(PC pc,const char prefix[])1528 PetscErrorCode PCSetOptionsPrefix(PC pc, const char prefix[])
1529 {
1530 PetscFunctionBegin;
1531 PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1532 PetscCall(PetscObjectSetOptionsPrefix((PetscObject)pc, prefix));
1533 PetscFunctionReturn(PETSC_SUCCESS);
1534 }
1535
1536 /*@
1537 PCAppendOptionsPrefix - Appends to the prefix used for searching for all
1538 `PC` options in the database.
1539
1540 Logically Collective
1541
1542 Input Parameters:
1543 + pc - the `PC` preconditioner context
1544 - prefix - the prefix string to prepend to all `PC` option requests
1545
1546 Note:
1547 A hyphen (-) must NOT be given at the beginning of the prefix name.
1548 The first character of all runtime options is AUTOMATICALLY the
1549 hyphen.
1550
1551 Level: advanced
1552
1553 .seealso: [](ch_ksp), `PC`, `PCSetFromOptions`, `PCSetOptionsPrefix()`, `PCGetOptionsPrefix()`
1554 @*/
PCAppendOptionsPrefix(PC pc,const char prefix[])1555 PetscErrorCode PCAppendOptionsPrefix(PC pc, const char prefix[])
1556 {
1557 PetscFunctionBegin;
1558 PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1559 PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)pc, prefix));
1560 PetscFunctionReturn(PETSC_SUCCESS);
1561 }
1562
1563 /*@
1564 PCGetOptionsPrefix - Gets the prefix used for searching for all
1565 PC options in the database.
1566
1567 Not Collective
1568
1569 Input Parameter:
1570 . pc - the `PC` preconditioner context
1571
1572 Output Parameter:
1573 . prefix - pointer to the prefix string used, is returned
1574
1575 Level: advanced
1576
1577 .seealso: [](ch_ksp), `PC`, `PCSetFromOptions`, `PCSetOptionsPrefix()`, `PCAppendOptionsPrefix()`
1578 @*/
PCGetOptionsPrefix(PC pc,const char * prefix[])1579 PetscErrorCode PCGetOptionsPrefix(PC pc, const char *prefix[])
1580 {
1581 PetscFunctionBegin;
1582 PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1583 PetscAssertPointer(prefix, 2);
1584 PetscCall(PetscObjectGetOptionsPrefix((PetscObject)pc, prefix));
1585 PetscFunctionReturn(PETSC_SUCCESS);
1586 }
1587
1588 /*
1589 Indicates the right-hand side will be changed by KSPSolve(), this occurs for a few
1590 preconditioners including BDDC and Eisentat that transform the equations before applying
1591 the Krylov methods
1592 */
PCPreSolveChangeRHS(PC pc,PetscBool * change)1593 PETSC_INTERN PetscErrorCode PCPreSolveChangeRHS(PC pc, PetscBool *change)
1594 {
1595 PetscFunctionBegin;
1596 PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1597 PetscAssertPointer(change, 2);
1598 *change = PETSC_FALSE;
1599 PetscTryMethod(pc, "PCPreSolveChangeRHS_C", (PC, PetscBool *), (pc, change));
1600 PetscFunctionReturn(PETSC_SUCCESS);
1601 }
1602
1603 /*@
1604 PCPreSolve - Optional pre-solve phase, intended for any preconditioner-specific actions that must be performed before
1605 the iterative solve itself. Used in conjunction with `PCPostSolve()`
1606
1607 Collective
1608
1609 Input Parameters:
1610 + pc - the `PC` preconditioner context
1611 - ksp - the Krylov subspace context
1612
1613 Level: developer
1614
1615 Notes:
1616 `KSPSolve()` calls this directly, so is rarely called by the user.
1617
1618 Certain preconditioners, such as the `PCType` of `PCEISENSTAT`, change the formulation of the linear system to be solved iteratively.
1619 This function performs that transformation. `PCPostSolve()` then transforms the system back to its original form after the solve.
1620 `PCPostSolve()` also transforms the resulting solution of the transformed system to the solution of the original problem.
1621
1622 `KSPSetPostSolve()` provides an alternative way to provide such transformations.
1623
1624 .seealso: [](ch_ksp), `PC`, `PCPostSolve()`, `KSP`, `PCSetPostSetUp()`, `KSPSetPreSolve()`, `KSPSetPostSolve()`
1625 @*/
PCPreSolve(PC pc,KSP ksp)1626 PetscErrorCode PCPreSolve(PC pc, KSP ksp)
1627 {
1628 Vec x, rhs;
1629
1630 PetscFunctionBegin;
1631 PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1632 PetscValidHeaderSpecific(ksp, KSP_CLASSID, 2);
1633 pc->presolvedone++;
1634 PetscCheck(pc->presolvedone <= 2, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Cannot embed PCPreSolve() more than twice");
1635 PetscCall(KSPGetSolution(ksp, &x));
1636 PetscCall(KSPGetRhs(ksp, &rhs));
1637 PetscTryTypeMethod(pc, presolve, ksp, rhs, x);
1638 PetscFunctionReturn(PETSC_SUCCESS);
1639 }
1640
1641 /*@C
1642 PCSetPostSetUp - Sets function called at the end of `PCSetUp()` to adjust the computed preconditioner
1643
1644 Logically Collective
1645
1646 Input Parameters:
1647 + pc - the preconditioner object
1648 - postsetup - the function to call after `PCSetUp()`
1649
1650 Calling sequence of `postsetup`:
1651 . pc - the `PC` context
1652
1653 Level: developer
1654
1655 .seealso: [](ch_ksp), `PC`, `PCSetUp()`
1656 @*/
PCSetPostSetUp(PC pc,PetscErrorCode (* postsetup)(PC pc))1657 PetscErrorCode PCSetPostSetUp(PC pc, PetscErrorCode (*postsetup)(PC pc))
1658 {
1659 PetscFunctionBegin;
1660 PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1661 pc->postsetup = postsetup;
1662 PetscFunctionReturn(PETSC_SUCCESS);
1663 }
1664
1665 /*@
1666 PCPostSolve - Optional post-solve phase, intended for any
1667 preconditioner-specific actions that must be performed after
1668 the iterative solve itself.
1669
1670 Collective
1671
1672 Input Parameters:
1673 + pc - the `PC` preconditioner context
1674 - ksp - the `KSP` Krylov subspace context
1675
1676 Example Usage:
1677 .vb
1678 PCPreSolve(pc,ksp);
1679 KSPSolve(ksp,b,x);
1680 PCPostSolve(pc,ksp);
1681 .ve
1682
1683 Level: developer
1684
1685 Note:
1686 `KSPSolve()` calls this routine directly, so it is rarely called by the user.
1687
1688 .seealso: [](ch_ksp), `PC`, `KSPSetPostSolve()`, `KSPSetPreSolve()`, `PCPreSolve()`, `KSPSolve()`
1689 @*/
PCPostSolve(PC pc,KSP ksp)1690 PetscErrorCode PCPostSolve(PC pc, KSP ksp)
1691 {
1692 Vec x, rhs;
1693
1694 PetscFunctionBegin;
1695 PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1696 PetscValidHeaderSpecific(ksp, KSP_CLASSID, 2);
1697 pc->presolvedone--;
1698 PetscCall(KSPGetSolution(ksp, &x));
1699 PetscCall(KSPGetRhs(ksp, &rhs));
1700 PetscTryTypeMethod(pc, postsolve, ksp, rhs, x);
1701 PetscFunctionReturn(PETSC_SUCCESS);
1702 }
1703
1704 /*@
1705 PCLoad - Loads a `PC` that has been stored in binary with `PCView()`.
1706
1707 Collective
1708
1709 Input Parameters:
1710 + newdm - the newly loaded `PC`, this needs to have been created with `PCCreate()` or
1711 some related function before a call to `PCLoad()`.
1712 - viewer - binary file viewer `PETSCVIEWERBINARY`, obtained from `PetscViewerBinaryOpen()`
1713
1714 Level: intermediate
1715
1716 Note:
1717 The type is determined by the data in the file, any `PCType` set into the `PC` before this call is ignored.
1718
1719 .seealso: [](ch_ksp), `PC`, `PetscViewerBinaryOpen()`, `PCView()`, `MatLoad()`, `VecLoad()`, `PETSCVIEWERBINARY`
1720 @*/
PCLoad(PC newdm,PetscViewer viewer)1721 PetscErrorCode PCLoad(PC newdm, PetscViewer viewer)
1722 {
1723 PetscBool isbinary;
1724 PetscInt classid;
1725 char type[256];
1726
1727 PetscFunctionBegin;
1728 PetscValidHeaderSpecific(newdm, PC_CLASSID, 1);
1729 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
1730 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
1731 PetscCheck(isbinary, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid viewer; open viewer with PetscViewerBinaryOpen()");
1732
1733 PetscCall(PetscViewerBinaryRead(viewer, &classid, 1, NULL, PETSC_INT));
1734 PetscCheck(classid == PC_FILE_CLASSID, PetscObjectComm((PetscObject)newdm), PETSC_ERR_ARG_WRONG, "Not PC next in file");
1735 PetscCall(PetscViewerBinaryRead(viewer, type, 256, NULL, PETSC_CHAR));
1736 PetscCall(PCSetType(newdm, type));
1737 PetscTryTypeMethod(newdm, load, viewer);
1738 PetscFunctionReturn(PETSC_SUCCESS);
1739 }
1740
1741 #include <petscdraw.h>
1742 #if defined(PETSC_HAVE_SAWS)
1743 #include <petscviewersaws.h>
1744 #endif
1745
1746 /*@
1747 PCViewFromOptions - View (print or provide information about) the `PC`, based on options in the options database
1748
1749 Collective
1750
1751 Input Parameters:
1752 + A - the `PC` context
1753 . obj - Optional object that provides the options prefix
1754 - name - command line option name
1755
1756 Level: developer
1757
1758 .seealso: [](ch_ksp), `PC`, `PCView`, `PetscObjectViewFromOptions()`, `PCCreate()`
1759 @*/
PCViewFromOptions(PC A,PetscObject obj,const char name[])1760 PetscErrorCode PCViewFromOptions(PC A, PetscObject obj, const char name[])
1761 {
1762 PetscFunctionBegin;
1763 PetscValidHeaderSpecific(A, PC_CLASSID, 1);
1764 PetscCall(PetscObjectViewFromOptions((PetscObject)A, obj, name));
1765 PetscFunctionReturn(PETSC_SUCCESS);
1766 }
1767
1768 /*@
1769 PCView - Prints information about the `PC`
1770
1771 Collective
1772
1773 Input Parameters:
1774 + pc - the `PC` preconditioner context
1775 - viewer - optional `PetscViewer` visualization context
1776
1777 Level: intermediate
1778
1779 Notes:
1780 The available visualization contexts include
1781 + `PETSC_VIEWER_STDOUT_SELF` - standard output (default)
1782 - `PETSC_VIEWER_STDOUT_WORLD` - synchronized standard
1783 output where only the first processor opens
1784 the file. All other processors send their
1785 data to the first processor to print.
1786
1787 The user can open an alternative visualization contexts with
1788 `PetscViewerASCIIOpen()` (output to a specified file).
1789
1790 .seealso: [](ch_ksp), `PC`, `PetscViewer`, `PetscViewerType`, `KSPView()`, `PetscViewerASCIIOpen()`
1791 @*/
PCView(PC pc,PetscViewer viewer)1792 PetscErrorCode PCView(PC pc, PetscViewer viewer)
1793 {
1794 PCType cstr;
1795 PetscViewerFormat format;
1796 PetscBool isascii, isstring, isbinary, isdraw, pop = PETSC_FALSE;
1797 #if defined(PETSC_HAVE_SAWS)
1798 PetscBool issaws;
1799 #endif
1800
1801 PetscFunctionBegin;
1802 PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1803 if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)pc), &viewer));
1804 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
1805 PetscCheckSameComm(pc, 1, viewer, 2);
1806
1807 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
1808 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSTRING, &isstring));
1809 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
1810 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
1811 #if defined(PETSC_HAVE_SAWS)
1812 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSAWS, &issaws));
1813 #endif
1814
1815 if (isascii) {
1816 PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)pc, viewer));
1817 if (!pc->setupcalled) PetscCall(PetscViewerASCIIPrintf(viewer, " PC has not been set up so information may be incomplete\n"));
1818 PetscCall(PetscViewerASCIIPushTab(viewer));
1819 PetscTryTypeMethod(pc, view, viewer);
1820 PetscCall(PetscViewerASCIIPopTab(viewer));
1821 if (pc->mat) {
1822 PetscCall(PetscViewerGetFormat(viewer, &format));
1823 if (format != PETSC_VIEWER_ASCII_INFO_DETAIL) {
1824 PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_INFO));
1825 pop = PETSC_TRUE;
1826 }
1827 if (pc->pmat == pc->mat) {
1828 PetscCall(PetscViewerASCIIPrintf(viewer, " linear system matrix, which is also used to construct the preconditioner:\n"));
1829 PetscCall(PetscViewerASCIIPushTab(viewer));
1830 PetscCall(MatView(pc->mat, viewer));
1831 PetscCall(PetscViewerASCIIPopTab(viewer));
1832 } else {
1833 if (pc->pmat) {
1834 PetscCall(PetscViewerASCIIPrintf(viewer, " linear system matrix, followed by the matrix used to construct the preconditioner:\n"));
1835 } else {
1836 PetscCall(PetscViewerASCIIPrintf(viewer, " linear system matrix:\n"));
1837 }
1838 PetscCall(PetscViewerASCIIPushTab(viewer));
1839 PetscCall(MatView(pc->mat, viewer));
1840 if (pc->pmat) PetscCall(MatView(pc->pmat, viewer));
1841 PetscCall(PetscViewerASCIIPopTab(viewer));
1842 }
1843 if (pop) PetscCall(PetscViewerPopFormat(viewer));
1844 }
1845 } else if (isstring) {
1846 PetscCall(PCGetType(pc, &cstr));
1847 PetscCall(PetscViewerStringSPrintf(viewer, " PCType: %-7.7s", cstr));
1848 PetscTryTypeMethod(pc, view, viewer);
1849 if (pc->mat) PetscCall(MatView(pc->mat, viewer));
1850 if (pc->pmat && pc->pmat != pc->mat) PetscCall(MatView(pc->pmat, viewer));
1851 } else if (isbinary) {
1852 PetscInt classid = PC_FILE_CLASSID;
1853 MPI_Comm comm;
1854 PetscMPIInt rank;
1855 char type[256];
1856
1857 PetscCall(PetscObjectGetComm((PetscObject)pc, &comm));
1858 PetscCallMPI(MPI_Comm_rank(comm, &rank));
1859 if (rank == 0) {
1860 PetscCall(PetscViewerBinaryWrite(viewer, &classid, 1, PETSC_INT));
1861 PetscCall(PetscStrncpy(type, ((PetscObject)pc)->type_name, 256));
1862 PetscCall(PetscViewerBinaryWrite(viewer, type, 256, PETSC_CHAR));
1863 }
1864 PetscTryTypeMethod(pc, view, viewer);
1865 } else if (isdraw) {
1866 PetscDraw draw;
1867 char str[25];
1868 PetscReal x, y, bottom, h;
1869 PetscInt n;
1870
1871 PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
1872 PetscCall(PetscDrawGetCurrentPoint(draw, &x, &y));
1873 if (pc->mat) {
1874 PetscCall(MatGetSize(pc->mat, &n, NULL));
1875 PetscCall(PetscSNPrintf(str, 25, "PC: %s (%" PetscInt_FMT ")", ((PetscObject)pc)->type_name, n));
1876 } else {
1877 PetscCall(PetscSNPrintf(str, 25, "PC: %s", ((PetscObject)pc)->type_name));
1878 }
1879 PetscCall(PetscDrawStringBoxed(draw, x, y, PETSC_DRAW_RED, PETSC_DRAW_BLACK, str, NULL, &h));
1880 bottom = y - h;
1881 PetscCall(PetscDrawPushCurrentPoint(draw, x, bottom));
1882 PetscTryTypeMethod(pc, view, viewer);
1883 PetscCall(PetscDrawPopCurrentPoint(draw));
1884 #if defined(PETSC_HAVE_SAWS)
1885 } else if (issaws) {
1886 PetscMPIInt rank;
1887
1888 PetscCall(PetscObjectName((PetscObject)pc));
1889 PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
1890 if (!((PetscObject)pc)->amsmem && rank == 0) PetscCall(PetscObjectViewSAWs((PetscObject)pc, viewer));
1891 if (pc->mat) PetscCall(MatView(pc->mat, viewer));
1892 if (pc->pmat && pc->pmat != pc->mat) PetscCall(MatView(pc->pmat, viewer));
1893 #endif
1894 }
1895 PetscFunctionReturn(PETSC_SUCCESS);
1896 }
1897
1898 /*@C
1899 PCRegister - Adds a method (`PCType`) to the PETSc preconditioner package.
1900
1901 Not collective. No Fortran Support
1902
1903 Input Parameters:
1904 + sname - name of a new user-defined solver
1905 - function - routine to create the method context which will be stored in a `PC` when `PCSetType()` is called
1906
1907 Example Usage:
1908 .vb
1909 PCRegister("my_solver", MySolverCreate);
1910 .ve
1911
1912 Then, your solver can be chosen with the procedural interface via
1913 .vb
1914 PCSetType(pc, "my_solver")
1915 .ve
1916 or at runtime via the option
1917 .vb
1918 -pc_type my_solver
1919 .ve
1920
1921 Level: advanced
1922
1923 Note:
1924 A simpler alternative to using `PCRegister()` for an application specific preconditioner is to use a `PC` of `PCType` `PCSHELL` and
1925 provide your customizations with `PCShellSetContext()` and `PCShellSetApply()`
1926
1927 `PCRegister()` may be called multiple times to add several user-defined preconditioners.
1928
1929 .seealso: [](ch_ksp), `PC`, `PCType`, `PCRegisterAll()`, `PCSetType()`, `PCShellSetContext()`, `PCShellSetApply()`, `PCSHELL`
1930 @*/
PCRegister(const char sname[],PetscErrorCode (* function)(PC))1931 PetscErrorCode PCRegister(const char sname[], PetscErrorCode (*function)(PC))
1932 {
1933 PetscFunctionBegin;
1934 PetscCall(PCInitializePackage());
1935 PetscCall(PetscFunctionListAdd(&PCList, sname, function));
1936 PetscFunctionReturn(PETSC_SUCCESS);
1937 }
1938
MatMult_PC(Mat A,Vec X,Vec Y)1939 static PetscErrorCode MatMult_PC(Mat A, Vec X, Vec Y)
1940 {
1941 PC pc;
1942
1943 PetscFunctionBegin;
1944 PetscCall(MatShellGetContext(A, &pc));
1945 PetscCall(PCApply(pc, X, Y));
1946 PetscFunctionReturn(PETSC_SUCCESS);
1947 }
1948
1949 /*@
1950 PCComputeOperator - Computes the explicit preconditioned operator as a matrix `Mat`.
1951
1952 Collective
1953
1954 Input Parameters:
1955 + pc - the `PC` preconditioner object
1956 - mattype - the `MatType` to be used for the operator
1957
1958 Output Parameter:
1959 . mat - the explicit preconditioned operator
1960
1961 Level: advanced
1962
1963 Note:
1964 This computation is done by applying the operators to columns of the identity matrix.
1965 This routine is costly in general, and is recommended for use only with relatively small systems.
1966 Currently, this routine uses a dense matrix format when `mattype` == `NULL`
1967
1968 Developer Note:
1969 This should be called `PCCreateExplicitOperator()`
1970
1971 .seealso: [](ch_ksp), `PC`, `KSPComputeOperator()`, `MatType`
1972 @*/
PCComputeOperator(PC pc,MatType mattype,Mat * mat)1973 PetscErrorCode PCComputeOperator(PC pc, MatType mattype, Mat *mat)
1974 {
1975 PetscInt N, M, m, n;
1976 Mat A, Apc;
1977
1978 PetscFunctionBegin;
1979 PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1980 PetscAssertPointer(mat, 3);
1981 PetscCall(PCGetOperators(pc, &A, NULL));
1982 PetscCall(MatGetLocalSize(A, &m, &n));
1983 PetscCall(MatGetSize(A, &M, &N));
1984 PetscCall(MatCreateShell(PetscObjectComm((PetscObject)pc), m, n, M, N, pc, &Apc));
1985 PetscCall(MatShellSetOperation(Apc, MATOP_MULT, (PetscErrorCodeFn *)MatMult_PC));
1986 PetscCall(MatComputeOperator(Apc, mattype, mat));
1987 PetscCall(MatDestroy(&Apc));
1988 PetscFunctionReturn(PETSC_SUCCESS);
1989 }
1990
1991 /*@
1992 PCSetCoordinates - sets the coordinates of all the nodes (degrees of freedom in the vector) on the local process
1993
1994 Collective
1995
1996 Input Parameters:
1997 + pc - the `PC` preconditioner context
1998 . dim - the dimension of the coordinates 1, 2, or 3
1999 . nloc - the blocked size of the coordinates array
2000 - coords - the coordinates array
2001
2002 Level: intermediate
2003
2004 Notes:
2005 `coords` is an array of the dim coordinates for the nodes on
2006 the local processor, of size `dim`*`nloc`.
2007 If there are 108 equations (dofs) on a processor
2008 for a 3d displacement finite element discretization of elasticity (so
2009 that there are nloc = 36 = 108/3 nodes) then the array must have 108
2010 double precision values (ie, 3 * 36). These x y z coordinates
2011 should be ordered for nodes 0 to N-1 like so: [ 0.x, 0.y, 0.z, 1.x,
2012 ... , N-1.z ].
2013
2014 The information provided here can be used by some preconditioners, such as `PCGAMG`, to produce a better preconditioner.
2015 See also `MatSetNearNullSpace()`.
2016
2017 .seealso: [](ch_ksp), `PC`, `MatSetNearNullSpace()`
2018 @*/
PCSetCoordinates(PC pc,PetscInt dim,PetscInt nloc,PetscReal coords[])2019 PetscErrorCode PCSetCoordinates(PC pc, PetscInt dim, PetscInt nloc, PetscReal coords[])
2020 {
2021 PetscFunctionBegin;
2022 PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
2023 PetscValidLogicalCollectiveInt(pc, dim, 2);
2024 PetscTryMethod(pc, "PCSetCoordinates_C", (PC, PetscInt, PetscInt, PetscReal[]), (pc, dim, nloc, coords));
2025 PetscFunctionReturn(PETSC_SUCCESS);
2026 }
2027
2028 /*@
2029 PCGetInterpolations - Gets interpolation matrices for all levels (except level 0)
2030
2031 Logically Collective
2032
2033 Input Parameter:
2034 . pc - the precondition context
2035
2036 Output Parameters:
2037 + num_levels - the number of levels
2038 - interpolations - the interpolation matrices (size of `num_levels`-1)
2039
2040 Level: advanced
2041
2042 Developer Note:
2043 Why is this here instead of in `PCMG` etc?
2044
2045 .seealso: [](ch_ksp), `PC`, `PCMG`, `PCMGGetRestriction()`, `PCMGSetInterpolation()`, `PCMGGetInterpolation()`, `PCGetCoarseOperators()`
2046 @*/
PCGetInterpolations(PC pc,PetscInt * num_levels,Mat * interpolations[])2047 PetscErrorCode PCGetInterpolations(PC pc, PetscInt *num_levels, Mat *interpolations[])
2048 {
2049 PetscFunctionBegin;
2050 PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
2051 PetscAssertPointer(num_levels, 2);
2052 PetscAssertPointer(interpolations, 3);
2053 PetscUseMethod(pc, "PCGetInterpolations_C", (PC, PetscInt *, Mat *[]), (pc, num_levels, interpolations));
2054 PetscFunctionReturn(PETSC_SUCCESS);
2055 }
2056
2057 /*@
2058 PCGetCoarseOperators - Gets coarse operator matrices for all levels (except the finest level)
2059
2060 Logically Collective
2061
2062 Input Parameter:
2063 . pc - the precondition context
2064
2065 Output Parameters:
2066 + num_levels - the number of levels
2067 - coarseOperators - the coarse operator matrices (size of `num_levels`-1)
2068
2069 Level: advanced
2070
2071 Developer Note:
2072 Why is this here instead of in `PCMG` etc?
2073
2074 .seealso: [](ch_ksp), `PC`, `PCMG`, `PCMGGetRestriction()`, `PCMGSetInterpolation()`, `PCMGGetRScale()`, `PCMGGetInterpolation()`, `PCGetInterpolations()`
2075 @*/
PCGetCoarseOperators(PC pc,PetscInt * num_levels,Mat * coarseOperators[])2076 PetscErrorCode PCGetCoarseOperators(PC pc, PetscInt *num_levels, Mat *coarseOperators[])
2077 {
2078 PetscFunctionBegin;
2079 PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
2080 PetscAssertPointer(num_levels, 2);
2081 PetscAssertPointer(coarseOperators, 3);
2082 PetscUseMethod(pc, "PCGetCoarseOperators_C", (PC, PetscInt *, Mat *[]), (pc, num_levels, coarseOperators));
2083 PetscFunctionReturn(PETSC_SUCCESS);
2084 }
2085