Actual source code: preonly.c
petsc-3.5.4 2015-05-23
2: #include <petsc-private/kspimpl.h>
6: static PetscErrorCode KSPSetUp_PREONLY(KSP ksp)
7: {
9: return(0);
10: }
14: static PetscErrorCode KSPSolve_PREONLY(KSP ksp)
15: {
17: PetscBool diagonalscale;
20: PCGetDiagonalScale(ksp->pc,&diagonalscale);
21: if (diagonalscale) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"Krylov method %s does not support diagonal scaling",((PetscObject)ksp)->type_name);
22: if (!ksp->guess_zero) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_USER,"Running KSP of preonly doesn't make sense with nonzero initial guess\n\
23: you probably want a KSP type of Richardson");
24: ksp->its = 0;
25: PCSetInitialGuessNonzero(ksp->pc,(PetscBool) !(int)ksp->guess_zero);
26: KSP_PCApply(ksp,ksp->vec_rhs,ksp->vec_sol);
27: ksp->its = 1;
28: ksp->reason = KSP_CONVERGED_ITS;
29: return(0);
30: }
32: /*MC
33: KSPPREONLY - This implements a stub method that applies ONLY the preconditioner.
34: This may be used in inner iterations, where it is desired to
35: allow multiple iterations as well as the "0-iteration" case. It is
36: commonly used with the direct solver preconditioners like PCLU and PCCHOLESKY
38: Options Database Keys:
39: . see KSPSolve()
41: Level: beginner
43: .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP
45: M*/
49: PETSC_EXTERN PetscErrorCode KSPCreate_PREONLY(KSP ksp)
50: {
54: KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_LEFT,3); /* LEFT/RIGHT is arbitrary, so "support" both */
55: KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_RIGHT,2);
57: ksp->data = (void*)0;
58: ksp->ops->setup = KSPSetUp_PREONLY;
59: ksp->ops->solve = KSPSolve_PREONLY;
60: ksp->ops->destroy = KSPDestroyDefault;
61: ksp->ops->buildsolution = KSPBuildSolutionDefault;
62: ksp->ops->buildresidual = KSPBuildResidualDefault;
63: ksp->ops->setfromoptions = 0;
64: ksp->ops->view = 0;
65: return(0);
66: }