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