xref: /petsc/src/snes/impls/ksponly/ksponly.c (revision 6cab3a1b8d4fca83e1d7b61474c608460de73d5f)
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 {
56*6cab3a1bSJed Brown   PetscErrorCode ierr;
57b79b07cfSJed Brown 
58b79b07cfSJed Brown   PetscFunctionBegin;
59*6cab3a1bSJed Brown   ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);
60b79b07cfSJed Brown   PetscFunctionReturn(0);
61b79b07cfSJed Brown }
62b79b07cfSJed Brown 
63b79b07cfSJed Brown #undef __FUNCT__
64b79b07cfSJed Brown #define __FUNCT__ "SNESDestroy_KSPONLY"
65b79b07cfSJed Brown static PetscErrorCode SNESDestroy_KSPONLY(SNES snes)
66b79b07cfSJed Brown {
67b79b07cfSJed Brown 
68b79b07cfSJed Brown   PetscFunctionBegin;
69b79b07cfSJed Brown   PetscFunctionReturn(0);
70b79b07cfSJed Brown }
71b79b07cfSJed Brown 
72b79b07cfSJed Brown /* -------------------------------------------------------------------------- */
73b79b07cfSJed Brown /*MC
74b79b07cfSJed Brown       SNESKSPONLY - Nonlinear solver that only performs one Newton step and does not compute any norms.
75b79b07cfSJed Brown       The main purpose of this solver is to solve linear problems using the SNES interface, without
76b79b07cfSJed Brown       any additional overhead in the form of vector operations.
77b79b07cfSJed Brown 
78b79b07cfSJed Brown    Level: beginner
79b79b07cfSJed Brown 
80b79b07cfSJed Brown .seealso:  SNESCreate(), SNES, SNESSetType(), SNESLS, SNESTR
81b79b07cfSJed Brown M*/
82b79b07cfSJed Brown EXTERN_C_BEGIN
83b79b07cfSJed Brown #undef __FUNCT__
84b79b07cfSJed Brown #define __FUNCT__ "SNESCreate_KSPONLY"
857087cfbeSBarry Smith PetscErrorCode  SNESCreate_KSPONLY(SNES snes)
86b79b07cfSJed Brown {
87b79b07cfSJed Brown 
88b79b07cfSJed Brown   PetscFunctionBegin;
89b79b07cfSJed Brown   snes->ops->setup           = SNESSetUp_KSPONLY;
90b79b07cfSJed Brown   snes->ops->solve           = SNESSolve_KSPONLY;
91b79b07cfSJed Brown   snes->ops->destroy         = SNESDestroy_KSPONLY;
9237596af1SLisandro Dalcin   snes->ops->setfromoptions  = 0;
9337596af1SLisandro Dalcin   snes->ops->view            = 0;
9437596af1SLisandro Dalcin   snes->ops->reset           = 0;
95b79b07cfSJed Brown 
9642f4f86dSBarry Smith   snes->usesksp         = PETSC_TRUE;
9742f4f86dSBarry Smith   snes->usespc          = PETSC_FALSE;
9842f4f86dSBarry Smith 
99b79b07cfSJed Brown   snes->data = 0;
100b79b07cfSJed Brown   PetscFunctionReturn(0);
101b79b07cfSJed Brown }
102b79b07cfSJed Brown EXTERN_C_END
103