Actual source code: petsc-kspimpl.h

petsc-3.4.2 2013-07-02
  2: #ifndef _KSPIMPL_H
  3: #define _KSPIMPL_H

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

  8: typedef struct _KSPOps *KSPOps;

 10: struct _KSPOps {
 11:   PetscErrorCode (*buildsolution)(KSP,Vec,Vec*);       /* Returns a pointer to the solution, or
 12:                                                           calculates the solution in a
 13:                                                           user-provided area. */
 14:   PetscErrorCode (*buildresidual)(KSP,Vec,Vec,Vec*);   /* Returns a pointer to the residual, or
 15:                                                           calculates the residual in a
 16:                                                           user-provided area.  */
 17:   PetscErrorCode (*solve)(KSP);                        /* actual solver */
 18:   PetscErrorCode (*setup)(KSP);
 19:   PetscErrorCode (*setfromoptions)(KSP);
 20:   PetscErrorCode (*publishoptions)(KSP);
 21:   PetscErrorCode (*computeextremesingularvalues)(KSP,PetscReal*,PetscReal*);
 22:   PetscErrorCode (*computeeigenvalues)(KSP,PetscInt,PetscReal*,PetscReal*,PetscInt *);
 23:   PetscErrorCode (*destroy)(KSP);
 24:   PetscErrorCode (*view)(KSP,PetscViewer);
 25:   PetscErrorCode (*reset)(KSP);
 26:   PetscErrorCode (*load)(KSP,PetscViewer);
 27: };

 29: typedef struct {PetscInt model,curl,maxl;Mat mat; KSP ksp;}* KSPGuessFischer;

 31: /*
 32:      Maximum number of monitors you can run with a single KSP
 33: */
 34: #define MAXKSPMONITORS 5
 35: typedef enum {KSP_SETUP_NEW, KSP_SETUP_NEWMATRIX, KSP_SETUP_NEWRHS} KSPSetUpStage;

 37: /*
 38:    Defines the KSP data structure.
 39: */
 40: struct _p_KSP {
 41:   PETSCHEADER(struct _KSPOps);
 42:   DM              dm;
 43:   PetscBool       dmAuto;       /* DM was created automatically by KSP */
 44:   PetscBool       dmActive;     /* KSP should use DM for computing operators */
 45:   /*------------------------- User parameters--------------------------*/
 46:   PetscInt        max_it;                     /* maximum number of iterations */
 47:   KSPFischerGuess guess;
 48:   PetscBool       guess_zero,                  /* flag for whether initial guess is 0 */
 49:                   calc_sings,                  /* calculate extreme Singular Values */
 50:                   guess_knoll;                /* use initial guess of PCApply(ksp->B,b */
 51:   PCSide          pc_side;                  /* flag for left, right, or symmetric preconditioning */
 52:   PetscInt        normsupporttable[KSP_NORM_MAX][PC_SIDE_MAX]; /* Table of supported norms and pc_side, see KSPSetSupportedNorm() */
 53:   PetscReal       rtol,                     /* relative tolerance */
 54:                   abstol,                     /* absolute tolerance */
 55:                   ttol,                     /* (not set by user)  */
 56:                   divtol;                   /* divergence tolerance */
 57:   PetscReal       rnorm0;                   /* initial residual norm (used for divergence testing) */
 58:   PetscReal       rnorm;                    /* current residual norm */
 59:   KSPConvergedReason reason;
 60:   PetscBool          printreason;     /* prints converged reason after solve */
 61:   PetscBool          errorifnotconverged;    /* create an error if the KSPSolve() does not converge */

 63:   Vec vec_sol,vec_rhs;            /* pointer to where user has stashed
 64:                                       the solution and rhs, these are
 65:                                       never touched by the code, only
 66:                                       passed back to the user */
 67:   PetscReal     *res_hist;            /* If !0 stores residual at iterations*/
 68:   PetscReal     *res_hist_alloc;      /* If !0 means user did not provide buffer, needs deallocation */
 69:   PetscInt      res_hist_len;         /* current size of residual history array */
 70:   PetscInt      res_hist_max;         /* actual amount of data in residual_history */
 71:   PetscBool     res_hist_reset;       /* reset history to size zero for each new solve */

 73:   PetscInt      chknorm;             /* only compute/check norm if iterations is great than this */
 74:   PetscBool     lagnorm;             /* Lag the residual norm calculation so that it is computed as part of the
 75:                                         MPI_Allreduce() for computing the inner products for the next iteration. */
 76:   /* --------User (or default) routines (most return -1 on error) --------*/
 77:   PetscErrorCode (*monitor[MAXKSPMONITORS])(KSP,PetscInt,PetscReal,void*); /* returns control to user after */
 78:   PetscErrorCode (*monitordestroy[MAXKSPMONITORS])(void**);         /* */
 79:   void *monitorcontext[MAXKSPMONITORS];                  /* residual calculation, allows user */
 80:   PetscInt  numbermonitors;                                   /* to, for instance, print residual norm, etc. */

 82:   PetscErrorCode (*converged)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*);
 83:   PetscErrorCode (*convergeddestroy)(void*);
 84:   void       *cnvP;

 86:   void       *user;             /* optional user-defined context */

 88:   PC         pc;

 90:   void       *data;                      /* holder for misc stuff associated
 91:                                    with a particular iterative solver */

 93:   /* ----------------Default work-area management -------------------- */
 94:   PetscInt       nwork;
 95:   Vec            *work;

 97:   KSPSetUpStage  setupstage;

 99:   PetscInt       its;       /* number of iterations so far computed */

