Actual source code: modpcf.c
petsc-3.7.5 2017-01-01
2: #include <petsc/private/kspimpl.h> /*I "petscksp.h" I*/
6: /*@C
7: KSPFGMRESSetModifyPC - Sets the routine used by FGMRES to modify the preconditioner.
9: Logically Collective on KSP
11: Input Parameters:
12: + ksp - iterative context obtained from KSPCreate
13: . fcn - modifypc function
14: . ctx - optional contex
15: - d - optional context destroy routine
17: Calling Sequence of function:
18: int fcn(KSP ksp,int total_its,int loc_its,PetscReal res_norm,void*ctx);
20: ksp - the ksp context being used.
21: total_its - the total number of FGMRES iterations that have occurred.
22: loc_its - the number of FGMRES iterations since last restart.
23: res_norm - the current residual norm.
24: ctx - optional context variable
26: Options Database Keys:
27: -ksp_fgmres_modifypcnochange
28: -ksp_fgmres_modifypcksp
30: Level: intermediate
32: Contributed by Allison Baker
34: Notes:
35: Several modifypc routines are predefined, including
36: KSPFGMRESModifyPCNoChange()
37: KSPFGMRESModifyPCKSP()
39: .seealso: KSPFGMRESModifyPCNoChange(), KSPFGMRESModifyPCKSP()
41: @*/
42: PetscErrorCode KSPFGMRESSetModifyPC(KSP ksp,PetscErrorCode (*fcn)(KSP,PetscInt,PetscInt,PetscReal,void*),void *ctx,PetscErrorCode (*d)(void*))
43: {
48: PetscTryMethod(ksp,"KSPFGMRESSetModifyPC_C",(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode (*)(void*)),(ksp,fcn,ctx,d));
49: return(0);
50: }
53: /* The following are different routines used to modify the preconditioner */
57: /*@
59: KSPFGMRESModifyPCNoChange - this is the default used by fgmres - it doesn't change the preconditioner.
61: Input Parameters:
62: + ksp - the ksp context being used.
63: . total_its - the total number of FGMRES iterations that have occurred.
64: . loc_its - the number of FGMRES iterations since last restart.
65: a restart (so number of Krylov directions to be computed)
66: . res_norm - the current residual norm.
67: - dummy - context variable, unused in this routine
69: Level: intermediate
71: Contributed by Allison Baker
73: You can use this as a template!
75: .seealso: KSPFGMRESSetModifyPC(), KSPFGMRESModifyPCKSP()
77: @*/
78: PetscErrorCode KSPFGMRESModifyPCNoChange(KSP ksp,PetscInt total_its,PetscInt loc_its,PetscReal res_norm,void *dummy)
79: {
81: return(0);
82: }
86: /*@
88: KSPFGMRESModifyPCKSP - modifies the attributes of the
89: GMRES preconditioner. It serves as an example (not as something
90: useful!)
92: Input Parameters:
93: + ksp - the ksp context being used.
94: . total_its - the total number of FGMRES iterations that have occurred.
95: . loc_its - the number of FGMRES iterations since last restart.
96: . res_norm - the current residual norm.
97: - dummy - context, not used here
99: Level: intermediate
101: Contributed by Allison Baker
103: This could be used as a template!
105: .seealso: KSPFGMRESSetModifyPC(), KSPFGMRESModifyPCKSP()
107: @*/
108: PetscErrorCode KSPFGMRESModifyPCKSP(KSP ksp,PetscInt total_its,PetscInt loc_its,PetscReal res_norm,void *dummy)
109: {
110: PC pc;
112: PetscInt maxits;
113: KSP sub_ksp;
114: PetscReal rtol,abstol,dtol;
115: PetscBool isksp;
118: KSPGetPC(ksp,&pc);
120: PetscObjectTypeCompare((PetscObject)pc,PCKSP,&isksp);
121: if (isksp) {
122: PCKSPGetKSP(pc,&sub_ksp);
124: /* note that at this point you could check the type of KSP with KSPGetType() */
126: /* Now we can use functions such as KSPGMRESSetRestart() or
127: KSPGMRESSetOrthogonalization() or KSPSetTolerances() */
129: KSPGetTolerances(sub_ksp,&rtol,&abstol,&dtol,&maxits);
130: if (!loc_its) rtol = .1;
131: else rtol *= .9;
132: KSPSetTolerances(sub_ksp,rtol,abstol,dtol,maxits);
133: }
134: return(0);
135: }