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