Actual source code: preonly.c

petsc-3.7.5 2017-01-01
Report Typos and Errors
  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;
 18:   PCFailedReason pcreason;

 21:   PCGetDiagonalScale(ksp->pc,&diagonalscale);
 22:   if (diagonalscale) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"Krylov method %s does not support diagonal scaling",((PetscObject)ksp)->type_name);
 23:   if (!ksp->guess_zero) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_USER,"Running KSP of preonly doesn't make sense with nonzero initial guess\n\
 24:                you probably want a KSP type of Richardson");
 25:   ksp->its = 0;
 26:   KSP_PCApply(ksp,ksp->vec_rhs,ksp->vec_sol);
 27:   PCGetSetUpFailedReason(ksp->pc,&pcreason);
 28:   if (pcreason) {
 29:     ksp->reason = KSP_DIVERGED_PCSETUP_FAILED;
 30:   } else {
 31:     ksp->its    = 1;
 32:     ksp->reason = KSP_CONVERGED_ITS;
 33:   }
 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: .   -ksp_type preonly

 46:    Level: beginner

 48: .seealso:  KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP

 50: M*/

 54: PETSC_EXTERN PetscErrorCode KSPCreate_PREONLY(KSP ksp)
 55: {

 59:   KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_LEFT,3); /* LEFT/RIGHT is arbitrary, so "support" both */
 60:   KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_RIGHT,2);

 62:   ksp->data                = (void*)0;
 63:   ksp->ops->setup          = KSPSetUp_PREONLY;
 64:   ksp->ops->solve          = KSPSolve_PREONLY;
 65:   ksp->ops->destroy        = KSPDestroyDefault;
 66:   ksp->ops->buildsolution  = KSPBuildSolutionDefault;
 67:   ksp->ops->buildresidual  = KSPBuildResidualDefault;
 68:   ksp->ops->setfromoptions = 0;
 69:   ksp->ops->view           = 0;
 70:   return(0);
 71: }