1 2 #include <petsc-private/snesimpl.h> 3 4 #undef __FUNCT__ 5 #define __FUNCT__ "SNESSolve_KSPONLY" 6 static PetscErrorCode SNESSolve_KSPONLY(SNES snes) 7 { 8 PetscErrorCode ierr; 9 PetscInt lits; 10 Vec Y,X,F; 11 KSPConvergedReason kspreason; 12 13 PetscFunctionBegin; 14 snes->numFailures = 0; 15 snes->numLinearSolveFailures = 0; 16 snes->reason = SNES_CONVERGED_ITERATING; 17 snes->iter = 0; 18 snes->norm = 0.0; 19 20 X = snes->vec_sol; 21 F = snes->vec_func; 22 Y = snes->vec_sol_update; 23 24 ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 25 if (snes->domainerror) { 26 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 27 PetscFunctionReturn(0); 28 } 29 if (snes->numbermonitors) { 30 PetscReal fnorm; 31 ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); 32 ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 33 } 34 35 /* Call general purpose update function */ 36 if (snes->ops->update) { 37 ierr = (*snes->ops->update)(snes, 0);CHKERRQ(ierr); 38 } 39 40 /* Solve J Y = F, where J is Jacobian matrix */ 41 ierr = SNESComputeJacobian(snes,X,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr); 42 ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr); 43 ierr = KSPSolve(snes->ksp,F,Y);CHKERRQ(ierr); 44 ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr); 45 if (kspreason < 0 && ++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 46 ierr = PetscInfo2(snes,"iter=%D, number linear solve failures %D greater than current SNES allowed, stopping solve\n",snes->iter,snes->numLinearSolveFailures);CHKERRQ(ierr); 47 snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 48 } else snes->reason = SNES_CONVERGED_ITS; 49 50 ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 51 snes->linear_its += lits; 52 ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr); 53 snes->iter++; 54 55 /* Take the computed step. */ 56 ierr = VecAXPY(X,-1.0,Y);CHKERRQ(ierr); 57 if (snes->numbermonitors) { 58 PetscReal fnorm; 59 ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 60 ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); 61 ierr = SNESMonitor(snes,1,fnorm);CHKERRQ(ierr); 62 } 63 PetscFunctionReturn(0); 64 } 65 66 #undef __FUNCT__ 67 #define __FUNCT__ "SNESSetUp_KSPONLY" 68 static PetscErrorCode SNESSetUp_KSPONLY(SNES snes) 69 { 70 PetscErrorCode ierr; 71 72 PetscFunctionBegin; 73 ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr); 74 PetscFunctionReturn(0); 75 } 76 77 #undef __FUNCT__ 78 #define __FUNCT__ "SNESDestroy_KSPONLY" 79 static PetscErrorCode SNESDestroy_KSPONLY(SNES snes) 80 { 81 82 PetscFunctionBegin; 83 PetscFunctionReturn(0); 84 } 85 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. 91 92 Level: beginner 93 94 .seealso: SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR 95 M*/ 96 #undef __FUNCT__ 97 #define __FUNCT__ "SNESCreate_KSPONLY" 98 PETSC_EXTERN PetscErrorCode SNESCreate_KSPONLY(SNES snes) 99 { 100 101 PetscFunctionBegin; 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; 108 109 snes->usesksp = PETSC_TRUE; 110 snes->usespc = PETSC_FALSE; 111 112 snes->data = 0; 113 PetscFunctionReturn(0); 114 } 115