Actual source code: preonly.c

petsc-main 2021-04-13
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:   if (ksp->numbermonitors) {
 33:     Vec       v;
 34:     PetscReal norm;
 35:     Mat       A;

 37:     VecNorm(ksp->vec_rhs,NORM_2,&norm);
 38:     KSPMonitor(ksp,0,norm);
 39:     VecDuplicate(ksp->vec_rhs,&v);
 40:     PCGetOperators(ksp->pc,&A,NULL);
 41:     KSP_MatMult(ksp,A,ksp->vec_sol,v);
 42:     VecAYPX(v,-1.0,ksp->vec_rhs);
 43:     VecNorm(v,NORM_2,&norm);
 44:     VecDestroy(&v);
 45:     KSPMonitor(ksp,1,norm);
 46:   }
 47:   return(0);
 48: }

 50: static PetscErrorCode KSPMatSolve_PREONLY(KSP ksp, Mat B, Mat X)
 51: {
 53:   PetscBool      diagonalscale;
 54:   PCFailedReason pcreason;

 57:   PCGetDiagonalScale(ksp->pc,&diagonalscale);
 58:   if (diagonalscale) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"Krylov method %s does not support diagonal scaling",((PetscObject)ksp)->type_name);
 59:   if (!ksp->guess_zero) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_USER,"Running KSP of preonly doesn't make sense with nonzero initial guess\n\
 60:                you probably want a KSP type of Richardson");
 61:   ksp->its = 0;
 62:   PCMatApply(ksp->pc,B,X);
 63:   PCGetFailedReason(ksp->pc,&pcreason);
 64:   /* 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 */
 65:   if (pcreason) {
 66:     MatSetInf(X);
 67:     ksp->reason = KSP_DIVERGED_PC_FAILED;
 68:   } else {
 69:     ksp->its    = 1;
 70:     ksp->reason = KSP_CONVERGED_ITS;
 71:   }
 72:   return(0);
 73: }

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

 81:    Options Database Keys:
 82: .   -ksp_type preonly

 84:    Level: beginner

 86:    Notes:
 87:     Since this does not involve an iteration the basic KSP parameters such as tolerances and iteration counts
 88:           do not apply

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

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

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

 98: M*/

100: PETSC_EXTERN PetscErrorCode KSPCreate_PREONLY(KSP ksp)
101: {

105:   KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_LEFT,3);
106:   KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_RIGHT,2);
107:   KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,2);
108:   KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_RIGHT,2);
109:   KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_LEFT,2);
110:   KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_RIGHT,2);
111:   KSPSetSupportedNorm(ksp,KSP_NORM_NATURAL,PC_LEFT,2);

113:   ksp->data                = NULL;
114:   ksp->ops->setup          = KSPSetUp_PREONLY;
115:   ksp->ops->solve          = KSPSolve_PREONLY;
116:   ksp->ops->matsolve       = KSPMatSolve_PREONLY;
117:   ksp->ops->destroy        = KSPDestroyDefault;
118:   ksp->ops->buildsolution  = KSPBuildSolutionDefault;
119:   ksp->ops->buildresidual  = KSPBuildResidualDefault;
120:   ksp->ops->setfromoptions = NULL;
121:   ksp->ops->view           = NULL;
122:   return(0);
123: }