101:   PetscBool      transpose_solve;    /* solve transpose system instead */

103:   KSPNormType    normtype;          /* type of norm used for convergence tests */

105:   /*   Allow diagonally scaling the matrix before computing the preconditioner or using
106:        the Krylov method. Note this is NOT just Jacobi preconditioning */

108:   PetscBool    dscale;       /* diagonal scale system; used with KSPSetDiagonalScale() */
109:   PetscBool    dscalefix;    /* unscale system after solve */
110:   PetscBool    dscalefix2;   /* system has been unscaled */
111:   Vec          diagonal;     /* 1/sqrt(diag of matrix) */
112:   Vec          truediagonal;

114:   MatNullSpace nullsp;      /* Null space of the operator, removed from Krylov space */

116:   PetscViewer  eigviewer;   /* Viewer where computed eigenvalues are displayed */

118:   PetscErrorCode (*presolve)(KSP,Vec,Vec,void*);
119:   PetscErrorCode (*postsolve)(KSP,Vec,Vec,void*);
120:   void           *prectx,*postctx;
121: };

123: typedef struct { /* dummy data structure used in KSPMonitorDynamicTolerance() */
124:   PetscReal coef;
125:   PetscReal bnrm;
126: } KSPDynTolCtx;

128: typedef struct {
129:   PetscBool  initialrtol;    /* default relative residual decrease is computing from initial residual, not rhs */
130:   PetscBool  mininitialrtol; /* default relative residual decrease is computing from min of initial residual and rhs */
131:   Vec        work;
132: } KSPDefaultConvergedCtx;

136: PETSC_STATIC_INLINE PetscErrorCode KSPLogResidualHistory(KSP ksp,PetscReal norm)
137: {

141:   PetscObjectAMSTakeAccess((PetscObject)ksp);
142:   if (ksp->res_hist && ksp->res_hist_max > ksp->res_hist_len) {
143:     ksp->res_hist[ksp->res_hist_len++] = norm;
144:   }
145:   PetscObjectAMSGrantAccess((PetscObject)ksp);
146:   return(0);
147: }

149: PETSC_INTERN PetscErrorCode KSPSetUpNorms_Private(KSP,KSPNormType*,PCSide*);

151: PETSC_INTERN PetscErrorCode KSPPlotEigenContours_Private(KSP,PetscInt,const PetscReal*,const PetscReal*);

153: typedef struct _p_DMKSP *DMKSP;
154: typedef struct _DMKSPOps *DMKSPOps;
155: struct _DMKSPOps {
156:   PetscErrorCode (*computeoperators)(KSP,Mat,Mat,MatStructure*,void*);
157:   PetscErrorCode (*computerhs)(KSP,Vec,void*);
158:   PetscErrorCode (*computeinitialguess)(KSP,Vec,void*);
159:   PetscErrorCode (*destroy)(DMKSP*);
160:   PetscErrorCode (*duplicate)(DMKSP,DMKSP);
161: };

