xref: /petsc/include/petsc/private/kspimpl.h (revision fbdf3bed89ca4550907b2a3c2c22a34a69e8be89)
1 
2 #ifndef _KSPIMPL_H
3 #define _KSPIMPL_H
4 
5 #include <petscksp.h>
6 #include <petsc/private/petscimpl.h>
7 
8 PETSC_EXTERN PetscBool KSPRegisterAllCalled;
9 PETSC_EXTERN PetscErrorCode KSPRegisterAll(void);
10 PETSC_EXTERN PetscErrorCode KSPMatRegisterAll(void);
11 
12 typedef struct _KSPOps *KSPOps;
13 
14 struct _KSPOps {
15   PetscErrorCode (*buildsolution)(KSP,Vec,Vec*);       /* Returns a pointer to the solution, or
16                                                           calculates the solution in a
17                                                           user-provided area. */
18   PetscErrorCode (*buildresidual)(KSP,Vec,Vec,Vec*);   /* Returns a pointer to the residual, or
19                                                           calculates the residual in a
20                                                           user-provided area.  */
21   PetscErrorCode (*solve)(KSP);                        /* actual solver */
22   PetscErrorCode (*setup)(KSP);
23   PetscErrorCode (*setfromoptions)(PetscOptions*,KSP);
24   PetscErrorCode (*publishoptions)(KSP);
25   PetscErrorCode (*computeextremesingularvalues)(KSP,PetscReal*,PetscReal*);
26   PetscErrorCode (*computeeigenvalues)(KSP,PetscInt,PetscReal*,PetscReal*,PetscInt *);
27   PetscErrorCode (*destroy)(KSP);
28   PetscErrorCode (*view)(KSP,PetscViewer);
29   PetscErrorCode (*reset)(KSP);
30   PetscErrorCode (*load)(KSP,PetscViewer);
31 };
32 
33 typedef struct {PetscInt model,curl,maxl;Mat mat; KSP ksp;}* KSPGuessFischer;
34 
35 /*
36      Maximum number of monitors you can run with a single KSP
37 */
38 #define MAXKSPMONITORS 5
39 typedef enum {KSP_SETUP_NEW, KSP_SETUP_NEWMATRIX, KSP_SETUP_NEWRHS} KSPSetUpStage;
40 
41 /*
42    Defines the KSP data structure.
43 */
44 struct _p_KSP {
45   PETSCHEADER(struct _KSPOps);
46   DM              dm;
47   PetscBool       dmAuto;       /* DM was created automatically by KSP */
48   PetscBool       dmActive;     /* KSP should use DM for computing operators */
49   /*------------------------- User parameters--------------------------*/
50   PetscInt        max_it;                     /* maximum number of iterations */
51   KSPFischerGuess guess;
52   PetscBool       guess_zero,                  /* flag for whether initial guess is 0 */
53                   calc_sings,                  /* calculate extreme Singular Values */
54                   guess_knoll;                /* use initial guess of PCApply(ksp->B,b */
55   PCSide          pc_side;                  /* flag for left, right, or symmetric preconditioning */
56   PetscInt        normsupporttable[KSP_NORM_MAX][PC_SIDE_MAX]; /* Table of supported norms and pc_side, see KSPSetSupportedNorm() */
57   PetscReal       rtol,                     /* relative tolerance */
58                   abstol,                     /* absolute tolerance */
59                   ttol,                     /* (not set by user)  */
60                   divtol;                   /* divergence tolerance */
61   PetscReal       rnorm0;                   /* initial residual norm (used for divergence testing) */
62   PetscReal       rnorm;                    /* current residual norm */
63   KSPConvergedReason reason;
64   PetscBool          errorifnotconverged;    /* create an error if the KSPSolve() does not converge */
65 
66   Vec vec_sol,vec_rhs;            /* pointer to where user has stashed
67                                       the solution and rhs, these are
68                                       never touched by the code, only
69                                       passed back to the user */
70   PetscReal     *res_hist;            /* If !0 stores residual at iterations*/
71   PetscReal     *res_hist_alloc;      /* If !0 means user did not provide buffer, needs deallocation */
72   PetscInt      res_hist_len;         /* current size of residual history array */
73   PetscInt      res_hist_max;         /* actual amount of data in residual_history */
74   PetscBool     res_hist_reset;       /* reset history to size zero for each new solve */
75 
76   PetscInt      chknorm;             /* only compute/check norm if iterations is great than this */
77   PetscBool     lagnorm;             /* Lag the residual norm calculation so that it is computed as part of the
78                                         MPI_Allreduce() for computing the inner products for the next iteration. */
79   /* --------User (or default) routines (most return -1 on error) --------*/
80   PetscErrorCode (*monitor[MAXKSPMONITORS])(KSP,PetscInt,PetscReal,void*); /* returns control to user after */
81   PetscErrorCode (*monitordestroy[MAXKSPMONITORS])(void**);         /* */
82   void *monitorcontext[MAXKSPMONITORS];                  /* residual calculation, allows user */
83   PetscInt  numbermonitors;                                   /* to, for instance, print residual norm, etc. */
84 
85   PetscErrorCode (*converged)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*);
86   PetscErrorCode (*convergeddestroy)(void*);
87   void       *cnvP;
88 
89   void       *user;             /* optional user-defined context */
90 
91   PC         pc;
92 
93   void       *data;                      /* holder for misc stuff associated
94                                    with a particular iterative solver */
95 
96   /* ----------------Default work-area management -------------------- */
97   PetscInt       nwork;
98   Vec            *work;
99 
100   KSPSetUpStage  setupstage;
101 
102   PetscInt       its;       /* number of iterations so far computed in THIS linear solve*/
103   PetscInt       totalits;   /* number of iterations used by this KSP object since it was created */
104 
105   PetscBool      transpose_solve;    /* solve transpose system instead */
106 
107   KSPNormType    normtype;          /* type of norm used for convergence tests */
108 
109   PCSide         pc_side_set;   /* PC type set explicitly by user */
110   KSPNormType    normtype_set;  /* Norm type set explicitly by user */
111 
112   /*   Allow diagonally scaling the matrix before computing the preconditioner or using
113        the Krylov method. Note this is NOT just Jacobi preconditioning */
114 
115   PetscBool    dscale;       /* diagonal scale system; used with KSPSetDiagonalScale() */
116   PetscBool    dscalefix;    /* unscale system after solve */
117   PetscBool    dscalefix2;   /* system has been unscaled */
118   Vec          diagonal;     /* 1/sqrt(diag of matrix) */
119   Vec          truediagonal;
120 
121   MatNullSpace nullsp;      /* Null space of the operator, removed from Krylov space */
122 
123   PetscBool    skippcsetfromoptions; /* if set then KSPSetFromOptions() does not call PCSetFromOptions() */
124 
125   PetscViewer  eigviewer;   /* Viewer where computed eigenvalues are displayed */
126 
127   PetscErrorCode (*presolve)(KSP,Vec,Vec,void*);
128   PetscErrorCode (*postsolve)(KSP,Vec,Vec,void*);
129   void           *prectx,*postctx;
130 };
131 
132 typedef struct { /* dummy data structure used in KSPMonitorDynamicTolerance() */
133   PetscReal coef;
134   PetscReal bnrm;
135 } KSPDynTolCtx;
136 
137 typedef struct {
138   PetscBool  initialrtol;    /* default relative residual decrease is computing from initial residual, not rhs */
139   PetscBool  mininitialrtol; /* default relative residual decrease is computing from min of initial residual and rhs */
140   Vec        work;
141 } KSPConvergedDefaultCtx;
142 
143 #undef __FUNCT__
144 #define __FUNCT__ "KSPLogResidualHistory"
145 PETSC_STATIC_INLINE PetscErrorCode KSPLogResidualHistory(KSP ksp,PetscReal norm)
146 {
147   PetscErrorCode ierr;
148 
149   PetscFunctionBegin;
150   ierr = PetscObjectSAWsTakeAccess((PetscObject)ksp);CHKERRQ(ierr);
151   if (ksp->res_hist && ksp->res_hist_max > ksp->res_hist_len) {
152     ksp->res_hist[ksp->res_hist_len++] = norm;
153   }
154   ierr = PetscObjectSAWsGrantAccess((PetscObject)ksp);CHKERRQ(ierr);
155   PetscFunctionReturn(0);
156 }
157 
158 PETSC_INTERN PetscErrorCode KSPSetUpNorms_Private(KSP,KSPNormType*,PCSide*);
159 
160 PETSC_INTERN PetscErrorCode KSPPlotEigenContours_Private(KSP,PetscInt,const PetscReal*,const PetscReal*);
161 
162 typedef struct _p_DMKSP *DMKSP;
163 typedef struct _DMKSPOps *DMKSPOps;
164 struct _DMKSPOps {
165   PetscErrorCode (*computeoperators)(KSP,Mat,Mat,void*);
166   PetscErrorCode (*computerhs)(KSP,Vec,void*);
167   PetscErrorCode (*computeinitialguess)(KSP,Vec,void*);
168   PetscErrorCode (*destroy)(DMKSP*);
169   PetscErrorCode (*duplicate)(DMKSP,DMKSP);
170 };
171 
172 struct _p_DMKSP {
173   PETSCHEADER(struct _DMKSPOps);
174   void *operatorsctx;
175   void *rhsctx;
176   void *initialguessctx;
177   void *data;
178 
179   /* This is NOT reference counted. The DM on which this context was first created is cached here to implement one-way
180    * copy-on-write. When DMGetDMKSPWrite() sees a request using a different DM, it makes a copy. Thus, if a user
181    * only interacts directly with one level, e.g., using KSPSetComputeOperators(), then coarse levels are constructed by
182    * PCMG, then the user changes the routine with another call to KSPSetComputeOperators(), it automatically propagates
183    * to all the levels. If instead, they get out a specific level and set the routines on that level, subsequent changes
184    * to the original level will no longer propagate to that level.
185    */
186   DM originaldm;
187 
188   void (*fortran_func_pointers[3])(void); /* Store our own function pointers so they are associated with the DMKSP instead of the DM */
189 };
190 PETSC_EXTERN PetscErrorCode DMGetDMKSP(DM,DMKSP*);
191 PETSC_EXTERN PetscErrorCode DMGetDMKSPWrite(DM,DMKSP*);
192 PETSC_EXTERN PetscErrorCode DMCopyDMKSP(DM,DM);
193 
194 /*
195        These allow the various Krylov methods to apply to either the linear system or its transpose.
196 */
197 #undef __FUNCT__
198 #define __FUNCT__ "KSP_RemoveNullSpace"
199 PETSC_STATIC_INLINE PetscErrorCode KSP_RemoveNullSpace(KSP ksp,Vec y)
200 {
201   PetscErrorCode ierr;
202   PetscFunctionBegin;
203   if (ksp->nullsp && ksp->pc_side == PC_LEFT) {ierr = MatNullSpaceRemove(ksp->nullsp,y);CHKERRQ(ierr);}
204   PetscFunctionReturn(0);
205 }
206 
207 #undef __FUNCT__
208 #define __FUNCT__ "KSP_MatMult"
209 PETSC_STATIC_INLINE PetscErrorCode KSP_MatMult(KSP ksp,Mat A,Vec x,Vec y)
210 {
211   PetscErrorCode ierr;
212   PetscFunctionBegin;
213   if (!ksp->transpose_solve) {ierr = MatMult(A,x,y);CHKERRQ(ierr);}
214   else                       {ierr = MatMultTranspose(A,x,y);CHKERRQ(ierr);}
215   PetscFunctionReturn(0);
216 }
217 
218 #undef __FUNCT__
219 #define __FUNCT__ "KSP_MatMultTranspose"
220 PETSC_STATIC_INLINE PetscErrorCode KSP_MatMultTranspose(KSP ksp,Mat A,Vec x,Vec y)
221 {
222   PetscErrorCode ierr;
223   PetscFunctionBegin;
224   if (!ksp->transpose_solve) {ierr = MatMultTranspose(A,x,y);CHKERRQ(ierr);}
225   else                       {ierr = MatMult(A,x,y);CHKERRQ(ierr);}
226   PetscFunctionReturn(0);
227 }
228 
229 #undef __FUNCT__
230 #define __FUNCT__ "KSP_PCApply"
231 PETSC_STATIC_INLINE PetscErrorCode KSP_PCApply(KSP ksp,Vec x,Vec y)
232 {
233   PetscErrorCode ierr;
234   PetscFunctionBegin;
235   if (!ksp->transpose_solve) {
236     ierr = PCApply(ksp->pc,x,y);CHKERRQ(ierr);
237     ierr = KSP_RemoveNullSpace(ksp,y);CHKERRQ(ierr);
238   } else {
239     ierr = PCApplyTranspose(ksp->pc,x,y);CHKERRQ(ierr);
240   }
241   PetscFunctionReturn(0);
242 }
243 
244 #undef __FUNCT__
245 #define __FUNCT__ "KSP_PCApplyTranspose"
246 PETSC_STATIC_INLINE PetscErrorCode KSP_PCApplyTranspose(KSP ksp,Vec x,Vec y)
247 {
248   PetscErrorCode ierr;
249   PetscFunctionBegin;
250   if (!ksp->transpose_solve) {
251     ierr = PCApplyTranspose(ksp->pc,x,y);CHKERRQ(ierr);
252   } else {
253     ierr = PCApply(ksp->pc,x,y);CHKERRQ(ierr);
254     ierr = KSP_RemoveNullSpace(ksp,y);CHKERRQ(ierr);
255   }
256   PetscFunctionReturn(0);
257 }
258 
259 #undef __FUNCT__
260 #define __FUNCT__ "KSP_PCApplyBAorAB"
261 PETSC_STATIC_INLINE PetscErrorCode KSP_PCApplyBAorAB(KSP ksp,Vec x,Vec y,Vec w)
262 {
263   PetscErrorCode ierr;
264   PetscFunctionBegin;
265   if (!ksp->transpose_solve) {
266     ierr = PCApplyBAorAB(ksp->pc,ksp->pc_side,x,y,w);CHKERRQ(ierr);
267     ierr = KSP_RemoveNullSpace(ksp,y);CHKERRQ(ierr);
268   } else {
269     ierr = PCApplyBAorABTranspose(ksp->pc,ksp->pc_side,x,y,w);CHKERRQ(ierr);
270   }
271   PetscFunctionReturn(0);
272 }
273 
274 #undef __FUNCT__
275 #define __FUNCT__ "KSP_PCApplyBAorABTranspose"
276 PETSC_STATIC_INLINE PetscErrorCode KSP_PCApplyBAorABTranspose(KSP ksp,Vec x,Vec y,Vec w)
277 {
278   PetscErrorCode ierr;
279   PetscFunctionBegin;
280   if (!ksp->transpose_solve) {
281     ierr = PCApplyBAorABTranspose(ksp->pc,ksp->pc_side,x,y,w);CHKERRQ(ierr);
282     ierr = KSP_RemoveNullSpace(ksp,y);CHKERRQ(ierr);
283   } else {
284     ierr = PCApplyBAorAB(ksp->pc,ksp->pc_side,x,y,w);CHKERRQ(ierr);
285   }
286   PetscFunctionReturn(0);
287 }
288 
289 PETSC_EXTERN PetscLogEvent KSP_GMRESOrthogonalization, KSP_SetUp, KSP_Solve;
290 
291 PETSC_INTERN PetscErrorCode MatGetSchurComplement_Basic(Mat,IS,IS,IS,IS,MatReuse,Mat*,MatSchurComplementAinvType,MatReuse,Mat*);
292 
293 #endif
294