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