xref: /petsc/src/snes/impls/ksponly/ksponly.c (revision 1ef274424f215a09ebb8fb75edff2501b0a2cb92)
1*1ef27442SStefano Zampini #include <petsc/private/snesimpl.h>             /*I   "petscsnes.h"   I*/
2b79b07cfSJed Brown 
3*1ef27442SStefano Zampini typedef struct {
4*1ef27442SStefano Zampini   PetscBool transpose_solve;
5*1ef27442SStefano Zampini } SNES_KSPONLY;
6b79b07cfSJed Brown 
7b79b07cfSJed Brown static PetscErrorCode SNESSolve_KSPONLY(SNES snes)
8b79b07cfSJed Brown {
9*1ef27442SStefano Zampini   SNES_KSPONLY   *ksponly = (SNES_KSPONLY*)snes->data;
10b79b07cfSJed Brown   PetscErrorCode ierr;
11b79b07cfSJed Brown   PetscInt       lits;
12b79b07cfSJed Brown   Vec            Y,X,F;
13b79b07cfSJed Brown 
14b79b07cfSJed Brown   PetscFunctionBegin;
156c4ed002SBarry 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);
16c579b300SPatrick Farrell 
17b79b07cfSJed Brown   snes->numFailures            = 0;
18b79b07cfSJed Brown   snes->numLinearSolveFailures = 0;
19b79b07cfSJed Brown   snes->reason                 = SNES_CONVERGED_ITERATING;
20b79b07cfSJed Brown   snes->iter                   = 0;
21b79b07cfSJed Brown   snes->norm                   = 0.0;
22b79b07cfSJed Brown 
23b79b07cfSJed Brown   X = snes->vec_sol;
24b79b07cfSJed Brown   F = snes->vec_func;
25b79b07cfSJed Brown   Y = snes->vec_sol_update;
26b79b07cfSJed Brown 
27b79b07cfSJed Brown   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
28a5eaddd9SBarry Smith   if (snes->numbermonitors) {
29a5eaddd9SBarry Smith     PetscReal fnorm;
30a5eaddd9SBarry Smith     ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);
31a5eaddd9SBarry Smith     ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
32a5eaddd9SBarry Smith   }
33b79b07cfSJed Brown 
345341784dSBarry Smith   /* Call general purpose update function */
355341784dSBarry Smith   if (snes->ops->update) {
365341784dSBarry Smith     ierr = (*snes->ops->update)(snes, 0);CHKERRQ(ierr);
375341784dSBarry Smith   }
385341784dSBarry Smith 
39b79b07cfSJed Brown   /* Solve J Y = F, where J is Jacobian matrix */
40d1e9a80fSBarry Smith   ierr = SNESComputeJacobian(snes,X,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr);
4123ee1639SBarry Smith   ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr);
42*1ef27442SStefano Zampini   if (ksponly->transpose_solve) {
43*1ef27442SStefano Zampini     ierr = KSPSolveTranspose(snes->ksp,F,Y);CHKERRQ(ierr);
44*1ef27442SStefano Zampini   } else {
45b79b07cfSJed Brown     ierr = KSPSolve(snes->ksp,F,Y);CHKERRQ(ierr);
46*1ef27442SStefano Zampini   }
47422a814eSBarry Smith   snes->reason = SNES_CONVERGED_ITS;
48422a814eSBarry Smith   SNESCheckKSPSolve(snes);
491aa26658SKarl Rupp 
50b79b07cfSJed Brown   ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
51b79b07cfSJed Brown   snes->linear_its += lits;
52b79b07cfSJed Brown   ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr);
53b79b07cfSJed Brown   snes->iter++;
54b79b07cfSJed Brown 
55b79b07cfSJed Brown   /* Take the computed step. */
56b79b07cfSJed Brown   ierr = VecAXPY(X,-1.0,Y);CHKERRQ(ierr);
57a5eaddd9SBarry Smith   if (snes->numbermonitors) {
58a5eaddd9SBarry Smith     PetscReal fnorm;
59fef909d3SLawrence Mitchell     ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
60a5eaddd9SBarry Smith     ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);
61a5eaddd9SBarry Smith     ierr = SNESMonitor(snes,1,fnorm);CHKERRQ(ierr);
62a5eaddd9SBarry Smith   }
63b79b07cfSJed Brown   PetscFunctionReturn(0);
64b79b07cfSJed Brown }
65b79b07cfSJed Brown 
66b79b07cfSJed Brown static PetscErrorCode SNESSetUp_KSPONLY(SNES snes)
67b79b07cfSJed Brown {
686cab3a1bSJed Brown   PetscErrorCode ierr;
69b79b07cfSJed Brown 
70b79b07cfSJed Brown   PetscFunctionBegin;
716cab3a1bSJed Brown   ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);
72b79b07cfSJed Brown   PetscFunctionReturn(0);
73b79b07cfSJed Brown }
74b79b07cfSJed Brown 
75b79b07cfSJed Brown static PetscErrorCode SNESDestroy_KSPONLY(SNES snes)
76b79b07cfSJed Brown {
77*1ef27442SStefano Zampini   PetscErrorCode ierr;
78b79b07cfSJed Brown 
79b79b07cfSJed Brown   PetscFunctionBegin;
80*1ef27442SStefano Zampini   ierr = PetscFree(snes->data);CHKERRQ(ierr);
81b79b07cfSJed Brown   PetscFunctionReturn(0);
82b79b07cfSJed Brown }
83b79b07cfSJed Brown 
84b79b07cfSJed Brown /* -------------------------------------------------------------------------- */
85b79b07cfSJed Brown /*MC
86b79b07cfSJed Brown       SNESKSPONLY - Nonlinear solver that only performs one Newton step and does not compute any norms.
87b79b07cfSJed Brown       The main purpose of this solver is to solve linear problems using the SNES interface, without
88b79b07cfSJed Brown       any additional overhead in the form of vector operations.
89b79b07cfSJed Brown 
90b79b07cfSJed Brown    Level: beginner
91b79b07cfSJed Brown 
9204d7464bSBarry Smith .seealso:  SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR
93b79b07cfSJed Brown M*/
948cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_KSPONLY(SNES snes)
95b79b07cfSJed Brown {
96*1ef27442SStefano Zampini   SNES_KSPONLY   *ksponly;
97*1ef27442SStefano Zampini   PetscErrorCode ierr;
98b79b07cfSJed Brown 
99b79b07cfSJed Brown   PetscFunctionBegin;
100b79b07cfSJed Brown   snes->ops->setup          = SNESSetUp_KSPONLY;
101b79b07cfSJed Brown   snes->ops->solve          = SNESSolve_KSPONLY;
102b79b07cfSJed Brown   snes->ops->destroy        = SNESDestroy_KSPONLY;
10337596af1SLisandro Dalcin   snes->ops->setfromoptions = 0;
10437596af1SLisandro Dalcin   snes->ops->view           = 0;
10537596af1SLisandro Dalcin   snes->ops->reset          = 0;
106b79b07cfSJed Brown 
10742f4f86dSBarry Smith   snes->usesksp = PETSC_TRUE;
108efd4aadfSBarry Smith   snes->usesnpc = PETSC_FALSE;
10942f4f86dSBarry Smith 
1104fc747eaSLawrence Mitchell   snes->alwayscomputesfinalresidual = PETSC_FALSE;
1114fc747eaSLawrence Mitchell 
112*1ef27442SStefano Zampini   ierr = PetscNewLog(snes,&ksponly);CHKERRQ(ierr);
113*1ef27442SStefano Zampini   snes->data = (void*)ksponly;
114*1ef27442SStefano Zampini   PetscFunctionReturn(0);
115*1ef27442SStefano Zampini }
116*1ef27442SStefano Zampini 
117*1ef27442SStefano Zampini /*MC
118*1ef27442SStefano Zampini       SNESKSPTRANSPOSEONLY - Nonlinear solver that only performs one Newton step and does not compute any norms.
119*1ef27442SStefano Zampini       The main purpose of this solver is to solve transposed linear problems using the SNES interface, without
120*1ef27442SStefano Zampini       any additional overhead in the form of vector operations within adjoint solvers.
121*1ef27442SStefano Zampini 
122*1ef27442SStefano Zampini    Level: beginner
123*1ef27442SStefano Zampini 
124*1ef27442SStefano Zampini .seealso:  SNESCreate(), SNES, SNESSetType(), SNESKSPTRANSPOSEONLY, SNESNEWTONLS, SNESNEWTONTR
125*1ef27442SStefano Zampini M*/
126*1ef27442SStefano Zampini PETSC_EXTERN PetscErrorCode SNESCreate_KSPTRANSPOSEONLY(SNES snes)
127*1ef27442SStefano Zampini {
128*1ef27442SStefano Zampini   SNES_KSPONLY   *kspo;
129*1ef27442SStefano Zampini   PetscErrorCode ierr;
130*1ef27442SStefano Zampini 
131*1ef27442SStefano Zampini   PetscFunctionBegin;
132*1ef27442SStefano Zampini   ierr = SNESCreate_KSPONLY(snes);CHKERRQ(ierr);
133*1ef27442SStefano Zampini   kspo = (SNES_KSPONLY*)snes->data;
134*1ef27442SStefano Zampini   kspo->transpose_solve = PETSC_TRUE;
135b79b07cfSJed Brown   PetscFunctionReturn(0);
136b79b07cfSJed Brown }
137