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