Actual source code: ksponly.c

petsc-3.7.5 2017-01-01
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;

 13:   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);

 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->numbermonitors) {
 27:     PetscReal fnorm;
 28:     VecNorm(F,NORM_2,&fnorm);
 29:     SNESMonitor(snes,0,fnorm);
 30:   }

 32:   /* Call general purpose update function */
 33:   if (snes->ops->update) {
 34:     (*snes->ops->update)(snes, 0);
 35:   }

 37:   /* Solve J Y = F, where J is Jacobian matrix */
 38:   SNESComputeJacobian(snes,X,snes->jacobian,snes->jacobian_pre);
 39:   KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre);
 40:   KSPSolve(snes->ksp,F,Y);
 41:   snes->reason = SNES_CONVERGED_ITS;
 42:   SNESCheckKSPSolve(snes);

 44:   KSPGetIterationNumber(snes->ksp,&lits);
 45:   snes->linear_its += lits;
 46:   PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);
 47:   snes->iter++;

 49:   /* Take the computed step. */
 50:   VecAXPY(X,-1.0,Y);
 51:   SNESComputeFunction(snes,X,F);
 52:   if (snes->numbermonitors) {
 53:     PetscReal fnorm;
 54:     VecNorm(F,NORM_2,&fnorm);
 55:     SNESMonitor(snes,1,fnorm);
 56:   }
 57:   return(0);
 58: }

 62: static PetscErrorCode SNESSetUp_KSPONLY(SNES snes)
 63: {

 67:   SNESSetUpMatrices(snes);
 68:   return(0);
 69: }

 73: static PetscErrorCode SNESDestroy_KSPONLY(SNES snes)
 74: {

 77:   return(0);
 78: }

 80: /* -------------------------------------------------------------------------- */
 81: /*MC
 82:       SNESKSPONLY - Nonlinear solver that only performs one Newton step and does not compute any norms.
 83:       The main purpose of this solver is to solve linear problems using the SNES interface, without
 84:       any additional overhead in the form of vector operations.

 86:    Level: beginner

 88: .seealso:  SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR
 89: M*/
 92: PETSC_EXTERN PetscErrorCode SNESCreate_KSPONLY(SNES snes)
 93: {

 96:   snes->ops->setup          = SNESSetUp_KSPONLY;
 97:   snes->ops->solve          = SNESSolve_KSPONLY;
 98:   snes->ops->destroy        = SNESDestroy_KSPONLY;
 99:   snes->ops->setfromoptions = 0;
100:   snes->ops->view           = 0;
101:   snes->ops->reset          = 0;

103:   snes->usesksp = PETSC_TRUE;
104:   snes->usespc  = PETSC_FALSE;

106:   snes->data = 0;
107:   return(0);
108: }