Actual source code: ksponly.c
petsc-3.12.5 2020-03-29
1: #include <petsc/private/snesimpl.h>
3: typedef struct {
4: PetscBool transpose_solve;
5: } SNES_KSPONLY;
7: static PetscErrorCode SNESSolve_KSPONLY(SNES snes)
8: {
9: SNES_KSPONLY *ksponly = (SNES_KSPONLY*)snes->data;
11: PetscInt lits;
12: Vec Y,X,F;
15: 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);
17: snes->numFailures = 0;
18: snes->numLinearSolveFailures = 0;
19: snes->reason = SNES_CONVERGED_ITERATING;
20: snes->iter = 0;
21: snes->norm = 0.0;
23: X = snes->vec_sol;
24: F = snes->vec_func;
25: Y = snes->vec_sol_update;
27: SNESComputeFunction(snes,X,F);
28: if (snes->numbermonitors) {
29: PetscReal fnorm;
30: VecNorm(F,NORM_2,&fnorm);
31: SNESMonitor(snes,0,fnorm);
32: }
34: /* Call general purpose update function */
35: if (snes->ops->update) {
36: (*snes->ops->update)(snes, 0);
37: }
39: /* Solve J Y = F, where J is Jacobian matrix */
40: SNESComputeJacobian(snes,X,snes->jacobian,snes->jacobian_pre);
42: SNESCheckJacobianDomainerror(snes);
44: KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre);
45: if (ksponly->transpose_solve) {
46: KSPSolveTranspose(snes->ksp,F,Y);
47: } else {
48: KSPSolve(snes->ksp,F,Y);
49: }
50: snes->reason = SNES_CONVERGED_ITS;
51: SNESCheckKSPSolve(snes);
53: KSPGetIterationNumber(snes->ksp,&lits);
54: PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);
55: snes->iter++;
57: /* Take the computed step. */
58: VecAXPY(X,-1.0,Y);
59: if (snes->numbermonitors) {
60: PetscReal fnorm;
61: SNESComputeFunction(snes,X,F);
62: VecNorm(F,NORM_2,&fnorm);
63: SNESMonitor(snes,1,fnorm);
64: }
65: return(0);
66: }
68: static PetscErrorCode SNESSetUp_KSPONLY(SNES snes)
69: {
73: SNESSetUpMatrices(snes);
74: return(0);
75: }
77: static PetscErrorCode SNESDestroy_KSPONLY(SNES snes)
78: {
82: PetscFree(snes->data);
83: return(0);
84: }
86: /* -------------------------------------------------------------------------- */
87: /*MC
88: SNESKSPONLY - Nonlinear solver that only performs one Newton step and does not compute any norms.
89: The main purpose of this solver is to solve linear problems using the SNES interface, without
90: any additional overhead in the form of vector operations.
92: Level: beginner
94: .seealso: SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR
95: M*/
96: PETSC_EXTERN PetscErrorCode SNESCreate_KSPONLY(SNES snes)
97: {
98: SNES_KSPONLY *ksponly;
102: snes->ops->setup = SNESSetUp_KSPONLY;
103: snes->ops->solve = SNESSolve_KSPONLY;
104: snes->ops->destroy = SNESDestroy_KSPONLY;
105: snes->ops->setfromoptions = 0;
106: snes->ops->view = 0;
107: snes->ops->reset = 0;
109: snes->usesksp = PETSC_TRUE;
110: snes->usesnpc = PETSC_FALSE;
112: snes->alwayscomputesfinalresidual = PETSC_FALSE;
114: PetscNewLog(snes,&ksponly);
115: snes->data = (void*)ksponly;
116: return(0);
117: }
119: /*MC
120: SNESKSPTRANSPOSEONLY - Nonlinear solver that only performs one Newton step and does not compute any norms.
121: The main purpose of this solver is to solve transposed linear problems using the SNES interface, without
122: any additional overhead in the form of vector operations within adjoint solvers.
124: Level: beginner
126: .seealso: SNESCreate(), SNES, SNESSetType(), SNESKSPTRANSPOSEONLY, SNESNEWTONLS, SNESNEWTONTR
127: M*/
128: PETSC_EXTERN PetscErrorCode SNESCreate_KSPTRANSPOSEONLY(SNES snes)
129: {
130: SNES_KSPONLY *kspo;
134: SNESCreate_KSPONLY(snes);
135: kspo = (SNES_KSPONLY*)snes->data;
136: kspo->transpose_solve = PETSC_TRUE;
137: return(0);
138: }