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) */ 201f109b39eSPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 202a312c225SMatthew G Knepley if (snes->domainerror) { 203a312c225SMatthew G Knepley snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 204a312c225SMatthew G Knepley PetscFunctionReturn(0); 205a312c225SMatthew G Knepley } 20619653cdaSPeter Brune 20798b3e84cSPeter Brune /* nu = (r, r) */ 208f109b39eSPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 209f109b39eSPeter Brune fminnorm = fnorm; 210f109b39eSPeter Brune nu = fnorm*fnorm; 211f109b39eSPeter Brune if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FP, "Infinite or not-a-number generated in function evaluation"); 21219653cdaSPeter Brune 21398b3e84cSPeter Brune /* q_{00} = nu */ 21419653cdaSPeter Brune Q(0,0) = nu; 215f109b39eSPeter Brune ngmres->fnorms[0] = fnorm; 216f109b39eSPeter Brune /* Fdot[0] = F */ 217f109b39eSPeter Brune ierr = VecCopy(X, Xdot[0]);CHKERRQ(ierr); 218f109b39eSPeter Brune ierr = VecCopy(F, Fdot[0]);CHKERRQ(ierr); 219087dfb9eSxuemin 220a312c225SMatthew G Knepley ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 221f109b39eSPeter Brune snes->norm = fnorm; 222a312c225SMatthew G Knepley ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 223f109b39eSPeter Brune SNESLogConvHistory(snes, fnorm, 0); 224f109b39eSPeter Brune ierr = SNESMonitor(snes, 0, fnorm);CHKERRQ(ierr); 225f109b39eSPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 226a312c225SMatthew G Knepley if (snes->reason) PetscFunctionReturn(0); 227a312c225SMatthew G Knepley 22819653cdaSPeter Brune k_restart = 1; 22919653cdaSPeter Brune l = 1; 23009c08436SPeter Brune for (k=1; k < snes->max_its+1; k++) { 23198b3e84cSPeter Brune /* select which vector of the stored subspace will be updated */ 23298b3e84cSPeter Brune ivec = k_restart % ngmres->msize; /* replace the last used part of the subspace */ 23319653cdaSPeter Brune 23498b3e84cSPeter Brune /* Computation of x^M */ 2358cc86e31SPeter Brune if (snes->pc) { 2369f425c49SPeter Brune ierr = VecCopy(X, XM);CHKERRQ(ierr); 2379f425c49SPeter Brune ierr = SNESSolve(snes->pc, B, XM);CHKERRQ(ierr); 2388cc86e31SPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 2398cc86e31SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 2408cc86e31SPeter Brune snes->reason = SNES_DIVERGED_INNER; 2418cc86e31SPeter Brune PetscFunctionReturn(0); 2428cc86e31SPeter Brune } 2434b32a720SPeter Brune ierr = SNESGetFunction(snes->pc, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr); 2444b32a720SPeter Brune ierr = VecCopy(FPC, FM);CHKERRQ(ierr); 2454b32a720SPeter Brune ierr = SNESGetFunctionNorm(snes->pc, &fMnorm);CHKERRQ(ierr); 2468cc86e31SPeter Brune } else { 247f109b39eSPeter Brune /* no preconditioner -- just take gradient descent with line search */ 248f109b39eSPeter Brune ierr = VecCopy(F, Y);CHKERRQ(ierr); 249e7058c64SPeter Brune ierr = VecCopy(F, FM);CHKERRQ(ierr); 250e7058c64SPeter Brune ierr = VecCopy(X, XM);CHKERRQ(ierr); 251e7058c64SPeter Brune fMnorm = fnorm; 252f1c6b773SPeter Brune ierr = SNESLineSearchApply(snes->linesearch,XM,FM,&fMnorm,Y);CHKERRQ(ierr); 253f1c6b773SPeter Brune ierr = SNESLineSearchGetSuccess(snes->linesearch, &lssucceed);CHKERRQ(ierr); 254f109b39eSPeter Brune if (!lssucceed) { 255f109b39eSPeter Brune if (++snes->numFailures >= snes->maxFailures) { 256f109b39eSPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 257f109b39eSPeter Brune PetscFunctionReturn(0); 258f109b39eSPeter Brune } 259f109b39eSPeter Brune } 2606634f59bSPeter Brune } 2611e633543SBarry Smith 26298b3e84cSPeter Brune /* r = F(x) */ 2639f425c49SPeter Brune nu = fMnorm*fMnorm; 2649f425c49SPeter Brune if (fminnorm > fMnorm) fminnorm = fMnorm; /* the minimum norm is now of F^M */ 26519653cdaSPeter Brune 26698b3e84cSPeter Brune /* construct the right hand side and xi factors */ 2678c09b6b7SPeter Brune ierr = VecMDot(FM, l, Fdot, xi);CHKERRQ(ierr); 2688c09b6b7SPeter Brune 26919653cdaSPeter Brune for (i = 0; i < l; i++) { 27019653cdaSPeter Brune beta[i] = nu - xi[i]; 27119653cdaSPeter Brune } 27219653cdaSPeter Brune 27398b3e84cSPeter Brune /* construct h */ 27419653cdaSPeter Brune for (j = 0; j < l; j++) { 27519653cdaSPeter Brune for (i = 0; i < l; i++) { 27619653cdaSPeter Brune H(i, j) = Q(i, j) - xi[i] - xi[j] + nu; 27719653cdaSPeter Brune } 27819653cdaSPeter Brune } 279f109b39eSPeter Brune 280f109b39eSPeter Brune if(l == 1) { 281f109b39eSPeter Brune /* simply set alpha[0] = beta[0] / H[0, 0] */ 282f109b39eSPeter Brune beta[0] = beta[0] / H(0, 0); 283f109b39eSPeter Brune } else { 28419653cdaSPeter Brune #ifdef PETSC_MISSING_LAPACK_GELSS 28519653cdaSPeter Brune SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "NGMRES with LS requires the LAPACK GELSS routine."); 2864a0c5b0cSMatthew G Knepley #else 28719653cdaSPeter Brune ngmres->m = PetscBLASIntCast(l); 28819653cdaSPeter Brune ngmres->n = PetscBLASIntCast(l); 28919653cdaSPeter Brune ngmres->info = PetscBLASIntCast(0); 29019653cdaSPeter Brune ngmres->rcond = -1.; 29119653cdaSPeter Brune #ifdef PETSC_USE_COMPLEX 29219653cdaSPeter Brune LAPACKgelss_(&ngmres->m, 29319653cdaSPeter Brune &ngmres->n, 29419653cdaSPeter Brune &ngmres->nrhs, 29519653cdaSPeter Brune ngmres->h, 29619653cdaSPeter Brune &ngmres->lda, 29719653cdaSPeter Brune ngmres->beta, 29819653cdaSPeter Brune &ngmres->ldb, 29919653cdaSPeter Brune ngmres->s, 30019653cdaSPeter Brune &ngmres->rcond, 30119653cdaSPeter Brune &ngmres->rank, 30219653cdaSPeter Brune ngmres->work, 30319653cdaSPeter Brune &ngmres->lwork, 30419653cdaSPeter Brune ngmres->rwork, 30519653cdaSPeter Brune &ngmres->info); 30619653cdaSPeter Brune #else 30719653cdaSPeter Brune LAPACKgelss_(&ngmres->m, 30819653cdaSPeter Brune &ngmres->n, 30919653cdaSPeter Brune &ngmres->nrhs, 31019653cdaSPeter Brune ngmres->h, 31119653cdaSPeter Brune &ngmres->lda, 31219653cdaSPeter Brune ngmres->beta, 31319653cdaSPeter Brune &ngmres->ldb, 31419653cdaSPeter Brune ngmres->s, 31519653cdaSPeter Brune &ngmres->rcond, 31619653cdaSPeter Brune &ngmres->rank, 31719653cdaSPeter Brune ngmres->work, 31819653cdaSPeter Brune &ngmres->lwork, 31919653cdaSPeter Brune &ngmres->info); 32019653cdaSPeter Brune #endif 32119653cdaSPeter Brune if (ngmres->info < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad argument to GELSS"); 32219653cdaSPeter Brune if (ngmres->info > 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"SVD failed to converge"); 323a312c225SMatthew G Knepley #endif 324f109b39eSPeter Brune } 325a312c225SMatthew G Knepley 32619653cdaSPeter Brune alph_total = 0.; 32719653cdaSPeter Brune for (i = 0; i < l; i++) { 32819653cdaSPeter Brune alph_total += beta[i]; 329c0bbabecSxuemin } 330f109b39eSPeter Brune 3319f425c49SPeter Brune ierr = VecCopy(XM, XA);CHKERRQ(ierr); 332f109b39eSPeter Brune ierr = VecScale(XA, 1. - alph_total);CHKERRQ(ierr); 333087dfb9eSxuemin 3348c09b6b7SPeter Brune ierr = VecMAXPY(XA, l, beta, Xdot);CHKERRQ(ierr); 335f109b39eSPeter Brune ierr = SNESComputeFunction(snes, XA, FA);CHKERRQ(ierr); 336f109b39eSPeter Brune ierr = VecNorm(FA, NORM_2, &fAnorm);CHKERRQ(ierr); 33719653cdaSPeter Brune 3389f425c49SPeter Brune /* differences for selection and restart */ 339f109b39eSPeter Brune ierr=VecCopy(XA, D);CHKERRQ(ierr); 340f109b39eSPeter Brune ierr=VecAXPY(D, -1.0, X);CHKERRQ(ierr); 341f109b39eSPeter Brune ierr=VecNorm(D, NORM_2, &dnorm);CHKERRQ(ierr); 342f109b39eSPeter Brune dminnorm = -1.0; 34319653cdaSPeter Brune for(i=0;i<l;i++) { 344f109b39eSPeter Brune ierr=VecCopy(XA, D);CHKERRQ(ierr); 345f109b39eSPeter Brune ierr=VecAXPY(D, -1.0, Xdot[i]);CHKERRQ(ierr); 346f109b39eSPeter Brune ierr=VecNorm(D, NORM_2, &dcurnorm);CHKERRQ(ierr); 347f109b39eSPeter Brune if((dcurnorm < dminnorm) || (dminnorm < 0.0)) dminnorm = dcurnorm; 34819653cdaSPeter Brune } 3499f425c49SPeter Brune 3509f425c49SPeter Brune /* combination (additive) or selection (multiplicative) of the N-GMRES solution */ 3519f425c49SPeter Brune if (ngmres->additive) { 3529f425c49SPeter Brune /* X = X + \lambda(XA - X) */ 3539f425c49SPeter Brune if (ngmres->monitor) { 3549f425c49SPeter Brune ierr = PetscViewerASCIIPrintf(ngmres->monitor, "||F_A||_2 = %e, ||F_M||_2 = %e\n", fAnorm, fMnorm);CHKERRQ(ierr); 3559f425c49SPeter Brune } 3569f425c49SPeter Brune ierr = VecCopy(XA, Y);CHKERRQ(ierr); 357eb1825c3SPeter Brune ierr = VecAYPX(Y, -1.0, X);CHKERRQ(ierr); 358f1c6b773SPeter Brune ierr = SNESLineSearchApply(ngmres->additive_linesearch,X,F,&fnorm,Y);CHKERRQ(ierr); 359f1c6b773SPeter Brune ierr = SNESLineSearchGetSuccess(ngmres->additive_linesearch, &lssucceed);CHKERRQ(ierr); 3609f425c49SPeter Brune if (!lssucceed) { 3619f425c49SPeter Brune if (++snes->numFailures >= snes->maxFailures) { 3629f425c49SPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 3639f425c49SPeter Brune PetscFunctionReturn(0); 3649f425c49SPeter Brune } 3659f425c49SPeter Brune } 366f1c6b773SPeter Brune ierr = SNESLineSearchGetNorms(ngmres->additive_linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr); 3679f425c49SPeter Brune if (ngmres->monitor) { 3689f425c49SPeter Brune ierr = PetscViewerASCIIPrintf(ngmres->monitor, "Additive solution: ||F||_2 = %e\n", fnorm);CHKERRQ(ierr); 3699f425c49SPeter Brune } 3709f425c49SPeter Brune ierr = VecCopy(FA, F);CHKERRQ(ierr); 3719f425c49SPeter Brune ierr = VecCopy(XA, X);CHKERRQ(ierr); 3729f425c49SPeter Brune } else { 3739f425c49SPeter Brune selectA = PETSC_TRUE; 3749f425c49SPeter Brune /* Conditions for choosing the accelerated answer */ 3759f425c49SPeter Brune /* Criterion A -- the norm of the function isn't increased above the minimum by too much */ 3769f425c49SPeter Brune if (fAnorm >= ngmres->gammaA*fminnorm) { 3779f425c49SPeter Brune selectA = PETSC_FALSE; 3789f425c49SPeter Brune } 3799f425c49SPeter Brune /* Criterion B -- the choice of x^A isn't too close to some other choice */ 380cf22de31SPeter Brune if (ngmres->epsilonB*dnorm<dminnorm || PetscSqrtReal(fnorm)<ngmres->deltaB*PetscSqrtReal(fminnorm)) { 38119653cdaSPeter Brune } else { 38219653cdaSPeter Brune selectA=PETSC_FALSE; 38319653cdaSPeter Brune } 38419653cdaSPeter Brune if (selectA) { 385dfbf837cSBarry Smith if (ngmres->monitor) { 3869e764e56SPeter Brune ierr = PetscViewerASCIIPrintf(ngmres->monitor, "picked X_A, ||F_A||_2 = %e, ||F_M||_2 = %e\n", fAnorm, fMnorm);CHKERRQ(ierr); 387dfbf837cSBarry Smith } 38898b3e84cSPeter Brune /* copy it over */ 389f109b39eSPeter Brune fnorm = fAnorm; 390f109b39eSPeter Brune nu = fnorm*fnorm; 391f109b39eSPeter Brune ierr = VecCopy(FA, F);CHKERRQ(ierr); 392f109b39eSPeter Brune ierr = VecCopy(XA, X);CHKERRQ(ierr); 39319653cdaSPeter Brune } else { 394dfbf837cSBarry Smith if (ngmres->monitor) { 395f109b39eSPeter Brune ierr = PetscViewerASCIIPrintf(ngmres->monitor, "picked X_M, ||F_A||_2 = %e, ||F_M||_2 = %e\n", fAnorm, fnorm);CHKERRQ(ierr); 396dfbf837cSBarry Smith } 3979f425c49SPeter Brune fnorm = fMnorm; 3989f425c49SPeter Brune nu = fnorm*fnorm; 3999f425c49SPeter Brune ierr = VecCopy(FM, F);CHKERRQ(ierr); 4009f425c49SPeter Brune ierr = VecCopy(XM, X);CHKERRQ(ierr); 4019f425c49SPeter Brune } 40219653cdaSPeter Brune } 403211b2d39SPeter Brune 40419653cdaSPeter Brune selectRestart = PETSC_FALSE; 40519653cdaSPeter Brune 40698b3e84cSPeter Brune /* difference stagnation restart */ 407cf22de31SPeter Brune if((ngmres->epsilonB*dnorm > dminnorm) && (PetscSqrtReal(fAnorm) > ngmres->deltaB*PetscSqrtReal(fminnorm))) { 408dfbf837cSBarry Smith if (ngmres->monitor) { 409f109b39eSPeter Brune ierr = PetscViewerASCIIPrintf(ngmres->monitor, "difference restart: %e > %e\n", ngmres->epsilonB*dnorm, dminnorm);CHKERRQ(ierr); 410dfbf837cSBarry Smith } 41119653cdaSPeter Brune selectRestart = PETSC_TRUE; 41219653cdaSPeter Brune } 41398b3e84cSPeter Brune /* residual stagnation restart */ 414cf22de31SPeter Brune if (PetscSqrtReal(fAnorm) > ngmres->gammaC*PetscSqrtReal(fminnorm)) { 415dfbf837cSBarry Smith if (ngmres->monitor) { 416cf22de31SPeter Brune ierr = PetscViewerASCIIPrintf(ngmres->monitor, "residual restart: %e > %e\n", PetscSqrtReal(fAnorm), ngmres->gammaC*PetscSqrtReal(fminnorm));CHKERRQ(ierr); 417dfbf837cSBarry Smith } 41819653cdaSPeter Brune selectRestart = PETSC_TRUE; 41919653cdaSPeter Brune } 42019653cdaSPeter Brune 42128ed4a04SPeter Brune /* if the restart conditions persist for more than restart_it iterations, restart. */ 42219653cdaSPeter Brune if (selectRestart) { 42328ed4a04SPeter Brune restart_count++; 42428ed4a04SPeter Brune } else { 42528ed4a04SPeter Brune restart_count = 0; 42628ed4a04SPeter Brune } 42728ed4a04SPeter Brune 42828ed4a04SPeter Brune /* restart after restart conditions have persisted for a fixed number of iterations */ 42928ed4a04SPeter Brune if (restart_count >= ngmres->restart_it) { 430dfbf837cSBarry Smith if (ngmres->monitor){ 431dfbf837cSBarry Smith ierr = PetscViewerASCIIPrintf(ngmres->monitor, "Restarted at iteration %d\n", k_restart);CHKERRQ(ierr); 432dfbf837cSBarry Smith } 43328ed4a04SPeter Brune restart_count = 0; 43419653cdaSPeter Brune k_restart = 1; 43519653cdaSPeter Brune l = 1; 43698b3e84cSPeter Brune /* q_{00} = nu */ 437d2e16ddcSPeter Brune if (ngmres->anderson) { 438d2e16ddcSPeter Brune ngmres->fnorms[0] = fMnorm; 439d2e16ddcSPeter Brune nu = fMnorm*fMnorm; 440d2e16ddcSPeter Brune Q(0,0) = nu; 441d2e16ddcSPeter Brune /* Fdot[0] = F */ 442d2e16ddcSPeter Brune ierr = VecCopy(XM, Xdot[0]);CHKERRQ(ierr); 443d2e16ddcSPeter Brune ierr = VecCopy(FM, Fdot[0]);CHKERRQ(ierr); 444d2e16ddcSPeter Brune } else { 445f109b39eSPeter Brune ngmres->fnorms[0] = fnorm; 446f109b39eSPeter Brune nu = fnorm*fnorm; 44719653cdaSPeter Brune Q(0,0) = nu; 448f109b39eSPeter Brune /* Fdot[0] = F */ 449f109b39eSPeter Brune ierr = VecCopy(X, Xdot[0]);CHKERRQ(ierr); 450f109b39eSPeter Brune ierr = VecCopy(F, Fdot[0]);CHKERRQ(ierr); 451d2e16ddcSPeter Brune } 45219653cdaSPeter Brune } else { 45398b3e84cSPeter Brune /* select the current size of the subspace */ 4541e633543SBarry Smith if (l < ngmres->msize) l++; 45519653cdaSPeter Brune k_restart++; 45698b3e84cSPeter Brune /* place the current entry in the list of previous entries */ 457d2e16ddcSPeter Brune if (ngmres->anderson) { 458d2e16ddcSPeter Brune ierr = VecCopy(FM, Fdot[ivec]);CHKERRQ(ierr); 459d2e16ddcSPeter Brune ierr = VecCopy(XM, Xdot[ivec]);CHKERRQ(ierr); 460d2e16ddcSPeter Brune ngmres->fnorms[ivec] = fMnorm; 461d2e16ddcSPeter Brune if (fminnorm > fMnorm) fminnorm = fMnorm; /* the minimum norm is now of FM */ 462d2e16ddcSPeter Brune for (i = 0; i < l; i++) { 463d2e16ddcSPeter Brune ierr = VecDot(FM, Fdot[i], &qentry);CHKERRQ(ierr); 464d2e16ddcSPeter Brune Q(i, ivec) = qentry; 465d2e16ddcSPeter Brune Q(ivec, i) = qentry; 466d2e16ddcSPeter Brune } 467d2e16ddcSPeter Brune } else { 468f109b39eSPeter Brune ierr = VecCopy(F, Fdot[ivec]);CHKERRQ(ierr); 469f109b39eSPeter Brune ierr = VecCopy(X, Xdot[ivec]);CHKERRQ(ierr); 470f109b39eSPeter Brune ngmres->fnorms[ivec] = fnorm; 471d2e16ddcSPeter Brune if (fminnorm > fnorm) fminnorm = fnorm; /* the minimum norm is now of FA */ 47219653cdaSPeter Brune for (i = 0; i < l; i++) { 473f109b39eSPeter Brune ierr = VecDot(F, Fdot[i], &qentry);CHKERRQ(ierr); 47419653cdaSPeter Brune Q(i, ivec) = qentry; 47519653cdaSPeter Brune Q(ivec, i) = qentry; 47619653cdaSPeter Brune } 47719653cdaSPeter Brune } 478d2e16ddcSPeter Brune } 47919653cdaSPeter Brune 4808409ca45SMatthew G Knepley ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 481087dfb9eSxuemin snes->iter = k; 482f109b39eSPeter Brune snes->norm = fnorm; 4838409ca45SMatthew G Knepley ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 4848409ca45SMatthew G Knepley SNESLogConvHistory(snes, snes->norm, snes->iter); 4858409ca45SMatthew G Knepley ierr = SNESMonitor(snes, snes->iter, snes->norm);CHKERRQ(ierr); 4868409ca45SMatthew G Knepley 487e7058c64SPeter Brune ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 488087dfb9eSxuemin if (snes->reason) PetscFunctionReturn(0); 489a312c225SMatthew G Knepley } 490a312c225SMatthew G Knepley snes->reason = SNES_DIVERGED_MAX_IT; 491a312c225SMatthew G Knepley PetscFunctionReturn(0); 492a312c225SMatthew G Knepley } 493a312c225SMatthew G Knepley 494dfbf837cSBarry Smith /*MC 4951867fe5bSPeter Brune SNESNGMRES - The Nonlinear Generalized Minimum Residual method. 496a312c225SMatthew G Knepley 497dfbf837cSBarry Smith Level: beginner 498dfbf837cSBarry Smith 4991867fe5bSPeter Brune Options Database: 5001867fe5bSPeter Brune + -snes_ngmres_additive - Use a variant that line-searches between the candidate solution and the combined solution. 5011867fe5bSPeter Brune . -snes_ngmres_anderson - Use Anderson mixing NGMRES variant which combines candidate solutions instead of actual solutions. 5021867fe5bSPeter Brune . -snes_ngmres_m - Number of stored previous solutions and residuals. 5031867fe5bSPeter Brune . -snes_ngmres_restart_it - Number of iterations the restart conditions hold before restart. 5041867fe5bSPeter Brune . -snes_ngmres_gammaA - Residual tolerance for solution selection between the candidate and combination. 5051867fe5bSPeter Brune . -snes_ngmres_gammaC - Residual tolerance for restart. 5061867fe5bSPeter Brune . -snes_ngmres_epsilonB - Difference tolerance between subsequent solutions triggering restart. 5071867fe5bSPeter Brune . -snes_ngmres_deltaB - Difference tolerance between residuals triggering restart. 5081867fe5bSPeter Brune . -snes_ngmres_monitor - Prints relevant information about the ngmres iteration. 509f3aaf736SPeter Brune - -snes_ls <basic,basicnonorms,quadratic,critical> - Line search type. 5101867fe5bSPeter Brune 5111867fe5bSPeter Brune Notes: 5121867fe5bSPeter Brune 5131867fe5bSPeter Brune The N-GMRES method combines m previous solutions into a minimum-residual solution by solving a small linearized 5141867fe5bSPeter Brune optimization problem at each iteration. 5151867fe5bSPeter Brune 5161867fe5bSPeter Brune References: 5171867fe5bSPeter Brune 518dfbf837cSBarry Smith "Krylov Subspace Acceleration of Nonlinear Multigrid with Application to Recirculating Flows", C. W. Oosterlee and T. Washio, 519dfbf837cSBarry Smith SIAM Journal on Scientific Computing, 21(5), 2000. 520dfbf837cSBarry Smith 521dfbf837cSBarry Smith This is also the same as the algorithm called Anderson acceleration introduced in "D. G. Anderson. Iterative procedures for nonlinear integral equations. 522dfbf837cSBarry Smith J. Assoc. Comput. Mach., 12:547–560, 1965." 523dfbf837cSBarry Smith 524dfbf837cSBarry Smith .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types) 525dfbf837cSBarry Smith M*/ 526a312c225SMatthew G Knepley 527a312c225SMatthew G Knepley EXTERN_C_BEGIN 528a312c225SMatthew G Knepley #undef __FUNCT__ 529a312c225SMatthew G Knepley #define __FUNCT__ "SNESCreate_NGMRES" 530a312c225SMatthew G Knepley PetscErrorCode SNESCreate_NGMRES(SNES snes) 531a312c225SMatthew G Knepley { 532a312c225SMatthew G Knepley SNES_NGMRES *ngmres; 533a312c225SMatthew G Knepley PetscErrorCode ierr; 534a312c225SMatthew G Knepley 535a312c225SMatthew G Knepley PetscFunctionBegin; 536a312c225SMatthew G Knepley snes->ops->destroy = SNESDestroy_NGMRES; 537a312c225SMatthew G Knepley snes->ops->setup = SNESSetUp_NGMRES; 538a312c225SMatthew G Knepley snes->ops->setfromoptions = SNESSetFromOptions_NGMRES; 539a312c225SMatthew G Knepley snes->ops->view = SNESView_NGMRES; 540a312c225SMatthew G Knepley snes->ops->solve = SNESSolve_NGMRES; 541a312c225SMatthew G Knepley snes->ops->reset = SNESReset_NGMRES; 542a312c225SMatthew G Knepley 54342f4f86dSBarry Smith snes->usespc = PETSC_TRUE; 5442c155ee1SBarry Smith snes->usesksp = PETSC_FALSE; 5452c155ee1SBarry Smith 546a312c225SMatthew G Knepley ierr = PetscNewLog(snes, SNES_NGMRES, &ngmres);CHKERRQ(ierr); 547a312c225SMatthew G Knepley snes->data = (void*) ngmres; 548d2e16ddcSPeter Brune ngmres->msize = 30; 54919653cdaSPeter Brune 550*88976e71SPeter Brune if (!snes->tolerancesset) { 5510e444f03SPeter Brune snes->max_funcs = 30000; 5520e444f03SPeter Brune snes->max_its = 10000; 553*88976e71SPeter Brune } 5540e444f03SPeter Brune 555d2e16ddcSPeter Brune ngmres->additive = PETSC_FALSE; 556d2e16ddcSPeter Brune ngmres->anderson = PETSC_FALSE; 557d2e16ddcSPeter Brune 558e7058c64SPeter Brune ngmres->additive_linesearch = PETSC_NULL; 559e7058c64SPeter Brune 56028ed4a04SPeter Brune ngmres->restart_it = 2; 561f109b39eSPeter Brune ngmres->gammaA = 2.0; 562f109b39eSPeter Brune ngmres->gammaC = 2.0; 563cac108bcSPeter Brune ngmres->deltaB = 0.9; 564cac108bcSPeter Brune ngmres->epsilonB = 0.1; 565e7058c64SPeter Brune 566a312c225SMatthew G Knepley PetscFunctionReturn(0); 567a312c225SMatthew G Knepley } 568a312c225SMatthew G Knepley EXTERN_C_END 569