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