xref: /petsc/src/snes/impls/ksponly/ksponly.c (revision f3fe499b4cc4d64bf04aa4f5e4963dcc4eb56541)
1 #define PETSCSNES_DLL
2 
3 #include "private/snesimpl.h"
4 
5 #undef __FUNCT__
6 #define __FUNCT__ "SNESSolve_KSPONLY"
7 static PetscErrorCode SNESSolve_KSPONLY(SNES snes)
8 {
9   PetscErrorCode     ierr;
10   PetscInt           lits;
11   MatStructure       flg = DIFFERENT_NONZERO_PATTERN;
12   Vec                Y,X,F;
13   KSPConvergedReason kspreason;
14 
15   PetscFunctionBegin;
16   snes->numFailures            = 0;
17   snes->numLinearSolveFailures = 0;
18   snes->reason                 = SNES_CONVERGED_ITERATING;
19   snes->iter                   = 0;
20   snes->norm                   = 0.0;
21 
22   X = snes->vec_sol;
23   F = snes->vec_func;
24   Y = snes->vec_sol_update;
25 
26   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
27   if (snes->domainerror) {
28     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
29     PetscFunctionReturn(0);
30   }
31 
32   /* Solve J Y = F, where J is Jacobian matrix */
33   ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
34   ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr);
35   ierr = KSPSolve(snes->ksp,F,Y);CHKERRQ(ierr);
36   ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr);
37   if (kspreason < 0 && ++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) {
38     ierr = PetscInfo2(snes,"iter=%D, number linear solve failures %D greater than current SNES allowed, stopping solve\n",snes->iter,snes->numLinearSolveFailures);CHKERRQ(ierr);
39     snes->reason = SNES_DIVERGED_LINEAR_SOLVE;
40   } else {
41     snes->reason = SNES_CONVERGED_ITS;
42   }
43   ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
44   snes->linear_its += lits;
45   ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr);
46   snes->iter++;
47 
48   /* Take the computed step. */
49   ierr = VecAXPY(X,-1.0,Y);CHKERRQ(ierr);
50   PetscFunctionReturn(0);
51 }
52 
53 #undef __FUNCT__
54 #define __FUNCT__ "SNESSetUp_KSPONLY"
55 static PetscErrorCode SNESSetUp_KSPONLY(SNES snes)
56 {
57   PetscErrorCode ierr;
58 
59   PetscFunctionBegin;
60   if (!snes->vec_sol_update) {
61     ierr = VecDuplicate(snes->vec_sol,&snes->vec_sol_update);CHKERRQ(ierr);
62     ierr = PetscLogObjectParent(snes,snes->vec_sol_update);CHKERRQ(ierr);
63   }
64   PetscFunctionReturn(0);
65 }
66 
67 #undef __FUNCT__
68 #define __FUNCT__ "SNESDestroy_KSPONLY"
69 static PetscErrorCode SNESDestroy_KSPONLY(SNES snes)
70 {
71   PetscErrorCode ierr;
72 
73   PetscFunctionBegin;
74   if (snes->vec_sol_update) {
75     ierr = VecDestroy(snes->vec_sol_update);CHKERRQ(ierr);
76     snes->vec_sol_update = PETSC_NULL;
77   }
78   PetscFunctionReturn(0);
79 }
80 
81 /* -------------------------------------------------------------------------- */
82 /*MC
83       SNESKSPONLY - Nonlinear solver that only performs one Newton step and does not compute any norms.
84       The main purpose of this solver is to solve linear problems using the SNES interface, without
85       any additional overhead in the form of vector operations.
86 
87    Level: beginner
88 
89 .seealso:  SNESCreate(), SNES, SNESSetType(), SNESLS, SNESTR
90 M*/
91 EXTERN_C_BEGIN
92 #undef __FUNCT__
93 #define __FUNCT__ "SNESCreate_KSPONLY"
94 PetscErrorCode PETSCSNES_DLLEXPORT SNESCreate_KSPONLY(SNES snes)
95 {
96 
97   PetscFunctionBegin;
98   snes->ops->setup   = SNESSetUp_KSPONLY;
99   snes->ops->solve   = SNESSolve_KSPONLY;
100   snes->ops->destroy = SNESDestroy_KSPONLY;
101 
102   snes->data = 0;
103   PetscFunctionReturn(0);
104 }
105 EXTERN_C_END
106