1 #define PETSCSNES_DLL 2 3 #include "src/snes/impls/tr/tr.h" /*I "petscsnes.h" I*/ 4 5 /* 6 This convergence test determines if the two norm of the 7 solution lies outside the trust region, if so it halts. 8 */ 9 #undef __FUNCT__ 10 #define __FUNCT__ "SNES_TR_KSPConverged_Private" 11 PetscErrorCode SNES_TR_KSPConverged_Private(KSP ksp,PetscInt n,PetscReal rnorm,KSPConvergedReason *reason,void *ctx) 12 { 13 SNES snes = (SNES) ctx; 14 SNES_TR *neP = (SNES_TR*)snes->data; 15 Vec x; 16 PetscReal nrm; 17 PetscErrorCode ierr; 18 19 PetscFunctionBegin; 20 ierr = KSPDefaultConverged(ksp,n,rnorm,reason,ctx);CHKERRQ(ierr); 21 if (*reason) { 22 ierr = PetscInfo2(snes,"default convergence test KSP iterations=%D, rnorm=%G\n",n,rnorm);CHKERRQ(ierr); 23 } 24 /* Determine norm of solution */ 25 ierr = KSPBuildSolution(ksp,0,&x);CHKERRQ(ierr); 26 ierr = VecNorm(x,NORM_2,&nrm);CHKERRQ(ierr); 27 if (nrm >= neP->delta) { 28 ierr = PetscInfo2(snes,"Ending linear iteration early, delta=%G, length=%G\n",neP->delta,nrm);CHKERRQ(ierr); 29 *reason = KSP_CONVERGED_STEP_LENGTH; 30 } 31 PetscFunctionReturn(0); 32 } 33 34 /* 35 SNESSolve_TR - Implements Newton's Method with a very simple trust 36 region approach for solving systems of nonlinear equations. 37 38 39 */ 40 #undef __FUNCT__ 41 #define __FUNCT__ "SNESSolve_TR" 42 static PetscErrorCode SNESSolve_TR(SNES snes) 43 { 44 SNES_TR *neP = (SNES_TR*)snes->data; 45 Vec X,F,Y,G,TMP,Ytmp; 46 PetscErrorCode ierr; 47 PetscInt maxits,i,lits; 48 MatStructure flg = DIFFERENT_NONZERO_PATTERN; 49 PetscReal rho,fnorm,gnorm,gpnorm,xnorm,delta,nrm,ynorm,norm1; 50 PetscScalar cnorm; 51 KSP ksp; 52 SNESConvergedReason reason; 53 PetscTruth conv,breakout = PETSC_FALSE; 54 55 PetscFunctionBegin; 56 maxits = snes->max_its; /* maximum number of iterations */ 57 X = snes->vec_sol; /* solution vector */ 58 F = snes->vec_func; /* residual vector */ 59 Y = snes->work[0]; /* work vectors */ 60 G = snes->work[1]; 61 Ytmp = snes->work[2]; 62 63 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 64 snes->iter = 0; 65 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 66 ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); /* xnorm = || X || */ 67 68 ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); /* F(X) */ 69 ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); /* fnorm <- || F || */ 70 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 71 snes->norm = fnorm; 72 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 73 delta = neP->delta0*fnorm; 74 neP->delta = delta; 75 SNESLogConvHistory(snes,fnorm,0); 76 SNESMonitor(snes,0,fnorm); 77 ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 78 79 if (fnorm < snes->abstol) {snes->reason = SNES_CONVERGED_FNORM_ABS; PetscFunctionReturn(0);} 80 81 /* set parameter for default relative tolerance convergence test */ 82 snes->ttol = fnorm*snes->rtol; 83 84 /* Set the stopping criteria to use the More' trick. */ 85 ierr = PetscOptionsHasName(PETSC_NULL,"-snes_tr_ksp_regular_convergence_test",&conv);CHKERRQ(ierr); 86 if (!conv) { 87 ierr = KSPSetConvergenceTest(ksp,SNES_TR_KSPConverged_Private,(void*)snes);CHKERRQ(ierr); 88 ierr = PetscInfo(snes,"Using Krylov convergence test SNES_TR_KSPConverged_Private\n");CHKERRQ(ierr); 89 } 90 91 for (i=0; i<maxits; i++) { 92 93 /* Call general purpose update function */ 94 if (snes->ops->update) { 95 ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 96 } 97 98 ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 99 ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr); 100 101 /* Solve J Y = F, where J is Jacobian matrix */ 102 ierr = SNES_KSPSolve(snes,ksp,F,Ytmp);CHKERRQ(ierr); 103 ierr = KSPGetIterationNumber(ksp,&lits);CHKERRQ(ierr); 104 snes->linear_its += lits; 105 ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr); 106 ierr = VecNorm(Ytmp,NORM_2,&nrm);CHKERRQ(ierr); 107 norm1 = nrm; 108 while(1) { 109 ierr = VecCopy(Ytmp,Y);CHKERRQ(ierr); 110 nrm = norm1; 111 112 /* Scale Y if need be and predict new value of F norm */ 113 if (nrm >= delta) { 114 nrm = delta/nrm; 115 gpnorm = (1.0 - nrm)*fnorm; 116 cnorm = nrm; 117 ierr = PetscInfo1(snes,"Scaling direction by %G\n",nrm);CHKERRQ(ierr); 118 ierr = VecScale(Y,cnorm);CHKERRQ(ierr); 119 nrm = gpnorm; 120 ynorm = delta; 121 } else { 122 gpnorm = 0.0; 123 ierr = PetscInfo(snes,"Direction is in Trust Region\n");CHKERRQ(ierr); 124 ynorm = nrm; 125 } 126 ierr = VecAYPX(Y,-1.0,X);CHKERRQ(ierr); /* Y <- X - Y */ 127 ierr = VecCopy(X,snes->vec_sol_update_always);CHKERRQ(ierr); 128 ierr = SNESComputeFunction(snes,Y,G);CHKERRQ(ierr); /* F(X) */ 129 ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr); /* gnorm <- || g || */ 130 if (fnorm == gpnorm) rho = 0.0; 131 else rho = (fnorm*fnorm - gnorm*gnorm)/(fnorm*fnorm - gpnorm*gpnorm); 132 133 /* Update size of trust region */ 134 if (rho < neP->mu) delta *= neP->delta1; 135 else if (rho < neP->eta) delta *= neP->delta2; 136 else delta *= neP->delta3; 137 ierr = PetscInfo3(snes,"fnorm=%G, gnorm=%G, ynorm=%G\n",fnorm,gnorm,ynorm);CHKERRQ(ierr); 138 ierr = PetscInfo3(snes,"gpred=%G, rho=%G, delta=%G\n",gpnorm,rho,delta);CHKERRQ(ierr); 139 neP->delta = delta; 140 if (rho > neP->sigma) break; 141 ierr = PetscInfo(snes,"Trying again in smaller region\n");CHKERRQ(ierr); 142 /* check to see if progress is hopeless */ 143 neP->itflag = PETSC_FALSE; 144 ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr); 145 if (reason) { 146 /* We're not progressing, so return with the current iterate */ 147 SNESMonitor(snes,i+1,fnorm); 148 breakout = PETSC_TRUE; 149 break; 150 } 151 snes->numFailures++; 152 } 153 if (!breakout) { 154 fnorm = gnorm; 155 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 156 snes->iter = i+1; 157 snes->norm = fnorm; 158 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 159 TMP = F; F = G; snes->vec_func_always = F; G = TMP; 160 TMP = X; X = Y; snes->vec_sol_always = X; Y = TMP; 161 ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); /* xnorm = || X || */ 162 SNESLogConvHistory(snes,fnorm,lits); 163 SNESMonitor(snes,i+1,fnorm); 164 165 /* Test for convergence */ 166 neP->itflag = PETSC_TRUE; 167 ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr); 168 if (reason) { 169 break; 170 } 171 } else { 172 break; 173 } 174 } 175 /* Verify solution is in corect location */ 176 if (X != snes->vec_sol) { 177 ierr = VecCopy(X,snes->vec_sol);CHKERRQ(ierr); 178 } 179 if (F != snes->vec_func) { 180 ierr = VecCopy(F,snes->vec_func);CHKERRQ(ierr); 181 } 182 snes->vec_sol_always = snes->vec_sol; 183 snes->vec_func_always = snes->vec_func; 184 if (i == maxits) { 185 ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr); 186 reason = SNES_DIVERGED_MAX_IT; 187 } 188 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 189 snes->reason = reason; 190 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 191 PetscFunctionReturn(0); 192 } 193 /*------------------------------------------------------------*/ 194 #undef __FUNCT__ 195 #define __FUNCT__ "SNESSetUp_TR" 196 static PetscErrorCode SNESSetUp_TR(SNES snes) 197 { 198 PetscErrorCode ierr; 199 200 PetscFunctionBegin; 201 if (!snes->work) { 202 snes->nwork = 4; 203 ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr); 204 ierr = PetscLogObjectParents(snes,snes->nwork,snes->work);CHKERRQ(ierr); 205 } 206 snes->vec_sol_update_always = snes->work[3]; 207 PetscFunctionReturn(0); 208 } 209 /*------------------------------------------------------------*/ 210 #undef __FUNCT__ 211 #define __FUNCT__ "SNESDestroy_TR" 212 static PetscErrorCode SNESDestroy_TR(SNES snes) 213 { 214 PetscErrorCode ierr; 215 216 PetscFunctionBegin; 217 if (snes->nwork) { 218 ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr); 219 } 220 ierr = PetscFree(snes->data);CHKERRQ(ierr); 221 PetscFunctionReturn(0); 222 } 223 /*------------------------------------------------------------*/ 224 225 #undef __FUNCT__ 226 #define __FUNCT__ "SNESSetFromOptions_TR" 227 static PetscErrorCode SNESSetFromOptions_TR(SNES snes) 228 { 229 SNES_TR *ctx = (SNES_TR *)snes->data; 230 PetscErrorCode ierr; 231 232 PetscFunctionBegin; 233 ierr = PetscOptionsHead("SNES trust region options for nonlinear equations");CHKERRQ(ierr); 234 ierr = PetscOptionsReal("-snes_trtol","Trust region tolerance","SNESSetTrustRegionTolerance",snes->deltatol,&snes->deltatol,0);CHKERRQ(ierr); 235 ierr = PetscOptionsReal("-snes_tr_mu","mu","None",ctx->mu,&ctx->mu,0);CHKERRQ(ierr); 236 ierr = PetscOptionsReal("-snes_tr_eta","eta","None",ctx->eta,&ctx->eta,0);CHKERRQ(ierr); 237 ierr = PetscOptionsReal("-snes_tr_sigma","sigma","None",ctx->sigma,&ctx->sigma,0);CHKERRQ(ierr); 238 ierr = PetscOptionsReal("-snes_tr_delta0","delta0","None",ctx->delta0,&ctx->delta0,0);CHKERRQ(ierr); 239 ierr = PetscOptionsReal("-snes_tr_delta1","delta1","None",ctx->delta1,&ctx->delta1,0);CHKERRQ(ierr); 240 ierr = PetscOptionsReal("-snes_tr_delta2","delta2","None",ctx->delta2,&ctx->delta2,0);CHKERRQ(ierr); 241 ierr = PetscOptionsReal("-snes_tr_delta3","delta3","None",ctx->delta3,&ctx->delta3,0);CHKERRQ(ierr); 242 ierr = PetscOptionsTail();CHKERRQ(ierr); 243 PetscFunctionReturn(0); 244 } 245 246 #undef __FUNCT__ 247 #define __FUNCT__ "SNESView_TR" 248 static PetscErrorCode SNESView_TR(SNES snes,PetscViewer viewer) 249 { 250 SNES_TR *tr = (SNES_TR *)snes->data; 251 PetscErrorCode ierr; 252 PetscTruth iascii; 253 254 PetscFunctionBegin; 255 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 256 if (iascii) { 257 ierr = PetscViewerASCIIPrintf(viewer," mu=%G, eta=%G, sigma=%G\n",tr->mu,tr->eta,tr->sigma);CHKERRQ(ierr); 258 ierr = PetscViewerASCIIPrintf(viewer," delta0=%G, delta1=%G, delta2=%G, delta3=%G\n",tr->delta0,tr->delta1,tr->delta2,tr->delta3);CHKERRQ(ierr); 259 } else { 260 SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for SNES EQ TR",((PetscObject)viewer)->type_name); 261 } 262 PetscFunctionReturn(0); 263 } 264 265 /* ---------------------------------------------------------------- */ 266 #undef __FUNCT__ 267 #define __FUNCT__ "SNESConverged_TR" 268 /*@C 269 SNESConverged_TR - Monitors the convergence of the trust region 270 method SNESTR for solving systems of nonlinear equations (default). 271 272 Collective on SNES 273 274 Input Parameters: 275 + snes - the SNES context 276 . xnorm - 2-norm of current iterate 277 . pnorm - 2-norm of current step 278 . fnorm - 2-norm of function 279 - dummy - unused context 280 281 Output Parameter: 282 . reason - one of 283 $ SNES_CONVERGED_FNORM_ABS - (fnorm < abstol), 284 $ SNES_CONVERGED_PNORM_RELATIVE - (pnorm < xtol*xnorm), 285 $ SNES_CONVERGED_FNORM_RELATIVE - (fnorm < rtol*fnorm0), 286 $ SNES_DIVERGED_FUNCTION_COUNT - (nfct > maxf), 287 $ SNES_DIVERGED_FNORM_NAN - (fnorm == NaN), 288 $ SNES_CONVERGED_TR_DELTA - (delta < xnorm*deltatol), 289 $ SNES_CONVERGED_ITERATING - (otherwise) 290 291 where 292 + delta - trust region paramenter 293 . deltatol - trust region size tolerance, 294 set with SNESSetTrustRegionTolerance() 295 . maxf - maximum number of function evaluations, 296 set with SNESSetTolerances() 297 . nfct - number of function evaluations, 298 . abstol - absolute function norm tolerance, 299 set with SNESSetTolerances() 300 - xtol - relative function norm tolerance, 301 set with SNESSetTolerances() 302 303 Level: intermediate 304 305 .keywords: SNES, nonlinear, default, converged, convergence 306 307 .seealso: SNESSetConvergenceTest(), SNESEisenstatWalkerConverged() 308 @*/ 309 PetscErrorCode PETSCSNES_DLLEXPORT SNESConverged_TR(SNES snes,PetscInt it,PetscReal xnorm,PetscReal pnorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy) 310 { 311 SNES_TR *neP = (SNES_TR *)snes->data; 312 PetscErrorCode ierr; 313 314 PetscFunctionBegin; 315 if (fnorm != fnorm) { 316 ierr = PetscInfo(snes,"Failed to converged, function norm is NaN\n");CHKERRQ(ierr); 317 *reason = SNES_DIVERGED_FNORM_NAN; 318 } else if (neP->delta < xnorm * snes->deltatol) { 319 ierr = PetscInfo3(snes,"Converged due to trust region param %G<%G*%G\n",neP->delta,xnorm,snes->deltatol);CHKERRQ(ierr); 320 *reason = SNES_CONVERGED_TR_DELTA; 321 } else if (neP->itflag) { 322 ierr = SNESConverged_LS(snes,it,xnorm,pnorm,fnorm,reason,dummy);CHKERRQ(ierr); 323 } else if (snes->nfuncs >= snes->max_funcs) { 324 ierr = PetscInfo2(snes,"Exceeded maximum number of function evaluations: %D > %D\n",snes->nfuncs,snes->max_funcs);CHKERRQ(ierr); 325 *reason = SNES_DIVERGED_FUNCTION_COUNT; 326 } else { 327 *reason = SNES_CONVERGED_ITERATING; 328 } 329 PetscFunctionReturn(0); 330 } 331 /* ------------------------------------------------------------ */ 332 /*MC 333 SNESTR - Newton based nonlinear solver that uses a trust region 334 335 Options Database: 336 + -snes_trtol <tol> Trust region tolerance 337 . -snes_tr_mu <mu> 338 . -snes_tr_eta <eta> 339 . -snes_tr_sigma <sigma> 340 . -snes_tr_delta0 <delta0> 341 . -snes_tr_delta1 <delta1> 342 . -snes_tr_delta2 <delta2> 343 - -snes_tr_delta3 <delta3> 344 345 The basic algorithm is taken from "The Minpack Project", by More', 346 Sorensen, Garbow, Hillstrom, pages 88-111 of "Sources and Development 347 of Mathematical Software", Wayne Cowell, editor. 348 349 This is intended as a model implementation, since it does not 350 necessarily have many of the bells and whistles of other 351 implementations. 352 353 Level: intermediate 354 355 .seealso: SNESCreate(), SNES, SNESSetType(), SNESLS, SNESSetTrustRegionTolerance() 356 357 M*/ 358 EXTERN_C_BEGIN 359 #undef __FUNCT__ 360 #define __FUNCT__ "SNESCreate_TR" 361 PetscErrorCode PETSCSNES_DLLEXPORT SNESCreate_TR(SNES snes) 362 { 363 SNES_TR *neP; 364 PetscErrorCode ierr; 365 366 PetscFunctionBegin; 367 snes->ops->setup = SNESSetUp_TR; 368 snes->ops->solve = SNESSolve_TR; 369 snes->ops->destroy = SNESDestroy_TR; 370 snes->ops->converged = SNESConverged_TR; 371 snes->ops->setfromoptions = SNESSetFromOptions_TR; 372 snes->ops->view = SNESView_TR; 373 snes->nwork = 0; 374 375 ierr = PetscNew(SNES_TR,&neP);CHKERRQ(ierr); 376 ierr = PetscLogObjectMemory(snes,sizeof(SNES_TR));CHKERRQ(ierr); 377 snes->data = (void*)neP; 378 neP->mu = 0.25; 379 neP->eta = 0.75; 380 neP->delta = 0.0; 381 neP->delta0 = 0.2; 382 neP->delta1 = 0.3; 383 neP->delta2 = 0.75; 384 neP->delta3 = 2.0; 385 neP->sigma = 0.0001; 386 neP->itflag = PETSC_FALSE; 387 neP->rnorm0 = 0; 388 neP->ttol = 0; 389 PetscFunctionReturn(0); 390 } 391 EXTERN_C_END 392 393