1 #include <petsc/private/snesimpl.h> /*I "petscsnes.h" I*/ 2 3 typedef struct { 4 PetscBool transpose_solve; 5 } SNES_KSPONLY; 6 7 static PetscErrorCode SNESSolve_KSPONLY(SNES snes) 8 { 9 SNES_KSPONLY *ksponly = (SNES_KSPONLY*)snes->data; 10 PetscErrorCode ierr; 11 PetscInt lits; 12 Vec Y,X,F; 13 14 PetscFunctionBegin; 15 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); 16 17 snes->numFailures = 0; 18 snes->numLinearSolveFailures = 0; 19 snes->reason = SNES_CONVERGED_ITERATING; 20 snes->iter = 0; 21 snes->norm = 0.0; 22 23 X = snes->vec_sol; 24 F = snes->vec_func; 25 Y = snes->vec_sol_update; 26 27 ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 28 if (snes->numbermonitors) { 29 PetscReal fnorm; 30 ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); 31 ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 32 } 33 34 /* Call general purpose update function */ 35 if (snes->ops->update) { 36 ierr = (*snes->ops->update)(snes, 0);CHKERRQ(ierr); 37 } 38 39 /* Solve J Y = F, where J is Jacobian matrix */ 40 ierr = SNESComputeJacobian(snes,X,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr); 41 ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr); 42 if (ksponly->transpose_solve) { 43 ierr = KSPSolveTranspose(snes->ksp,F,Y);CHKERRQ(ierr); 44 } else { 45 ierr = KSPSolve(snes->ksp,F,Y);CHKERRQ(ierr); 46 } 47 snes->reason = SNES_CONVERGED_ITS; 48 SNESCheckKSPSolve(snes); 49 50 ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 51 ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr); 52 snes->iter++; 53 54 /* Take the computed step. */ 55 ierr = VecAXPY(X,-1.0,Y);CHKERRQ(ierr); 56 if (snes->numbermonitors) { 57 PetscReal fnorm; 58 ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 59 ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); 60 ierr = SNESMonitor(snes,1,fnorm);CHKERRQ(ierr); 61 } 62 PetscFunctionReturn(0); 63 } 64 65 static PetscErrorCode SNESSetUp_KSPONLY(SNES snes) 66 { 67 PetscErrorCode ierr; 68 69 PetscFunctionBegin; 70 ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr); 71 PetscFunctionReturn(0); 72 } 73 74 static PetscErrorCode SNESDestroy_KSPONLY(SNES snes) 75 { 76 PetscErrorCode ierr; 77 78 PetscFunctionBegin; 79 ierr = PetscFree(snes->data);CHKERRQ(ierr); 80 PetscFunctionReturn(0); 81 } 82 83 /* -------------------------------------------------------------------------- */ 84 /*MC 85 SNESKSPONLY - Nonlinear solver that only performs one Newton step and does not compute any norms. 86 The main purpose of this solver is to solve linear problems using the SNES interface, without 87 any additional overhead in the form of vector operations. 88 89 Level: beginner 90 91 .seealso: SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR 92 M*/ 93 PETSC_EXTERN PetscErrorCode SNESCreate_KSPONLY(SNES snes) 94 { 95 SNES_KSPONLY *ksponly; 96 PetscErrorCode ierr; 97 98 PetscFunctionBegin; 99 snes->ops->setup = SNESSetUp_KSPONLY; 100 snes->ops->solve = SNESSolve_KSPONLY; 101 snes->ops->destroy = SNESDestroy_KSPONLY; 102 snes->ops->setfromoptions = 0; 103 snes->ops->view = 0; 104 snes->ops->reset = 0; 105 106 snes->usesksp = PETSC_TRUE; 107 snes->usesnpc = PETSC_FALSE; 108 109 snes->alwayscomputesfinalresidual = PETSC_FALSE; 110 111 ierr = PetscNewLog(snes,&ksponly);CHKERRQ(ierr); 112 snes->data = (void*)ksponly; 113 PetscFunctionReturn(0); 114 } 115 116 /*MC 117 SNESKSPTRANSPOSEONLY - Nonlinear solver that only performs one Newton step and does not compute any norms. 118 The main purpose of this solver is to solve transposed linear problems using the SNES interface, without 119 any additional overhead in the form of vector operations within adjoint solvers. 120 121 Level: beginner 122 123 .seealso: SNESCreate(), SNES, SNESSetType(), SNESKSPTRANSPOSEONLY, SNESNEWTONLS, SNESNEWTONTR 124 M*/ 125 PETSC_EXTERN PetscErrorCode SNESCreate_KSPTRANSPOSEONLY(SNES snes) 126 { 127 SNES_KSPONLY *kspo; 128 PetscErrorCode ierr; 129 130 PetscFunctionBegin; 131 ierr = SNESCreate_KSPONLY(snes);CHKERRQ(ierr); 132 kspo = (SNES_KSPONLY*)snes->data; 133 kspo->transpose_solve = PETSC_TRUE; 134 PetscFunctionReturn(0); 135 } 136