Actual source code: preonly.c

petsc-3.14.6 2021-03-30
Report Typos and Errors

  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:   PCGetFailedReasonRank(ksp->pc,&pcreason);
 24:   /* Note: only some ranks may have this set; this may lead to problems if the caller assumes ksp->reason is set on all processes or just uses the result */
 25:   if (pcreason) {
 26:     VecSetInf(ksp->vec_sol);
 27:     ksp->reason = KSP_DIVERGED_PC_FAILED;
 28:   } else {
 29:     ksp->its    = 1;
 30:     ksp->reason = KSP_CONVERGED_ITS;
 31:   }
 32:   return(0);
 33: }

 35: static PetscErrorCode KSPMatSolve_PREONLY(KSP ksp, Mat B, Mat X)
 36: {
 38:   PetscBool      diagonalscale;
 39:   PCFailedReason pcreason;

 42:   PCGetDiagonalScale(ksp->pc,&diagonalscale);
 43:   if (diagonalscale) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"Krylov method %s does not support diagonal scaling",((PetscObject)ksp)->type_name);
 44:   if (!ksp->guess_zero) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_USER,"Running KSP of preonly doesn't make sense with nonzero initial guess\n\
 45:                you probably want a KSP type of Richardson");
 46:   ksp->its = 0;
 47:   PCMatApply(ksp->pc,B,X);
 48:   PCGetFailedReason(ksp->pc,&pcreason);
 49:   /* Note: only some ranks may have this set; this may lead to problems if the caller assumes ksp->reason is set on all processes or just uses the result */
 50:   if (pcreason) {
 51:     MatSetInf(X);
 52:     ksp->reason = KSP_DIVERGED_PC_FAILED;
 53:   } else {
 54:     ksp->its    = 1;
 55:     ksp->reason = KSP_CONVERGED_ITS;
 56:   }
 57:   return(0);
 58: }

 60: /*MC
 61:      KSPPREONLY - This implements a method that applies ONLY the preconditioner exactly once.
 62:                   This may be used in inner iterations, where it is desired to
 63:                   allow multiple iterations as well as the "0-iteration" case. It is
 64:                   commonly used with the direct solver preconditioners like PCLU and PCCHOLESKY

 66:    Options Database Keys:
 67: .   -ksp_type preonly

 69:    Level: beginner

 71:    Notes:
 72:     Since this does not involve an iteration the basic KSP parameters such as tolerances and iteration counts
 73:           do not apply

 75:     To apply multiple preconditioners in a simple iteration use KSPRICHARDSON

 77:    Developer Notes:
 78:     Even though this method does not use any norms, the user is allowed to set the KSPNormType to any value.
 79:     This is so the users does not have to change KSPNormType options when they switch from other KSP methods to this one.

 81: .seealso:  KSPCreate(), KSPSetType(), KSPType, KSP, KSPRICHARDSON, KSPCHEBYSHEV

 83: M*/

 85: PETSC_EXTERN PetscErrorCode KSPCreate_PREONLY(KSP ksp)
 86: {

 90:   KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_LEFT,3);
 91:   KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_RIGHT,2);
 92:   KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,2);
 93:   KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_RIGHT,2);
 94:   KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_LEFT,2);
 95:   KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_RIGHT,2);
 96:   KSPSetSupportedNorm(ksp,KSP_NORM_NATURAL,PC_LEFT,2);

 98:   ksp->data                = NULL;
 99:   ksp->ops->setup          = KSPSetUp_PREONLY;
100:   ksp->ops->solve          = KSPSolve_PREONLY;
101:   ksp->ops->matsolve       = KSPMatSolve_PREONLY;
102:   ksp->ops->destroy        = KSPDestroyDefault;
103:   ksp->ops->buildsolution  = KSPBuildSolutionDefault;
104:   ksp->ops->buildresidual  = KSPBuildResidualDefault;
105:   ksp->ops->setfromoptions = NULL;
106:   ksp->ops->view           = NULL;
107:   return(0);
108: }