119653cdaSPeter Brune #include <../src/snes/impls/ngmres/snesngmres.h> 219653cdaSPeter Brune #include <petscblaslapack.h> 3a312c225SMatthew G Knepley 4a312c225SMatthew G Knepley #undef __FUNCT__ 5a312c225SMatthew G Knepley #define __FUNCT__ "SNESReset_NGMRES" 6a312c225SMatthew G Knepley PetscErrorCode SNESReset_NGMRES(SNES snes) 7a312c225SMatthew G Knepley { 8a312c225SMatthew G Knepley SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data; 9a312c225SMatthew G Knepley PetscErrorCode ierr; 10a312c225SMatthew G Knepley 11a312c225SMatthew G Knepley PetscFunctionBegin; 12f109b39eSPeter Brune ierr = VecDestroyVecs(ngmres->msize, &ngmres->Fdot);CHKERRQ(ierr); 13f109b39eSPeter Brune ierr = VecDestroyVecs(ngmres->msize, &ngmres->Xdot);CHKERRQ(ierr); 14f1c6b773SPeter Brune ierr = SNESLineSearchDestroy(&ngmres->additive_linesearch);CHKERRQ(ierr); 159e764e56SPeter Brune 16a312c225SMatthew G Knepley PetscFunctionReturn(0); 17a312c225SMatthew G Knepley } 18a312c225SMatthew G Knepley 19a312c225SMatthew G Knepley #undef __FUNCT__ 20a312c225SMatthew G Knepley #define __FUNCT__ "SNESDestroy_NGMRES" 21a312c225SMatthew G Knepley PetscErrorCode SNESDestroy_NGMRES(SNES snes) 22a312c225SMatthew G Knepley { 23a312c225SMatthew G Knepley PetscErrorCode ierr; 2478440776SJed Brown SNES_NGMRES *ngmres = (SNES_NGMRES*)snes->data; 25a312c225SMatthew G Knepley 26a312c225SMatthew G Knepley PetscFunctionBegin; 27a312c225SMatthew G Knepley ierr = SNESReset_NGMRES(snes);CHKERRQ(ierr); 28f109b39eSPeter Brune ierr = PetscFree5(ngmres->h, ngmres->beta, ngmres->xi, ngmres->fnorms, ngmres->q);CHKERRQ(ierr); 2919653cdaSPeter Brune ierr = PetscFree(ngmres->s);CHKERRQ(ierr); 3019653cdaSPeter Brune #if PETSC_USE_COMPLEX 3119653cdaSPeter Brune ierr = PetscFree(ngmres->rwork); 3219653cdaSPeter Brune #endif 3319653cdaSPeter Brune ierr = PetscFree(ngmres->work); 3419653cdaSPeter Brune ierr = PetscFree(snes->data); 35a312c225SMatthew G Knepley PetscFunctionReturn(0); 36a312c225SMatthew G Knepley } 37a312c225SMatthew G Knepley 38a312c225SMatthew G Knepley #undef __FUNCT__ 39a312c225SMatthew G Knepley #define __FUNCT__ "SNESSetUp_NGMRES" 40a312c225SMatthew G Knepley PetscErrorCode SNESSetUp_NGMRES(SNES snes) 41a312c225SMatthew G Knepley { 42a312c225SMatthew G Knepley SNES_NGMRES *ngmres = (SNES_NGMRES *) snes->data; 43e7058c64SPeter Brune const char *optionsprefix; 4419653cdaSPeter Brune PetscInt msize,hsize; 45a312c225SMatthew G Knepley PetscErrorCode ierr; 46a312c225SMatthew G Knepley 47a312c225SMatthew G Knepley PetscFunctionBegin; 4878440776SJed Brown ierr = SNESDefaultGetWork(snes,5);CHKERRQ(ierr); 4978440776SJed Brown if (!ngmres->Xdot) {ierr = VecDuplicateVecs(snes->vec_sol,ngmres->msize,&ngmres->Xdot);CHKERRQ(ierr);} 5078440776SJed Brown if (!ngmres->Fdot) {ierr = VecDuplicateVecs(snes->vec_sol,ngmres->msize,&ngmres->Fdot);CHKERRQ(ierr);} 5178440776SJed Brown if (!ngmres->setup_called) { 52087dfb9eSxuemin msize = ngmres->msize; /* restart size */ 5319653cdaSPeter Brune hsize = msize * msize; 54087dfb9eSxuemin 5598b3e84cSPeter Brune /* explicit least squares minimization solve */ 5619653cdaSPeter Brune ierr = PetscMalloc5(hsize,PetscScalar,&ngmres->h, 5719653cdaSPeter Brune msize,PetscScalar,&ngmres->beta, 5819653cdaSPeter Brune msize,PetscScalar,&ngmres->xi, 59f109b39eSPeter Brune msize,PetscReal, &ngmres->fnorms, 6019653cdaSPeter Brune hsize,PetscScalar,&ngmres->q);CHKERRQ(ierr); 6119653cdaSPeter Brune ngmres->nrhs = 1; 6219653cdaSPeter Brune ngmres->lda = msize; 6319653cdaSPeter Brune ngmres->ldb = msize; 6419653cdaSPeter Brune ierr = PetscMalloc(msize*sizeof(PetscScalar),&ngmres->s);CHKERRQ(ierr); 6519653cdaSPeter Brune ierr = PetscMemzero(ngmres->h, hsize*sizeof(PetscScalar));CHKERRQ(ierr); 6619653cdaSPeter Brune ierr = PetscMemzero(ngmres->q, hsize*sizeof(PetscScalar));CHKERRQ(ierr); 6719653cdaSPeter Brune ierr = PetscMemzero(ngmres->xi, msize*sizeof(PetscScalar));CHKERRQ(ierr); 6819653cdaSPeter Brune ierr = PetscMemzero(ngmres->beta,msize*sizeof(PetscScalar));CHKERRQ(ierr); 6919653cdaSPeter Brune ngmres->lwork = 12*msize; 7019653cdaSPeter Brune #if PETSC_USE_COMPLEX 7119653cdaSPeter Brune ierr = PetscMalloc(sizeof(PetscReal)*ngmres->lwork,&ngmres->rwork); 7219653cdaSPeter Brune #endif 7319653cdaSPeter Brune ierr = PetscMalloc(sizeof(PetscScalar)*ngmres->lwork,&ngmres->work); 7478440776SJed Brown } 75e7058c64SPeter Brune 76e7058c64SPeter Brune /* linesearch setup */ 77e7058c64SPeter Brune ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr); 78e7058c64SPeter Brune 79e7058c64SPeter Brune if (ngmres->additive) { 80f1c6b773SPeter Brune ierr = SNESLineSearchCreate(((PetscObject)snes)->comm, &ngmres->additive_linesearch);CHKERRQ(ierr); 81f1c6b773SPeter Brune ierr = SNESLineSearchSetSNES(ngmres->additive_linesearch, snes);CHKERRQ(ierr); 82f1c6b773SPeter Brune ierr = SNESLineSearchSetType(ngmres->additive_linesearch, SNES_LINESEARCH_L2);CHKERRQ(ierr); 83f1c6b773SPeter Brune ierr = SNESLineSearchAppendOptionsPrefix(ngmres->additive_linesearch, "additive_");CHKERRQ(ierr); 84f1c6b773SPeter Brune ierr = SNESLineSearchAppendOptionsPrefix(ngmres->additive_linesearch, optionsprefix);CHKERRQ(ierr); 85f1c6b773SPeter Brune ierr = SNESLineSearchSetFromOptions(ngmres->additive_linesearch);CHKERRQ(ierr); 86e7058c64SPeter Brune } 87e7058c64SPeter Brune 8878440776SJed Brown ngmres->setup_called = PETSC_TRUE; 89a312c225SMatthew G Knepley PetscFunctionReturn(0); 90a312c225SMatthew G Knepley } 91a312c225SMatthew G Knepley 92a312c225SMatthew G Knepley #undef __FUNCT__ 93a312c225SMatthew G Knepley #define __FUNCT__ "SNESSetFromOptions_NGMRES" 94a312c225SMatthew G Knepley PetscErrorCode SNESSetFromOptions_NGMRES(SNES snes) 95a312c225SMatthew G Knepley { 96a312c225SMatthew G Knepley SNES_NGMRES *ngmres = (SNES_NGMRES *) snes->data; 97a312c225SMatthew G Knepley PetscErrorCode ierr; 98dfbf837cSBarry Smith PetscBool debug; 99f1c6b773SPeter Brune SNESLineSearch linesearch; 100a312c225SMatthew G Knepley PetscFunctionBegin; 101a312c225SMatthew G Knepley ierr = PetscOptionsHead("SNES NGMRES options");CHKERRQ(ierr); 1029f425c49SPeter Brune ierr = PetscOptionsBool("-snes_ngmres_additive", "Use additive variant vs. choice", "SNES", ngmres->additive, &ngmres->additive, PETSC_NULL);CHKERRQ(ierr); 103d2e16ddcSPeter Brune ierr = PetscOptionsBool("-snes_ngmres_anderson", "Use Anderson mixing storage", "SNES", ngmres->anderson, &ngmres->anderson, PETSC_NULL);CHKERRQ(ierr); 10419653cdaSPeter Brune ierr = PetscOptionsInt("-snes_ngmres_m", "Number of directions", "SNES", ngmres->msize, &ngmres->msize, PETSC_NULL);CHKERRQ(ierr); 10528ed4a04SPeter Brune ierr = PetscOptionsInt("-snes_ngmres_restart_it","Tolerance iterations before restart","SNES", ngmres->restart_it, &ngmres->restart_it, PETSC_NULL);CHKERRQ(ierr); 106dfbf837cSBarry Smith ierr = PetscOptionsBool("-snes_ngmres_monitor", "Monitor actions of NGMRES", "SNES", ngmres->monitor ? PETSC_TRUE: PETSC_FALSE, &debug, PETSC_NULL);CHKERRQ(ierr); 107dfbf837cSBarry Smith if (debug) { 108dfbf837cSBarry Smith ngmres->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr); 109dfbf837cSBarry Smith } 1106a7cf640SPeter Brune ierr = PetscOptionsReal("-snes_ngmres_gammaA", "Residual selection constant", "SNES", ngmres->gammaA, &ngmres->gammaA, PETSC_NULL);CHKERRQ(ierr); 1116a7cf640SPeter Brune ierr = PetscOptionsReal("-snes_ngmres_gammaC", " Residual restart constant", "SNES", ngmres->gammaC, &ngmres->gammaC, PETSC_NULL);CHKERRQ(ierr); 1126a7cf640SPeter Brune ierr = PetscOptionsReal("-snes_ngmres_epsilonB", "Difference selection constant", "SNES", ngmres->epsilonB, &ngmres->epsilonB, PETSC_NULL);CHKERRQ(ierr); 1136a7cf640SPeter Brune ierr = PetscOptionsReal("-snes_ngmres_deltaB", "Difference residual selection constant", "SNES", ngmres->deltaB, &ngmres->deltaB, PETSC_NULL);CHKERRQ(ierr); 114a312c225SMatthew G Knepley ierr = PetscOptionsTail();CHKERRQ(ierr); 1156a7cf640SPeter Brune if ((ngmres->gammaA > ngmres->gammaC) && (ngmres->gammaC > 2.)) ngmres->gammaC = ngmres->gammaA; 1169e764e56SPeter Brune 1179e764e56SPeter Brune /* set the default type of the line search if the user hasn't already. */ 1189e764e56SPeter Brune if (!snes->linesearch) { 119f1c6b773SPeter Brune ierr = SNESGetSNESLineSearch(snes, &linesearch);CHKERRQ(ierr); 120f1c6b773SPeter Brune ierr = SNESLineSearchSetType(linesearch, SNES_LINESEARCH_BASIC);CHKERRQ(ierr); 1219e764e56SPeter Brune } 122a312c225SMatthew G Knepley PetscFunctionReturn(0); 123a312c225SMatthew G Knepley } 124a312c225SMatthew G Knepley 125a312c225SMatthew G Knepley #undef __FUNCT__ 126a312c225SMatthew G Knepley #define __FUNCT__ "SNESView_NGMRES" 127a312c225SMatthew G Knepley PetscErrorCode SNESView_NGMRES(SNES snes, PetscViewer viewer) 128a312c225SMatthew G Knepley { 129a312c225SMatthew G Knepley SNES_NGMRES *ngmres = (SNES_NGMRES *) snes->data; 130a312c225SMatthew G Knepley PetscBool iascii; 131a312c225SMatthew G Knepley PetscErrorCode ierr; 132a312c225SMatthew G Knepley 133a312c225SMatthew G Knepley PetscFunctionBegin; 134a312c225SMatthew G Knepley ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 135a312c225SMatthew G Knepley if (iascii) { 1369e764e56SPeter Brune 137f109b39eSPeter Brune ierr = PetscViewerASCIIPrintf(viewer, " Number of stored past updates: %d\n", ngmres->msize);CHKERRQ(ierr); 138f109b39eSPeter Brune ierr = PetscViewerASCIIPrintf(viewer, " Residual selection: gammaA=%1.0e, gammaC=%1.0e\n", ngmres->gammaA, ngmres->gammaC);CHKERRQ(ierr); 139f109b39eSPeter Brune ierr = PetscViewerASCIIPrintf(viewer, " Difference restart: epsilonB=%1.0e, deltaB=%1.0e\n", ngmres->epsilonB, ngmres->deltaB);CHKERRQ(ierr); 140a312c225SMatthew G Knepley } 141a312c225SMatthew G Knepley PetscFunctionReturn(0); 142a312c225SMatthew G Knepley } 143a312c225SMatthew G Knepley 144a312c225SMatthew G Knepley #undef __FUNCT__ 145a312c225SMatthew G Knepley #define __FUNCT__ "SNESSolve_NGMRES" 146211b2d39SPeter Brune 147a312c225SMatthew G Knepley PetscErrorCode SNESSolve_NGMRES(SNES snes) 148a312c225SMatthew G Knepley { 149087dfb9eSxuemin SNES_NGMRES *ngmres = (SNES_NGMRES *) snes->data; 15098b3e84cSPeter Brune /* present solution, residual, and preconditioned residual */ 1519f425c49SPeter Brune Vec X, F, B, D, Y; 152f109b39eSPeter Brune 153f109b39eSPeter Brune /* candidate linear combination answers */ 1544b32a720SPeter Brune Vec XA, FA, XM, FM, FPC; 15519653cdaSPeter Brune 15698b3e84cSPeter Brune /* previous iterations to construct the subspace */ 157f109b39eSPeter Brune Vec *Fdot = ngmres->Fdot; 158f109b39eSPeter Brune Vec *Xdot = ngmres->Xdot; 15919653cdaSPeter Brune 16098b3e84cSPeter Brune /* coefficients and RHS to the minimization problem */ 16119653cdaSPeter Brune PetscScalar *beta = ngmres->beta; 16219653cdaSPeter Brune PetscScalar *xi = ngmres->xi; 1639f425c49SPeter Brune PetscReal fnorm, fMnorm, fAnorm, ynorm, xnorm = 0.0; 16419653cdaSPeter Brune PetscReal nu; 16519653cdaSPeter Brune PetscScalar alph_total = 0.; 16619653cdaSPeter Brune PetscScalar qentry; 16728ed4a04SPeter Brune PetscInt i, j, k, k_restart, l, ivec, restart_count = 0; 16819653cdaSPeter Brune 16998b3e84cSPeter Brune /* solution selection data */ 17019653cdaSPeter Brune PetscBool selectA, selectRestart; 171f109b39eSPeter Brune PetscReal dnorm, dminnorm, dcurnorm; 172f109b39eSPeter Brune PetscReal fminnorm; 17319653cdaSPeter Brune 1741e633543SBarry Smith SNESConvergedReason reason; 175f109b39eSPeter Brune PetscBool lssucceed; 176a312c225SMatthew G Knepley PetscErrorCode ierr; 177a312c225SMatthew G Knepley 178a312c225SMatthew G Knepley PetscFunctionBegin; 17998b3e84cSPeter Brune /* variable initialization */ 180a312c225SMatthew G Knepley snes->reason = SNES_CONVERGED_ITERATING; 181f109b39eSPeter Brune X = snes->vec_sol; 182f109b39eSPeter Brune F = snes->vec_func; 183f109b39eSPeter Brune B = snes->vec_rhs; 184f109b39eSPeter Brune XA = snes->vec_sol_update; 185f109b39eSPeter Brune FA = snes->work[0]; 186f109b39eSPeter Brune D = snes->work[1]; 187f109b39eSPeter Brune 188f109b39eSPeter Brune /* work for the line search */ 189f109b39eSPeter Brune Y = snes->work[2]; 1909f425c49SPeter Brune XM = snes->work[3]; 1919f425c49SPeter Brune FM = snes->work[4]; 192a312c225SMatthew G Knepley 193a312c225SMatthew G Knepley ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 194a312c225SMatthew G Knepley snes->iter = 0; 195a312c225SMatthew G Knepley snes->norm = 0.; 196a312c225SMatthew G Knepley ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 19719653cdaSPeter Brune 19898b3e84cSPeter Brune /* initialization */ 19919653cdaSPeter Brune 20098b3e84cSPeter Brune /* r = F(x) */ 201*e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 202f109b39eSPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 203a312c225SMatthew G Knepley if (snes->domainerror) { 204a312c225SMatthew G Knepley snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 205a312c225SMatthew G Knepley PetscFunctionReturn(0); 206a312c225SMatthew G Knepley } 207*e4ed7901SPeter Brune } else { 208*e4ed7901SPeter Brune snes->vec_func_init_set = PETSC_FALSE; 209*e4ed7901SPeter Brune } 21019653cdaSPeter Brune 211*e4ed7901SPeter Brune if (!snes->norm_init_set) { 212f109b39eSPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 213f109b39eSPeter Brune if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FP, "Infinite or not-a-number generated in function evaluation"); 214*e4ed7901SPeter Brune } else { 215*e4ed7901SPeter Brune fnorm = snes->norm_init; 216*e4ed7901SPeter Brune snes->norm_init_set = PETSC_FALSE; 217*e4ed7901SPeter Brune } 218*e4ed7901SPeter Brune 219*e4ed7901SPeter Brune fminnorm = fnorm; 220*e4ed7901SPeter Brune /* nu = (r, r) */ 221*e4ed7901SPeter Brune nu = fnorm*fnorm; 22219653cdaSPeter Brune 22398b3e84cSPeter Brune /* q_{00} = nu */ 22419653cdaSPeter Brune Q(0,0) = nu; 225f109b39eSPeter Brune ngmres->fnorms[0] = fnorm; 226f109b39eSPeter Brune /* Fdot[0] = F */ 227f109b39eSPeter Brune ierr = VecCopy(X, Xdot[0]);CHKERRQ(ierr); 228f109b39eSPeter Brune ierr = VecCopy(F, Fdot[0]);CHKERRQ(ierr); 229087dfb9eSxuemin 230a312c225SMatthew G Knepley ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 231f109b39eSPeter Brune snes->norm = fnorm; 232a312c225SMatthew G Knepley ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 233f109b39eSPeter Brune SNESLogConvHistory(snes, fnorm, 0); 234f109b39eSPeter Brune ierr = SNESMonitor(snes, 0, fnorm);CHKERRQ(ierr); 235f109b39eSPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 236a312c225SMatthew G Knepley if (snes->reason) PetscFunctionReturn(0); 237a312c225SMatthew G Knepley 23819653cdaSPeter Brune k_restart = 1; 23919653cdaSPeter Brune l = 1; 24009c08436SPeter Brune for (k=1; k < snes->max_its+1; k++) { 24198b3e84cSPeter Brune /* select which vector of the stored subspace will be updated */ 24298b3e84cSPeter Brune ivec = k_restart % ngmres->msize; /* replace the last used part of the subspace */ 24319653cdaSPeter Brune 24498b3e84cSPeter Brune /* Computation of x^M */ 2458cc86e31SPeter Brune if (snes->pc) { 2469f425c49SPeter Brune ierr = VecCopy(X, XM);CHKERRQ(ierr); 247*e4ed7901SPeter Brune ierr = SNESSetInitialFunction(snes->pc, F);CHKERRQ(ierr); 248*e4ed7901SPeter Brune ierr = SNESSetInitialFunctionNorm(snes->pc, fnorm);CHKERRQ(ierr); 2499f425c49SPeter Brune ierr = SNESSolve(snes->pc, B, XM);CHKERRQ(ierr); 2508cc86e31SPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 2518cc86e31SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 2528cc86e31SPeter Brune snes->reason = SNES_DIVERGED_INNER; 2538cc86e31SPeter Brune PetscFunctionReturn(0); 2548cc86e31SPeter Brune } 2554b32a720SPeter Brune ierr = SNESGetFunction(snes->pc, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr); 2564b32a720SPeter Brune ierr = VecCopy(FPC, FM);CHKERRQ(ierr); 2574b32a720SPeter Brune ierr = SNESGetFunctionNorm(snes->pc, &fMnorm);CHKERRQ(ierr); 2588cc86e31SPeter Brune } else { 259f109b39eSPeter Brune /* no preconditioner -- just take gradient descent with line search */ 260f109b39eSPeter Brune ierr = VecCopy(F, Y);CHKERRQ(ierr); 261e7058c64SPeter Brune ierr = VecCopy(F, FM);CHKERRQ(ierr); 262e7058c64SPeter Brune ierr = VecCopy(X, XM);CHKERRQ(ierr); 263e7058c64SPeter Brune fMnorm = fnorm; 264f1c6b773SPeter Brune ierr = SNESLineSearchApply(snes->linesearch,XM,FM,&fMnorm,Y);CHKERRQ(ierr); 265f1c6b773SPeter Brune ierr = SNESLineSearchGetSuccess(snes->linesearch, &lssucceed);CHKERRQ(ierr); 266f109b39eSPeter Brune if (!lssucceed) { 267f109b39eSPeter Brune if (++snes->numFailures >= snes->maxFailures) { 268f109b39eSPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 269f109b39eSPeter Brune PetscFunctionReturn(0); 270f109b39eSPeter Brune } 271f109b39eSPeter Brune } 2726634f59bSPeter Brune } 2731e633543SBarry Smith 27498b3e84cSPeter Brune /* r = F(x) */ 2759f425c49SPeter Brune nu = fMnorm*fMnorm; 2769f425c49SPeter Brune if (fminnorm > fMnorm) fminnorm = fMnorm; /* the minimum norm is now of F^M */ 27719653cdaSPeter Brune 27898b3e84cSPeter Brune /* construct the right hand side and xi factors */ 2798c09b6b7SPeter Brune ierr = VecMDot(FM, l, Fdot, xi);CHKERRQ(ierr); 2808c09b6b7SPeter Brune 28119653cdaSPeter Brune for (i = 0; i < l; i++) { 28219653cdaSPeter Brune beta[i] = nu - xi[i]; 28319653cdaSPeter Brune } 28419653cdaSPeter Brune 28598b3e84cSPeter Brune /* construct h */ 28619653cdaSPeter Brune for (j = 0; j < l; j++) { 28719653cdaSPeter Brune for (i = 0; i < l; i++) { 28819653cdaSPeter Brune H(i, j) = Q(i, j) - xi[i] - xi[j] + nu; 28919653cdaSPeter Brune } 29019653cdaSPeter Brune } 291f109b39eSPeter Brune 292f109b39eSPeter Brune if(l == 1) { 293f109b39eSPeter Brune /* simply set alpha[0] = beta[0] / H[0, 0] */ 294f109b39eSPeter Brune beta[0] = beta[0] / H(0, 0); 295f109b39eSPeter Brune } else { 29619653cdaSPeter Brune #ifdef PETSC_MISSING_LAPACK_GELSS 29719653cdaSPeter Brune SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "NGMRES with LS requires the LAPACK GELSS routine."); 2984a0c5b0cSMatthew G Knepley #else 29919653cdaSPeter Brune ngmres->m = PetscBLASIntCast(l); 30019653cdaSPeter Brune ngmres->n = PetscBLASIntCast(l); 30119653cdaSPeter Brune ngmres->info = PetscBLASIntCast(0); 30219653cdaSPeter Brune ngmres->rcond = -1.; 3038e57ea43SSatish Balay ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 30419653cdaSPeter Brune #ifdef PETSC_USE_COMPLEX 30519653cdaSPeter Brune LAPACKgelss_(&ngmres->m, 30619653cdaSPeter Brune &ngmres->n, 30719653cdaSPeter Brune &ngmres->nrhs, 30819653cdaSPeter Brune ngmres->h, 30919653cdaSPeter Brune &ngmres->lda, 31019653cdaSPeter Brune ngmres->beta, 31119653cdaSPeter Brune &ngmres->ldb, 31219653cdaSPeter Brune ngmres->s, 31319653cdaSPeter Brune &ngmres->rcond, 31419653cdaSPeter Brune &ngmres->rank, 31519653cdaSPeter Brune ngmres->work, 31619653cdaSPeter Brune &ngmres->lwork, 31719653cdaSPeter Brune ngmres->rwork, 31819653cdaSPeter Brune &ngmres->info); 31919653cdaSPeter Brune #else 32019653cdaSPeter Brune LAPACKgelss_(&ngmres->m, 32119653cdaSPeter Brune &ngmres->n, 32219653cdaSPeter Brune &ngmres->nrhs, 32319653cdaSPeter Brune ngmres->h, 32419653cdaSPeter Brune &ngmres->lda, 32519653cdaSPeter Brune ngmres->beta, 32619653cdaSPeter Brune &ngmres->ldb, 32719653cdaSPeter Brune ngmres->s, 32819653cdaSPeter Brune &ngmres->rcond, 32919653cdaSPeter Brune &ngmres->rank, 33019653cdaSPeter Brune ngmres->work, 33119653cdaSPeter Brune &ngmres->lwork, 33219653cdaSPeter Brune &ngmres->info); 33319653cdaSPeter Brune #endif 3348e57ea43SSatish Balay ierr = PetscFPTrapPop();CHKERRQ(ierr); 33519653cdaSPeter Brune if (ngmres->info < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad argument to GELSS"); 33619653cdaSPeter Brune if (ngmres->info > 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"SVD failed to converge"); 337a312c225SMatthew G Knepley #endif 338f109b39eSPeter Brune } 339a312c225SMatthew G Knepley 34019653cdaSPeter Brune alph_total = 0.; 34119653cdaSPeter Brune for (i = 0; i < l; i++) { 34219653cdaSPeter Brune alph_total += beta[i]; 343c0bbabecSxuemin } 344f109b39eSPeter Brune 3459f425c49SPeter Brune ierr = VecCopy(XM, XA);CHKERRQ(ierr); 346f109b39eSPeter Brune ierr = VecScale(XA, 1. - alph_total);CHKERRQ(ierr); 347087dfb9eSxuemin 3488c09b6b7SPeter Brune ierr = VecMAXPY(XA, l, beta, Xdot);CHKERRQ(ierr); 349f109b39eSPeter Brune ierr = SNESComputeFunction(snes, XA, FA);CHKERRQ(ierr); 350f109b39eSPeter Brune ierr = VecNorm(FA, NORM_2, &fAnorm);CHKERRQ(ierr); 35119653cdaSPeter Brune 3529f425c49SPeter Brune /* differences for selection and restart */ 353f109b39eSPeter Brune ierr=VecCopy(XA, D);CHKERRQ(ierr); 354f109b39eSPeter Brune ierr=VecAXPY(D, -1.0, X);CHKERRQ(ierr); 355f109b39eSPeter Brune ierr=VecNorm(D, NORM_2, &dnorm);CHKERRQ(ierr); 356f109b39eSPeter Brune dminnorm = -1.0; 35719653cdaSPeter Brune for(i=0;i<l;i++) { 358f109b39eSPeter Brune ierr=VecCopy(XA, D);CHKERRQ(ierr); 359f109b39eSPeter Brune ierr=VecAXPY(D, -1.0, Xdot[i]);CHKERRQ(ierr); 360f109b39eSPeter Brune ierr=VecNorm(D, NORM_2, &dcurnorm);CHKERRQ(ierr); 361f109b39eSPeter Brune if((dcurnorm < dminnorm) || (dminnorm < 0.0)) dminnorm = dcurnorm; 36219653cdaSPeter Brune } 3639f425c49SPeter Brune 3649f425c49SPeter Brune /* combination (additive) or selection (multiplicative) of the N-GMRES solution */ 3659f425c49SPeter Brune if (ngmres->additive) { 3669f425c49SPeter Brune /* X = X + \lambda(XA - X) */ 3679f425c49SPeter Brune if (ngmres->monitor) { 3689f425c49SPeter Brune ierr = PetscViewerASCIIPrintf(ngmres->monitor, "||F_A||_2 = %e, ||F_M||_2 = %e\n", fAnorm, fMnorm);CHKERRQ(ierr); 3699f425c49SPeter Brune } 3709f425c49SPeter Brune ierr = VecCopy(XA, Y);CHKERRQ(ierr); 371eb1825c3SPeter Brune ierr = VecAYPX(Y, -1.0, X);CHKERRQ(ierr); 372f1c6b773SPeter Brune ierr = SNESLineSearchApply(ngmres->additive_linesearch,X,F,&fnorm,Y);CHKERRQ(ierr); 373f1c6b773SPeter Brune ierr = SNESLineSearchGetSuccess(ngmres->additive_linesearch, &lssucceed);CHKERRQ(ierr); 3749f425c49SPeter Brune if (!lssucceed) { 3759f425c49SPeter Brune if (++snes->numFailures >= snes->maxFailures) { 3769f425c49SPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 3779f425c49SPeter Brune PetscFunctionReturn(0); 3789f425c49SPeter Brune } 3799f425c49SPeter Brune } 380f1c6b773SPeter Brune ierr = SNESLineSearchGetNorms(ngmres->additive_linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr); 3819f425c49SPeter Brune if (ngmres->monitor) { 3829f425c49SPeter Brune ierr = PetscViewerASCIIPrintf(ngmres->monitor, "Additive solution: ||F||_2 = %e\n", fnorm);CHKERRQ(ierr); 3839f425c49SPeter Brune } 3849f425c49SPeter Brune ierr = VecCopy(FA, F);CHKERRQ(ierr); 3859f425c49SPeter Brune ierr = VecCopy(XA, X);CHKERRQ(ierr); 3869f425c49SPeter Brune } else { 3879f425c49SPeter Brune selectA = PETSC_TRUE; 3889f425c49SPeter Brune /* Conditions for choosing the accelerated answer */ 3899f425c49SPeter Brune /* Criterion A -- the norm of the function isn't increased above the minimum by too much */ 3909f425c49SPeter Brune if (fAnorm >= ngmres->gammaA*fminnorm) { 3919f425c49SPeter Brune selectA = PETSC_FALSE; 3929f425c49SPeter Brune } 3939f425c49SPeter Brune /* Criterion B -- the choice of x^A isn't too close to some other choice */ 394cf22de31SPeter Brune if (ngmres->epsilonB*dnorm<dminnorm || PetscSqrtReal(fnorm)<ngmres->deltaB*PetscSqrtReal(fminnorm)) { 39519653cdaSPeter Brune } else { 39619653cdaSPeter Brune selectA=PETSC_FALSE; 39719653cdaSPeter Brune } 39819653cdaSPeter Brune if (selectA) { 399dfbf837cSBarry Smith if (ngmres->monitor) { 4009e764e56SPeter Brune ierr = PetscViewerASCIIPrintf(ngmres->monitor, "picked X_A, ||F_A||_2 = %e, ||F_M||_2 = %e\n", fAnorm, fMnorm);CHKERRQ(ierr); 401dfbf837cSBarry Smith } 40298b3e84cSPeter Brune /* copy it over */ 403f109b39eSPeter Brune fnorm = fAnorm; 404f109b39eSPeter Brune nu = fnorm*fnorm; 405f109b39eSPeter Brune ierr = VecCopy(FA, F);CHKERRQ(ierr); 406f109b39eSPeter Brune ierr = VecCopy(XA, X);CHKERRQ(ierr); 40719653cdaSPeter Brune } else { 408dfbf837cSBarry Smith if (ngmres->monitor) { 409f109b39eSPeter Brune ierr = PetscViewerASCIIPrintf(ngmres->monitor, "picked X_M, ||F_A||_2 = %e, ||F_M||_2 = %e\n", fAnorm, fnorm);CHKERRQ(ierr); 410dfbf837cSBarry Smith } 4119f425c49SPeter Brune fnorm = fMnorm; 4129f425c49SPeter Brune nu = fnorm*fnorm; 4139f425c49SPeter Brune ierr = VecCopy(FM, F);CHKERRQ(ierr); 4149f425c49SPeter Brune ierr = VecCopy(XM, X);CHKERRQ(ierr); 4159f425c49SPeter Brune } 41619653cdaSPeter Brune } 417211b2d39SPeter Brune 41819653cdaSPeter Brune selectRestart = PETSC_FALSE; 41919653cdaSPeter Brune 42098b3e84cSPeter Brune /* difference stagnation restart */ 421cf22de31SPeter Brune if((ngmres->epsilonB*dnorm > dminnorm) && (PetscSqrtReal(fAnorm) > ngmres->deltaB*PetscSqrtReal(fminnorm))) { 422dfbf837cSBarry Smith if (ngmres->monitor) { 423f109b39eSPeter Brune ierr = PetscViewerASCIIPrintf(ngmres->monitor, "difference restart: %e > %e\n", ngmres->epsilonB*dnorm, dminnorm);CHKERRQ(ierr); 424dfbf837cSBarry Smith } 42519653cdaSPeter Brune selectRestart = PETSC_TRUE; 42619653cdaSPeter Brune } 42798b3e84cSPeter Brune /* residual stagnation restart */ 428cf22de31SPeter Brune if (PetscSqrtReal(fAnorm) > ngmres->gammaC*PetscSqrtReal(fminnorm)) { 429dfbf837cSBarry Smith if (ngmres->monitor) { 430cf22de31SPeter Brune ierr = PetscViewerASCIIPrintf(ngmres->monitor, "residual restart: %e > %e\n", PetscSqrtReal(fAnorm), ngmres->gammaC*PetscSqrtReal(fminnorm));CHKERRQ(ierr); 431dfbf837cSBarry Smith } 43219653cdaSPeter Brune selectRestart = PETSC_TRUE; 43319653cdaSPeter Brune } 43419653cdaSPeter Brune 43528ed4a04SPeter Brune /* if the restart conditions persist for more than restart_it iterations, restart. */ 43619653cdaSPeter Brune if (selectRestart) { 43728ed4a04SPeter Brune restart_count++; 43828ed4a04SPeter Brune } else { 43928ed4a04SPeter Brune restart_count = 0; 44028ed4a04SPeter Brune } 44128ed4a04SPeter Brune 44228ed4a04SPeter Brune /* restart after restart conditions have persisted for a fixed number of iterations */ 44328ed4a04SPeter Brune if (restart_count >= ngmres->restart_it) { 444dfbf837cSBarry Smith if (ngmres->monitor){ 445dfbf837cSBarry Smith ierr = PetscViewerASCIIPrintf(ngmres->monitor, "Restarted at iteration %d\n", k_restart);CHKERRQ(ierr); 446dfbf837cSBarry Smith } 44728ed4a04SPeter Brune restart_count = 0; 44819653cdaSPeter Brune k_restart = 1; 44919653cdaSPeter Brune l = 1; 45098b3e84cSPeter Brune /* q_{00} = nu */ 451d2e16ddcSPeter Brune if (ngmres->anderson) { 452d2e16ddcSPeter Brune ngmres->fnorms[0] = fMnorm; 453d2e16ddcSPeter Brune nu = fMnorm*fMnorm; 454d2e16ddcSPeter Brune Q(0,0) = nu; 455d2e16ddcSPeter Brune /* Fdot[0] = F */ 456d2e16ddcSPeter Brune ierr = VecCopy(XM, Xdot[0]);CHKERRQ(ierr); 457d2e16ddcSPeter Brune ierr = VecCopy(FM, Fdot[0]);CHKERRQ(ierr); 458d2e16ddcSPeter Brune } else { 459f109b39eSPeter Brune ngmres->fnorms[0] = fnorm; 460f109b39eSPeter Brune nu = fnorm*fnorm; 46119653cdaSPeter Brune Q(0,0) = nu; 462f109b39eSPeter Brune /* Fdot[0] = F */ 463f109b39eSPeter Brune ierr = VecCopy(X, Xdot[0]);CHKERRQ(ierr); 464f109b39eSPeter Brune ierr = VecCopy(F, Fdot[0]);CHKERRQ(ierr); 465d2e16ddcSPeter Brune } 46619653cdaSPeter Brune } else { 46798b3e84cSPeter Brune /* select the current size of the subspace */ 4681e633543SBarry Smith if (l < ngmres->msize) l++; 46919653cdaSPeter Brune k_restart++; 47098b3e84cSPeter Brune /* place the current entry in the list of previous entries */ 471d2e16ddcSPeter Brune if (ngmres->anderson) { 472d2e16ddcSPeter Brune ierr = VecCopy(FM, Fdot[ivec]);CHKERRQ(ierr); 473d2e16ddcSPeter Brune ierr = VecCopy(XM, Xdot[ivec]);CHKERRQ(ierr); 474d2e16ddcSPeter Brune ngmres->fnorms[ivec] = fMnorm; 475d2e16ddcSPeter Brune if (fminnorm > fMnorm) fminnorm = fMnorm; /* the minimum norm is now of FM */ 476d2e16ddcSPeter Brune for (i = 0; i < l; i++) { 477d2e16ddcSPeter Brune ierr = VecDot(FM, Fdot[i], &qentry);CHKERRQ(ierr); 478d2e16ddcSPeter Brune Q(i, ivec) = qentry; 479d2e16ddcSPeter Brune Q(ivec, i) = qentry; 480d2e16ddcSPeter Brune } 481d2e16ddcSPeter Brune } else { 482f109b39eSPeter Brune ierr = VecCopy(F, Fdot[ivec]);CHKERRQ(ierr); 483f109b39eSPeter Brune ierr = VecCopy(X, Xdot[ivec]);CHKERRQ(ierr); 484f109b39eSPeter Brune ngmres->fnorms[ivec] = fnorm; 485d2e16ddcSPeter Brune if (fminnorm > fnorm) fminnorm = fnorm; /* the minimum norm is now of FA */ 48619653cdaSPeter Brune for (i = 0; i < l; i++) { 487f109b39eSPeter Brune ierr = VecDot(F, Fdot[i], &qentry);CHKERRQ(ierr); 48819653cdaSPeter Brune Q(i, ivec) = qentry; 48919653cdaSPeter Brune Q(ivec, i) = qentry; 49019653cdaSPeter Brune } 49119653cdaSPeter Brune } 492d2e16ddcSPeter Brune } 49319653cdaSPeter Brune 4948409ca45SMatthew G Knepley ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 495087dfb9eSxuemin snes->iter = k; 496f109b39eSPeter Brune snes->norm = fnorm; 4978409ca45SMatthew G Knepley ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 4988409ca45SMatthew G Knepley SNESLogConvHistory(snes, snes->norm, snes->iter); 4998409ca45SMatthew G Knepley ierr = SNESMonitor(snes, snes->iter, snes->norm);CHKERRQ(ierr); 5008409ca45SMatthew G Knepley 501e7058c64SPeter Brune ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 502087dfb9eSxuemin if (snes->reason) PetscFunctionReturn(0); 503a312c225SMatthew G Knepley } 504a312c225SMatthew G Knepley snes->reason = SNES_DIVERGED_MAX_IT; 505a312c225SMatthew G Knepley PetscFunctionReturn(0); 506a312c225SMatthew G Knepley } 507a312c225SMatthew G Knepley 508dfbf837cSBarry Smith /*MC 5091867fe5bSPeter Brune SNESNGMRES - The Nonlinear Generalized Minimum Residual method. 510a312c225SMatthew G Knepley 511dfbf837cSBarry Smith Level: beginner 512dfbf837cSBarry Smith 5131867fe5bSPeter Brune Options Database: 5141867fe5bSPeter Brune + -snes_ngmres_additive - Use a variant that line-searches between the candidate solution and the combined solution. 5151867fe5bSPeter Brune . -snes_ngmres_anderson - Use Anderson mixing NGMRES variant which combines candidate solutions instead of actual solutions. 5161867fe5bSPeter Brune . -snes_ngmres_m - Number of stored previous solutions and residuals. 5171867fe5bSPeter Brune . -snes_ngmres_restart_it - Number of iterations the restart conditions hold before restart. 5181867fe5bSPeter Brune . -snes_ngmres_gammaA - Residual tolerance for solution selection between the candidate and combination. 5191867fe5bSPeter Brune . -snes_ngmres_gammaC - Residual tolerance for restart. 5201867fe5bSPeter Brune . -snes_ngmres_epsilonB - Difference tolerance between subsequent solutions triggering restart. 5211867fe5bSPeter Brune . -snes_ngmres_deltaB - Difference tolerance between residuals triggering restart. 5221867fe5bSPeter Brune . -snes_ngmres_monitor - Prints relevant information about the ngmres iteration. 523f3aaf736SPeter Brune - -snes_ls <basic,basicnonorms,quadratic,critical> - Line search type. 5241867fe5bSPeter Brune 5251867fe5bSPeter Brune Notes: 5261867fe5bSPeter Brune 5271867fe5bSPeter Brune The N-GMRES method combines m previous solutions into a minimum-residual solution by solving a small linearized 5281867fe5bSPeter Brune optimization problem at each iteration. 5291867fe5bSPeter Brune 5301867fe5bSPeter Brune References: 5311867fe5bSPeter Brune 532dfbf837cSBarry Smith "Krylov Subspace Acceleration of Nonlinear Multigrid with Application to Recirculating Flows", C. W. Oosterlee and T. Washio, 533dfbf837cSBarry Smith SIAM Journal on Scientific Computing, 21(5), 2000. 534dfbf837cSBarry Smith 535dfbf837cSBarry Smith This is also the same as the algorithm called Anderson acceleration introduced in "D. G. Anderson. Iterative procedures for nonlinear integral equations. 536dfbf837cSBarry Smith J. Assoc. Comput. Mach., 12:547–560, 1965." 537dfbf837cSBarry Smith 538dfbf837cSBarry Smith .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types) 539dfbf837cSBarry Smith M*/ 540a312c225SMatthew G Knepley 541a312c225SMatthew G Knepley EXTERN_C_BEGIN 542a312c225SMatthew G Knepley #undef __FUNCT__ 543a312c225SMatthew G Knepley #define __FUNCT__ "SNESCreate_NGMRES" 544a312c225SMatthew G Knepley PetscErrorCode SNESCreate_NGMRES(SNES snes) 545a312c225SMatthew G Knepley { 546a312c225SMatthew G Knepley SNES_NGMRES *ngmres; 547a312c225SMatthew G Knepley PetscErrorCode ierr; 548a312c225SMatthew G Knepley 549a312c225SMatthew G Knepley PetscFunctionBegin; 550a312c225SMatthew G Knepley snes->ops->destroy = SNESDestroy_NGMRES; 551a312c225SMatthew G Knepley snes->ops->setup = SNESSetUp_NGMRES; 552a312c225SMatthew G Knepley snes->ops->setfromoptions = SNESSetFromOptions_NGMRES; 553a312c225SMatthew G Knepley snes->ops->view = SNESView_NGMRES; 554a312c225SMatthew G Knepley snes->ops->solve = SNESSolve_NGMRES; 555a312c225SMatthew G Knepley snes->ops->reset = SNESReset_NGMRES; 556a312c225SMatthew G Knepley 55742f4f86dSBarry Smith snes->usespc = PETSC_TRUE; 5582c155ee1SBarry Smith snes->usesksp = PETSC_FALSE; 5592c155ee1SBarry Smith 560a312c225SMatthew G Knepley ierr = PetscNewLog(snes, SNES_NGMRES, &ngmres);CHKERRQ(ierr); 561a312c225SMatthew G Knepley snes->data = (void*) ngmres; 562d2e16ddcSPeter Brune ngmres->msize = 30; 56319653cdaSPeter Brune 56488976e71SPeter Brune if (!snes->tolerancesset) { 5650e444f03SPeter Brune snes->max_funcs = 30000; 5660e444f03SPeter Brune snes->max_its = 10000; 56788976e71SPeter Brune } 5680e444f03SPeter Brune 569d2e16ddcSPeter Brune ngmres->additive = PETSC_FALSE; 570d2e16ddcSPeter Brune ngmres->anderson = PETSC_FALSE; 571d2e16ddcSPeter Brune 572e7058c64SPeter Brune ngmres->additive_linesearch = PETSC_NULL; 573e7058c64SPeter Brune 57428ed4a04SPeter Brune ngmres->restart_it = 2; 575f109b39eSPeter Brune ngmres->gammaA = 2.0; 576f109b39eSPeter Brune ngmres->gammaC = 2.0; 577cac108bcSPeter Brune ngmres->deltaB = 0.9; 578cac108bcSPeter Brune ngmres->epsilonB = 0.1; 579e7058c64SPeter Brune 580a312c225SMatthew G Knepley PetscFunctionReturn(0); 581a312c225SMatthew G Knepley } 582a312c225SMatthew G Knepley EXTERN_C_END 583