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 = PetscLineSearchDestroy(&ngmres->linesearch);CHKERRQ(ierr); 15 ierr = PetscLineSearchDestroy(&ngmres->additive_linesearch);CHKERRQ(ierr); 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 ierr = PetscLineSearchCreate(((PetscObject)snes)->comm, &ngmres->linesearch);CHKERRQ(ierr); 80 ierr = PetscLineSearchSetSNES(ngmres->linesearch, snes);CHKERRQ(ierr); 81 ierr = PetscLineSearchSetType(ngmres->linesearch, PETSCLINESEARCHBASIC);CHKERRQ(ierr); 82 ierr = PetscLineSearchAppendOptionsPrefix(ngmres->linesearch, optionsprefix);CHKERRQ(ierr); 83 ierr = PetscLineSearchSetFromOptions(ngmres->linesearch);CHKERRQ(ierr); 84 85 if (ngmres->additive) { 86 ierr = PetscLineSearchCreate(((PetscObject)snes)->comm, &ngmres->additive_linesearch);CHKERRQ(ierr); 87 ierr = PetscLineSearchSetSNES(ngmres->additive_linesearch, snes);CHKERRQ(ierr); 88 ierr = PetscLineSearchSetType(ngmres->additive_linesearch, PETSCLINESEARCHL2);CHKERRQ(ierr); 89 ierr = PetscLineSearchAppendOptionsPrefix(ngmres->additive_linesearch, "additive_");CHKERRQ(ierr); 90 ierr = PetscLineSearchAppendOptionsPrefix(ngmres->additive_linesearch, optionsprefix);CHKERRQ(ierr); 91 ierr = PetscLineSearchSetFromOptions(ngmres->additive_linesearch);CHKERRQ(ierr); 92 } 93 94 ngmres->setup_called = PETSC_TRUE; 95 PetscFunctionReturn(0); 96 } 97 98 #undef __FUNCT__ 99 #define __FUNCT__ "SNESSetFromOptions_NGMRES" 100 PetscErrorCode SNESSetFromOptions_NGMRES(SNES snes) 101 { 102 SNES_NGMRES *ngmres = (SNES_NGMRES *) snes->data; 103 PetscErrorCode ierr; 104 PetscBool debug; 105 106 PetscFunctionBegin; 107 ierr = PetscOptionsHead("SNES NGMRES options");CHKERRQ(ierr); 108 ierr = PetscOptionsBool("-snes_ngmres_additive", "Use additive variant vs. choice", "SNES", ngmres->additive, &ngmres->additive, PETSC_NULL);CHKERRQ(ierr); 109 ierr = PetscOptionsBool("-snes_ngmres_anderson", "Use Anderson mixing storage", "SNES", ngmres->anderson, &ngmres->anderson, PETSC_NULL);CHKERRQ(ierr); 110 ierr = PetscOptionsInt("-snes_ngmres_m", "Number of directions", "SNES", ngmres->msize, &ngmres->msize, PETSC_NULL);CHKERRQ(ierr); 111 ierr = PetscOptionsInt("-snes_ngmres_restart_it","Tolerance iterations before restart","SNES", ngmres->restart_it, &ngmres->restart_it, PETSC_NULL);CHKERRQ(ierr); 112 ierr = PetscOptionsBool("-snes_ngmres_monitor", "Monitor actions of NGMRES", "SNES", ngmres->monitor ? PETSC_TRUE: PETSC_FALSE, &debug, PETSC_NULL);CHKERRQ(ierr); 113 if (debug) { 114 ngmres->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr); 115 } 116 ierr = PetscOptionsReal("-snes_ngmres_gammaA", "Residual selection constant", "SNES", ngmres->gammaA, &ngmres->gammaA, PETSC_NULL);CHKERRQ(ierr); 117 ierr = PetscOptionsReal("-snes_ngmres_gammaC", " Residual restart constant", "SNES", ngmres->gammaC, &ngmres->gammaC, PETSC_NULL);CHKERRQ(ierr); 118 ierr = PetscOptionsReal("-snes_ngmres_epsilonB", "Difference selection constant", "SNES", ngmres->epsilonB, &ngmres->epsilonB, PETSC_NULL);CHKERRQ(ierr); 119 ierr = PetscOptionsReal("-snes_ngmres_deltaB", "Difference residual selection constant", "SNES", ngmres->deltaB, &ngmres->deltaB, PETSC_NULL);CHKERRQ(ierr); 120 ierr = PetscOptionsTail();CHKERRQ(ierr); 121 if ((ngmres->gammaA > ngmres->gammaC) && (ngmres->gammaC > 2.)) ngmres->gammaC = ngmres->gammaA; 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 ierr = PetscViewerASCIIPrintf(viewer, " line search variant: %s\n",((PetscObject)ngmres->linesearch)->type_name);CHKERRQ(ierr); 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 ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 202 if (snes->domainerror) { 203 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 204 PetscFunctionReturn(0); 205 } 206 207 /* nu = (r, r) */ 208 ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 209 fminnorm = fnorm; 210 nu = fnorm*fnorm; 211 if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FP, "Infinite or not-a-number generated in function evaluation"); 212 213 /* q_{00} = nu */ 214 Q(0,0) = nu; 215 ngmres->fnorms[0] = fnorm; 216 /* Fdot[0] = F */ 217 ierr = VecCopy(X, Xdot[0]);CHKERRQ(ierr); 218 ierr = VecCopy(F, Fdot[0]);CHKERRQ(ierr); 219 220 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 221 snes->norm = fnorm; 222 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 223 SNESLogConvHistory(snes, fnorm, 0); 224 ierr = SNESMonitor(snes, 0, fnorm);CHKERRQ(ierr); 225 ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 226 if (snes->reason) PetscFunctionReturn(0); 227 228 k_restart = 1; 229 l = 1; 230 for (k=1; k < snes->max_its+1; k++) { 231 /* select which vector of the stored subspace will be updated */ 232 ivec = k_restart % ngmres->msize; /* replace the last used part of the subspace */ 233 234 /* Computation of x^M */ 235 if (snes->pc) { 236 ierr = VecCopy(X, XM);CHKERRQ(ierr); 237 ierr = SNESSolve(snes->pc, B, XM);CHKERRQ(ierr); 238 ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 239 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 240 snes->reason = SNES_DIVERGED_INNER; 241 PetscFunctionReturn(0); 242 } 243 ierr = SNESGetFunction(snes->pc, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr); 244 ierr = VecCopy(FPC, FM);CHKERRQ(ierr); 245 ierr = SNESGetFunctionNorm(snes->pc, &fMnorm);CHKERRQ(ierr); 246 } else { 247 /* no preconditioner -- just take gradient descent with line search */ 248 ierr = VecCopy(F, Y);CHKERRQ(ierr); 249 ierr = VecCopy(F, FM);CHKERRQ(ierr); 250 ierr = VecCopy(X, XM);CHKERRQ(ierr); 251 fMnorm = fnorm; 252 ierr = PetscLineSearchApply(ngmres->linesearch,XM,FM,&fMnorm,Y);CHKERRQ(ierr); 253 ierr = PetscLineSearchGetSuccess(ngmres->linesearch, &lssucceed);CHKERRQ(ierr); 254 if (!lssucceed) { 255 if (++snes->numFailures >= snes->maxFailures) { 256 snes->reason = SNES_DIVERGED_LINE_SEARCH; 257 PetscFunctionReturn(0); 258 } 259 } 260 } 261 262 /* r = F(x) */ 263 nu = fMnorm*fMnorm; 264 if (fminnorm > fMnorm) fminnorm = fMnorm; /* the minimum norm is now of F^M */ 265 266 /* construct the right hand side and xi factors */ 267 ierr = VecMDot(FM, l, Fdot, xi);CHKERRQ(ierr); 268 269 for (i = 0; i < l; i++) { 270 beta[i] = nu - xi[i]; 271 } 272 273 /* construct h */ 274 for (j = 0; j < l; j++) { 275 for (i = 0; i < l; i++) { 276 H(i, j) = Q(i, j) - xi[i] - xi[j] + nu; 277 } 278 } 279 280 if(l == 1) { 281 /* simply set alpha[0] = beta[0] / H[0, 0] */ 282 beta[0] = beta[0] / H(0, 0); 283 } else { 284 #ifdef PETSC_MISSING_LAPACK_GELSS 285 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "NGMRES with LS requires the LAPACK GELSS routine."); 286 #else 287 ngmres->m = PetscBLASIntCast(l); 288 ngmres->n = PetscBLASIntCast(l); 289 ngmres->info = PetscBLASIntCast(0); 290 ngmres->rcond = -1.; 291 #ifdef PETSC_USE_COMPLEX 292 LAPACKgelss_(&ngmres->m, 293 &ngmres->n, 294 &ngmres->nrhs, 295 ngmres->h, 296 &ngmres->lda, 297 ngmres->beta, 298 &ngmres->ldb, 299 ngmres->s, 300 &ngmres->rcond, 301 &ngmres->rank, 302 ngmres->work, 303 &ngmres->lwork, 304 ngmres->rwork, 305 &ngmres->info); 306 #else 307 LAPACKgelss_(&ngmres->m, 308 &ngmres->n, 309 &ngmres->nrhs, 310 ngmres->h, 311 &ngmres->lda, 312 ngmres->beta, 313 &ngmres->ldb, 314 ngmres->s, 315 &ngmres->rcond, 316 &ngmres->rank, 317 ngmres->work, 318 &ngmres->lwork, 319 &ngmres->info); 320 #endif 321 if (ngmres->info < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad argument to GELSS"); 322 if (ngmres->info > 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"SVD failed to converge"); 323 #endif 324 } 325 326 alph_total = 0.; 327 for (i = 0; i < l; i++) { 328 alph_total += beta[i]; 329 } 330 331 ierr = VecCopy(XM, XA);CHKERRQ(ierr); 332 ierr = VecScale(XA, 1. - alph_total);CHKERRQ(ierr); 333 334 ierr = VecMAXPY(XA, l, beta, Xdot);CHKERRQ(ierr); 335 ierr = SNESComputeFunction(snes, XA, FA);CHKERRQ(ierr); 336 ierr = VecNorm(FA, NORM_2, &fAnorm);CHKERRQ(ierr); 337 338 /* differences for selection and restart */ 339 ierr=VecCopy(XA, D);CHKERRQ(ierr); 340 ierr=VecAXPY(D, -1.0, X);CHKERRQ(ierr); 341 ierr=VecNorm(D, NORM_2, &dnorm);CHKERRQ(ierr); 342 dminnorm = -1.0; 343 for(i=0;i<l;i++) { 344 ierr=VecCopy(XA, D);CHKERRQ(ierr); 345 ierr=VecAXPY(D, -1.0, Xdot[i]);CHKERRQ(ierr); 346 ierr=VecNorm(D, NORM_2, &dcurnorm);CHKERRQ(ierr); 347 if((dcurnorm < dminnorm) || (dminnorm < 0.0)) dminnorm = dcurnorm; 348 } 349 350 /* combination (additive) or selection (multiplicative) of the N-GMRES solution */ 351 if (ngmres->additive) { 352 /* X = X + \lambda(XA - X) */ 353 if (ngmres->monitor) { 354 ierr = PetscViewerASCIIPrintf(ngmres->monitor, "||F_A||_2 = %e, ||F_M||_2 = %e\n", fAnorm, fMnorm);CHKERRQ(ierr); 355 } 356 ierr = VecCopy(XA, Y);CHKERRQ(ierr); 357 ierr = VecAYPX(Y, -1.0, X);CHKERRQ(ierr); 358 ierr = PetscLineSearchApply(ngmres->additive_linesearch,X,F,&fnorm,Y);CHKERRQ(ierr); 359 ierr = PetscLineSearchGetSuccess(ngmres->additive_linesearch, &lssucceed);CHKERRQ(ierr); 360 if (!lssucceed) { 361 if (++snes->numFailures >= snes->maxFailures) { 362 snes->reason = SNES_DIVERGED_LINE_SEARCH; 363 PetscFunctionReturn(0); 364 } 365 } 366 ierr = PetscLineSearchGetNorms(ngmres->additive_linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr); 367 if (ngmres->monitor) { 368 ierr = PetscViewerASCIIPrintf(ngmres->monitor, "Additive solution: ||F||_2 = %e\n", fnorm);CHKERRQ(ierr); 369 } 370 ierr = VecCopy(FA, F);CHKERRQ(ierr); 371 ierr = VecCopy(XA, X);CHKERRQ(ierr); 372 } else { 373 selectA = PETSC_TRUE; 374 /* Conditions for choosing the accelerated answer */ 375 /* Criterion A -- the norm of the function isn't increased above the minimum by too much */ 376 if (fAnorm >= ngmres->gammaA*fminnorm) { 377 selectA = PETSC_FALSE; 378 } 379 /* Criterion B -- the choice of x^A isn't too close to some other choice */ 380 if (ngmres->epsilonB*dnorm<dminnorm || PetscSqrtReal(fnorm)<ngmres->deltaB*PetscSqrtReal(fminnorm)) { 381 } else { 382 selectA=PETSC_FALSE; 383 } 384 if (selectA) { 385 if (ngmres->monitor) { 386 ierr = PetscViewerASCIIPrintf(ngmres->monitor, "picked X_A, ||F_A||_2 = %e, ||F_M||_2 = %e\n", fAnorm, fnorm);CHKERRQ(ierr); 387 } 388 /* copy it over */ 389 fnorm = fAnorm; 390 nu = fnorm*fnorm; 391 ierr = VecCopy(FA, F);CHKERRQ(ierr); 392 ierr = VecCopy(XA, X);CHKERRQ(ierr); 393 } else { 394 if (ngmres->monitor) { 395 ierr = PetscViewerASCIIPrintf(ngmres->monitor, "picked X_M, ||F_A||_2 = %e, ||F_M||_2 = %e\n", fAnorm, fnorm);CHKERRQ(ierr); 396 } 397 fnorm = fMnorm; 398 nu = fnorm*fnorm; 399 ierr = VecCopy(FM, F);CHKERRQ(ierr); 400 ierr = VecCopy(XM, X);CHKERRQ(ierr); 401 } 402 } 403 404 selectRestart = PETSC_FALSE; 405 406 /* difference stagnation restart */ 407 if((ngmres->epsilonB*dnorm > dminnorm) && (PetscSqrtReal(fAnorm) > ngmres->deltaB*PetscSqrtReal(fminnorm))) { 408 if (ngmres->monitor) { 409 ierr = PetscViewerASCIIPrintf(ngmres->monitor, "difference restart: %e > %e\n", ngmres->epsilonB*dnorm, dminnorm);CHKERRQ(ierr); 410 } 411 selectRestart = PETSC_TRUE; 412 } 413 /* residual stagnation restart */ 414 if (PetscSqrtReal(fAnorm) > ngmres->gammaC*PetscSqrtReal(fminnorm)) { 415 if (ngmres->monitor) { 416 ierr = PetscViewerASCIIPrintf(ngmres->monitor, "residual restart: %e > %e\n", PetscSqrtReal(fAnorm), ngmres->gammaC*PetscSqrtReal(fminnorm));CHKERRQ(ierr); 417 } 418 selectRestart = PETSC_TRUE; 419 } 420 421 /* if the restart conditions persist for more than restart_it iterations, restart. */ 422 if (selectRestart) { 423 restart_count++; 424 } else { 425 restart_count = 0; 426 } 427 428 /* restart after restart conditions have persisted for a fixed number of iterations */ 429 if (restart_count >= ngmres->restart_it) { 430 if (ngmres->monitor){ 431 ierr = PetscViewerASCIIPrintf(ngmres->monitor, "Restarted at iteration %d\n", k_restart);CHKERRQ(ierr); 432 } 433 restart_count = 0; 434 k_restart = 1; 435 l = 1; 436 /* q_{00} = nu */ 437 if (ngmres->anderson) { 438 ngmres->fnorms[0] = fMnorm; 439 nu = fMnorm*fMnorm; 440 Q(0,0) = nu; 441 /* Fdot[0] = F */ 442 ierr = VecCopy(XM, Xdot[0]);CHKERRQ(ierr); 443 ierr = VecCopy(FM, Fdot[0]);CHKERRQ(ierr); 444 } else { 445 ngmres->fnorms[0] = fnorm; 446 nu = fnorm*fnorm; 447 Q(0,0) = nu; 448 /* Fdot[0] = F */ 449 ierr = VecCopy(X, Xdot[0]);CHKERRQ(ierr); 450 ierr = VecCopy(F, Fdot[0]);CHKERRQ(ierr); 451 } 452 } else { 453 /* select the current size of the subspace */ 454 if (l < ngmres->msize) l++; 455 k_restart++; 456 /* place the current entry in the list of previous entries */ 457 if (ngmres->anderson) { 458 ierr = VecCopy(FM, Fdot[ivec]);CHKERRQ(ierr); 459 ierr = VecCopy(XM, Xdot[ivec]);CHKERRQ(ierr); 460 ngmres->fnorms[ivec] = fMnorm; 461 if (fminnorm > fMnorm) fminnorm = fMnorm; /* the minimum norm is now of FM */ 462 for (i = 0; i < l; i++) { 463 ierr = VecDot(FM, Fdot[i], &qentry);CHKERRQ(ierr); 464 Q(i, ivec) = qentry; 465 Q(ivec, i) = qentry; 466 } 467 } else { 468 ierr = VecCopy(F, Fdot[ivec]);CHKERRQ(ierr); 469 ierr = VecCopy(X, Xdot[ivec]);CHKERRQ(ierr); 470 ngmres->fnorms[ivec] = fnorm; 471 if (fminnorm > fnorm) fminnorm = fnorm; /* the minimum norm is now of FA */ 472 for (i = 0; i < l; i++) { 473 ierr = VecDot(F, Fdot[i], &qentry);CHKERRQ(ierr); 474 Q(i, ivec) = qentry; 475 Q(ivec, i) = qentry; 476 } 477 } 478 } 479 480 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 481 snes->iter = k; 482 snes->norm = fnorm; 483 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 484 SNESLogConvHistory(snes, snes->norm, snes->iter); 485 ierr = SNESMonitor(snes, snes->iter, snes->norm);CHKERRQ(ierr); 486 487 ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 488 if (snes->reason) PetscFunctionReturn(0); 489 } 490 snes->reason = SNES_DIVERGED_MAX_IT; 491 PetscFunctionReturn(0); 492 } 493 494 /*MC 495 SNESNGMRES - The Nonlinear Generalized Minimum Residual method. 496 497 Level: beginner 498 499 Options Database: 500 + -snes_ngmres_additive - Use a variant that line-searches between the candidate solution and the combined solution. 501 . -snes_ngmres_anderson - Use Anderson mixing NGMRES variant which combines candidate solutions instead of actual solutions. 502 . -snes_ngmres_m - Number of stored previous solutions and residuals. 503 . -snes_ngmres_restart_it - Number of iterations the restart conditions hold before restart. 504 . -snes_ngmres_gammaA - Residual tolerance for solution selection between the candidate and combination. 505 . -snes_ngmres_gammaC - Residual tolerance for restart. 506 . -snes_ngmres_epsilonB - Difference tolerance between subsequent solutions triggering restart. 507 . -snes_ngmres_deltaB - Difference tolerance between residuals triggering restart. 508 . -snes_ngmres_monitor - Prints relevant information about the ngmres iteration. 509 - -snes_ls <basic,basicnonorms,quadratic,critical> - Line search type. 510 511 Notes: 512 513 The N-GMRES method combines m previous solutions into a minimum-residual solution by solving a small linearized 514 optimization problem at each iteration. 515 516 References: 517 518 "Krylov Subspace Acceleration of Nonlinear Multigrid with Application to Recirculating Flows", C. W. Oosterlee and T. Washio, 519 SIAM Journal on Scientific Computing, 21(5), 2000. 520 521 This is also the same as the algorithm called Anderson acceleration introduced in "D. G. Anderson. Iterative procedures for nonlinear integral equations. 522 J. Assoc. Comput. Mach., 12:547–560, 1965." 523 524 .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types) 525 M*/ 526 527 EXTERN_C_BEGIN 528 #undef __FUNCT__ 529 #define __FUNCT__ "SNESCreate_NGMRES" 530 PetscErrorCode SNESCreate_NGMRES(SNES snes) 531 { 532 SNES_NGMRES *ngmres; 533 PetscErrorCode ierr; 534 535 PetscFunctionBegin; 536 snes->ops->destroy = SNESDestroy_NGMRES; 537 snes->ops->setup = SNESSetUp_NGMRES; 538 snes->ops->setfromoptions = SNESSetFromOptions_NGMRES; 539 snes->ops->view = SNESView_NGMRES; 540 snes->ops->solve = SNESSolve_NGMRES; 541 snes->ops->reset = SNESReset_NGMRES; 542 543 snes->usespc = PETSC_TRUE; 544 snes->usesksp = PETSC_FALSE; 545 546 ierr = PetscNewLog(snes, SNES_NGMRES, &ngmres);CHKERRQ(ierr); 547 snes->data = (void*) ngmres; 548 ngmres->msize = 30; 549 550 snes->max_funcs = 30000; 551 snes->max_its = 10000; 552 553 ngmres->additive = PETSC_FALSE; 554 ngmres->anderson = PETSC_FALSE; 555 556 ngmres->linesearch = PETSC_NULL; 557 ngmres->additive_linesearch = PETSC_NULL; 558 559 ngmres->restart_it = 2; 560 ngmres->gammaA = 2.0; 561 ngmres->gammaC = 2.0; 562 ngmres->deltaB = 0.9; 563 ngmres->epsilonB = 0.1; 564 565 PetscFunctionReturn(0); 566 } 567 EXTERN_C_END 568