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