xref: /petsc/src/snes/impls/ksponly/ksponly.c (revision 96d2d4e583fe3011d58d6b0fe81dc6dfaa0a9d13)
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   PetscErrorCode ierr;
11   PetscInt       lits;
12   Vec            Y,X,F;
13 
14   PetscFunctionBegin;
15   if (snes->xl || snes->xu || snes->ops->computevariablebounds) SETERRQ1(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   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
28   if (snes->numbermonitors) {
29     PetscReal fnorm;
30     ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);
31     ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
32   }
33 
34   /* Call general purpose update function */
35   if (snes->ops->update) {
36     ierr = (*snes->ops->update)(snes, 0);CHKERRQ(ierr);
37   }
38 
39   /* Solve J Y = F, where J is Jacobian matrix */
40   ierr = SNESComputeJacobian(snes,X,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr);
41   ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr);
42   if (ksponly->transpose_solve) {
43     ierr = KSPSolveTranspose(snes->ksp,F,Y);CHKERRQ(ierr);
44   } else {
45     ierr = KSPSolve(snes->ksp,F,Y);CHKERRQ(ierr);
46   }
47   snes->reason = SNES_CONVERGED_ITS;
48   SNESCheckKSPSolve(snes);
49 
50   ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
51   snes->linear_its += lits;
52   ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr);
53   snes->iter++;
54 
55   /* Take the computed step. */
56   ierr = VecAXPY(X,-1.0,Y);CHKERRQ(ierr);
57   if (snes->numbermonitors) {
58     PetscReal fnorm;
59     ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
60     ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);
61     ierr = SNESMonitor(snes,1,fnorm);CHKERRQ(ierr);
62   }
63   PetscFunctionReturn(0);
64 }
65 
66 static PetscErrorCode SNESSetUp_KSPONLY(SNES snes)
67 {
68   PetscErrorCode ierr;
69 
70   PetscFunctionBegin;
71   ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);
72   PetscFunctionReturn(0);
73 }
74 
75 static PetscErrorCode SNESDestroy_KSPONLY(SNES snes)
76 {
77   PetscErrorCode ierr;
78 
79   PetscFunctionBegin;
80   ierr = PetscFree(snes->data);CHKERRQ(ierr);
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   PetscErrorCode ierr;
98 
99   PetscFunctionBegin;
100   snes->ops->setup          = SNESSetUp_KSPONLY;
101   snes->ops->solve          = SNESSolve_KSPONLY;
102   snes->ops->destroy        = SNESDestroy_KSPONLY;
103   snes->ops->setfromoptions = 0;
104   snes->ops->view           = 0;
105   snes->ops->reset          = 0;
106 
107   snes->usesksp = PETSC_TRUE;
108   snes->usesnpc = PETSC_FALSE;
109 
110   snes->alwayscomputesfinalresidual = PETSC_FALSE;
111 
112   ierr = PetscNewLog(snes,&ksponly);CHKERRQ(ierr);
113   snes->data = (void*)ksponly;
114   PetscFunctionReturn(0);
115 }
116 
117 /*MC
118       SNESKSPTRANSPOSEONLY - Nonlinear solver that only performs one Newton step and does not compute any norms.
119       The main purpose of this solver is to solve transposed linear problems using the SNES interface, without
120       any additional overhead in the form of vector operations within adjoint solvers.
121 
122    Level: beginner
123 
124 .seealso:  SNESCreate(), SNES, SNESSetType(), SNESKSPTRANSPOSEONLY, SNESNEWTONLS, SNESNEWTONTR
125 M*/
126 PETSC_EXTERN PetscErrorCode SNESCreate_KSPTRANSPOSEONLY(SNES snes)
127 {
128   SNES_KSPONLY   *kspo;
129   PetscErrorCode ierr;
130 
131   PetscFunctionBegin;
132   ierr = SNESCreate_KSPONLY(snes);CHKERRQ(ierr);
133   kspo = (SNES_KSPONLY*)snes->data;
134   kspo->transpose_solve = PETSC_TRUE;
135   PetscFunctionReturn(0);
136 }
137