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