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