1 2 #include <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 MatStructure flg = DIFFERENT_NONZERO_PATTERN; 11 Vec Y,X,F; 12 KSPConvergedReason kspreason; 13 14 PetscFunctionBegin; 15 snes->numFailures = 0; 16 snes->numLinearSolveFailures = 0; 17 snes->reason = SNES_CONVERGED_ITERATING; 18 snes->iter = 0; 19 snes->norm = 0.0; 20 21 X = snes->vec_sol; 22 F = snes->vec_func; 23 Y = snes->vec_sol_update; 24 25 ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 26 if (snes->domainerror) { 27 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 28 PetscFunctionReturn(0); 29 } 30 31 /* Solve J Y = F, where J is Jacobian matrix */ 32 ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 33 ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr); 34 ierr = KSPSolve(snes->ksp,F,Y);CHKERRQ(ierr); 35 ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr); 36 if (kspreason < 0 && ++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 37 ierr = PetscInfo2(snes,"iter=%D, number linear solve failures %D greater than current SNES allowed, stopping solve\n",snes->iter,snes->numLinearSolveFailures);CHKERRQ(ierr); 38 snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 39 } else { 40 snes->reason = SNES_CONVERGED_ITS; 41 } 42 ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 43 snes->linear_its += lits; 44 ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr); 45 snes->iter++; 46 47 /* Take the computed step. */ 48 ierr = VecAXPY(X,-1.0,Y);CHKERRQ(ierr); 49 PetscFunctionReturn(0); 50 } 51 52 #undef __FUNCT__ 53 #define __FUNCT__ "SNESSetUp_KSPONLY" 54 static PetscErrorCode SNESSetUp_KSPONLY(SNES snes) 55 { 56 57 PetscFunctionBegin; 58 PetscFunctionReturn(0); 59 } 60 61 #undef __FUNCT__ 62 #define __FUNCT__ "SNESDestroy_KSPONLY" 63 static PetscErrorCode SNESDestroy_KSPONLY(SNES snes) 64 { 65 66 PetscFunctionBegin; 67 PetscFunctionReturn(0); 68 } 69 70 /* -------------------------------------------------------------------------- */ 71 /*MC 72 SNESKSPONLY - Nonlinear solver that only performs one Newton step and does not compute any norms. 73 The main purpose of this solver is to solve linear problems using the SNES interface, without 74 any additional overhead in the form of vector operations. 75 76 Level: beginner 77 78 .seealso: SNESCreate(), SNES, SNESSetType(), SNESLS, SNESTR 79 M*/ 80 EXTERN_C_BEGIN 81 #undef __FUNCT__ 82 #define __FUNCT__ "SNESCreate_KSPONLY" 83 PetscErrorCode SNESCreate_KSPONLY(SNES snes) 84 { 85 86 PetscFunctionBegin; 87 snes->ops->setup = SNESSetUp_KSPONLY; 88 snes->ops->solve = SNESSolve_KSPONLY; 89 snes->ops->destroy = SNESDestroy_KSPONLY; 90 snes->ops->setfromoptions = 0; 91 snes->ops->view = 0; 92 snes->ops->reset = 0; 93 94 snes->usesksp = PETSC_TRUE; 95 snes->usespc = PETSC_FALSE; 96 97 snes->data = 0; 98 PetscFunctionReturn(0); 99 } 100 EXTERN_C_END 101