xref: /petsc/src/snes/impls/ksponly/ksponly.c (revision 8cc058d9cd56c1ccb3be12a47760ddfc446aaffc)
1b79b07cfSJed Brown 
2b45d2f2cSJed Brown #include <petsc-private/snesimpl.h>
3b79b07cfSJed Brown 
4b79b07cfSJed Brown #undef __FUNCT__
5b79b07cfSJed Brown #define __FUNCT__ "SNESSolve_KSPONLY"
6b79b07cfSJed Brown static PetscErrorCode SNESSolve_KSPONLY(SNES snes)
7b79b07cfSJed Brown {
8b79b07cfSJed Brown   PetscErrorCode     ierr;
9b79b07cfSJed Brown   PetscInt           lits;
10b79b07cfSJed Brown   MatStructure       flg = DIFFERENT_NONZERO_PATTERN;
11b79b07cfSJed Brown   Vec                Y,X,F;
12b79b07cfSJed Brown   KSPConvergedReason kspreason;
13b79b07cfSJed Brown 
14b79b07cfSJed Brown   PetscFunctionBegin;
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 
25b79b07cfSJed Brown   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
26b79b07cfSJed Brown   if (snes->domainerror) {
27b79b07cfSJed Brown     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
28b79b07cfSJed Brown     PetscFunctionReturn(0);
29b79b07cfSJed Brown   }
30a5eaddd9SBarry Smith   if (snes->numbermonitors) {
31a5eaddd9SBarry Smith     PetscReal fnorm;
32a5eaddd9SBarry Smith     ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);
33a5eaddd9SBarry Smith     ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
34a5eaddd9SBarry Smith   }
35b79b07cfSJed Brown 
365341784dSBarry Smith   /* Call general purpose update function */
375341784dSBarry Smith   if (snes->ops->update) {
385341784dSBarry Smith     ierr = (*snes->ops->update)(snes, 0);CHKERRQ(ierr);
395341784dSBarry Smith   }
405341784dSBarry Smith 
41b79b07cfSJed Brown   /* Solve J Y = F, where J is Jacobian matrix */
42b79b07cfSJed Brown   ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
43b79b07cfSJed Brown   ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr);
44b79b07cfSJed Brown   ierr = KSPSolve(snes->ksp,F,Y);CHKERRQ(ierr);
45b79b07cfSJed Brown   ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr);
46b79b07cfSJed Brown   if (kspreason < 0 && ++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) {
47b79b07cfSJed Brown     ierr         = PetscInfo2(snes,"iter=%D, number linear solve failures %D greater than current SNES allowed, stopping solve\n",snes->iter,snes->numLinearSolveFailures);CHKERRQ(ierr);
48b79b07cfSJed Brown     snes->reason = SNES_DIVERGED_LINEAR_SOLVE;
491aa26658SKarl Rupp   } else snes->reason = SNES_CONVERGED_ITS;
501aa26658SKarl Rupp 
51b79b07cfSJed Brown   ierr              = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
52b79b07cfSJed Brown   snes->linear_its += lits;
53b79b07cfSJed Brown   ierr              = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr);
54b79b07cfSJed Brown   snes->iter++;
55b79b07cfSJed Brown 
56b79b07cfSJed Brown   /* Take the computed step. */
57b79b07cfSJed Brown   ierr = VecAXPY(X,-1.0,Y);CHKERRQ(ierr);
58a5eaddd9SBarry Smith   if (snes->numbermonitors) {
59a5eaddd9SBarry Smith     PetscReal fnorm;
60a5eaddd9SBarry Smith     ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
61a5eaddd9SBarry Smith     ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);
62a5eaddd9SBarry Smith     ierr = SNESMonitor(snes,1,fnorm);CHKERRQ(ierr);
63a5eaddd9SBarry Smith   }
64b79b07cfSJed Brown   PetscFunctionReturn(0);
65b79b07cfSJed Brown }
66b79b07cfSJed Brown 
67b79b07cfSJed Brown #undef __FUNCT__
68b79b07cfSJed Brown #define __FUNCT__ "SNESSetUp_KSPONLY"
69b79b07cfSJed Brown static PetscErrorCode SNESSetUp_KSPONLY(SNES snes)
70b79b07cfSJed Brown {
716cab3a1bSJed Brown   PetscErrorCode ierr;
72b79b07cfSJed Brown 
73b79b07cfSJed Brown   PetscFunctionBegin;
746cab3a1bSJed Brown   ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);
75b79b07cfSJed Brown   PetscFunctionReturn(0);
76b79b07cfSJed Brown }
77b79b07cfSJed Brown 
78b79b07cfSJed Brown #undef __FUNCT__
79b79b07cfSJed Brown #define __FUNCT__ "SNESDestroy_KSPONLY"
80b79b07cfSJed Brown static PetscErrorCode SNESDestroy_KSPONLY(SNES snes)
81b79b07cfSJed Brown {
82b79b07cfSJed Brown 
83b79b07cfSJed Brown   PetscFunctionBegin;
84b79b07cfSJed Brown   PetscFunctionReturn(0);
85b79b07cfSJed Brown }
86b79b07cfSJed Brown 
87b79b07cfSJed Brown /* -------------------------------------------------------------------------- */
88b79b07cfSJed Brown /*MC
89b79b07cfSJed Brown       SNESKSPONLY - Nonlinear solver that only performs one Newton step and does not compute any norms.
90b79b07cfSJed Brown       The main purpose of this solver is to solve linear problems using the SNES interface, without
91b79b07cfSJed Brown       any additional overhead in the form of vector operations.
92b79b07cfSJed Brown 
93b79b07cfSJed Brown    Level: beginner
94b79b07cfSJed Brown 
9504d7464bSBarry Smith .seealso:  SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR
96b79b07cfSJed Brown M*/
97b79b07cfSJed Brown #undef __FUNCT__
98b79b07cfSJed Brown #define __FUNCT__ "SNESCreate_KSPONLY"
99*8cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_KSPONLY(SNES snes)
100b79b07cfSJed Brown {
101b79b07cfSJed Brown 
102b79b07cfSJed Brown   PetscFunctionBegin;
103b79b07cfSJed Brown   snes->ops->setup          = SNESSetUp_KSPONLY;
104b79b07cfSJed Brown   snes->ops->solve          = SNESSolve_KSPONLY;
105b79b07cfSJed Brown   snes->ops->destroy        = SNESDestroy_KSPONLY;
10637596af1SLisandro Dalcin   snes->ops->setfromoptions = 0;
10737596af1SLisandro Dalcin   snes->ops->view           = 0;
10837596af1SLisandro Dalcin   snes->ops->reset          = 0;
109b79b07cfSJed Brown 
11042f4f86dSBarry Smith   snes->usesksp = PETSC_TRUE;
11142f4f86dSBarry Smith   snes->usespc  = PETSC_FALSE;
11242f4f86dSBarry Smith 
113b79b07cfSJed Brown   snes->data = 0;
114b79b07cfSJed Brown   PetscFunctionReturn(0);
115b79b07cfSJed Brown }
116