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 (*applyBAtranspose)(PC, PetscInt, Vec, Vec, Vec);
 19:   PetscErrorCode (*setfromoptions)(PC, PetscOptionItems *);
 20:   PetscErrorCode (*presolve)(PC, KSP, Vec, Vec);
 21:   PetscErrorCode (*postsolve)(PC, KSP, Vec, Vec);
 22:   PetscErrorCode (*getfactoredmatrix)(PC, Mat *);
 23:   PetscErrorCode (*applysymmetricleft)(PC, Vec, Vec);
 24:   PetscErrorCode (*applysymmetricright)(PC, Vec, Vec);
 25:   PetscErrorCode (*setuponblocks)(PC);
 26:   PetscErrorCode (*destroy)(PC);
 27:   PetscErrorCode (*view)(PC, PetscViewer);
 28:   PetscErrorCode (*reset)(PC);
 29:   PetscErrorCode (*load)(PC, PetscViewer);
 30: };

 32: /*
 33:    Preconditioner context
 34: */
 35: struct _p_PC {
 36:   PETSCHEADER(struct _PCOps);
 37:   DM               dm;
 38:   PetscInt         setupcalled;
 39:   PetscObjectState matstate, matnonzerostate; /* last known nonzero state of the pmat associated with this PC */
 40:   PetscBool        reusepreconditioner;
 41:   MatStructure     flag; /* reset each PCSetUp() to indicate to PC implementations if nonzero structure has changed */

 43:   PetscInt  setfromoptionscalled;
 44:   PetscBool erroriffailure; /* Generate an error if FPE detected (for example a zero pivot) instead of returning*/
 45:   Mat       mat, pmat;
 46:   Vec       diagonalscaleright, diagonalscaleleft; /* used for time integration scaling */
 47:   PetscBool diagonalscale;
 48:   PetscBool useAmat;                                                                        /* used by several PC that including applying the operator inside the preconditioner */
 49:   PetscErrorCode (*modifysubmatrices)(PC, PetscInt, const IS[], const IS[], Mat[], void *); /* user provided routine */
 50:   void          *modifysubmatricesP;                                                        /* context for user routine */
 51:   void          *data;
 52:   PetscInt       presolvedone;     /* has PCPreSolve() already been run */
 53:   void          *user;             /* 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 */

 57:   PetscErrorCode (*presolve)(PC, KSP);

 59:   PetscInt kspnestlevel; /* how many levels of nesting does the KSP have that contains the PC */
 60: };

 62: PETSC_EXTERN PetscLogEvent PC_SetUp;
 63: PETSC_EXTERN PetscLogEvent PC_SetUpOnBlocks;
 64: PETSC_EXTERN PetscLogEvent PC_Apply;
 65: PETSC_EXTERN PetscLogEvent PC_MatApply;
 66: PETSC_EXTERN PetscLogEvent PC_ApplyCoarse;
 67: PETSC_EXTERN PetscLogEvent PC_ApplyMultiple;
 68: PETSC_EXTERN PetscLogEvent PC_ApplySymmetricLeft;
 69: PETSC_EXTERN PetscLogEvent PC_ApplySymmetricRight;
 70: PETSC_EXTERN PetscLogEvent PC_ModifySubMatrices;
 71: PETSC_EXTERN PetscLogEvent PC_ApplyOnBlocks;
 72: PETSC_EXTERN PetscLogEvent PC_ApplyTransposeOnBlocks;
 73: PETSC_EXTERN PetscLogStage PCMPIStage;