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