Actual source code: preonly.c
petsc-3.8.4 2018-03-24
2: #include <petsc/private/kspimpl.h>
4: static PetscErrorCode KSPSetUp_PREONLY(KSP ksp)
5: {
7: return(0);
8: }
10: static PetscErrorCode KSPSolve_PREONLY(KSP ksp)
11: {
13: PetscBool diagonalscale;
14: PCFailedReason pcreason;
17: PCGetDiagonalScale(ksp->pc,&diagonalscale);
18: if (diagonalscale) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"Krylov method %s does not support diagonal scaling",((PetscObject)ksp)->type_name);
19: if (!ksp->guess_zero) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_USER,"Running KSP of preonly doesn't make sense with nonzero initial guess\n\
20: you probably want a KSP type of Richardson");
21: ksp->its = 0;
22: KSP_PCApply(ksp,ksp->vec_rhs,ksp->vec_sol);
23: PCGetSetUpFailedReason(ksp->pc,&pcreason);
24: if (pcreason) {
25: ksp->reason = KSP_DIVERGED_PCSETUP_FAILED;
26: } else {
27: ksp->its = 1;
28: ksp->reason = KSP_CONVERGED_ITS;
29: }
30: return(0);
31: }
33: /*MC
34: KSPPREONLY - This implements a stub method that applies ONLY the preconditioner.
35: This may be used in inner iterations, where it is desired to
36: allow multiple iterations as well as the "0-iteration" case. It is
37: commonly used with the direct solver preconditioners like PCLU and PCCHOLESKY
39: Options Database Keys:
40: . -ksp_type preonly
42: Level: beginner
44: Notes: Since this does not involve an iteration the basic KSP parameters such as tolerances and iteration counts
45: do not apply
47: Developer Notes: Even though this method does not use any norms, the user is allowed to set the KSPNormType to any value.
48: This is so the users does not have to change KSPNormType options when they switch from other KSP methods to this one.
50: .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP
52: M*/
54: PETSC_EXTERN PetscErrorCode KSPCreate_PREONLY(KSP ksp)
55: {
59: KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_LEFT,3);
60: KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_RIGHT,2);
61: KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,2);
62: KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_RIGHT,2);
63: KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_LEFT,2);
64: KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_RIGHT,2);
65: KSPSetSupportedNorm(ksp,KSP_NORM_NATURAL,PC_LEFT,2);
67: ksp->data = NULL;
68: ksp->ops->setup = KSPSetUp_PREONLY;
69: ksp->ops->solve = KSPSolve_PREONLY;
70: ksp->ops->destroy = KSPDestroyDefault;
71: ksp->ops->buildsolution = KSPBuildSolutionDefault;
72: ksp->ops->buildresidual = KSPBuildResidualDefault;
73: ksp->ops->setfromoptions = 0;
74: ksp->ops->view = 0;
75: return(0);
76: }