xref: /petsc/src/snes/impls/ksponly/ksponly.c (revision efd4aadf157bf1ba2d80c2be092fcf4247860003)
1b79b07cfSJed Brown 
2af0996ceSBarry Smith #include <petsc/private/snesimpl.h>
3b79b07cfSJed Brown 
4b79b07cfSJed Brown static PetscErrorCode SNESSolve_KSPONLY(SNES snes)
5b79b07cfSJed Brown {
6b79b07cfSJed Brown   PetscErrorCode     ierr;
7b79b07cfSJed Brown   PetscInt           lits;
8b79b07cfSJed Brown   Vec                Y,X,F;
9b79b07cfSJed Brown 
10b79b07cfSJed Brown   PetscFunctionBegin;
116c4ed002SBarry Smith   if (snes->xl || snes->xu || snes->ops->computevariablebounds) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name);
12c579b300SPatrick Farrell 
13b79b07cfSJed Brown   snes->numFailures            = 0;
14b79b07cfSJed Brown   snes->numLinearSolveFailures = 0;
15b79b07cfSJed Brown   snes->reason                 = SNES_CONVERGED_ITERATING;
16b79b07cfSJed Brown   snes->iter                   = 0;
17b79b07cfSJed Brown   snes->norm                   = 0.0;
18b79b07cfSJed Brown 
19b79b07cfSJed Brown   X = snes->vec_sol;
20b79b07cfSJed Brown   F = snes->vec_func;
21b79b07cfSJed Brown   Y = snes->vec_sol_update;
22b79b07cfSJed Brown 
23b79b07cfSJed Brown   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
24a5eaddd9SBarry Smith   if (snes->numbermonitors) {
25a5eaddd9SBarry Smith     PetscReal fnorm;
26a5eaddd9SBarry Smith     ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);
27a5eaddd9SBarry Smith     ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
28a5eaddd9SBarry Smith   }
29b79b07cfSJed Brown 
305341784dSBarry Smith   /* Call general purpose update function */
315341784dSBarry Smith   if (snes->ops->update) {
325341784dSBarry Smith     ierr = (*snes->ops->update)(snes, 0);CHKERRQ(ierr);
335341784dSBarry Smith   }
345341784dSBarry Smith 
35b79b07cfSJed Brown   /* Solve J Y = F, where J is Jacobian matrix */
36d1e9a80fSBarry Smith   ierr = SNESComputeJacobian(snes,X,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr);
3723ee1639SBarry Smith   ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr);
38b79b07cfSJed Brown   ierr = KSPSolve(snes->ksp,F,Y);CHKERRQ(ierr);
39422a814eSBarry Smith   snes->reason = SNES_CONVERGED_ITS;
40422a814eSBarry Smith   SNESCheckKSPSolve(snes);
411aa26658SKarl Rupp 
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);
49a5eaddd9SBarry Smith   if (snes->numbermonitors) {
50a5eaddd9SBarry Smith     PetscReal fnorm;
51fef909d3SLawrence Mitchell     ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
52a5eaddd9SBarry Smith     ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);
53a5eaddd9SBarry Smith     ierr = SNESMonitor(snes,1,fnorm);CHKERRQ(ierr);
54a5eaddd9SBarry Smith   }
55b79b07cfSJed Brown   PetscFunctionReturn(0);
56b79b07cfSJed Brown }
57b79b07cfSJed Brown 
58b79b07cfSJed Brown static PetscErrorCode SNESSetUp_KSPONLY(SNES snes)
59b79b07cfSJed Brown {
606cab3a1bSJed Brown   PetscErrorCode ierr;
61b79b07cfSJed Brown 
62b79b07cfSJed Brown   PetscFunctionBegin;
636cab3a1bSJed Brown   ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);
64b79b07cfSJed Brown   PetscFunctionReturn(0);
65b79b07cfSJed Brown }
66b79b07cfSJed Brown 
67b79b07cfSJed Brown static PetscErrorCode SNESDestroy_KSPONLY(SNES snes)
68b79b07cfSJed Brown {
69b79b07cfSJed Brown 
70b79b07cfSJed Brown   PetscFunctionBegin;
71b79b07cfSJed Brown   PetscFunctionReturn(0);
72b79b07cfSJed Brown }
73b79b07cfSJed Brown 
74b79b07cfSJed Brown /* -------------------------------------------------------------------------- */
75b79b07cfSJed Brown /*MC
76b79b07cfSJed Brown       SNESKSPONLY - Nonlinear solver that only performs one Newton step and does not compute any norms.
77b79b07cfSJed Brown       The main purpose of this solver is to solve linear problems using the SNES interface, without
78b79b07cfSJed Brown       any additional overhead in the form of vector operations.
79b79b07cfSJed Brown 
80b79b07cfSJed Brown    Level: beginner
81b79b07cfSJed Brown 
8204d7464bSBarry Smith .seealso:  SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR
83b79b07cfSJed Brown M*/
848cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_KSPONLY(SNES snes)
85b79b07cfSJed Brown {
86b79b07cfSJed Brown 
87b79b07cfSJed Brown   PetscFunctionBegin;
88b79b07cfSJed Brown   snes->ops->setup          = SNESSetUp_KSPONLY;
89b79b07cfSJed Brown   snes->ops->solve          = SNESSolve_KSPONLY;
90b79b07cfSJed Brown   snes->ops->destroy        = SNESDestroy_KSPONLY;
9137596af1SLisandro Dalcin   snes->ops->setfromoptions = 0;
9237596af1SLisandro Dalcin   snes->ops->view           = 0;
9337596af1SLisandro Dalcin   snes->ops->reset          = 0;
94b79b07cfSJed Brown 
9542f4f86dSBarry Smith   snes->usesksp = PETSC_TRUE;
96*efd4aadfSBarry Smith   snes->usesnpc = PETSC_FALSE;
9742f4f86dSBarry Smith 
984fc747eaSLawrence Mitchell   snes->alwayscomputesfinalresidual = PETSC_FALSE;
994fc747eaSLawrence Mitchell 
100b79b07cfSJed Brown   snes->data = 0;
101b79b07cfSJed Brown   PetscFunctionReturn(0);
102b79b07cfSJed Brown }
103