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