xref: /petsc/src/snes/impls/ksponly/ksponly.c (revision 03436ad9ecd159e250ca4f04ee48a77d5bbe8d4f)
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 
42   SNESCheckJacobianDomainerror(snes);
43 
44   ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr);
45   if (ksponly->transpose_solve) {
46     ierr = KSPSolveTranspose(snes->ksp,F,Y);CHKERRQ(ierr);
47   } else {
48     ierr = KSPSolve(snes->ksp,F,Y);CHKERRQ(ierr);
49   }
50   snes->reason = SNES_CONVERGED_ITS;
51   SNESCheckKSPSolve(snes);
52 
53   ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
54   ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr);
55   snes->iter++;
56 
57   /* Take the computed step. */
58   ierr = VecAXPY(X,-1.0,Y);CHKERRQ(ierr);
59   if (snes->numbermonitors) {
60     PetscReal fnorm;
61     ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
62     ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);
63     ierr = SNESMonitor(snes,1,fnorm);CHKERRQ(ierr);
64   }
65   PetscFunctionReturn(0);
66 }
67 
68 static PetscErrorCode SNESSetUp_KSPONLY(SNES snes)
69 {
70   PetscErrorCode ierr;
71 
72   PetscFunctionBegin;
73   ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);
74   PetscFunctionReturn(0);
75 }
76 
77 static PetscErrorCode SNESDestroy_KSPONLY(SNES snes)
78 {
79   PetscErrorCode ierr;
80 
81   PetscFunctionBegin;
82   ierr = PetscFree(snes->data);CHKERRQ(ierr);
83   PetscFunctionReturn(0);
84 }
85 
86 /* -------------------------------------------------------------------------- */
87 /*MC
88       SNESKSPONLY - Nonlinear solver that only performs one Newton step and does not compute any norms.
89       The main purpose of this solver is to solve linear problems using the SNES interface, without
90       any additional overhead in the form of vector operations.
91 
92    Level: beginner
93 
94 .seealso:  SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR
95 M*/
96 PETSC_EXTERN PetscErrorCode SNESCreate_KSPONLY(SNES snes)
97 {
98   SNES_KSPONLY   *ksponly;
99   PetscErrorCode ierr;
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 = 0;
106   snes->ops->view           = 0;
107   snes->ops->reset          = 0;
108 
109   snes->usesksp = PETSC_TRUE;
110   snes->usesnpc = PETSC_FALSE;
111 
112   snes->alwayscomputesfinalresidual = PETSC_FALSE;
113 
114   ierr = PetscNewLog(snes,&ksponly);CHKERRQ(ierr);
115   snes->data = (void*)ksponly;
116   PetscFunctionReturn(0);
117 }
118 
119 /*MC
120       SNESKSPTRANSPOSEONLY - Nonlinear solver that only performs one Newton step and does not compute any norms.
121       The main purpose of this solver is to solve transposed linear problems using the SNES interface, without
122       any additional overhead in the form of vector operations within adjoint solvers.
123 
124    Level: beginner
125 
126 .seealso:  SNESCreate(), SNES, SNESSetType(), SNESKSPTRANSPOSEONLY, SNESNEWTONLS, SNESNEWTONTR
127 M*/
128 PETSC_EXTERN PetscErrorCode SNESCreate_KSPTRANSPOSEONLY(SNES snes)
129 {
130   SNES_KSPONLY   *kspo;
131   PetscErrorCode ierr;
132 
133   PetscFunctionBegin;
134   ierr = SNESCreate_KSPONLY(snes);CHKERRQ(ierr);
135   ierr = PetscObjectChangeTypeName((PetscObject)snes,SNESKSPTRANSPOSEONLY);CHKERRQ(ierr);
136   kspo = (SNES_KSPONLY*)snes->data;
137   kspo->transpose_solve = PETSC_TRUE;
138   PetscFunctionReturn(0);
139 }
140