Actual source code: ksponly.c
1: #include <petsc/private/snesimpl.h>
3: typedef struct {
4: PetscBool transpose_solve;
5: } SNES_KSPONLY;
7: static PetscErrorCode SNESSolve_KSPONLY(SNES snes)
8: {
9: SNES_KSPONLY *ksponly = (SNES_KSPONLY *)snes->data;
10: PetscInt lits;
11: Vec Y, X, F;
13: PetscFunctionBegin;
14: PetscCheck(!snes->xl && !snes->xu && !snes->ops->computevariablebounds, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name);
16: snes->numFailures = 0;
17: snes->numLinearSolveFailures = 0;
18: snes->reason = SNES_CONVERGED_ITERATING;
19: snes->iter = 0;
20: snes->norm = 0.0;
22: X = snes->vec_sol;
23: F = snes->vec_func;
24: Y = snes->vec_sol_update;
26: if (!snes->vec_func_init_set) {
27: PetscCall(SNESComputeFunction(snes, X, F));
28: } else snes->vec_func_init_set = PETSC_FALSE;
30: if (snes->numbermonitors) {
31: PetscReal fnorm;
32: PetscCall(VecNorm(F, NORM_2, &fnorm));
33: SNESCheckFunctionNorm(snes, fnorm);
34: PetscCall(SNESMonitor(snes, 0, fnorm));
35: }
37: /* Call general purpose update function */
38: PetscTryTypeMethod(snes, update, 0);
40: /* Solve J Y = F, where J is Jacobian matrix */
41: PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre));
43: SNESCheckJacobianDomainerror(snes);
45: PetscCall(KSPSetOperators(snes->ksp, snes->jacobian, snes->jacobian_pre));
46: if (ksponly->transpose_solve) {
47: PetscCall(KSPSolveTranspose(snes->ksp, F, Y));
48: } else {
49: PetscCall(KSPSolve(snes->ksp, F, Y));
50: }
51: snes->reason = SNES_CONVERGED_ITS;
52: SNESCheckKSPSolve(snes);
54: PetscCall(KSPGetIterationNumber(snes->ksp, &lits));
55: PetscCall(PetscInfo(snes, "iter=%" PetscInt_FMT ", linear solve iterations=%" PetscInt_FMT "\n", snes->iter, lits));
56: snes->iter++;
58: /* Take the computed step. */
59: PetscCall(VecAXPY(X, -1.0, Y));
60: if (snes->numbermonitors) {
61: PetscReal fnorm;
62: PetscCall(SNESComputeFunction(snes, X, F));
63: PetscCall(VecNorm(F, NORM_2, &fnorm));
64: SNESCheckFunctionNorm(snes, fnorm);
65: PetscCall(SNESMonitor(snes, 1, fnorm));
66: }
67: PetscFunctionReturn(PETSC_SUCCESS);
68: }
70: static PetscErrorCode SNESSetUp_KSPONLY(SNES snes)
71: {
72: PetscFunctionBegin;
73: PetscCall(SNESSetUpMatrices(snes));
74: PetscFunctionReturn(PETSC_SUCCESS);
75: }
77: static PetscErrorCode SNESDestroy_KSPONLY(SNES snes)
78: {
79: PetscFunctionBegin;
80: PetscCall(PetscFree(snes->data));
81: PetscFunctionReturn(PETSC_SUCCESS);
82: }
84: /*MC
85: SNESKSPONLY - Nonlinear solver that performs one Newton step with `KSPSolve()` and does not compute any norms.
87: Level: beginner
89: Note:
90: The main purpose of this solver is to solve linear problems using the `SNES` interface, without
91: any additional overhead in the form of vector norm operations.
93: .seealso: [](ch_snes), `SNES`, `SNESType`, `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONLS`, `SNESNEWTONTR`, `SNESKSPTRANSPOSEONLY`
94: M*/
95: PETSC_EXTERN PetscErrorCode SNESCreate_KSPONLY(SNES snes)
96: {
97: SNES_KSPONLY *ksponly;
99: PetscFunctionBegin;
100: snes->ops->setup = SNESSetUp_KSPONLY;
101: snes->ops->solve = SNESSolve_KSPONLY;
102: snes->ops->destroy = SNESDestroy_KSPONLY;
103: snes->ops->setfromoptions = NULL;
104: snes->ops->view = NULL;
105: snes->ops->reset = NULL;
107: snes->usesksp = PETSC_TRUE;
108: snes->usesnpc = PETSC_FALSE;
110: snes->alwayscomputesfinalresidual = PETSC_FALSE;
112: PetscCall(SNESParametersInitialize(snes));
114: PetscCall(PetscNew(&ksponly));
115: snes->data = (void *)ksponly;
116: PetscFunctionReturn(PETSC_SUCCESS);
117: }
119: /*MC
120: SNESKSPTRANSPOSEONLY - Nonlinear solver that performs one Newton step with `KSPSolveTranspose()` and does not compute any norms.
122: Level: beginner
124: Note:
125: The main purpose of this solver is to solve transposed linear problems using the `SNES` interface, without
126: any additional overhead in the form of vector operations within adjoint solvers.
128: .seealso: [](ch_snes), `SNES`, `SNESType`, `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESKS`, `SNESNEWTONLS`, `SNESNEWTONTR`
129: M*/
130: PETSC_EXTERN PetscErrorCode SNESCreate_KSPTRANSPOSEONLY(SNES snes)
131: {
132: SNES_KSPONLY *kspo;
134: PetscFunctionBegin;
135: PetscCall(SNESCreate_KSPONLY(snes));
136: PetscCall(PetscObjectChangeTypeName((PetscObject)snes, SNESKSPTRANSPOSEONLY));
137: kspo = (SNES_KSPONLY *)snes->data;
138: kspo->transpose_solve = PETSC_TRUE;
139: PetscFunctionReturn(PETSC_SUCCESS);
140: }