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