xref: /petsc/src/snes/impls/ksponly/ksponly.c (revision 3fc23ce8277e00a65c4b1be6e2b4535b56ba7978)
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) {
38     PetscCall((*snes->ops->update)(snes, 0));
39   }
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   if (snes->numbermonitors) {
62     PetscReal fnorm;
63     PetscCall(SNESComputeFunction(snes,X,F));
64     PetscCall(VecNorm(F,NORM_2,&fnorm));
65     PetscCall(SNESMonitor(snes,1,fnorm));
66   }
67   PetscFunctionReturn(0);
68 }
69 
70 static PetscErrorCode SNESSetUp_KSPONLY(SNES snes)
71 {
72   PetscFunctionBegin;
73   PetscCall(SNESSetUpMatrices(snes));
74   PetscFunctionReturn(0);
75 }
76 
77 static PetscErrorCode SNESDestroy_KSPONLY(SNES snes)
78 {
79   PetscFunctionBegin;
80   PetscCall(PetscFree(snes->data));
81   PetscFunctionReturn(0);
82 }
83 
84 /* -------------------------------------------------------------------------- */
85 /*MC
86       SNESKSPONLY - Nonlinear solver that only performs one Newton step and does not compute any norms.
87       The main purpose of this solver is to solve linear problems using the SNES interface, without
88       any additional overhead in the form of vector operations.
89 
90    Level: beginner
91 
92 .seealso:  SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR
93 M*/
94 PETSC_EXTERN PetscErrorCode SNESCreate_KSPONLY(SNES snes)
95 {
96   SNES_KSPONLY   *ksponly;
97 
98   PetscFunctionBegin;
99   snes->ops->setup          = SNESSetUp_KSPONLY;
100   snes->ops->solve          = SNESSolve_KSPONLY;
101   snes->ops->destroy        = SNESDestroy_KSPONLY;
102   snes->ops->setfromoptions = NULL;
103   snes->ops->view           = NULL;
104   snes->ops->reset          = NULL;
105 
106   snes->usesksp = PETSC_TRUE;
107   snes->usesnpc = PETSC_FALSE;
108 
109   snes->alwayscomputesfinalresidual = PETSC_FALSE;
110 
111   PetscCall(PetscNewLog(snes,&ksponly));
112   snes->data = (void*)ksponly;
113   PetscFunctionReturn(0);
114 }
115 
116 /*MC
117       SNESKSPTRANSPOSEONLY - Nonlinear solver that only performs one Newton step and does not compute any norms.
118       The main purpose of this solver is to solve transposed linear problems using the SNES interface, without
119       any additional overhead in the form of vector operations within adjoint solvers.
120 
121    Level: beginner
122 
123 .seealso:  SNESCreate(), SNES, SNESSetType(), SNESKSPTRANSPOSEONLY, SNESNEWTONLS, SNESNEWTONTR
124 M*/
125 PETSC_EXTERN PetscErrorCode SNESCreate_KSPTRANSPOSEONLY(SNES snes)
126 {
127   SNES_KSPONLY   *kspo;
128 
129   PetscFunctionBegin;
130   PetscCall(SNESCreate_KSPONLY(snes));
131   PetscCall(PetscObjectChangeTypeName((PetscObject)snes,SNESKSPTRANSPOSEONLY));
132   kspo = (SNES_KSPONLY*)snes->data;
133   kspo->transpose_solve = PETSC_TRUE;
134   PetscFunctionReturn(0);
135 }
136