xref: /petsc/src/snes/impls/ksponly/ksponly.c (revision 2f75a4aaab0e3f692f7125f85b5f5852e7ceb6cb)
1 
2 #include <petsc/private/snesimpl.h>
3 
4 static PetscErrorCode SNESSolve_KSPONLY(SNES snes)
5 {
6   PetscErrorCode     ierr;
7   PetscInt           lits;
8   Vec                Y,X,F;
9 
10   PetscFunctionBegin;
11   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);
12 
13   snes->numFailures            = 0;
14   snes->numLinearSolveFailures = 0;
15   snes->reason                 = SNES_CONVERGED_ITERATING;
16   snes->iter                   = 0;
17   snes->norm                   = 0.0;
18 
19   X = snes->vec_sol;
20   F = snes->vec_func;
21   Y = snes->vec_sol_update;
22 
23   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
24   if (snes->numbermonitors) {
25     PetscReal fnorm;
26     ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);
27     ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
28   }
29 
30   /* Call general purpose update function */
31   if (snes->ops->update) {
32     ierr = (*snes->ops->update)(snes, 0);CHKERRQ(ierr);
33   }
34 
35   /* Solve J Y = F, where J is Jacobian matrix */
36   ierr = SNESComputeJacobian(snes,X,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr);
37   ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr);
38   ierr = KSPSolve(snes->ksp,F,Y);CHKERRQ(ierr);
39   snes->reason = SNES_CONVERGED_ITS;
40   SNESCheckKSPSolve(snes);
41 
42   ierr              = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
43   snes->linear_its += lits;
44   ierr              = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr);
45   snes->iter++;
46 
47   /* Take the computed step. */
48   ierr = VecAXPY(X,-1.0,Y);CHKERRQ(ierr);
49   if (snes->numbermonitors) {
50     PetscReal fnorm;
51     ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
52     ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);
53     ierr = SNESMonitor(snes,1,fnorm);CHKERRQ(ierr);
54   }
55   PetscFunctionReturn(0);
56 }
57 
58 static PetscErrorCode SNESSetUp_KSPONLY(SNES snes)
59 {
60   PetscErrorCode ierr;
61 
62   PetscFunctionBegin;
63   ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);
64   PetscFunctionReturn(0);
65 }
66 
67 static PetscErrorCode SNESDestroy_KSPONLY(SNES snes)
68 {
69 
70   PetscFunctionBegin;
71   PetscFunctionReturn(0);
72 }
73 
74 /* -------------------------------------------------------------------------- */
75 /*MC
76       SNESKSPONLY - Nonlinear solver that only performs one Newton step and does not compute any norms.
77       The main purpose of this solver is to solve linear problems using the SNES interface, without
78       any additional overhead in the form of vector operations.
79 
80    Level: beginner
81 
82 .seealso:  SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR
83 M*/
84 PETSC_EXTERN PetscErrorCode SNESCreate_KSPONLY(SNES snes)
85 {
86 
87   PetscFunctionBegin;
88   snes->ops->setup          = SNESSetUp_KSPONLY;
89   snes->ops->solve          = SNESSolve_KSPONLY;
90   snes->ops->destroy        = SNESDestroy_KSPONLY;
91   snes->ops->setfromoptions = 0;
92   snes->ops->view           = 0;
93   snes->ops->reset          = 0;
94 
95   snes->usesksp = PETSC_TRUE;
96   snes->usesnpc = PETSC_FALSE;
97 
98   snes->alwayscomputesfinalresidual = PETSC_FALSE;
99 
100   snes->data = 0;
101   PetscFunctionReturn(0);
102 }
103