Actual source code: pcimpl.h

petsc-3.6.4 2016-04-12
Report Typos and Errors
  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 (*applyrichardson)(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool ,PetscInt*,PCRichardsonConvergedReason*);
 17:   PetscErrorCode (*applyBA)(PC,PCSide,Vec,Vec,Vec);
 18:   PetscErrorCode (*applytranspose)(PC,Vec,Vec);
 19:   PetscErrorCode (*applyBAtranspose)(PC,PetscInt,Vec,Vec,Vec);
 20:   PetscErrorCode (*setfromoptions)(PetscOptions*,PC);
 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:   PetscInt         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        nonzero_guess; /* used by PCKSP, PCREDUNDANT */
 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: };

 58: PETSC_EXTERN PetscLogEvent PC_SetUp, PC_SetUpOnBlocks, PC_Apply, PC_ApplyCoarse, PC_ApplyMultiple, PC_ApplySymmetricLeft;
 59: PETSC_EXTERN PetscLogEvent PC_ApplySymmetricRight, PC_ModifySubMatrices, PC_ApplyOnBlocks, PC_ApplyTransposeOnBlocks, PC_ApplyOnMproc;

 61: #endif