Actual source code: preonly.c
2: #include src/ksp/ksp/kspimpl.h
6: static PetscErrorCode KSPSetUp_PREONLY(KSP ksp)
7: {
9: return(0);
10: }
14: static PetscErrorCode KSPSolve_PREONLY(KSP ksp)
15: {
17: Vec X,B;
18: PetscTruth diagonalscale;
21: PCDiagonalScale(ksp->pc,&diagonalscale);
22: if (diagonalscale) SETERRQ1(PETSC_ERR_SUP,"Krylov method %s does not support diagonal scaling",ksp->type_name);
23: if (!ksp->guess_zero) {
24: SETERRQ(PETSC_ERR_USER,"Running KSP of preonly doesn't make sense with nonzero initial guess\n\
25: you probably want a KSP type of Richardson");
26: }
27: ksp->its = 0;
28: X = ksp->vec_sol;
29: B = ksp->vec_rhs;
30: KSP_PCApply(ksp,B,X);
31: ksp->its = 1;
32: ksp->reason = KSP_CONVERGED_ITS;
33: return(0);
34: }
36: /*MC
37: KSPPREONLY - This implements a stub method that applies ONLY the preconditioner.
38: This may be used in inner iterations, where it is desired to
39: allow multiple iterations as well as the "0-iteration" case. It is
40: commonly used with the direct solver preconditioners like PCLU and PCCHOLESKY
42: Options Database Keys:
43: . see KSPSolve()
45: Level: beginner
47: .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP
49: M*/
54: PetscErrorCode KSPCreate_PREONLY(KSP ksp)
55: {
57: ksp->data = (void*)0;
58: ksp->ops->setup = KSPSetUp_PREONLY;
59: ksp->ops->solve = KSPSolve_PREONLY;
60: ksp->ops->destroy = KSPDefaultDestroy;
61: ksp->ops->buildsolution = KSPDefaultBuildSolution;
62: ksp->ops->buildresidual = KSPDefaultBuildResidual;
63: ksp->ops->setfromoptions = 0;
64: ksp->ops->view = 0;
65: return(0);
66: }