Actual source code: ksponly.c

petsc-3.5.4 2015-05-23
Report Typos and Errors
  2: #include <petsc-private/snesimpl.h>

  6: static PetscErrorCode SNESSolve_KSPONLY(SNES snes)
  7: {
  8:   PetscErrorCode     ierr;
  9:   PetscInt           lits;
 10:   Vec                Y,X,F;
 11:   KSPConvergedReason kspreason;

 14:   snes->numFailures            = 0;
 15:   snes->numLinearSolveFailures = 0;
 16:   snes->reason                 = SNES_CONVERGED_ITERATING;
 17:   snes->iter                   = 0;
 18:   snes->norm                   = 0.0;

 20:   X = snes->vec_sol;
 21:   F = snes->vec_func;
 22:   Y = snes->vec_sol_update;

 24:   SNESComputeFunction(snes,X,F);
 25:   if (snes->domainerror) {
 26:     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
 27:     return(0);
 28:   }
 29:   if (snes->numbermonitors) {
 30:     PetscReal fnorm;
 31:     VecNorm(F,NORM_2,&fnorm);
 32:     SNESMonitor(snes,0,fnorm);
 33:   }

 35:   /* Call general purpose update function */
 36:   if (snes->ops->update) {
 37:     (*snes->ops->update)(snes, 0);
 38:   }

 40:   /* Solve J Y = F, where J is Jacobian matrix */
 41:   SNESComputeJacobian(snes,X,snes->jacobian,snes->jacobian_pre);
 42:   KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre);
 43:   KSPSolve(snes->ksp,F,Y);
 44:   KSPGetConvergedReason(snes->ksp,&kspreason);
 45:   if (kspreason < 0 && ++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) {
 46:     PetscInfo2(snes,"iter=%D, number linear solve failures %D greater than current SNES allowed, stopping solve\n",snes->iter,snes->numLinearSolveFailures);
 47:     snes->reason = SNES_DIVERGED_LINEAR_SOLVE;
 48:   } else snes->reason = SNES_CONVERGED_ITS;

 50:   KSPGetIterationNumber(snes->ksp,&lits);
 51:   snes->linear_its += lits;
 52:   PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);
 53:   snes->iter++;

 55:   /* Take the computed step. */
 56:   VecAXPY(X,-1.0,Y);
 57:   if (snes->numbermonitors) {
 58:     PetscReal fnorm;
 59:     SNESComputeFunction(snes,X,F);
 60:     VecNorm(F,NORM_2,&fnorm);
 61:     SNESMonitor(snes,1,fnorm);
 62:   }
 63:   return(0);
 64: }

 68: static PetscErrorCode SNESSetUp_KSPONLY(SNES snes)
 69: {

 73:   SNESSetUpMatrices(snes);
 74:   return(0);
 75: }

 79: static PetscErrorCode SNESDestroy_KSPONLY(SNES snes)
 80: {

 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*/
 98: PETSC_EXTERN PetscErrorCode SNESCreate_KSPONLY(SNES snes)
 99: {

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->usespc  = PETSC_FALSE;

112:   snes->data = 0;
113:   return(0);
114: }