xref: /petsc/src/snes/impls/ksponly/ksponly.c (revision 07b62357f41014a58be1c7cccdd6f30e8aca3e32)
11ef27442SStefano Zampini #include <petsc/private/snesimpl.h>             /*I   "petscsnes.h"   I*/
2b79b07cfSJed Brown 
31ef27442SStefano Zampini typedef struct {
41ef27442SStefano Zampini   PetscBool transpose_solve;
51ef27442SStefano Zampini } SNES_KSPONLY;
6b79b07cfSJed Brown 
7b79b07cfSJed Brown static PetscErrorCode SNESSolve_KSPONLY(SNES snes)
8b79b07cfSJed Brown {
91ef27442SStefano 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);
41*07b62357SFande Kong 
42*07b62357SFande Kong   SNESCheckJacobianDomainerror(snes);
43*07b62357SFande Kong 
4423ee1639SBarry Smith   ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr);
451ef27442SStefano Zampini   if (ksponly->transpose_solve) {
461ef27442SStefano Zampini     ierr = KSPSolveTranspose(snes->ksp,F,Y);CHKERRQ(ierr);
471ef27442SStefano Zampini   } else {
48b79b07cfSJed Brown     ierr = KSPSolve(snes->ksp,F,Y);CHKERRQ(ierr);
491ef27442SStefano Zampini   }
50422a814eSBarry Smith   snes->reason = SNES_CONVERGED_ITS;
51422a814eSBarry Smith   SNESCheckKSPSolve(snes);
521aa26658SKarl Rupp 
53b79b07cfSJed Brown   ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
54b79b07cfSJed Brown   ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr);
55b79b07cfSJed Brown   snes->iter++;
56b79b07cfSJed Brown 
57b79b07cfSJed Brown   /* Take the computed step. */
58b79b07cfSJed Brown   ierr = VecAXPY(X,-1.0,Y);CHKERRQ(ierr);
59a5eaddd9SBarry Smith   if (snes->numbermonitors) {
60a5eaddd9SBarry Smith     PetscReal fnorm;
61fef909d3SLawrence Mitchell     ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
62a5eaddd9SBarry Smith     ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);
63a5eaddd9SBarry Smith     ierr = SNESMonitor(snes,1,fnorm);CHKERRQ(ierr);
64a5eaddd9SBarry Smith   }
65b79b07cfSJed Brown   PetscFunctionReturn(0);
66b79b07cfSJed Brown }
67b79b07cfSJed Brown 
68b79b07cfSJed Brown static PetscErrorCode SNESSetUp_KSPONLY(SNES snes)
69b79b07cfSJed Brown {
706cab3a1bSJed Brown   PetscErrorCode ierr;
71b79b07cfSJed Brown 
72b79b07cfSJed Brown   PetscFunctionBegin;
736cab3a1bSJed Brown   ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);
74b79b07cfSJed Brown   PetscFunctionReturn(0);
75b79b07cfSJed Brown }
76b79b07cfSJed Brown 
77b79b07cfSJed Brown static PetscErrorCode SNESDestroy_KSPONLY(SNES snes)
78b79b07cfSJed Brown {
791ef27442SStefano Zampini   PetscErrorCode ierr;
80b79b07cfSJed Brown 
81b79b07cfSJed Brown   PetscFunctionBegin;
821ef27442SStefano Zampini   ierr = PetscFree(snes->data);CHKERRQ(ierr);
83b79b07cfSJed Brown   PetscFunctionReturn(0);
84b79b07cfSJed Brown }
85b79b07cfSJed Brown 
86b79b07cfSJed Brown /* -------------------------------------------------------------------------- */
87b79b07cfSJed Brown /*MC
88b79b07cfSJed Brown       SNESKSPONLY - Nonlinear solver that only performs one Newton step and does not compute any norms.
89b79b07cfSJed Brown       The main purpose of this solver is to solve linear problems using the SNES interface, without
90b79b07cfSJed Brown       any additional overhead in the form of vector operations.
91b79b07cfSJed Brown 
92b79b07cfSJed Brown    Level: beginner
93b79b07cfSJed Brown 
9404d7464bSBarry Smith .seealso:  SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR
95b79b07cfSJed Brown M*/
968cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_KSPONLY(SNES snes)
97b79b07cfSJed Brown {
981ef27442SStefano Zampini   SNES_KSPONLY   *ksponly;
991ef27442SStefano Zampini   PetscErrorCode ierr;
100b79b07cfSJed Brown 
101b79b07cfSJed Brown   PetscFunctionBegin;
102b79b07cfSJed Brown   snes->ops->setup          = SNESSetUp_KSPONLY;
103b79b07cfSJed Brown   snes->ops->solve          = SNESSolve_KSPONLY;
104b79b07cfSJed Brown   snes->ops->destroy        = SNESDestroy_KSPONLY;
10537596af1SLisandro Dalcin   snes->ops->setfromoptions = 0;
10637596af1SLisandro Dalcin   snes->ops->view           = 0;
10737596af1SLisandro Dalcin   snes->ops->reset          = 0;
108b79b07cfSJed Brown 
10942f4f86dSBarry Smith   snes->usesksp = PETSC_TRUE;
110efd4aadfSBarry Smith   snes->usesnpc = PETSC_FALSE;
11142f4f86dSBarry Smith 
1124fc747eaSLawrence Mitchell   snes->alwayscomputesfinalresidual = PETSC_FALSE;
1134fc747eaSLawrence Mitchell 
1141ef27442SStefano Zampini   ierr = PetscNewLog(snes,&ksponly);CHKERRQ(ierr);
1151ef27442SStefano Zampini   snes->data = (void*)ksponly;
1161ef27442SStefano Zampini   PetscFunctionReturn(0);
1171ef27442SStefano Zampini }
1181ef27442SStefano Zampini 
1191ef27442SStefano Zampini /*MC
1201ef27442SStefano Zampini       SNESKSPTRANSPOSEONLY - Nonlinear solver that only performs one Newton step and does not compute any norms.
1211ef27442SStefano Zampini       The main purpose of this solver is to solve transposed linear problems using the SNES interface, without
1221ef27442SStefano Zampini       any additional overhead in the form of vector operations within adjoint solvers.
1231ef27442SStefano Zampini 
1241ef27442SStefano Zampini    Level: beginner
1251ef27442SStefano Zampini 
1261ef27442SStefano Zampini .seealso:  SNESCreate(), SNES, SNESSetType(), SNESKSPTRANSPOSEONLY, SNESNEWTONLS, SNESNEWTONTR
1271ef27442SStefano Zampini M*/
1281ef27442SStefano Zampini PETSC_EXTERN PetscErrorCode SNESCreate_KSPTRANSPOSEONLY(SNES snes)
1291ef27442SStefano Zampini {
1301ef27442SStefano Zampini   SNES_KSPONLY   *kspo;
1311ef27442SStefano Zampini   PetscErrorCode ierr;
1321ef27442SStefano Zampini 
1331ef27442SStefano Zampini   PetscFunctionBegin;
1341ef27442SStefano Zampini   ierr = SNESCreate_KSPONLY(snes);CHKERRQ(ierr);
1351ef27442SStefano Zampini   kspo = (SNES_KSPONLY*)snes->data;
1361ef27442SStefano Zampini   kspo->transpose_solve = PETSC_TRUE;
137b79b07cfSJed Brown   PetscFunctionReturn(0);
138b79b07cfSJed Brown }
139