Actual source code: kspimpl.h
1: /* $Id: kspimpl.h,v 1.50 2001/08/06 18:04:29 bsmith Exp $ */
3: #ifndef _KSPIMPL
4: #define _KSPIMPL
6: #include petscksp.h
8: typedef struct _KSPOps *KSPOps;
10: struct _KSPOps {
11: int (*buildsolution)(KSP,Vec,Vec*); /* Returns a pointer to the solution, or
12: calculates the solution in a
13: user-provided area. */
14: int (*buildresidual)(KSP,Vec,Vec,Vec*); /* Returns a pointer to the residual, or
15: calculates the residual in a
16: user-provided area. */
17: int (*solve)(KSP,int*); /* actual solver */
18: int (*setup)(KSP);
19: int (*setfromoptions)(KSP);
20: int (*publishoptions)(KSP);
21: int (*computeextremesingularvalues)(KSP,PetscReal*,PetscReal*);
22: int (*computeeigenvalues)(KSP,int,PetscReal*,PetscReal*,int *);
23: int (*destroy)(KSP);
24: int (*view)(KSP,PetscViewer);
25: };
27: /*
28: Maximum number of monitors you can run with a single KSP
29: */
30: #define MAXKSPMONITORS 5
32: /*
33: Defines the KSP data structure.
34: */
35: struct _p_KSP {
36: PETSCHEADER(struct _KSPOps)
37: /*------------------------- User parameters--------------------------*/
38: int max_it; /* maximum number of iterations */
39: PetscTruth guess_zero, /* flag for whether initial guess is 0 */
40: calc_sings; /* calculate extreme Singular Values */
41: PCSide pc_side; /* flag for left, right, or symmetric
42: preconditioning */
43: PetscReal rtol, /* relative tolerance */
44: atol, /* absolute tolerance */
45: ttol, /* (not set by user) */
46: divtol; /* divergence tolerance */
47: PetscReal rnorm0; /* initial residual norm (used for divergence testing) */
48: PetscReal rnorm; /* current residual norm */
49: KSPConvergedReason reason;
51: Vec vec_sol,vec_rhs; /* pointer to where user has stashed
52: the solution and rhs, these are
53: never touched by the code, only
54: passed back to the user */
55: PetscReal *res_hist; /* If !0 stores residual at iterations*/
56: int res_hist_len; /* current size of residual history array */
57: int res_hist_max; /* actual amount of data in residual_history */
58: PetscTruth res_hist_reset; /* reset history to size zero for each new solve */
60: /* --------User (or default) routines (most return -1 on error) --------*/
61: int (*monitor[MAXKSPMONITORS])(KSP,int,PetscReal,void*); /* returns control to user after */
62: int (*monitordestroy[MAXKSPMONITORS])(void*); /* */
63: void *monitorcontext[MAXKSPMONITORS]; /* residual calculation, allows user */
64: int numbermonitors; /* to, for instance, print residual norm, etc. */
66: int (*converged)(KSP,int,PetscReal,KSPConvergedReason*,void*);
67: void *cnvP;
69: PC B;
71: void *data; /* holder for misc stuff associated
72: with a particular iterative solver */
74: /* ----------------Default work-area management -------------------- */
75: int nwork;
76: Vec *work;
78: int setupcalled;
80: int its; /* number of iterations so far computed */
82: PetscTruth transpose_solve; /* solve transpose system instead */
84: KSPNormType normtype; /* type of norm used for convergence tests */
85: };
87: #define KSPLogResidualHistory(ksp,norm)
88: {if (ksp->res_hist && ksp->res_hist_max > ksp->res_hist_len)
89: ksp->res_hist[ksp->res_hist_len++] = norm;}
91: #define KSPMonitor(ksp,it,rnorm)
92: { int _ierr,_i,_im = ksp->numbermonitors;
93: for (_i=0; _i<_im; _i++) {
94: _(*ksp->monitor[_i])(ksp,it,rnorm,ksp->monitorcontext[_i]);CHKERRQ(_ierr);
95: }
96: }
98: EXTERN int KSPDefaultBuildSolution(KSP,Vec,Vec*);
99: EXTERN int KSPDefaultBuildResidual(KSP,Vec,Vec,Vec *);
100: EXTERN int KSPDefaultDestroy(KSP);
101: EXTERN int KSPDefaultGetWork(KSP,int);
102: EXTERN int KSPDefaultFreeWork(KSP);
103: EXTERN int KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec,Vec);
104: EXTERN int KSPUnwindPreconditioner(KSP,Vec,Vec);
106: /*
107: These allow the various Krylov methods to apply to either the linear system
108: or its transpose.
109: */
110: #define KSP_MatMult(ksp,A,x,y) (!ksp->transpose_solve) ? MatMult(A,x,y) : MatMultTranspose(A,x,y)
111: #define KSP_MatMultTranspose(ksp,A,x,y) (!ksp->transpose_solve) ? MatMultTranspose(A,x,y) : MatMult(A,x,y)
112: #define KSP_PCApply(ksp,B,x,y) (!ksp->transpose_solve) ? PCApply(B,x,y) : PCApplyTranspose(B,x,y)
113: #define KSP_PCApplyTranspose(ksp,B,x,y) (!ksp->transpose_solve) ? PCApplyTranspose(B,x,y) : PCApply(B,x,y)
114: #define KSP_PCApplyBAorAB(ksp,pc,side,x,y,w) (!ksp->transpose_solve) ? PCApplyBAorAB(pc,side,x,y,w) : PCApplyBAorABTranspose(pc,side,x,y,w)
116: #endif