xref: /petsc/src/ksp/ksp/impls/preonly/preonly.c (revision 76d6960897ba55d8238485170da43545084300c6)
1 #include <petsc/private/kspimpl.h>
2 
KSPSetUp_PREONLY(KSP ksp)3 static PetscErrorCode KSPSetUp_PREONLY(KSP ksp)
4 {
5   PetscFunctionBegin;
6   PetscFunctionReturn(PETSC_SUCCESS);
7 }
8 
KSPSolve_PREONLY(KSP ksp)9 static PetscErrorCode KSPSolve_PREONLY(KSP ksp)
10 {
11   PetscReal      norm;
12   PetscBool      flg;
13   PCFailedReason pcreason;
14 
15   PetscFunctionBegin;
16   PetscCall(PCGetDiagonalScale(ksp->pc, &flg));
17   PetscCheck(!flg, PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "Krylov method %s does not support diagonal scaling", ((PetscObject)ksp)->type_name);
18   if (!ksp->guess_zero) {
19     PetscCall(PetscObjectTypeCompareAny((PetscObject)ksp->pc, &flg, PCREDISTRIBUTE, PCMPI, ""));
20     PetscCheck(flg, PetscObjectComm((PetscObject)ksp), PETSC_ERR_USER, "KSP of type preonly (application of preconditioner only) doesn't make sense with nonzero initial guess you probably want a KSP of type Richardson");
21   }
22   ksp->its = 0;
23   if (ksp->numbermonitors) {
24     PetscCall(VecNorm(ksp->vec_rhs, NORM_2, &norm));
25     PetscCall(KSPMonitor(ksp, 0, norm));
26   }
27   PetscCall(KSP_PCApply(ksp, ksp->vec_rhs, ksp->vec_sol));
28   PetscCall(PCReduceFailedReason(ksp->pc));
29   PetscCall(PCGetFailedReason(ksp->pc, &pcreason));
30   PetscCall(VecFlag(ksp->vec_sol, pcreason));
31   if (pcreason) {
32     PetscCheck(!ksp->errorifnotconverged, PetscObjectComm((PetscObject)ksp), PETSC_ERR_NOT_CONVERGED, "KSPSolve has not converged with PCFailedReason %s", PCFailedReasons[pcreason]);
33     ksp->reason = KSP_DIVERGED_PC_FAILED;
34   } else {
35     ksp->its    = 1;
36     ksp->reason = KSP_CONVERGED_ITS;
37   }
38   if (ksp->numbermonitors) {
39     Vec v;
40     Mat A;
41 
42     PetscCall(VecDuplicate(ksp->vec_rhs, &v));
43     PetscCall(PCGetOperators(ksp->pc, &A, NULL));
44     PetscCall(KSP_MatMult(ksp, A, ksp->vec_sol, v));
45     PetscCall(VecAYPX(v, -1.0, ksp->vec_rhs));
46     PetscCall(VecNorm(v, NORM_2, &norm));
47     PetscCall(VecDestroy(&v));
48     PetscCall(KSPMonitor(ksp, 1, norm));
49   }
50   PetscFunctionReturn(PETSC_SUCCESS);
51 }
52 
KSPMatSolve_PREONLY(KSP ksp,Mat B,Mat X)53 static PetscErrorCode KSPMatSolve_PREONLY(KSP ksp, Mat B, Mat X)
54 {
55   PetscBool      diagonalscale;
56   PCFailedReason pcreason;
57 
58   PetscFunctionBegin;
59   PetscCall(PCGetDiagonalScale(ksp->pc, &diagonalscale));
60   PetscCheck(!diagonalscale, PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "Krylov method %s does not support diagonal scaling", ((PetscObject)ksp)->type_name);
61   PetscCheck(ksp->guess_zero, PetscObjectComm((PetscObject)ksp), PETSC_ERR_USER, "Running KSP of preonly doesn't make sense with nonzero initial guess you probably want a KSP type of Richardson");
62   ksp->its = 0;
63   PetscCall(KSP_PCMatApply(ksp, B, X));
64   PetscCall(PCGetFailedReason(ksp->pc, &pcreason));
65   /* 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 */
66   if (pcreason) {
67     PetscCall(MatSetInf(X));
68     ksp->reason = KSP_DIVERGED_PC_FAILED;
69   } else {
70     ksp->its    = 1;
71     ksp->reason = KSP_CONVERGED_ITS;
72   }
73   PetscFunctionReturn(PETSC_SUCCESS);
74 }
75 
76 /*MC
77    KSPNONE - An alias for `KSPPREONLY`
78 
79    Options Database Key:
80 .   -ksp_type none - use a single application of the preconditioner only
81 
82    Level: beginner
83 
84    Note:
85    See `KSPPREONLY` for more details
86 
87 .seealso: [](ch_ksp), `KSPCreate()`, `KSPSetType()`, `KSPType`, `KSPPREONLY`, `KSP`, `KSPRICHARDSON`, `KSPCHEBYSHEV`, `KSPGetPC()`, `KSPSetInitialGuessNonzero()`,
88           `PCREDISTRIBUTE`, `PCRedistributeGetKSP()`, `KSPPREONLY`
89 M*/
90 
91 /*MC
92    KSPPREONLY - This implements a method that applies ONLY the preconditioner exactly once.
93 
94    It is commonly used with the direct solver preconditioners like `PCLU` and `PCCHOLESKY`, but it may also be used when a single iteration of the
95    preconditioner is needed for smoothing in multigrid, `PCMG` or `PCGAMG` or within some other nested linear solve such as `PCFIELDSPLIT` or `PCBJACOBI`.
96 
97    There is an alias of this with the name `KSPNONE`.
98 
99    Options Database Key:
100 .   -ksp_type preonly - use a single application of the preconditioner only
101 
102    Level: beginner
103 
104    Notes:
105    Since this does not involve an iteration the basic `KSP` parameters such as tolerances and maximum iteration counts
106    do not apply
107 
108    To apply the preconditioner multiple times in a simple iteration use `KSPRICHARDSON`
109 
110    This `KSPType` cannot be used with the flag `-ksp_initial_guess_nonzero` or the call `KSPSetInitialGuessNonzero()` since it simply applies
111    the preconditioner to the given right-hand side during `KSPSolve()`. Except when the
112    `PCType` is `PCREDISTRIBUTE`; in that situation pass the nonzero initial guess flag with `-ksp_initial_guess_nonzero` or `KSPSetInitialGuessNonzero()`
113    both to the outer `KSP` (which is `KSPPREONLY`) and the inner `KSP` object obtained with `KSPGetPC()` followed by `PCRedistributedGetKSP()`
114    followed by `KSPSetInitialGuessNonzero()` or the option  `-redistribute_ksp_initial_guess_nonzero`.
115 
116    Developer Note:
117    Even though this method does not use any norms, the user is allowed to set the `KSPNormType` to any value.
118    This is so the users does not have to change `KSPNormType` options when they switch from other `KSP` methods to this one.
119 
120 .seealso: [](ch_ksp), `KSPCreate()`, `KSPSetType()`, `KSPType`, `KSP`, `KSPRICHARDSON`, `KSPCHEBYSHEV`, `KSPGetPC()`, `KSPSetInitialGuessNonzero()`,
121           `PCREDISTRIBUTE`, `PCRedistributeGetKSP()`, `KSPNONE`
122 M*/
123 
KSPCreate_PREONLY(KSP ksp)124 PETSC_EXTERN PetscErrorCode KSPCreate_PREONLY(KSP ksp)
125 {
126   PetscFunctionBegin;
127   PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_NONE, PC_LEFT, 3));
128   PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_NONE, PC_RIGHT, 2));
129   PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_PRECONDITIONED, PC_LEFT, 2));
130   PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_PRECONDITIONED, PC_RIGHT, 2));
131   PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_UNPRECONDITIONED, PC_LEFT, 2));
132   PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_UNPRECONDITIONED, PC_RIGHT, 2));
133   PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_NATURAL, PC_LEFT, 2));
134 
135   ksp->data                = NULL;
136   ksp->ops->setup          = KSPSetUp_PREONLY;
137   ksp->ops->solve          = KSPSolve_PREONLY;
138   ksp->ops->matsolve       = KSPMatSolve_PREONLY;
139   ksp->ops->destroy        = KSPDestroyDefault;
140   ksp->ops->buildsolution  = KSPBuildSolutionDefault;
141   ksp->ops->buildresidual  = KSPBuildResidualDefault;
142   ksp->ops->setfromoptions = NULL;
143   ksp->ops->view           = NULL;
144   ksp->guess_not_read      = PETSC_TRUE; // A KSP of preonly never needs to zero the input x since PC do not use an initial guess
145   PetscFunctionReturn(PETSC_SUCCESS);
146 }
147