Actual source code: pcimpl.h


  2: #ifndef _PCIMPL_H
  3: #define _PCIMPL_H

  5: #include <petscksp.h>
  6: #include <petscpc.h>
  7: #include <petsc/private/petscimpl.h>

  9: PETSC_EXTERN PetscBool      PCRegisterAllCalled;
 10: PETSC_EXTERN PetscErrorCode PCRegisterAll(void);

 12: typedef struct _PCOps *PCOps;
 13: struct _PCOps {
 14:   PetscErrorCode (*setup)(PC);
 15:   PetscErrorCode (*apply)(PC, Vec, Vec);
 16:   PetscErrorCode (*matapply)(PC, Mat, Mat);
 17:   PetscErrorCode (*applyrichardson)(PC, Vec, Vec, Vec, PetscReal, PetscReal, PetscReal, PetscInt, PetscBool, PetscInt *, PCRichardsonConvergedReason *);
 18:   PetscErrorCode (*applyBA)(PC, PCSide, Vec, Vec, Vec);
 19:   PetscErrorCode (*applytranspose)(PC, Vec, Vec);
 20:   PetscErrorCode (*applyBAtranspose)(PC, PetscInt, Vec, Vec, Vec);
 21:   PetscErrorCode (*setfromoptions)(PC, PetscOptionItems *);
 22:   PetscErrorCode (*presolve)(PC, KSP, Vec, Vec);
 23:   PetscErrorCode (*postsolve)(PC, KSP, Vec, Vec);
 24:   PetscErrorCode (*getfactoredmatrix)(PC, Mat *);
 25:   PetscErrorCode (*applysymmetricleft)(PC, Vec, Vec);
 26:   PetscErrorCode (*applysymmetricright)(PC, Vec, Vec);
 27:   PetscErrorCode (*setuponblocks)(PC);
 28:   PetscErrorCode (*destroy)(PC);
 29:   PetscErrorCode (*view)(PC, PetscViewer);
 30:   PetscErrorCode (*reset)(PC);
 31:   PetscErrorCode (*load)(PC, PetscViewer);
 32: };

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

 45:   PetscInt  setfromoptionscalled;
 46:   PetscBool erroriffailure; /* Generate an error if FPE detected (for example a zero pivot) instead of returning*/
 47:   Mat       mat, pmat;
 48:   Vec       diagonalscaleright, diagonalscaleleft; /* used for time integration scaling */
 49:   PetscBool diagonalscale;
 50:   PetscBool useAmat;                                                                        /* used by several PC that including applying the operator inside the preconditioner */
 51:   PetscErrorCode (*modifysubmatrices)(PC, PetscInt, const IS[], const IS[], Mat[], void *); /* user provided routine */
 52:   void          *modifysubmatricesP;                                                        /* context for user routine */
 53:   void          *data;
 54:   PetscInt       presolvedone;     /* has PCPreSolve() already been run */
 55:   void          *user;             /* optional user-defined context */
 56:   PCFailedReason failedreason;     /* after VecNorm or VecDot contains maximum of all rank failed reasons */
 57:   PCFailedReason failedreasonrank; /* failed reason on this rank */

 59:   PetscErrorCode (*presolve)(PC, KSP);
 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;

 75: #endif