Actual source code: preonly.c
2: #include <petsc/private/kspimpl.h>
4: static PetscErrorCode KSPSetUp_PREONLY(KSP ksp)
5: {
6: return 0;
7: }
9: static PetscErrorCode KSPSolve_PREONLY(KSP ksp)
10: {
11: PetscBool diagonalscale;
12: PCFailedReason pcreason;
14: PCGetDiagonalScale(ksp->pc,&diagonalscale);
17: you probably want a KSP type of Richardson");
18: ksp->its = 0;
19: KSP_PCApply(ksp,ksp->vec_rhs,ksp->vec_sol);
20: PCGetFailedReasonRank(ksp->pc,&pcreason);
21: /* Note: only some ranks may have this set; this may lead to problems if the caller assumes ksp->reason is set on all processes or just uses the result */
22: if (pcreason) {
23: VecSetInf(ksp->vec_sol);
24: ksp->reason = KSP_DIVERGED_PC_FAILED;
25: } else {
26: ksp->its = 1;
27: ksp->reason = KSP_CONVERGED_ITS;
28: }
29: if (ksp->numbermonitors) {
30: Vec v;
31: PetscReal norm;
32: Mat A;
34: VecNorm(ksp->vec_rhs,NORM_2,&norm);
35: KSPMonitor(ksp,0,norm);
36: VecDuplicate(ksp->vec_rhs,&v);
37: PCGetOperators(ksp->pc,&A,NULL);
38: KSP_MatMult(ksp,A,ksp->vec_sol,v);
39: VecAYPX(v,-1.0,ksp->vec_rhs);
40: VecNorm(v,NORM_2,&norm);
41: VecDestroy(&v);
42: KSPMonitor(ksp,1,norm);
43: }
44: return 0;
45: }
47: static PetscErrorCode KSPMatSolve_PREONLY(KSP ksp, Mat B, Mat X)
48: {
49: PetscBool diagonalscale;
50: PCFailedReason pcreason;
52: PCGetDiagonalScale(ksp->pc,&diagonalscale);
55: you probably want a KSP type of Richardson");
56: ksp->its = 0;
57: PCMatApply(ksp->pc,B,X);
58: PCGetFailedReason(ksp->pc,&pcreason);
59: /* Note: only some ranks may have this set; this may lead to problems if the caller assumes ksp->reason is set on all processes or just uses the result */
60: if (pcreason) {
61: MatSetInf(X);
62: ksp->reason = KSP_DIVERGED_PC_FAILED;
63: } else {
64: ksp->its = 1;
65: ksp->reason = KSP_CONVERGED_ITS;
66: }
67: return 0;
68: }
70: /*MC
71: KSPPREONLY - This implements a method that applies ONLY the preconditioner exactly once.
72: This may be used in inner iterations, where it is desired to
73: allow multiple iterations as well as the "0-iteration" case. It is
74: commonly used with the direct solver preconditioners like PCLU and PCCHOLESKY
76: Options Database Keys:
77: . -ksp_type preonly - use preconditioner only
79: Level: beginner
81: Notes:
82: Since this does not involve an iteration the basic KSP parameters such as tolerances and iteration counts
83: do not apply
85: To apply multiple preconditioners in a simple iteration use KSPRICHARDSON
87: Developer Notes:
88: Even though this method does not use any norms, the user is allowed to set the KSPNormType to any value.
89: This is so the users does not have to change KSPNormType options when they switch from other KSP methods to this one.
91: .seealso: KSPCreate(), KSPSetType(), KSPType, KSP, KSPRICHARDSON, KSPCHEBYSHEV
93: M*/
95: PETSC_EXTERN PetscErrorCode KSPCreate_PREONLY(KSP ksp)
96: {
97: KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_LEFT,3);
98: KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_RIGHT,2);
99: KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,2);
100: KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_RIGHT,2);
101: KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_LEFT,2);
102: KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_RIGHT,2);
103: KSPSetSupportedNorm(ksp,KSP_NORM_NATURAL,PC_LEFT,2);
105: ksp->data = NULL;
106: ksp->ops->setup = KSPSetUp_PREONLY;
107: ksp->ops->solve = KSPSolve_PREONLY;
108: ksp->ops->matsolve = KSPMatSolve_PREONLY;
109: ksp->ops->destroy = KSPDestroyDefault;
110: ksp->ops->buildsolution = KSPBuildSolutionDefault;
111: ksp->ops->buildresidual = KSPBuildResidualDefault;
112: ksp->ops->setfromoptions = NULL;
113: ksp->ops->view = NULL;
114: return 0;
115: }