Actual source code: ksponly.c
petsc-3.4.5 2014-06-29
2: #include <petsc-private/snesimpl.h>
6: static PetscErrorCode SNESSolve_KSPONLY(SNES snes)
7: {
8: PetscErrorCode ierr;
9: PetscInt lits;
10: MatStructure flg = DIFFERENT_NONZERO_PATTERN;
11: Vec Y,X,F;
12: KSPConvergedReason kspreason;
15: snes->numFailures = 0;
16: snes->numLinearSolveFailures = 0;
17: snes->reason = SNES_CONVERGED_ITERATING;
18: snes->iter = 0;
19: snes->norm = 0.0;
21: X = snes->vec_sol;
22: F = snes->vec_func;
23: Y = snes->vec_sol_update;
25: SNESComputeFunction(snes,X,F);
26: if (snes->domainerror) {
27: snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
28: return(0);
29: }
30: if (snes->numbermonitors) {
31: PetscReal fnorm;
32: VecNorm(F,NORM_2,&fnorm);
33: SNESMonitor(snes,0,fnorm);
34: }
36: /* Call general purpose update function */
37: if (snes->ops->update) {
38: (*snes->ops->update)(snes, 0);
39: }
41: /* Solve J Y = F, where J is Jacobian matrix */
42: SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);
43: KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);
44: KSPSolve(snes->ksp,F,Y);
45: KSPGetConvergedReason(snes->ksp,&kspreason);
46: if (kspreason < 0 && ++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) {
47: PetscInfo2(snes,"iter=%D, number linear solve failures %D greater than current SNES allowed, stopping solve\n",snes->iter,snes->numLinearSolveFailures);
48: snes->reason = SNES_DIVERGED_LINEAR_SOLVE;
49: } else snes->reason = SNES_CONVERGED_ITS;
51: KSPGetIterationNumber(snes->ksp,&lits);
52: snes->linear_its += lits;
53: PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);
54: snes->iter++;
56: /* Take the computed step. */
57: VecAXPY(X,-1.0,Y);
58: if (snes->numbermonitors) {
59: PetscReal fnorm;
60: SNESComputeFunction(snes,X,F);
61: VecNorm(F,NORM_2,&fnorm);
62: SNESMonitor(snes,1,fnorm);
63: }
64: return(0);
65: }
69: static PetscErrorCode SNESSetUp_KSPONLY(SNES snes)
70: {
74: SNESSetUpMatrices(snes);
75: return(0);
76: }
80: static PetscErrorCode SNESDestroy_KSPONLY(SNES snes)
81: {
84: return(0);
85: }
87: /* -------------------------------------------------------------------------- */
88: /*MC
89: SNESKSPONLY - Nonlinear solver that only performs one Newton step and does not compute any norms.
90: The main purpose of this solver is to solve linear problems using the SNES interface, without
91: any additional overhead in the form of vector operations.
93: Level: beginner
95: .seealso: SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR
96: M*/
99: PETSC_EXTERN PetscErrorCode SNESCreate_KSPONLY(SNES snes)
100: {
103: snes->ops->setup = SNESSetUp_KSPONLY;
104: snes->ops->solve = SNESSolve_KSPONLY;
105: snes->ops->destroy = SNESDestroy_KSPONLY;
106: snes->ops->setfromoptions = 0;
107: snes->ops->view = 0;
108: snes->ops->reset = 0;
110: snes->usesksp = PETSC_TRUE;
111: snes->usespc = PETSC_FALSE;
113: snes->data = 0;
114: return(0);
115: }