Actual source code: preonly.c

  1: #include <petsc/private/kspimpl.h>

  3: static PetscErrorCode KSPSetUp_PREONLY(KSP ksp)
  4: {
  5:   PetscFunctionBegin;
  6:   PetscFunctionReturn(PETSC_SUCCESS);
  7: }

  9: static PetscErrorCode KSPSolve_PREONLY(KSP ksp)
 10: {
 11:   PetscBool      diagonalscale;
 12:   PCFailedReason pcreason;

 14:   PetscFunctionBegin;
 15:   PetscCall(PCGetDiagonalScale(ksp->pc, &diagonalscale));
 16:   PetscCheck(!diagonalscale, PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "Krylov method %s does not support diagonal scaling", ((PetscObject)ksp)->type_name);
 17:   if (!ksp->guess_zero) {
 18:     PetscBool red;
 19:     PetscCall(PetscObjectTypeCompare((PetscObject)ksp->pc, PCREDISTRIBUTE, &red));
 20:     PetscCheck(red, PetscObjectComm((PetscObject)ksp), PETSC_ERR_USER, "KSP of type preonly (application of preconditioner only) doesn't make sense with nonzero initial guess\n\
 21:                 you probably want a KSP of type Richardson");
 22:   }
 23:   ksp->its = 0;
 24:   PetscCall(KSP_PCApply(ksp, ksp->vec_rhs, ksp->vec_sol));

 26:   PetscCall(PCReduceFailedReason(ksp->pc));
 27:   PetscCall(PCGetFailedReason(ksp->pc, &pcreason));
 28:   if (pcreason) {
 29:     PetscCheck(!ksp->errorifnotconverged, PetscObjectComm((PetscObject)ksp), PETSC_ERR_NOT_CONVERGED, "KSPSolve has not converged with PCFailedReason %s", PCFailedReasons[pcreason]);
 30:     PetscCall(VecSetInf(ksp->vec_sol));
 31:     ksp->reason = KSP_DIVERGED_PC_FAILED;
 32:   } else {
 33:     ksp->its    = 1;
 34:     ksp->reason = KSP_CONVERGED_ITS;
 35:   }

 37:   if (ksp->numbermonitors) {
 38:     Vec       v;
 39:     PetscReal norm;
 40:     Mat       A;

 42:     PetscCall(VecNorm(ksp->vec_rhs, NORM_2, &norm));
 43:     PetscCall(KSPMonitor(ksp, 0, norm));
 44:     PetscCall(VecDuplicate(ksp->vec_rhs, &v));
 45:     PetscCall(PCGetOperators(ksp->pc, &A, NULL));
 46:     PetscCall(KSP_MatMult(ksp, A, ksp->vec_sol, v));
 47:     PetscCall(VecAYPX(v, -1.0, ksp->vec_rhs));
 48:     PetscCall(VecNorm(v, NORM_2, &norm));
 49:     PetscCall(VecDestroy(&v));
 50:     PetscCall(KSPMonitor(ksp, 1, norm));
 51:   }
 52:   PetscFunctionReturn(PETSC_SUCCESS);
 53: }

 55: static PetscErrorCode KSPMatSolve_PREONLY(KSP ksp, Mat B, Mat X)
 56: {
 57:   PetscBool      diagonalscale;
 58:   PCFailedReason pcreason;

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

 79: /*MC
 80:    KSPNONE - An alias for `KSPPREONLY`

 82:    Options Database Key:
 83: .   -ksp_type none - use a single application of the preconditioner only

 85:    Level: beginner

 87:    Note:
 88:    See `KSPPREONLY` for more details

 90: .seealso: [](ch_ksp), `KSPCreate()`, `KSPSetType()`, `KSPType`, `KSPPREONLY`, `KSP`, `KSPRICHARDSON`, `KSPCHEBYSHEV`, `KSPGetPC()`, `KSPSetInitialGuessNonzero()`,
 91:           `PCREDISTRIBUTE`, `PCRedistributeGetKSP()`, `KSPPREONLY`
 92: M*/

 94: /*MC
 95:    KSPPREONLY - This implements a method that applies ONLY the preconditioner exactly once.

 97:    This may be used in inner iterations, where it is desired to
 98:    allow multiple iterations as well as the "0-iteration" case. It is
 99:    commonly used with the direct solver preconditioners like `PCLU` and `PCCHOLESKY`.
100:    There is an alias of this with the name `KSPNONE`.

102:    Options Database Key:
103: .   -ksp_type preonly - use a single application of the preconditioner only

105:    Level: beginner

107:    Notes:
108:    Since this does not involve an iteration the basic `KSP` parameters such as tolerances and iteration counts
109:    do not apply

111:    To apply multiple preconditioners in a simple iteration use `KSPRICHARDSON`

113:    This `KSPType` cannot be used with the flag `-ksp_initial_guess_nonzero` or the call `KSPSetInitialGuessNonzero()` since it simply applies
114:    the preconditioner to the given right-hand side during `KSPSolve()`. Except when the
115:    `PCType` is `PCREDISTRIBUTE`; in that situation pass the nonzero initial guess flag with `-ksp_initial_guess_nonzero` or `KSPSetInitialGuessNonzero()`
116:    both to the outer `KSP` (which is `KSPPREONLY`) and the inner `KSP` object obtained with `KSPGetPC()` followed by `PCRedistributedGetKSP()`
117:    followed by `KSPSetInitialGuessNonzero()` or the option  `-redistribute_ksp_initial_guess_nonzero`.

119:    Developer Note:
120:    Even though this method does not use any norms, the user is allowed to set the `KSPNormType` to any value.
121:    This is so the users does not have to change `KSPNormType` options when they switch from other `KSP` methods to this one.

123: .seealso: [](ch_ksp), `KSPCreate()`, `KSPSetType()`, `KSPType`, `KSP`, `KSPRICHARDSON`, `KSPCHEBYSHEV`, `KSPGetPC()`, `KSPSetInitialGuessNonzero()`,
124:           `PCREDISTRIBUTE`, `PCRedistributeGetKSP()`, `KSPNONE`
125: M*/

127: PETSC_EXTERN PetscErrorCode KSPCreate_PREONLY(KSP ksp)
128: {
129:   PetscFunctionBegin;
130:   PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_NONE, PC_LEFT, 3));
131:   PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_NONE, PC_RIGHT, 2));
132:   PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_PRECONDITIONED, PC_LEFT, 2));
133:   PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_PRECONDITIONED, PC_RIGHT, 2));
134:   PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_UNPRECONDITIONED, PC_LEFT, 2));
135:   PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_UNPRECONDITIONED, PC_RIGHT, 2));
136:   PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_NATURAL, PC_LEFT, 2));

138:   ksp->data                = NULL;
139:   ksp->ops->setup          = KSPSetUp_PREONLY;
140:   ksp->ops->solve          = KSPSolve_PREONLY;
141:   ksp->ops->matsolve       = KSPMatSolve_PREONLY;
142:   ksp->ops->destroy        = KSPDestroyDefault;
143:   ksp->ops->buildsolution  = KSPBuildSolutionDefault;
144:   ksp->ops->buildresidual  = KSPBuildResidualDefault;
145:   ksp->ops->setfromoptions = NULL;
146:   ksp->ops->view           = NULL;
147:   ksp->guess_not_read      = PETSC_TRUE; // A KSP of preonly never needs to zero the input x since PC do not use an initial guess
148:   PetscFunctionReturn(PETSC_SUCCESS);
149: }