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