1 #include <../src/snes/impls/ngmres/snesngmres.h> 2 #include <petscblaslapack.h> 3 4 #undef __FUNCT__ 5 #define __FUNCT__ "SNESReset_NGMRES" 6 PetscErrorCode SNESReset_NGMRES(SNES snes) 7 { 8 SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data; 9 PetscErrorCode ierr; 10 11 PetscFunctionBegin; 12 ierr = VecDestroyVecs(ngmres->msize, &ngmres->Fdot);CHKERRQ(ierr); 13 ierr = VecDestroyVecs(ngmres->msize, &ngmres->Xdot);CHKERRQ(ierr); 14 ierr = SNESLineSearchDestroy(&ngmres->additive_linesearch);CHKERRQ(ierr); 15 16 PetscFunctionReturn(0); 17 } 18 19 #undef __FUNCT__ 20 #define __FUNCT__ "SNESDestroy_NGMRES" 21 PetscErrorCode SNESDestroy_NGMRES(SNES snes) 22 { 23 PetscErrorCode ierr; 24 SNES_NGMRES *ngmres = (SNES_NGMRES*)snes->data; 25 26 PetscFunctionBegin; 27 ierr = SNESReset_NGMRES(snes);CHKERRQ(ierr); 28 ierr = PetscFree5(ngmres->h, ngmres->beta, ngmres->xi, ngmres->fnorms, ngmres->q);CHKERRQ(ierr); 29 ierr = PetscFree(ngmres->s);CHKERRQ(ierr); 30 ierr = PetscFree(ngmres->xnorms);CHKERRQ(ierr); 31 #if PETSC_USE_COMPLEX 32 ierr = PetscFree(ngmres->rwork); 33 #endif 34 ierr = PetscFree(ngmres->work); 35 ierr = PetscFree(snes->data); 36 PetscFunctionReturn(0); 37 } 38 39 #undef __FUNCT__ 40 #define __FUNCT__ "SNESSetUp_NGMRES" 41 PetscErrorCode SNESSetUp_NGMRES(SNES snes) 42 { 43 SNES_NGMRES *ngmres = (SNES_NGMRES *) snes->data; 44 const char *optionsprefix; 45 PetscInt msize,hsize; 46 PetscErrorCode ierr; 47 48 PetscFunctionBegin; 49 ierr = SNESDefaultGetWork(snes,5);CHKERRQ(ierr); 50 if (!ngmres->Xdot) {ierr = VecDuplicateVecs(snes->vec_sol,ngmres->msize,&ngmres->Xdot);CHKERRQ(ierr);} 51 if (!ngmres->Fdot) {ierr = VecDuplicateVecs(snes->vec_sol,ngmres->msize,&ngmres->Fdot);CHKERRQ(ierr);} 52 if (!ngmres->setup_called) { 53 msize = ngmres->msize; /* restart size */ 54 hsize = msize * msize; 55 56 /* explicit least squares minimization solve */ 57 ierr = PetscMalloc5(hsize,PetscScalar,&ngmres->h, 58 msize,PetscScalar,&ngmres->beta, 59 msize,PetscScalar,&ngmres->xi, 60 msize,PetscReal, &ngmres->fnorms, 61 hsize,PetscScalar,&ngmres->q);CHKERRQ(ierr); 62 if (ngmres->singlereduction) { 63 ierr = PetscMalloc(msize*sizeof(PetscReal),&ngmres->xnorms);CHKERRQ(ierr); 64 } 65 ngmres->nrhs = 1; 66 ngmres->lda = msize; 67 ngmres->ldb = msize; 68 ierr = PetscMalloc(msize*sizeof(PetscScalar),&ngmres->s);CHKERRQ(ierr); 69 ierr = PetscMemzero(ngmres->h, hsize*sizeof(PetscScalar));CHKERRQ(ierr); 70 ierr = PetscMemzero(ngmres->q, hsize*sizeof(PetscScalar));CHKERRQ(ierr); 71 ierr = PetscMemzero(ngmres->xi, msize*sizeof(PetscScalar));CHKERRQ(ierr); 72 ierr = PetscMemzero(ngmres->beta,msize*sizeof(PetscScalar));CHKERRQ(ierr); 73 ngmres->lwork = 12*msize; 74 #if PETSC_USE_COMPLEX 75 ierr = PetscMalloc(sizeof(PetscReal)*ngmres->lwork,&ngmres->rwork); 76 #endif 77 ierr = PetscMalloc(sizeof(PetscScalar)*ngmres->lwork,&ngmres->work); 78 } 79 80 /* linesearch setup */ 81 ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr); 82 83 if (ngmres->additive) { 84 ierr = SNESLineSearchCreate(((PetscObject)snes)->comm, &ngmres->additive_linesearch);CHKERRQ(ierr); 85 ierr = SNESLineSearchSetSNES(ngmres->additive_linesearch, snes);CHKERRQ(ierr); 86 ierr = SNESLineSearchSetType(ngmres->additive_linesearch, SNESLINESEARCHL2);CHKERRQ(ierr); 87 ierr = SNESLineSearchAppendOptionsPrefix(ngmres->additive_linesearch, "additive_");CHKERRQ(ierr); 88 ierr = SNESLineSearchAppendOptionsPrefix(ngmres->additive_linesearch, optionsprefix);CHKERRQ(ierr); 89 ierr = SNESLineSearchSetFromOptions(ngmres->additive_linesearch);CHKERRQ(ierr); 90 } 91 92 ngmres->setup_called = PETSC_TRUE; 93 PetscFunctionReturn(0); 94 } 95 96 #undef __FUNCT__ 97 #define __FUNCT__ "SNESSetFromOptions_NGMRES" 98 PetscErrorCode SNESSetFromOptions_NGMRES(SNES snes) 99 { 100 SNES_NGMRES *ngmres = (SNES_NGMRES *) snes->data; 101 PetscErrorCode ierr; 102 PetscBool debug; 103 SNESLineSearch linesearch; 104 PetscFunctionBegin; 105 ierr = PetscOptionsHead("SNES NGMRES options");CHKERRQ(ierr); 106 ierr = PetscOptionsBool("-snes_ngmres_additive", "Use additive variant vs. choice", "SNES", ngmres->additive, &ngmres->additive, PETSC_NULL);CHKERRQ(ierr); 107 ierr = PetscOptionsBool("-snes_ngmres_anderson", "Use Anderson mixing storage", "SNES", ngmres->anderson, &ngmres->anderson, PETSC_NULL);CHKERRQ(ierr); 108 ierr = PetscOptionsInt("-snes_ngmres_m", "Number of directions", "SNES", ngmres->msize, &ngmres->msize, PETSC_NULL);CHKERRQ(ierr); 109 ierr = PetscOptionsInt("-snes_ngmres_restart_it","Tolerance iterations before restart","SNES", ngmres->restart_it, &ngmres->restart_it, PETSC_NULL);CHKERRQ(ierr); 110 ierr = PetscOptionsBool("-snes_ngmres_monitor", "Monitor actions of NGMRES", "SNES", ngmres->monitor ? PETSC_TRUE: PETSC_FALSE, &debug, PETSC_NULL);CHKERRQ(ierr); 111 if (debug) { 112 ngmres->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr); 113 } 114 ierr = PetscOptionsReal("-snes_ngmres_gammaA", "Residual selection constant", "SNES", ngmres->gammaA, &ngmres->gammaA, PETSC_NULL);CHKERRQ(ierr); 115 ierr = PetscOptionsReal("-snes_ngmres_gammaC", " Residual restart constant", "SNES", ngmres->gammaC, &ngmres->gammaC, PETSC_NULL);CHKERRQ(ierr); 116 ierr = PetscOptionsReal("-snes_ngmres_epsilonB", "Difference selection constant", "SNES", ngmres->epsilonB, &ngmres->epsilonB, PETSC_NULL);CHKERRQ(ierr); 117 ierr = PetscOptionsReal("-snes_ngmres_deltaB", "Difference residual selection constant", "SNES", ngmres->deltaB, &ngmres->deltaB, PETSC_NULL);CHKERRQ(ierr); 118 ierr = PetscOptionsBool("-snes_ngmres_single_reduction", "Aggregate reductions", "SNES", ngmres->singlereduction, &ngmres->singlereduction, PETSC_NULL);CHKERRQ(ierr); 119 ierr = PetscOptionsTail();CHKERRQ(ierr); 120 if ((ngmres->gammaA > ngmres->gammaC) && (ngmres->gammaC > 2.)) ngmres->gammaC = ngmres->gammaA; 121 122 /* set the default type of the line search if the user hasn't already. */ 123 if (!snes->linesearch) { 124 ierr = SNESGetSNESLineSearch(snes, &linesearch);CHKERRQ(ierr); 125 ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHBASIC);CHKERRQ(ierr); 126 } 127 PetscFunctionReturn(0); 128 } 129 130 #undef __FUNCT__ 131 #define __FUNCT__ "SNESView_NGMRES" 132 PetscErrorCode SNESView_NGMRES(SNES snes, PetscViewer viewer) 133 { 134 SNES_NGMRES *ngmres = (SNES_NGMRES *) snes->data; 135 PetscBool iascii; 136 PetscErrorCode ierr; 137 138 PetscFunctionBegin; 139 ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 140 if (iascii) { 141 142 ierr = PetscViewerASCIIPrintf(viewer, " Number of stored past updates: %d\n", ngmres->msize);CHKERRQ(ierr); 143 ierr = PetscViewerASCIIPrintf(viewer, " Residual selection: gammaA=%1.0e, gammaC=%1.0e\n", ngmres->gammaA, ngmres->gammaC);CHKERRQ(ierr); 144 ierr = PetscViewerASCIIPrintf(viewer, " Difference restart: epsilonB=%1.0e, deltaB=%1.0e\n", ngmres->epsilonB, ngmres->deltaB);CHKERRQ(ierr); 145 } 146 PetscFunctionReturn(0); 147 } 148 149 #undef __FUNCT__ 150 #define __FUNCT__ "SNESSolve_NGMRES" 151 152 PetscErrorCode SNESSolve_NGMRES(SNES snes) 153 { 154 SNES_NGMRES *ngmres = (SNES_NGMRES *) snes->data; 155 /* present solution, residual, and preconditioned residual */ 156 Vec X, F, B, D, Y; 157 158 /* candidate linear combination answers */ 159 Vec XA, FA, XM, FM, FPC; 160 161 /* previous iterations to construct the subspace */ 162 Vec *Fdot = ngmres->Fdot; 163 Vec *Xdot = ngmres->Xdot; 164 165 /* coefficients and RHS to the minimization problem */ 166 PetscScalar *beta = ngmres->beta; 167 PetscScalar *xi = ngmres->xi; 168 PetscReal fnorm, fMnorm, fAnorm; 169 PetscReal nu; 170 PetscScalar alph_total = 0.; 171 PetscInt i, j, k, k_restart, l, ivec, restart_count = 0; 172 173 /* solution selection data */ 174 PetscBool selectA, selectRestart; 175 PetscReal dnorm, dminnorm, dcurnorm; 176 PetscReal fminnorm; 177 178 SNESConvergedReason reason; 179 PetscBool lssucceed; 180 PetscErrorCode ierr; 181 182 PetscFunctionBegin; 183 /* variable initialization */ 184 snes->reason = SNES_CONVERGED_ITERATING; 185 X = snes->vec_sol; 186 F = snes->vec_func; 187 B = snes->vec_rhs; 188 XA = snes->vec_sol_update; 189 FA = snes->work[0]; 190 D = snes->work[1]; 191 192 /* work for the line search */ 193 Y = snes->work[2]; 194 XM = snes->work[3]; 195 FM = snes->work[4]; 196 197 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 198 snes->iter = 0; 199 snes->norm = 0.; 200 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 201 202 /* initialization */ 203 204 /* r = F(x) */ 205 if (!snes->vec_func_init_set) { 206 ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 207 if (snes->domainerror) { 208 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 209 PetscFunctionReturn(0); 210 } 211 } else { 212 snes->vec_func_init_set = PETSC_FALSE; 213 } 214 215 if (!snes->norm_init_set) { 216 ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 217 if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FP, "Infinite or not-a-number generated in function evaluation"); 218 } else { 219 fnorm = snes->norm_init; 220 snes->norm_init_set = PETSC_FALSE; 221 } 222 fminnorm = fnorm; 223 /* nu = (r, r) */ 224 nu = fnorm*fnorm; 225 226 /* q_{00} = nu */ 227 Q(0,0) = nu; 228 ngmres->fnorms[0] = fnorm; 229 /* Fdot[0] = F */ 230 ierr = VecCopy(X, Xdot[0]);CHKERRQ(ierr); 231 ierr = VecCopy(F, Fdot[0]);CHKERRQ(ierr); 232 233 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 234 snes->norm = fnorm; 235 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 236 SNESLogConvHistory(snes, fnorm, 0); 237 ierr = SNESMonitor(snes, 0, fnorm);CHKERRQ(ierr); 238 ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 239 if (snes->reason) PetscFunctionReturn(0); 240 241 k_restart = 1; 242 l = 1; 243 for (k=1; k < snes->max_its+1; k++) { 244 /* select which vector of the stored subspace will be updated */ 245 ivec = k_restart % ngmres->msize; /* replace the last used part of the subspace */ 246 247 /* Computation of x^M */ 248 if (snes->pc) { 249 ierr = VecCopy(X, XM);CHKERRQ(ierr); 250 ierr = SNESSetInitialFunction(snes->pc, F);CHKERRQ(ierr); 251 ierr = SNESSetInitialFunctionNorm(snes->pc, fnorm);CHKERRQ(ierr); 252 ierr = SNESSolve(snes->pc, B, XM);CHKERRQ(ierr); 253 ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 254 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 255 snes->reason = SNES_DIVERGED_INNER; 256 PetscFunctionReturn(0); 257 } 258 ierr = SNESGetFunction(snes->pc, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr); 259 ierr = VecCopy(FPC, FM);CHKERRQ(ierr); 260 ierr = SNESGetFunctionNorm(snes->pc, &fMnorm);CHKERRQ(ierr); 261 } else { 262 /* no preconditioner -- just take gradient descent with line search */ 263 ierr = VecCopy(F, Y);CHKERRQ(ierr); 264 ierr = VecCopy(F, FM);CHKERRQ(ierr); 265 ierr = VecCopy(X, XM);CHKERRQ(ierr); 266 fMnorm = fnorm; 267 ierr = SNESLineSearchApply(snes->linesearch,XM,FM,&fMnorm,Y);CHKERRQ(ierr); 268 ierr = SNESLineSearchGetSuccess(snes->linesearch, &lssucceed);CHKERRQ(ierr); 269 if (!lssucceed) { 270 if (++snes->numFailures >= snes->maxFailures) { 271 snes->reason = SNES_DIVERGED_LINE_SEARCH; 272 PetscFunctionReturn(0); 273 } 274 } 275 } 276 277 /* r = F(x) */ 278 nu = fMnorm*fMnorm; 279 if (fminnorm > fMnorm) fminnorm = fMnorm; /* the minimum norm is now of F^M */ 280 281 /* construct the right hand side and xi factors */ 282 ierr = VecMDot(FM, l, Fdot, xi);CHKERRQ(ierr); 283 284 for (i = 0; i < l; i++) { 285 beta[i] = nu - xi[i]; 286 } 287 288 /* construct h */ 289 for (j = 0; j < l; j++) { 290 for (i = 0; i < l; i++) { 291 H(i, j) = Q(i, j) - xi[i] - xi[j] + nu; 292 } 293 } 294 295 if(l == 1) { 296 /* simply set alpha[0] = beta[0] / H[0, 0] */ 297 beta[0] = beta[0] / H(0, 0); 298 } else { 299 #ifdef PETSC_MISSING_LAPACK_GELSS 300 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "NGMRES with LS requires the LAPACK GELSS routine."); 301 #else 302 ngmres->m = PetscBLASIntCast(l); 303 ngmres->n = PetscBLASIntCast(l); 304 ngmres->info = PetscBLASIntCast(0); 305 ngmres->rcond = -1.; 306 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 307 #ifdef PETSC_USE_COMPLEX 308 LAPACKgelss_(&ngmres->m, 309 &ngmres->n, 310 &ngmres->nrhs, 311 ngmres->h, 312 &ngmres->lda, 313 ngmres->beta, 314 &ngmres->ldb, 315 ngmres->s, 316 &ngmres->rcond, 317 &ngmres->rank, 318 ngmres->work, 319 &ngmres->lwork, 320 ngmres->rwork, 321 &ngmres->info); 322 #else 323 LAPACKgelss_(&ngmres->m, 324 &ngmres->n, 325 &ngmres->nrhs, 326 ngmres->h, 327 &ngmres->lda, 328 ngmres->beta, 329 &ngmres->ldb, 330 ngmres->s, 331 &ngmres->rcond, 332 &ngmres->rank, 333 ngmres->work, 334 &ngmres->lwork, 335 &ngmres->info); 336 #endif 337 ierr = PetscFPTrapPop();CHKERRQ(ierr); 338 if (ngmres->info < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad argument to GELSS"); 339 if (ngmres->info > 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"SVD failed to converge"); 340 #endif 341 } 342 343 alph_total = 0.; 344 for (i = 0; i < l; i++) { 345 alph_total += beta[i]; 346 } 347 348 ierr = VecCopy(XM, XA);CHKERRQ(ierr); 349 ierr = VecScale(XA, 1. - alph_total);CHKERRQ(ierr); 350 351 ierr = VecMAXPY(XA, l, beta, Xdot);CHKERRQ(ierr); 352 ierr = SNESComputeFunction(snes, XA, FA);CHKERRQ(ierr); 353 354 /* differences for selection and restart */ 355 356 if (ngmres->singlereduction) { 357 dminnorm = -1.0; 358 ierr=VecCopy(XA, D);CHKERRQ(ierr); 359 ierr=VecAXPY(D, -1.0, XM);CHKERRQ(ierr); 360 for(i=0;i<l;i++) { 361 ierr=VecAXPY(Xdot[i],-1.0,XA);CHKERRQ(ierr); 362 } 363 ierr = VecNormBegin(FA, NORM_2, &fAnorm);CHKERRQ(ierr); 364 ierr = VecNormBegin(D, NORM_2, &dnorm);CHKERRQ(ierr); 365 for (i=0;i<l;i++) { 366 ierr = VecNormBegin(Xdot[i], NORM_2, &ngmres->xnorms[i]);CHKERRQ(ierr); 367 } 368 ierr = VecNormEnd(FA, NORM_2, &fAnorm);CHKERRQ(ierr); 369 ierr = VecNormEnd(D, NORM_2, &dnorm);CHKERRQ(ierr); 370 for (i=0;i<l;i++) { 371 ierr = VecNormEnd(Xdot[i], NORM_2, &ngmres->xnorms[i]);CHKERRQ(ierr); 372 } 373 for(i=0;i<l;i++) { 374 dcurnorm = ngmres->xnorms[i]; 375 if ((dcurnorm < dminnorm) || (dminnorm < 0.0)) dminnorm = dcurnorm; 376 ierr=VecAXPY(Xdot[i],1.0,XA);CHKERRQ(ierr); 377 } 378 } else { 379 ierr=VecCopy(XA, D);CHKERRQ(ierr); 380 ierr=VecAXPY(D, -1.0, XM);CHKERRQ(ierr); 381 ierr=VecNormBegin(D, NORM_2, &dnorm);CHKERRQ(ierr); 382 ierr=VecNormBegin(FA, NORM_2, &fAnorm);CHKERRQ(ierr); 383 ierr=VecNormEnd(D, NORM_2, &dnorm);CHKERRQ(ierr); 384 ierr=VecNormEnd(FA, NORM_2, &fAnorm);CHKERRQ(ierr); 385 dminnorm = -1.0; 386 for(i=0;i<l;i++) { 387 ierr=VecCopy(XA, D);CHKERRQ(ierr); 388 ierr=VecAXPY(D, -1.0, Xdot[i]);CHKERRQ(ierr); 389 ierr=VecNorm(D, NORM_2, &dcurnorm);CHKERRQ(ierr); 390 if((dcurnorm < dminnorm) || (dminnorm < 0.0)) dminnorm = dcurnorm; 391 } 392 } 393 394 /* combination (additive) or selection (multiplicative) of the N-GMRES solution */ 395 if (ngmres->additive) { 396 /* X = X + \lambda(XA - X) */ 397 if (ngmres->monitor) { 398 ierr = PetscViewerASCIIPrintf(ngmres->monitor, "||F_A||_2 = %e, ||F_M||_2 = %e\n", fAnorm, fMnorm);CHKERRQ(ierr); 399 } 400 ierr = VecCopy(XA, Y);CHKERRQ(ierr); 401 ierr = VecAYPX(Y, -1.0, X);CHKERRQ(ierr); 402 ierr = SNESLineSearchApply(ngmres->additive_linesearch,X,F,&fnorm,Y);CHKERRQ(ierr); 403 ierr = SNESLineSearchGetSuccess(ngmres->additive_linesearch, &lssucceed);CHKERRQ(ierr); 404 if (!lssucceed) { 405 if (++snes->numFailures >= snes->maxFailures) { 406 snes->reason = SNES_DIVERGED_LINE_SEARCH; 407 PetscFunctionReturn(0); 408 } 409 } 410 if (ngmres->monitor) { 411 ierr = PetscViewerASCIIPrintf(ngmres->monitor, "Additive solution: ||F||_2 = %e\n", fnorm);CHKERRQ(ierr); 412 } 413 } else { 414 selectA = PETSC_TRUE; 415 /* Conditions for choosing the accelerated answer */ 416 /* Criterion A -- the norm of the function isn't increased above the minimum by too much */ 417 if (fAnorm >= ngmres->gammaA*fminnorm) { 418 selectA = PETSC_FALSE; 419 } 420 /* Criterion B -- the choice of x^A isn't too close to some other choice */ 421 if (ngmres->epsilonB*dnorm<dminnorm || PetscSqrtReal(fnorm)<ngmres->deltaB*PetscSqrtReal(fminnorm)) { 422 } else { 423 selectA=PETSC_FALSE; 424 } 425 if (selectA) { 426 if (ngmres->monitor) { 427 ierr = PetscViewerASCIIPrintf(ngmres->monitor, "picked X_A, ||F_A||_2 = %e, ||F_M||_2 = %e\n", fAnorm, fMnorm);CHKERRQ(ierr); 428 } 429 /* copy it over */ 430 fnorm = fAnorm; 431 nu = fnorm*fnorm; 432 ierr = VecCopy(FA, F);CHKERRQ(ierr); 433 ierr = VecCopy(XA, X);CHKERRQ(ierr); 434 } else { 435 if (ngmres->monitor) { 436 ierr = PetscViewerASCIIPrintf(ngmres->monitor, "picked X_M, ||F_A||_2 = %e, ||F_M||_2 = %e\n", fAnorm, fnorm);CHKERRQ(ierr); 437 } 438 fnorm = fMnorm; 439 nu = fnorm*fnorm; 440 ierr = VecCopy(FM, F);CHKERRQ(ierr); 441 ierr = VecCopy(XM, X);CHKERRQ(ierr); 442 } 443 } 444 445 selectRestart = PETSC_FALSE; 446 447 /* difference stagnation restart */ 448 if((ngmres->epsilonB*dnorm > dminnorm) && (PetscSqrtReal(fAnorm) > ngmres->deltaB*PetscSqrtReal(fminnorm))) { 449 if (ngmres->monitor) { 450 ierr = PetscViewerASCIIPrintf(ngmres->monitor, "difference restart: %e > %e\n", ngmres->epsilonB*dnorm, dminnorm);CHKERRQ(ierr); 451 } 452 selectRestart = PETSC_TRUE; 453 } 454 /* residual stagnation restart */ 455 if (PetscSqrtReal(fAnorm) > ngmres->gammaC*PetscSqrtReal(fminnorm)) { 456 if (ngmres->monitor) { 457 ierr = PetscViewerASCIIPrintf(ngmres->monitor, "residual restart: %e > %e\n", PetscSqrtReal(fAnorm), ngmres->gammaC*PetscSqrtReal(fminnorm));CHKERRQ(ierr); 458 } 459 selectRestart = PETSC_TRUE; 460 } 461 462 /* if the restart conditions persist for more than restart_it iterations, restart. */ 463 if (selectRestart) { 464 restart_count++; 465 } else { 466 restart_count = 0; 467 } 468 469 /* restart after restart conditions have persisted for a fixed number of iterations */ 470 if (restart_count >= ngmres->restart_it) { 471 if (ngmres->monitor){ 472 ierr = PetscViewerASCIIPrintf(ngmres->monitor, "Restarted at iteration %d\n", k_restart);CHKERRQ(ierr); 473 } 474 restart_count = 0; 475 k_restart = 1; 476 l = 1; 477 /* q_{00} = nu */ 478 if (ngmres->anderson) { 479 ngmres->fnorms[0] = fMnorm; 480 nu = fMnorm*fMnorm; 481 Q(0,0) = nu; 482 /* Fdot[0] = F */ 483 ierr = VecCopy(XM, Xdot[0]);CHKERRQ(ierr); 484 ierr = VecCopy(FM, Fdot[0]);CHKERRQ(ierr); 485 } else { 486 ngmres->fnorms[0] = fnorm; 487 nu = fnorm*fnorm; 488 Q(0,0) = nu; 489 /* Fdot[0] = F */ 490 ierr = VecCopy(X, Xdot[0]);CHKERRQ(ierr); 491 ierr = VecCopy(F, Fdot[0]);CHKERRQ(ierr); 492 } 493 } else { 494 /* select the current size of the subspace */ 495 if (l < ngmres->msize) l++; 496 k_restart++; 497 /* place the current entry in the list of previous entries */ 498 if (ngmres->anderson) { 499 ierr = VecCopy(FM, Fdot[ivec]);CHKERRQ(ierr); 500 ierr = VecCopy(XM, Xdot[ivec]);CHKERRQ(ierr); 501 ngmres->fnorms[ivec] = fMnorm; 502 if (fminnorm > fMnorm) fminnorm = fMnorm; /* the minimum norm is now of FM */ 503 xi[ivec] = fMnorm*fMnorm; 504 for (i = 0; i < l; i++) { 505 Q(i, ivec) = xi[i]; 506 Q(ivec, i) = xi[i]; 507 } 508 } else { 509 ierr = VecCopy(F, Fdot[ivec]);CHKERRQ(ierr); 510 ierr = VecCopy(X, Xdot[ivec]);CHKERRQ(ierr); 511 ngmres->fnorms[ivec] = fnorm; 512 if (fminnorm > fnorm) fminnorm = fnorm; /* the minimum norm is now of FA */ 513 ierr = VecMDot(F, l, Fdot, xi);CHKERRQ(ierr); 514 for (i = 0; i < l; i++) { 515 Q(i, ivec) = xi[i]; 516 Q(ivec, i) = xi[i]; 517 } 518 } 519 } 520 521 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 522 snes->iter = k; 523 snes->norm = fnorm; 524 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 525 SNESLogConvHistory(snes, snes->norm, snes->iter); 526 ierr = SNESMonitor(snes, snes->iter, snes->norm);CHKERRQ(ierr); 527 528 ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 529 if (snes->reason) PetscFunctionReturn(0); 530 } 531 snes->reason = SNES_DIVERGED_MAX_IT; 532 PetscFunctionReturn(0); 533 } 534 535 /*MC 536 SNESNGMRES - The Nonlinear Generalized Minimum Residual method. 537 538 Level: beginner 539 540 Options Database: 541 + -snes_ngmres_additive - Use a variant that line-searches between the candidate solution and the combined solution. 542 . -snes_ngmres_anderson - Use Anderson mixing NGMRES variant which combines candidate solutions instead of actual solutions. 543 . -snes_ngmres_m - Number of stored previous solutions and residuals. 544 . -snes_ngmres_restart_it - Number of iterations the restart conditions hold before restart. 545 . -snes_ngmres_gammaA - Residual tolerance for solution selection between the candidate and combination. 546 . -snes_ngmres_gammaC - Residual tolerance for restart. 547 . -snes_ngmres_epsilonB - Difference tolerance between subsequent solutions triggering restart. 548 . -snes_ngmres_deltaB - Difference tolerance between residuals triggering restart. 549 . -snes_ngmres_monitor - Prints relevant information about the ngmres iteration. 550 - -snes_ls <basic,basicnonorms,quadratic,critical> - Line search type. 551 552 Notes: 553 554 The N-GMRES method combines m previous solutions into a minimum-residual solution by solving a small linearized 555 optimization problem at each iteration. 556 557 References: 558 559 "Krylov Subspace Acceleration of Nonlinear Multigrid with Application to Recirculating Flows", C. W. Oosterlee and T. Washio, 560 SIAM Journal on Scientific Computing, 21(5), 2000. 561 562 This is also the same as the algorithm called Anderson acceleration introduced in "D. G. Anderson. Iterative procedures for nonlinear integral equations. 563 J. Assoc. Comput. Mach., 12:547–560, 1965." 564 565 .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types) 566 M*/ 567 568 EXTERN_C_BEGIN 569 #undef __FUNCT__ 570 #define __FUNCT__ "SNESCreate_NGMRES" 571 PetscErrorCode SNESCreate_NGMRES(SNES snes) 572 { 573 SNES_NGMRES *ngmres; 574 PetscErrorCode ierr; 575 576 PetscFunctionBegin; 577 snes->ops->destroy = SNESDestroy_NGMRES; 578 snes->ops->setup = SNESSetUp_NGMRES; 579 snes->ops->setfromoptions = SNESSetFromOptions_NGMRES; 580 snes->ops->view = SNESView_NGMRES; 581 snes->ops->solve = SNESSolve_NGMRES; 582 snes->ops->reset = SNESReset_NGMRES; 583 584 snes->usespc = PETSC_TRUE; 585 snes->usesksp = PETSC_FALSE; 586 587 ierr = PetscNewLog(snes, SNES_NGMRES, &ngmres);CHKERRQ(ierr); 588 snes->data = (void*) ngmres; 589 ngmres->msize = 30; 590 591 if (!snes->tolerancesset) { 592 snes->max_funcs = 30000; 593 snes->max_its = 10000; 594 } 595 596 ngmres->additive = PETSC_FALSE; 597 ngmres->anderson = PETSC_FALSE; 598 599 ngmres->additive_linesearch = PETSC_NULL; 600 601 ngmres->restart_it = 2; 602 ngmres->gammaA = 2.0; 603 ngmres->gammaC = 2.0; 604 ngmres->deltaB = 0.9; 605 ngmres->epsilonB = 0.1; 606 607 PetscFunctionReturn(0); 608 } 609 EXTERN_C_END 610