1 #include <../src/snes/impls/richardson/snesrichardsonimpl.h> 2 3 4 #undef __FUNCT__ 5 #define __FUNCT__ "SNESReset_NRichardson" 6 PetscErrorCode SNESReset_NRichardson(SNES snes) 7 { 8 PetscFunctionBegin; 9 PetscFunctionReturn(0); 10 } 11 12 /* 13 SNESDestroy_NRichardson - Destroys the private SNES_NRichardson context that was created with SNESCreate_NRichardson(). 14 15 Input Parameter: 16 . snes - the SNES context 17 18 Application Interface Routine: SNESDestroy() 19 */ 20 #undef __FUNCT__ 21 #define __FUNCT__ "SNESDestroy_NRichardson" 22 PetscErrorCode SNESDestroy_NRichardson(SNES snes) 23 { 24 PetscErrorCode ierr; 25 26 PetscFunctionBegin; 27 ierr = SNESReset_NRichardson(snes);CHKERRQ(ierr); 28 ierr = PetscFree(snes->data);CHKERRQ(ierr); 29 PetscFunctionReturn(0); 30 } 31 32 /* 33 SNESSetUp_NRichardson - Sets up the internal data structures for the later use 34 of the SNESNRICHARDSON nonlinear solver. 35 36 Input Parameters: 37 + snes - the SNES context 38 - x - the solution vector 39 40 Application Interface Routine: SNESSetUp() 41 */ 42 #undef __FUNCT__ 43 #define __FUNCT__ "SNESSetUp_NRichardson" 44 PetscErrorCode SNESSetUp_NRichardson(SNES snes) 45 { 46 PetscFunctionBegin; 47 if (snes->pcside == PC_RIGHT) {SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"NRichardson only supports left preconditioning");} 48 PetscFunctionReturn(0); 49 } 50 51 /* 52 SNESSetFromOptions_NRichardson - Sets various parameters for the SNESNEWTONLS method. 53 54 Input Parameter: 55 . snes - the SNES context 56 57 Application Interface Routine: SNESSetFromOptions() 58 */ 59 #undef __FUNCT__ 60 #define __FUNCT__ "SNESSetFromOptions_NRichardson" 61 static PetscErrorCode SNESSetFromOptions_NRichardson(SNES snes) 62 { 63 PetscErrorCode ierr; 64 SNESLineSearch linesearch; 65 66 PetscFunctionBegin; 67 ierr = PetscOptionsHead("SNES Richardson options");CHKERRQ(ierr); 68 ierr = PetscOptionsTail();CHKERRQ(ierr); 69 if (!snes->linesearch) { 70 ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr); 71 ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHL2);CHKERRQ(ierr); 72 } 73 PetscFunctionReturn(0); 74 } 75 76 /* 77 SNESView_NRichardson - Prints info from the SNESRichardson data structure. 78 79 Input Parameters: 80 + SNES - the SNES context 81 - viewer - visualization context 82 83 Application Interface Routine: SNESView() 84 */ 85 #undef __FUNCT__ 86 #define __FUNCT__ "SNESView_NRichardson" 87 static PetscErrorCode SNESView_NRichardson(SNES snes, PetscViewer viewer) 88 { 89 PetscBool iascii; 90 PetscErrorCode ierr; 91 92 PetscFunctionBegin; 93 ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 94 if (iascii) { 95 } 96 PetscFunctionReturn(0); 97 } 98 99 /* 100 SNESSolve_NRichardson - Solves a nonlinear system with the Richardson method. 101 102 Input Parameters: 103 . snes - the SNES context 104 105 Output Parameter: 106 . outits - number of iterations until termination 107 108 Application Interface Routine: SNESSolve() 109 */ 110 #undef __FUNCT__ 111 #define __FUNCT__ "SNESSolve_NRichardson" 112 PetscErrorCode SNESSolve_NRichardson(SNES snes) 113 { 114 Vec X, Y, F; 115 PetscReal xnorm, fnorm, ynorm; 116 PetscInt maxits, i; 117 PetscErrorCode ierr; 118 PetscBool lsSuccess; 119 SNESConvergedReason reason; 120 121 PetscFunctionBegin; 122 snes->reason = SNES_CONVERGED_ITERATING; 123 124 maxits = snes->max_its; /* maximum number of iterations */ 125 X = snes->vec_sol; /* X^n */ 126 Y = snes->vec_sol_update; /* \tilde X */ 127 F = snes->vec_func; /* residual vector */ 128 129 ierr = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr); 130 snes->iter = 0; 131 snes->norm = 0.; 132 ierr = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr); 133 if (!snes->vec_func_init_set) { 134 ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 135 if (snes->domainerror) { 136 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 137 PetscFunctionReturn(0); 138 } 139 } else snes->vec_func_init_set = PETSC_FALSE; 140 141 if (!snes->norm_init_set) { 142 ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 143 if (PetscIsInfOrNanReal(fnorm)) { 144 snes->reason = SNES_DIVERGED_FNORM_NAN; 145 PetscFunctionReturn(0); 146 } 147 } else { 148 fnorm = snes->norm_init; 149 snes->norm_init_set = PETSC_FALSE; 150 } 151 152 ierr = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr); 153 snes->norm = fnorm; 154 ierr = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr); 155 ierr = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr); 156 ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 157 158 /* set parameter for default relative tolerance convergence test */ 159 snes->ttol = fnorm*snes->rtol; 160 /* test convergence */ 161 ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 162 if (snes->reason) PetscFunctionReturn(0); 163 164 for (i = 0; i < maxits; i++) { 165 lsSuccess = PETSC_TRUE; 166 /* Call general purpose update function */ 167 if (snes->ops->update) { 168 ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 169 } 170 if (snes->pc) { 171 ierr = VecCopy(X,Y);CHKERRQ(ierr); 172 if (snes->functype == SNES_FUNCTION_UNPRECONDITIONED) { 173 ierr = SNESSetInitialFunction(snes->pc, F);CHKERRQ(ierr); 174 ierr = SNESSetInitialFunctionNorm(snes->pc, fnorm);CHKERRQ(ierr); 175 } 176 ierr = PetscLogEventBegin(SNES_NPCSolve,snes->pc,Y,snes->vec_rhs,0);CHKERRQ(ierr); 177 ierr = SNESSolve(snes->pc, snes->vec_rhs, Y);CHKERRQ(ierr); 178 ierr = PetscLogEventEnd(SNES_NPCSolve,snes->pc,Y,snes->vec_rhs,0);CHKERRQ(ierr); 179 ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 180 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 181 snes->reason = SNES_DIVERGED_INNER; 182 PetscFunctionReturn(0); 183 } 184 ierr = VecAYPX(Y,-1.0,X);CHKERRQ(ierr); 185 } else { 186 ierr = VecCopy(F,Y);CHKERRQ(ierr); 187 } 188 ierr = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Y);CHKERRQ(ierr); 189 ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr); 190 ierr = SNESLineSearchGetSuccess(snes->linesearch, &lsSuccess);CHKERRQ(ierr); 191 if (!lsSuccess) { 192 if (++snes->numFailures >= snes->maxFailures) { 193 snes->reason = SNES_DIVERGED_LINE_SEARCH; 194 break; 195 } 196 } 197 if (snes->nfuncs >= snes->max_funcs) { 198 snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 199 break; 200 } 201 if (snes->domainerror) { 202 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 203 PetscFunctionReturn(0); 204 } 205 /* Monitor convergence */ 206 ierr = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr); 207 snes->iter = i+1; 208 snes->norm = fnorm; 209 ierr = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr); 210 ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 211 ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 212 /* Test for convergence */ 213 ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 214 if (snes->reason) break; 215 } 216 if (i == maxits) { 217 ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr); 218 if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 219 } 220 PetscFunctionReturn(0); 221 } 222 223 /*MC 224 SNESNRICHARDSON - Richardson nonlinear solver that uses successive substitutions, also sometimes known as Picard iteration. 225 226 Level: beginner 227 228 Options Database: 229 + -snes_linesearch_type <l2,cp,basic> Line search type. 230 - -snes_linesearch_damping<1.0> Damping for the line search. 231 232 Notes: If no inner nonlinear preconditioner is provided then solves F(x) - b = 0 using x^{n+1} = x^{n} - lambda 233 (F(x^n) - b) where lambda is obtained either SNESLineSearchSetDamping(), -snes_damping or a line search. If 234 an inner nonlinear preconditioner is provided (either with -npc_snes_type or SNESSetPC()) then the inner 235 solver is called an initial solution x^n and the nonlinear Richardson uses x^{n+1} = x^{n} + lambda d^{n} 236 where d^{n} = \hat{x}^{n} - x^{n} where \hat{x}^{n} is the solution returned from the inner solver. 237 238 The update, especially without inner nonlinear preconditioner, may be ill-scaled. If using the basic 239 linesearch, one may have to scale the update with -snes_linesearch_damping 240 241 This uses no derivative information thus will be much slower then Newton's method obtained with -snes_type ls 242 243 .seealso: SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR, SNESNGMRES, SNESQN, SNESNCG 244 M*/ 245 #undef __FUNCT__ 246 #define __FUNCT__ "SNESCreate_NRichardson" 247 PETSC_EXTERN PetscErrorCode SNESCreate_NRichardson(SNES snes) 248 { 249 PetscErrorCode ierr; 250 SNES_NRichardson *neP; 251 252 PetscFunctionBegin; 253 snes->ops->destroy = SNESDestroy_NRichardson; 254 snes->ops->setup = SNESSetUp_NRichardson; 255 snes->ops->setfromoptions = SNESSetFromOptions_NRichardson; 256 snes->ops->view = SNESView_NRichardson; 257 snes->ops->solve = SNESSolve_NRichardson; 258 snes->ops->reset = SNESReset_NRichardson; 259 260 snes->usesksp = PETSC_FALSE; 261 snes->usespc = PETSC_TRUE; 262 263 snes->functype = SNES_FUNCTION_UNPRECONDITIONED; 264 snes->pcside = PC_LEFT; 265 266 ierr = PetscNewLog(snes, SNES_NRichardson, &neP);CHKERRQ(ierr); 267 snes->data = (void*) neP; 268 269 if (!snes->tolerancesset) { 270 snes->max_funcs = 30000; 271 snes->max_its = 10000; 272 snes->stol = 1e-20; 273 } 274 PetscFunctionReturn(0); 275 } 276