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