xref: /petsc/src/snes/impls/ksponly/ksponly.c (revision d6f11fddcbc5503a4d8fb26ab62e20dc044b7e4b)
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   ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr);
44   snes->iter++;
45 
46   /* Take the computed step. */
47   ierr = VecAXPY(X,-1.0,Y);CHKERRQ(ierr);
48   if (snes->numbermonitors) {
49     PetscReal fnorm;
50     ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
51     ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);
52     ierr = SNESMonitor(snes,1,fnorm);CHKERRQ(ierr);
53   }
54   PetscFunctionReturn(0);
55 }
56 
57 static PetscErrorCode SNESSetUp_KSPONLY(SNES snes)
58 {
59   PetscErrorCode ierr;
60 
61   PetscFunctionBegin;
62   ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);
63   PetscFunctionReturn(0);
64 }
65 
66 static PetscErrorCode SNESDestroy_KSPONLY(SNES snes)
67 {
68 
69   PetscFunctionBegin;
70   PetscFunctionReturn(0);
71 }
72 
73 /* -------------------------------------------------------------------------- */
74 /*MC
75       SNESKSPONLY - Nonlinear solver that only performs one Newton step and does not compute any norms.
76       The main purpose of this solver is to solve linear problems using the SNES interface, without
77       any additional overhead in the form of vector operations.
78 
79    Level: beginner
80 
81 .seealso:  SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR
82 M*/
83 PETSC_EXTERN PetscErrorCode SNESCreate_KSPONLY(SNES snes)
84 {
85 
86   PetscFunctionBegin;
87   snes->ops->setup          = SNESSetUp_KSPONLY;
88   snes->ops->solve          = SNESSolve_KSPONLY;
89   snes->ops->destroy        = SNESDestroy_KSPONLY;
90   snes->ops->setfromoptions = 0;
91   snes->ops->view           = 0;
92   snes->ops->reset          = 0;
93 
94   snes->usesksp = PETSC_TRUE;
95   snes->usesnpc = PETSC_FALSE;
96 
97   snes->alwayscomputesfinalresidual = PETSC_FALSE;
98 
99   snes->data = 0;
100   PetscFunctionReturn(0);
101 }
102