Actual source code: pcimpl.h
1: #pragma once
3: #include <petscksp.h>
4: #include <petscpc.h>
5: #include <petsc/private/petscimpl.h>
7: PETSC_EXTERN PetscBool PCRegisterAllCalled;
8: PETSC_EXTERN PetscErrorCode PCRegisterAll(void);
10: typedef struct _PCOps *PCOps;
11: struct _PCOps {
12: PetscErrorCode (*setup)(PC);
13: PetscErrorCode (*apply)(PC, Vec, Vec);
14: PetscErrorCode (*matapply)(PC, Mat, Mat);
15: PetscErrorCode (*applyrichardson)(PC, Vec, Vec, Vec, PetscReal, PetscReal, PetscReal, PetscInt, PetscBool, PetscInt *, PCRichardsonConvergedReason *);
16: PetscErrorCode (*applyBA)(PC, PCSide, Vec, Vec, Vec);
17: PetscErrorCode (*applytranspose)(PC, Vec, Vec);
18: PetscErrorCode (*matapplytranspose)(PC, Mat, Mat);
19: PetscErrorCode (*applyBAtranspose)(PC, PetscInt, Vec, Vec, Vec);
20: PetscErrorCode (*setfromoptions)(PC, PetscOptionItems);
21: PetscErrorCode (*presolve)(PC, KSP, Vec, Vec);
22: PetscErrorCode (*postsolve)(PC, KSP, Vec, Vec);
23: PetscErrorCode (*getfactoredmatrix)(PC, Mat *);
24: PetscErrorCode (*applysymmetricleft)(PC, Vec, Vec);
25: PetscErrorCode (*applysymmetricright)(PC, Vec, Vec);
26: PetscErrorCode (*setuponblocks)(PC);
27: PetscErrorCode (*destroy)(PC);
28: PetscErrorCode (*view)(PC, PetscViewer);
29: PetscErrorCode (*reset)(PC);
30: PetscErrorCode (*load)(PC, PetscViewer);
31: };
33: /*
34: Preconditioner context
35: */
36: struct _p_PC {
37: PETSCHEADER(struct _PCOps);
38: DM dm;
39: PetscBool setupcalled;
40: PetscObjectState matstate, matnonzerostate; /* last known nonzero state of the pmat associated with this PC */
41: PetscBool reusepreconditioner;
42: MatStructure flag; /* reset each PCSetUp() to indicate to PC implementations if nonzero structure has changed */
44: PetscInt setfromoptionscalled;
45: PetscBool erroriffailure; /* Generate an error if FPE detected (for example a zero pivot) instead of returning*/
46: Mat mat, pmat;
47: Vec diagonalscaleright, diagonalscaleleft; /* used for time integration scaling */
48: PetscBool diagonalscale;
49: PetscBool useAmat; /* used by several PC that including applying the operator inside the preconditioner */
50: PetscErrorCode (*modifysubmatrices)(PC, PetscInt, const IS[], const IS[], Mat[], void *); /* user provided routine */
51: void *modifysubmatricesP; /* context for user routine */
52: void *data;
53: void *ctx; /* optional user-defined context */
54: PCFailedReason failedreason; /* after VecNorm or VecDot contains maximum of all rank failed reasons */
55: PCFailedReason failedreasonrank; /* failed reason on this rank */
56: PetscInt presolvedone;
57: PetscErrorCode (*postsetup)(PC);
58: PetscInt kspnestlevel; /* how many levels of nesting does the KSP have that contains the PC */
59: };
61: PETSC_EXTERN PetscLogEvent PC_SetUp;
62: PETSC_EXTERN PetscLogEvent PC_SetUpOnBlocks;
63: PETSC_EXTERN PetscLogEvent PC_Apply;
64: PETSC_EXTERN PetscLogEvent PC_MatApply;
65: PETSC_EXTERN PetscLogEvent PC_ApplyCoarse;
66: PETSC_EXTERN PetscLogEvent PC_ApplySymmetricLeft;
67: PETSC_EXTERN PetscLogEvent PC_ApplySymmetricRight;
68: PETSC_EXTERN PetscLogEvent PC_ModifySubMatrices;
69: PETSC_EXTERN PetscLogEvent PC_ApplyOnBlocks;
70: PETSC_EXTERN PetscLogEvent PC_ApplyTransposeOnBlocks;
71: PETSC_EXTERN PetscLogStage PCMPIStage;