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