163: struct _p_DMKSP {
164:   PETSCHEADER(struct _DMKSPOps);
165:   void *operatorsctx;
166:   void *rhsctx;
167:   void *initialguessctx;
168:   void *data;

170:   /* This is NOT reference counted. The DM on which this context was first created is cached here to implement one-way
171:    * copy-on-write. When DMGetDMKSPWrite() sees a request using a different DM, it makes a copy. Thus, if a user
172:    * only interacts directly with one level, e.g., using KSPSetComputeOperators(), then coarse levels are constructed by
173:    * PCMG, then the user changes the routine with another call to KSPSetComputeOperators(), it automatically propagates
174:    * to all the levels. If instead, they get out a specific level and set the routines on that level, subsequent changes
175:    * to the original level will no longer propagate to that level.
176:    */
177:   DM originaldm;

179:   void (*fortran_func_pointers[3])(void); /* Store our own function pointers so they are associated with the DMKSP instead of the DM */
180: };
181: PETSC_EXTERN PetscErrorCode DMGetDMKSP(DM,DMKSP*);
182: PETSC_EXTERN PetscErrorCode DMGetDMKSPWrite(DM,DMKSP*);
183: PETSC_EXTERN PetscErrorCode DMCopyDMKSP(DM,DM);

185: /*
186:        These allow the various Krylov methods to apply to either the linear system or its transpose.
187: */
190: PETSC_STATIC_INLINE PetscErrorCode KSP_RemoveNullSpace(KSP ksp,Vec y)
191: {
194:   if (ksp->nullsp && ksp->pc_side == PC_LEFT) {MatNullSpaceRemove(ksp->nullsp,y,NULL);}
195:   return(0);
196: }

200: PETSC_STATIC_INLINE PetscErrorCode KSP_MatMult(KSP ksp,Mat A,Vec x,Vec y)
201: {
204:   if (!ksp->transpose_solve) {MatMult(A,x,y);}
205:   else                       {MatMultTranspose(A,x,y);}
206:   return(0);
207: }

211: PETSC_STATIC_INLINE PetscErrorCode KSP_MatMultTranspose(KSP ksp,Mat A,Vec x,Vec y)
212: {
215:   if (!ksp->transpose_solve) {MatMultTranspose(A,x,y);}
216:   else                       {MatMult(A,x,y);}
217:   return(0);
218: }

222: PETSC_STATIC_INLINE PetscErrorCode KSP_PCApply(KSP ksp,Vec x,Vec y)
223: {
226:   if (!ksp->transpose_solve) {
227:     PCApply(ksp->pc,x,y);
228:     KSP_RemoveNullSpace(ksp,y);
229:   } else {
230:     PCApplyTranspose(ksp->pc,x,y);
231:   }
232:   return(0);
233: }

237: PETSC_STATIC_INLINE PetscErrorCode KSP_PCApplyTranspose(KSP ksp,Vec x,Vec y)
238: {
241:   if (!ksp->transpose_solve) {
242:     PCApplyTranspose(ksp->pc,x,y);
243:   } else {
244:     PCApply(ksp->pc,x,y);
245:     KSP_RemoveNullSpace(ksp,y);
246:   }
247:   return(0);
248: }

252: PETSC_STATIC_INLINE PetscErrorCode KSP_PCApplyBAorAB(KSP ksp,Vec x,Vec y,Vec w)
253: {
256:   if (!ksp->transpose_solve) {
257:     PCApplyBAorAB(ksp->pc,ksp->pc_side,x,y,w);
258:     KSP_RemoveNullSpace(ksp,y);
259:   } else {
260:     PCApplyBAorABTranspose(ksp->pc,ksp->pc_side,x,y,w);
261:   }
262:   return(0);
263: }

267: PETSC_STATIC_INLINE PetscErrorCode KSP_PCApplyBAorABTranspose(KSP ksp,Vec x,Vec y,Vec w)
268: {
271:   if (!ksp->transpose_solve) {
272:     PCApplyBAorABTranspose(ksp->pc,ksp->pc_side,x,y,w);
273:     KSP_RemoveNullSpace(ksp,y);
274:   } else {
275:     PCApplyBAorAB(ksp->pc,ksp->pc_side,x,y,w);
276:   }
277:   return(0);
278: }

280: PETSC_EXTERN PetscLogEvent KSP_GMRESOrthogonalization, KSP_SetUp, KSP_Solve;

282: PETSC_INTERN PetscErrorCode MatGetSchurComplement_Basic(Mat,IS,IS,IS,IS,MatReuse,Mat*,MatReuse,Mat*);

284: #endif