Actual source code: ksponly.c

petsc-3.3-p7 2013-05-11
  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:   }

 31:   /* Solve J Y = F, where J is Jacobian matrix */
 32:   SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);
 33:   KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);
 34:   KSPSolve(snes->ksp,F,Y);
 35:   KSPGetConvergedReason(snes->ksp,&kspreason);
 36:   if (kspreason < 0 && ++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) {
 37:     PetscInfo2(snes,"iter=%D, number linear solve failures %D greater than current SNES allowed, stopping solve\n",snes->iter,snes->numLinearSolveFailures);
 38:     snes->reason = SNES_DIVERGED_LINEAR_SOLVE;
 39:   } else {
 40:     snes->reason = SNES_CONVERGED_ITS;
 41:   }
 42:   KSPGetIterationNumber(snes->ksp,&lits);
 43:   snes->linear_its += lits;
 44:   PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);
 45:   snes->iter++;

 47:   /* Take the computed step. */
 48:   VecAXPY(X,-1.0,Y);
 49:   return(0);
 50: }

 54: static PetscErrorCode SNESSetUp_KSPONLY(SNES snes)
 55: {

 59:   SNESSetUpMatrices(snes);
 60:   return(0);
 61: }

 65: static PetscErrorCode SNESDestroy_KSPONLY(SNES snes)
 66: {

 69:   return(0);
 70: }

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

 78:    Level: beginner

 80: .seealso:  SNESCreate(), SNES, SNESSetType(), SNESLS, SNESTR
 81: M*/
 82: EXTERN_C_BEGIN
 85: PetscErrorCode  SNESCreate_KSPONLY(SNES snes)
 86: {

 89:   snes->ops->setup           = SNESSetUp_KSPONLY;
 90:   snes->ops->solve           = SNESSolve_KSPONLY;
 91:   snes->ops->destroy         = SNESDestroy_KSPONLY;
 92:   snes->ops->setfromoptions  = 0;
 93:   snes->ops->view            = 0;
 94:   snes->ops->reset           = 0;

 96:   snes->usesksp         = PETSC_TRUE;
 97:   snes->usespc          = PETSC_FALSE;

 99:   snes->data = 0;
100:   return(0);
101: }
102: EXTERN_C_END