xref: /petsc/src/snes/impls/ksponly/ksponly.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
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 
7*9371c9d4SSatish Balay static PetscErrorCode SNESSolve_KSPONLY(SNES snes) {
81ef27442SStefano Zampini   SNES_KSPONLY *ksponly = (SNES_KSPONLY *)snes->data;
9b79b07cfSJed Brown   PetscInt      lits;
10b79b07cfSJed Brown   Vec           Y, X, F;
11b79b07cfSJed Brown 
12b79b07cfSJed Brown   PetscFunctionBegin;
130b121fc5SBarry Smith   PetscCheck(!snes->xl && !snes->xu && !snes->ops->computevariablebounds, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name);
14c579b300SPatrick Farrell 
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 
25451891c6SStefano Zampini   if (!snes->vec_func_init_set) {
269566063dSJacob Faibussowitsch     PetscCall(SNESComputeFunction(snes, X, F));
27451891c6SStefano Zampini   } else snes->vec_func_init_set = PETSC_FALSE;
28451891c6SStefano Zampini 
29a5eaddd9SBarry Smith   if (snes->numbermonitors) {
30a5eaddd9SBarry Smith     PetscReal fnorm;
319566063dSJacob Faibussowitsch     PetscCall(VecNorm(F, NORM_2, &fnorm));
329566063dSJacob Faibussowitsch     PetscCall(SNESMonitor(snes, 0, fnorm));
33a5eaddd9SBarry Smith   }
34b79b07cfSJed Brown 
355341784dSBarry Smith   /* Call general purpose update function */
36dbbe0bcdSBarry Smith   PetscTryTypeMethod(snes, update, 0);
375341784dSBarry Smith 
38b79b07cfSJed Brown   /* Solve J Y = F, where J is Jacobian matrix */
399566063dSJacob Faibussowitsch   PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre));
4007b62357SFande Kong 
4107b62357SFande Kong   SNESCheckJacobianDomainerror(snes);
4207b62357SFande Kong 
439566063dSJacob Faibussowitsch   PetscCall(KSPSetOperators(snes->ksp, snes->jacobian, snes->jacobian_pre));
441ef27442SStefano Zampini   if (ksponly->transpose_solve) {
459566063dSJacob Faibussowitsch     PetscCall(KSPSolveTranspose(snes->ksp, F, Y));
461ef27442SStefano Zampini   } else {
479566063dSJacob Faibussowitsch     PetscCall(KSPSolve(snes->ksp, F, Y));
481ef27442SStefano Zampini   }
49422a814eSBarry Smith   snes->reason = SNES_CONVERGED_ITS;
50422a814eSBarry Smith   SNESCheckKSPSolve(snes);
511aa26658SKarl Rupp 
529566063dSJacob Faibussowitsch   PetscCall(KSPGetIterationNumber(snes->ksp, &lits));
5363a3b9bcSJacob Faibussowitsch   PetscCall(PetscInfo(snes, "iter=%" PetscInt_FMT ", linear solve iterations=%" PetscInt_FMT "\n", snes->iter, lits));
54b79b07cfSJed Brown   snes->iter++;
55b79b07cfSJed Brown 
56b79b07cfSJed Brown   /* Take the computed step. */
579566063dSJacob Faibussowitsch   PetscCall(VecAXPY(X, -1.0, Y));
58a5eaddd9SBarry Smith   if (snes->numbermonitors) {
59a5eaddd9SBarry Smith     PetscReal fnorm;
609566063dSJacob Faibussowitsch     PetscCall(SNESComputeFunction(snes, X, F));
619566063dSJacob Faibussowitsch     PetscCall(VecNorm(F, NORM_2, &fnorm));
629566063dSJacob Faibussowitsch     PetscCall(SNESMonitor(snes, 1, fnorm));
63a5eaddd9SBarry Smith   }
64b79b07cfSJed Brown   PetscFunctionReturn(0);
65b79b07cfSJed Brown }
66b79b07cfSJed Brown 
67*9371c9d4SSatish Balay static PetscErrorCode SNESSetUp_KSPONLY(SNES snes) {
68b79b07cfSJed Brown   PetscFunctionBegin;
699566063dSJacob Faibussowitsch   PetscCall(SNESSetUpMatrices(snes));
70b79b07cfSJed Brown   PetscFunctionReturn(0);
71b79b07cfSJed Brown }
72b79b07cfSJed Brown 
73*9371c9d4SSatish Balay static PetscErrorCode SNESDestroy_KSPONLY(SNES snes) {
74b79b07cfSJed Brown   PetscFunctionBegin;
759566063dSJacob Faibussowitsch   PetscCall(PetscFree(snes->data));
76b79b07cfSJed Brown   PetscFunctionReturn(0);
77b79b07cfSJed Brown }
78b79b07cfSJed Brown 
79b79b07cfSJed Brown /* -------------------------------------------------------------------------- */
80b79b07cfSJed Brown /*MC
81b79b07cfSJed Brown       SNESKSPONLY - Nonlinear solver that only performs one Newton step and does not compute any norms.
82b79b07cfSJed Brown       The main purpose of this solver is to solve linear problems using the SNES interface, without
83b79b07cfSJed Brown       any additional overhead in the form of vector operations.
84b79b07cfSJed Brown 
85b79b07cfSJed Brown    Level: beginner
86b79b07cfSJed Brown 
87db781477SPatrick Sanan .seealso: `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONLS`, `SNESNEWTONTR`
88b79b07cfSJed Brown M*/
89*9371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode SNESCreate_KSPONLY(SNES snes) {
901ef27442SStefano Zampini   SNES_KSPONLY *ksponly;
91b79b07cfSJed Brown 
92b79b07cfSJed Brown   PetscFunctionBegin;
93b79b07cfSJed Brown   snes->ops->setup          = SNESSetUp_KSPONLY;
94b79b07cfSJed Brown   snes->ops->solve          = SNESSolve_KSPONLY;
95b79b07cfSJed Brown   snes->ops->destroy        = SNESDestroy_KSPONLY;
969e5d0892SLisandro Dalcin   snes->ops->setfromoptions = NULL;
979e5d0892SLisandro Dalcin   snes->ops->view           = NULL;
989e5d0892SLisandro Dalcin   snes->ops->reset          = NULL;
99b79b07cfSJed Brown 
10042f4f86dSBarry Smith   snes->usesksp = PETSC_TRUE;
101efd4aadfSBarry Smith   snes->usesnpc = PETSC_FALSE;
10242f4f86dSBarry Smith 
1034fc747eaSLawrence Mitchell   snes->alwayscomputesfinalresidual = PETSC_FALSE;
1044fc747eaSLawrence Mitchell 
1059566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(snes, &ksponly));
1061ef27442SStefano Zampini   snes->data = (void *)ksponly;
1071ef27442SStefano Zampini   PetscFunctionReturn(0);
1081ef27442SStefano Zampini }
1091ef27442SStefano Zampini 
1101ef27442SStefano Zampini /*MC
1111ef27442SStefano Zampini       SNESKSPTRANSPOSEONLY - Nonlinear solver that only performs one Newton step and does not compute any norms.
1121ef27442SStefano Zampini       The main purpose of this solver is to solve transposed linear problems using the SNES interface, without
1131ef27442SStefano Zampini       any additional overhead in the form of vector operations within adjoint solvers.
1141ef27442SStefano Zampini 
1151ef27442SStefano Zampini    Level: beginner
1161ef27442SStefano Zampini 
117db781477SPatrick Sanan .seealso: `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESKSPTRANSPOSEONLY`, `SNESNEWTONLS`, `SNESNEWTONTR`
1181ef27442SStefano Zampini M*/
119*9371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode SNESCreate_KSPTRANSPOSEONLY(SNES snes) {
1201ef27442SStefano Zampini   SNES_KSPONLY *kspo;
1211ef27442SStefano Zampini 
1221ef27442SStefano Zampini   PetscFunctionBegin;
1239566063dSJacob Faibussowitsch   PetscCall(SNESCreate_KSPONLY(snes));
1249566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)snes, SNESKSPTRANSPOSEONLY));
1251ef27442SStefano Zampini   kspo                  = (SNES_KSPONLY *)snes->data;
1261ef27442SStefano Zampini   kspo->transpose_solve = PETSC_TRUE;
127b79b07cfSJed Brown   PetscFunctionReturn(0);
128b79b07cfSJed Brown }
129