1a312c225SMatthew G Knepley /* Defines the basic SNES object */ 219653cdaSPeter Brune #include <../src/snes/impls/ngmres/snesngmres.h> 319653cdaSPeter Brune #include <petscblaslapack.h> 4a312c225SMatthew G Knepley 5a312c225SMatthew G Knepley #undef __FUNCT__ 6a312c225SMatthew G Knepley #define __FUNCT__ "SNESReset_NGMRES" 7a312c225SMatthew G Knepley PetscErrorCode SNESReset_NGMRES(SNES snes) 8a312c225SMatthew G Knepley { 9a312c225SMatthew G Knepley SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data; 10a312c225SMatthew G Knepley PetscErrorCode ierr; 11a312c225SMatthew G Knepley 12a312c225SMatthew G Knepley PetscFunctionBegin; 13f109b39eSPeter Brune ierr = VecDestroyVecs(ngmres->msize, &ngmres->Fdot);CHKERRQ(ierr); 14f109b39eSPeter Brune ierr = VecDestroyVecs(ngmres->msize, &ngmres->Xdot);CHKERRQ(ierr); 15732c72c5SBarry Smith if (snes->work) {ierr = VecDestroyVecs(snes->nwork, &snes->work);CHKERRQ(ierr);} 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; 24a312c225SMatthew G Knepley 25a312c225SMatthew G Knepley PetscFunctionBegin; 26a312c225SMatthew G Knepley ierr = SNESReset_NGMRES(snes);CHKERRQ(ierr); 2719653cdaSPeter Brune if (snes->data) { 2819653cdaSPeter Brune SNES_NGMRES * ngmres = (SNES_NGMRES *)snes->data; 29f109b39eSPeter Brune ierr = PetscFree5(ngmres->h, ngmres->beta, ngmres->xi, ngmres->fnorms, ngmres->q);CHKERRQ(ierr); 3019653cdaSPeter Brune ierr = PetscFree(ngmres->s);CHKERRQ(ierr); 3119653cdaSPeter Brune #if PETSC_USE_COMPLEX 3219653cdaSPeter Brune ierr = PetscFree(ngmres->rwork); 3319653cdaSPeter Brune #endif 3419653cdaSPeter Brune ierr = PetscFree(ngmres->work); 3519653cdaSPeter Brune } 3619653cdaSPeter Brune ierr = PetscFree(snes->data); 37a312c225SMatthew G Knepley PetscFunctionReturn(0); 38a312c225SMatthew G Knepley } 39a312c225SMatthew G Knepley 40a312c225SMatthew G Knepley #undef __FUNCT__ 41a312c225SMatthew G Knepley #define __FUNCT__ "SNESSetUp_NGMRES" 42a312c225SMatthew G Knepley PetscErrorCode SNESSetUp_NGMRES(SNES snes) 43a312c225SMatthew G Knepley { 44a312c225SMatthew G Knepley SNES_NGMRES *ngmres = (SNES_NGMRES *) snes->data; 4519653cdaSPeter Brune PetscInt msize,hsize; 46a312c225SMatthew G Knepley PetscErrorCode ierr; 47a312c225SMatthew G Knepley 48a312c225SMatthew G Knepley PetscFunctionBegin; 49087dfb9eSxuemin msize = ngmres->msize; /* restart size */ 5019653cdaSPeter Brune hsize = msize * msize; 51087dfb9eSxuemin 52087dfb9eSxuemin 5398b3e84cSPeter Brune /* explicit least squares minimization solve */ 5419653cdaSPeter Brune ierr = PetscMalloc5(hsize,PetscScalar,&ngmres->h, 5519653cdaSPeter Brune msize,PetscScalar,&ngmres->beta, 5619653cdaSPeter Brune msize,PetscScalar,&ngmres->xi, 57f109b39eSPeter Brune msize,PetscReal, &ngmres->fnorms, 5819653cdaSPeter Brune hsize,PetscScalar,&ngmres->q);CHKERRQ(ierr); 5919653cdaSPeter Brune ngmres->nrhs = 1; 6019653cdaSPeter Brune ngmres->lda = msize; 6119653cdaSPeter Brune ngmres->ldb = msize; 6219653cdaSPeter Brune ierr = PetscMalloc(msize*sizeof(PetscScalar),&ngmres->s);CHKERRQ(ierr); 6319653cdaSPeter Brune ierr = PetscMemzero(ngmres->h, hsize*sizeof(PetscScalar));CHKERRQ(ierr); 6419653cdaSPeter Brune ierr = PetscMemzero(ngmres->q, hsize*sizeof(PetscScalar));CHKERRQ(ierr); 6519653cdaSPeter Brune ierr = PetscMemzero(ngmres->xi, msize*sizeof(PetscScalar));CHKERRQ(ierr); 6619653cdaSPeter Brune ierr = PetscMemzero(ngmres->beta,msize*sizeof(PetscScalar));CHKERRQ(ierr); 6719653cdaSPeter Brune ngmres->lwork = 12*msize; 6819653cdaSPeter Brune #if PETSC_USE_COMPLEX 6919653cdaSPeter Brune ierr = PetscMalloc(sizeof(PetscReal)*ngmres->lwork,&ngmres->rwork); 7019653cdaSPeter Brune #endif 7119653cdaSPeter Brune ierr = PetscMalloc(sizeof(PetscScalar)*ngmres->lwork,&ngmres->work); 72211b2d39SPeter Brune 73f109b39eSPeter Brune ierr = VecDuplicateVecs(snes->vec_sol, ngmres->msize, &ngmres->Xdot);CHKERRQ(ierr); 74f109b39eSPeter Brune ierr = VecDuplicateVecs(snes->vec_sol, ngmres->msize, &ngmres->Fdot);CHKERRQ(ierr); 75f109b39eSPeter Brune ierr = SNESDefaultGetWork(snes, 5);CHKERRQ(ierr); 76a312c225SMatthew G Knepley PetscFunctionReturn(0); 77a312c225SMatthew G Knepley } 78a312c225SMatthew G Knepley 79a312c225SMatthew G Knepley #undef __FUNCT__ 80a312c225SMatthew G Knepley #define __FUNCT__ "SNESSetFromOptions_NGMRES" 81a312c225SMatthew G Knepley PetscErrorCode SNESSetFromOptions_NGMRES(SNES snes) 82a312c225SMatthew G Knepley { 83a312c225SMatthew G Knepley SNES_NGMRES *ngmres = (SNES_NGMRES *) snes->data; 84a312c225SMatthew G Knepley PetscErrorCode ierr; 85dfbf837cSBarry Smith PetscBool debug; 86a312c225SMatthew G Knepley 87a312c225SMatthew G Knepley PetscFunctionBegin; 88a312c225SMatthew G Knepley ierr = PetscOptionsHead("SNES NGMRES options");CHKERRQ(ierr); 899f425c49SPeter Brune ierr = PetscOptionsBool("-snes_ngmres_additive", "Use additive variant vs. choice", "SNES", ngmres->additive, &ngmres->additive, PETSC_NULL);CHKERRQ(ierr); 9019653cdaSPeter Brune ierr = PetscOptionsInt("-snes_ngmres_m", "Number of directions", "SNES", ngmres->msize, &ngmres->msize, PETSC_NULL);CHKERRQ(ierr); 9128ed4a04SPeter Brune ierr = PetscOptionsInt("-snes_ngmres_restart_it","Tolerance iterations before restart","SNES", ngmres->restart_it, &ngmres->restart_it, PETSC_NULL);CHKERRQ(ierr); 92dfbf837cSBarry Smith ierr = PetscOptionsBool("-snes_ngmres_monitor", "Monitor actions of NGMRES", "SNES", ngmres->monitor ? PETSC_TRUE: PETSC_FALSE, &debug, PETSC_NULL);CHKERRQ(ierr); 93dfbf837cSBarry Smith if (debug) { 94dfbf837cSBarry Smith ngmres->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr); 95dfbf837cSBarry Smith } 966a7cf640SPeter Brune ierr = PetscOptionsReal("-snes_ngmres_gammaA", "Residual selection constant", "SNES", ngmres->gammaA, &ngmres->gammaA, PETSC_NULL);CHKERRQ(ierr); 976a7cf640SPeter Brune ierr = PetscOptionsReal("-snes_ngmres_gammaC", " Residual restart constant", "SNES", ngmres->gammaC, &ngmres->gammaC, PETSC_NULL);CHKERRQ(ierr); 986a7cf640SPeter Brune ierr = PetscOptionsReal("-snes_ngmres_epsilonB", "Difference selection constant", "SNES", ngmres->epsilonB, &ngmres->epsilonB, PETSC_NULL);CHKERRQ(ierr); 996a7cf640SPeter Brune ierr = PetscOptionsReal("-snes_ngmres_deltaB", "Difference residual selection constant", "SNES", ngmres->deltaB, &ngmres->deltaB, PETSC_NULL);CHKERRQ(ierr); 100a312c225SMatthew G Knepley ierr = PetscOptionsTail();CHKERRQ(ierr); 1016a7cf640SPeter Brune if ((ngmres->gammaA > ngmres->gammaC) && (ngmres->gammaC > 2.)) ngmres->gammaC = ngmres->gammaA; 102a312c225SMatthew G Knepley PetscFunctionReturn(0); 103a312c225SMatthew G Knepley } 104a312c225SMatthew G Knepley 105a312c225SMatthew G Knepley #undef __FUNCT__ 106a312c225SMatthew G Knepley #define __FUNCT__ "SNESView_NGMRES" 107a312c225SMatthew G Knepley PetscErrorCode SNESView_NGMRES(SNES snes, PetscViewer viewer) 108a312c225SMatthew G Knepley { 109a312c225SMatthew G Knepley SNES_NGMRES *ngmres = (SNES_NGMRES *) snes->data; 110a312c225SMatthew G Knepley PetscBool iascii; 111f109b39eSPeter Brune const char *cstr; 112a312c225SMatthew G Knepley PetscErrorCode ierr; 113a312c225SMatthew G Knepley 114a312c225SMatthew G Knepley PetscFunctionBegin; 115a312c225SMatthew G Knepley ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 116a312c225SMatthew G Knepley if (iascii) { 117f109b39eSPeter Brune cstr = SNESLineSearchTypeName(snes->ls_type); 118f109b39eSPeter Brune ierr = PetscViewerASCIIPrintf(viewer, " line search variant: %s\n",cstr);CHKERRQ(ierr); 119f109b39eSPeter Brune ierr = PetscViewerASCIIPrintf(viewer, " Number of stored past updates: %d\n", ngmres->msize);CHKERRQ(ierr); 120f109b39eSPeter Brune ierr = PetscViewerASCIIPrintf(viewer, " Residual selection: gammaA=%1.0e, gammaC=%1.0e\n", ngmres->gammaA, ngmres->gammaC);CHKERRQ(ierr); 121f109b39eSPeter Brune ierr = PetscViewerASCIIPrintf(viewer, " Difference restart: epsilonB=%1.0e, deltaB=%1.0e\n", ngmres->epsilonB, ngmres->deltaB);CHKERRQ(ierr); 122a312c225SMatthew G Knepley } 123a312c225SMatthew G Knepley PetscFunctionReturn(0); 124a312c225SMatthew G Knepley } 125a312c225SMatthew G Knepley 126f109b39eSPeter Brune EXTERN_C_BEGIN 127f109b39eSPeter Brune #undef __FUNCT__ 128f109b39eSPeter Brune #define __FUNCT__ "SNESLineSearchSetType_NGMRES" 129f109b39eSPeter Brune PetscErrorCode SNESLineSearchSetType_NGMRES(SNES snes, SNESLineSearchType type) 130f109b39eSPeter Brune { 131f109b39eSPeter Brune PetscErrorCode ierr; 132f109b39eSPeter Brune PetscFunctionBegin; 133f109b39eSPeter Brune 134f109b39eSPeter Brune switch (type) { 135f109b39eSPeter Brune case SNES_LS_BASIC: 136f109b39eSPeter Brune ierr = SNESLineSearchSet(snes,SNESLineSearchNo,PETSC_NULL);CHKERRQ(ierr); 137f109b39eSPeter Brune break; 138f109b39eSPeter Brune case SNES_LS_BASIC_NONORMS: 139f109b39eSPeter Brune ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr); 140f109b39eSPeter Brune break; 141f109b39eSPeter Brune case SNES_LS_QUADRATIC: 142f109b39eSPeter Brune ierr = SNESLineSearchSet(snes,SNESLineSearchQuadraticSecant,PETSC_NULL);CHKERRQ(ierr); 143f109b39eSPeter Brune break; 144f109b39eSPeter Brune default: 145f109b39eSPeter Brune SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP,"Unknown line search type."); 146f109b39eSPeter Brune break; 147f109b39eSPeter Brune } 148f109b39eSPeter Brune snes->ls_type = type; 149f109b39eSPeter Brune PetscFunctionReturn(0); 150f109b39eSPeter Brune } 151f109b39eSPeter Brune EXTERN_C_END 152f109b39eSPeter Brune 15319653cdaSPeter Brune 154a312c225SMatthew G Knepley #undef __FUNCT__ 155a312c225SMatthew G Knepley #define __FUNCT__ "SNESSolve_NGMRES" 156211b2d39SPeter Brune 157a312c225SMatthew G Knepley PetscErrorCode SNESSolve_NGMRES(SNES snes) 158a312c225SMatthew G Knepley { 159087dfb9eSxuemin SNES_NGMRES *ngmres = (SNES_NGMRES *) snes->data; 16098b3e84cSPeter Brune /* present solution, residual, and preconditioned residual */ 1619f425c49SPeter Brune Vec X, F, B, D, Y; 162f109b39eSPeter Brune 163f109b39eSPeter Brune /* candidate linear combination answers */ 1649f425c49SPeter Brune Vec XA, FA, XM, FM; 16519653cdaSPeter Brune 16698b3e84cSPeter Brune /* previous iterations to construct the subspace */ 167f109b39eSPeter Brune Vec *Fdot = ngmres->Fdot; 168f109b39eSPeter Brune Vec *Xdot = ngmres->Xdot; 16919653cdaSPeter Brune 17098b3e84cSPeter Brune /* coefficients and RHS to the minimization problem */ 17119653cdaSPeter Brune PetscScalar *beta = ngmres->beta; 17219653cdaSPeter Brune PetscScalar *xi = ngmres->xi; 1739f425c49SPeter Brune PetscReal fnorm, fMnorm, fAnorm, ynorm, xnorm = 0.0; 17419653cdaSPeter Brune PetscReal nu; 17519653cdaSPeter Brune PetscScalar alph_total = 0.; 17619653cdaSPeter Brune PetscScalar qentry; 17728ed4a04SPeter Brune PetscInt i, j, k, k_restart, l, ivec, restart_count = 0; 17819653cdaSPeter Brune 17998b3e84cSPeter Brune /* solution selection data */ 18019653cdaSPeter Brune PetscBool selectA, selectRestart; 181f109b39eSPeter Brune PetscReal dnorm, dminnorm, dcurnorm; 182f109b39eSPeter Brune PetscReal fminnorm; 18319653cdaSPeter Brune 1841e633543SBarry Smith SNESConvergedReason reason; 185f109b39eSPeter Brune PetscBool lssucceed; 186a312c225SMatthew G Knepley PetscErrorCode ierr; 187a312c225SMatthew G Knepley 188a312c225SMatthew G Knepley PetscFunctionBegin; 18998b3e84cSPeter Brune /* variable initialization */ 190a312c225SMatthew G Knepley snes->reason = SNES_CONVERGED_ITERATING; 191f109b39eSPeter Brune X = snes->vec_sol; 192f109b39eSPeter Brune F = snes->vec_func; 193f109b39eSPeter Brune B = snes->vec_rhs; 194f109b39eSPeter Brune XA = snes->vec_sol_update; 195f109b39eSPeter Brune FA = snes->work[0]; 196f109b39eSPeter Brune D = snes->work[1]; 197f109b39eSPeter Brune 198f109b39eSPeter Brune /* work for the line search */ 199f109b39eSPeter Brune Y = snes->work[2]; 2009f425c49SPeter Brune XM = snes->work[3]; 2019f425c49SPeter Brune FM = snes->work[4]; 202a312c225SMatthew G Knepley 203a312c225SMatthew G Knepley ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 204a312c225SMatthew G Knepley snes->iter = 0; 205a312c225SMatthew G Knepley snes->norm = 0.; 206a312c225SMatthew G Knepley ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 20719653cdaSPeter Brune 20898b3e84cSPeter Brune /* initialization */ 20919653cdaSPeter Brune 21098b3e84cSPeter Brune /* r = F(x) */ 211f109b39eSPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 212a312c225SMatthew G Knepley if (snes->domainerror) { 213a312c225SMatthew G Knepley snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 214a312c225SMatthew G Knepley PetscFunctionReturn(0); 215a312c225SMatthew G Knepley } 21619653cdaSPeter Brune 21798b3e84cSPeter Brune /* nu = (r, r) */ 218f109b39eSPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 219f109b39eSPeter Brune fminnorm = fnorm; 220f109b39eSPeter Brune nu = fnorm*fnorm; 221f109b39eSPeter Brune if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FP, "Infinite or not-a-number generated in function evaluation"); 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); 2479f425c49SPeter Brune ierr = SNESSolve(snes->pc, B, XM);CHKERRQ(ierr); 2488cc86e31SPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 2498cc86e31SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 2508cc86e31SPeter Brune snes->reason = SNES_DIVERGED_INNER; 2518cc86e31SPeter Brune PetscFunctionReturn(0); 2528cc86e31SPeter Brune } 2539f425c49SPeter Brune ierr = SNESComputeFunction(snes, XM, FM);CHKERRQ(ierr); 2549f425c49SPeter Brune ierr = VecNorm(FM, NORM_2, &fMnorm);CHKERRQ(ierr); 2559f425c49SPeter Brune 256d28543b3SPeter Brune } else if (snes->usegs && snes->ops->computegs) { 2578cc86e31SPeter Brune /* compute the update using the supplied Gauss-Seidel routine */ 2589f425c49SPeter Brune ierr = VecCopy(X, XM);CHKERRQ(ierr); 2599f425c49SPeter Brune ierr = SNESComputeGS(snes, B, XM);CHKERRQ(ierr); 2609f425c49SPeter Brune ierr = SNESComputeFunction(snes, XM, FM);CHKERRQ(ierr); 2619f425c49SPeter Brune ierr = VecNorm(FM, NORM_2, &fMnorm);CHKERRQ(ierr); 2629f425c49SPeter Brune 2638cc86e31SPeter Brune } else { 264f109b39eSPeter Brune /* no preconditioner -- just take gradient descent with line search */ 265f109b39eSPeter Brune ierr = VecCopy(F, Y);CHKERRQ(ierr); 266f109b39eSPeter Brune ierr = VecScale(Y, -1.0);CHKERRQ(ierr); 2679f425c49SPeter Brune ierr = (*snes->ops->linesearch)(snes,snes->lsP,X,F,Y,fnorm,xnorm,FM,XM,&ynorm,&fMnorm,&lssucceed);CHKERRQ(ierr); 268f109b39eSPeter Brune if (!lssucceed) { 269f109b39eSPeter Brune if (++snes->numFailures >= snes->maxFailures) { 270f109b39eSPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 271f109b39eSPeter Brune PetscFunctionReturn(0); 272f109b39eSPeter Brune } 273f109b39eSPeter Brune } 2746634f59bSPeter Brune } 2751e633543SBarry Smith 27698b3e84cSPeter Brune /* r = F(x) */ 2779f425c49SPeter Brune nu = fMnorm*fMnorm; 2789f425c49SPeter Brune if (fminnorm > fMnorm) fminnorm = fMnorm; /* the minimum norm is now of F^M */ 27919653cdaSPeter Brune 28098b3e84cSPeter Brune /* construct the right hand side and xi factors */ 28119653cdaSPeter Brune for (i = 0; i < l; i++) { 2829f425c49SPeter Brune ierr = VecDot(Fdot[i], FM, &xi[i]);CHKERRQ(ierr); 28319653cdaSPeter Brune beta[i] = nu - xi[i]; 28419653cdaSPeter Brune } 28519653cdaSPeter Brune 28698b3e84cSPeter Brune /* construct h */ 28719653cdaSPeter Brune for (j = 0; j < l; j++) { 28819653cdaSPeter Brune for (i = 0; i < l; i++) { 28919653cdaSPeter Brune H(i, j) = Q(i, j) - xi[i] - xi[j] + nu; 29019653cdaSPeter Brune } 29119653cdaSPeter Brune } 292f109b39eSPeter Brune 293f109b39eSPeter Brune if(l == 1) { 294f109b39eSPeter Brune /* simply set alpha[0] = beta[0] / H[0, 0] */ 295f109b39eSPeter Brune beta[0] = beta[0] / H(0, 0); 296f109b39eSPeter Brune } else { 29719653cdaSPeter Brune #ifdef PETSC_MISSING_LAPACK_GELSS 29819653cdaSPeter Brune SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "NGMRES with LS requires the LAPACK GELSS routine."); 2994a0c5b0cSMatthew G Knepley #else 30019653cdaSPeter Brune ngmres->m = PetscBLASIntCast(l); 30119653cdaSPeter Brune ngmres->n = PetscBLASIntCast(l); 30219653cdaSPeter Brune ngmres->info = PetscBLASIntCast(0); 30319653cdaSPeter Brune ngmres->rcond = -1.; 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 33419653cdaSPeter Brune if (ngmres->info < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad argument to GELSS"); 33519653cdaSPeter Brune if (ngmres->info > 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"SVD failed to converge"); 336a312c225SMatthew G Knepley #endif 337f109b39eSPeter Brune } 338a312c225SMatthew G Knepley 33919653cdaSPeter Brune alph_total = 0.; 34019653cdaSPeter Brune for (i = 0; i < l; i++) { 34119653cdaSPeter Brune alph_total += beta[i]; 342c0bbabecSxuemin } 343f109b39eSPeter Brune 3449f425c49SPeter Brune ierr = VecCopy(XM, XA);CHKERRQ(ierr); 345f109b39eSPeter Brune ierr = VecScale(XA, 1. - alph_total);CHKERRQ(ierr); 346087dfb9eSxuemin 34719653cdaSPeter Brune for(i=0;i<l;i++){ 348f109b39eSPeter Brune ierr= VecAXPY(XA, beta[i], Xdot[i]);CHKERRQ(ierr); 34919653cdaSPeter Brune } 350f109b39eSPeter Brune ierr = SNESComputeFunction(snes, XA, FA);CHKERRQ(ierr); 351f109b39eSPeter Brune ierr = VecNorm(FA, NORM_2, &fAnorm);CHKERRQ(ierr); 35219653cdaSPeter Brune 3539f425c49SPeter Brune /* differences for selection and restart */ 354f109b39eSPeter Brune ierr=VecCopy(XA, D);CHKERRQ(ierr); 355f109b39eSPeter Brune ierr=VecAXPY(D, -1.0, X);CHKERRQ(ierr); 356f109b39eSPeter Brune ierr=VecNorm(D, NORM_2, &dnorm);CHKERRQ(ierr); 357f109b39eSPeter Brune dminnorm = -1.0; 35819653cdaSPeter Brune for(i=0;i<l;i++) { 359f109b39eSPeter Brune ierr=VecCopy(XA, D);CHKERRQ(ierr); 360f109b39eSPeter Brune ierr=VecAXPY(D, -1.0, Xdot[i]);CHKERRQ(ierr); 361f109b39eSPeter Brune ierr=VecNorm(D, NORM_2, &dcurnorm);CHKERRQ(ierr); 362f109b39eSPeter Brune if((dcurnorm < dminnorm) || (dminnorm < 0.0)) dminnorm = dcurnorm; 36319653cdaSPeter Brune } 3649f425c49SPeter Brune 3659f425c49SPeter Brune /* combination (additive) or selection (multiplicative) of the N-GMRES solution */ 3669f425c49SPeter Brune if (ngmres->additive) { 3679f425c49SPeter Brune /* X = X + \lambda(XA - X) */ 3689f425c49SPeter Brune if (ngmres->monitor) { 3699f425c49SPeter Brune ierr = PetscViewerASCIIPrintf(ngmres->monitor, "||F_A||_2 = %e, ||F_M||_2 = %e\n", fAnorm, fMnorm);CHKERRQ(ierr); 3709f425c49SPeter Brune } 3719f425c49SPeter Brune ierr = VecCopy(XA, Y);CHKERRQ(ierr); 3729f425c49SPeter Brune ierr = VecAXPY(Y, -1.0, X);CHKERRQ(ierr); 3739f425c49SPeter Brune ierr = (*snes->ops->linesearch)(snes,snes->lsP,X,F,Y,fnorm,xnorm,FA,XA,&ynorm,&fnorm,&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 } 3809f425c49SPeter Brune if (ngmres->monitor) { 3819f425c49SPeter Brune ierr = PetscViewerASCIIPrintf(ngmres->monitor, "Additive solution: ||F||_2 = %e\n", fnorm);CHKERRQ(ierr); 3829f425c49SPeter Brune } 3839f425c49SPeter Brune ierr = VecCopy(FA, F);CHKERRQ(ierr); 3849f425c49SPeter Brune ierr = VecCopy(XA, X);CHKERRQ(ierr); 3859f425c49SPeter Brune } else { 3869f425c49SPeter Brune selectA = PETSC_TRUE; 3879f425c49SPeter Brune /* Conditions for choosing the accelerated answer */ 3889f425c49SPeter Brune /* Criterion A -- the norm of the function isn't increased above the minimum by too much */ 3899f425c49SPeter Brune if (fAnorm >= ngmres->gammaA*fminnorm) { 3909f425c49SPeter Brune selectA = PETSC_FALSE; 3919f425c49SPeter Brune } 3929f425c49SPeter Brune /* Criterion B -- the choice of x^A isn't too close to some other choice */ 393*cf22de31SPeter Brune if (ngmres->epsilonB*dnorm<dminnorm || PetscSqrtReal(fnorm)<ngmres->deltaB*PetscSqrtReal(fminnorm)) { 39419653cdaSPeter Brune } else { 39519653cdaSPeter Brune selectA=PETSC_FALSE; 39619653cdaSPeter Brune } 39719653cdaSPeter Brune if (selectA) { 398dfbf837cSBarry Smith if (ngmres->monitor) { 399f109b39eSPeter Brune ierr = PetscViewerASCIIPrintf(ngmres->monitor, "picked X_A, ||F_A||_2 = %e, ||F_M||_2 = %e\n", fAnorm, fnorm);CHKERRQ(ierr); 400dfbf837cSBarry Smith } 40198b3e84cSPeter Brune /* copy it over */ 402f109b39eSPeter Brune fnorm = fAnorm; 403f109b39eSPeter Brune nu = fnorm*fnorm; 404f109b39eSPeter Brune ierr = VecCopy(FA, F);CHKERRQ(ierr); 405f109b39eSPeter Brune ierr = VecCopy(XA, X);CHKERRQ(ierr); 40619653cdaSPeter Brune } else { 407dfbf837cSBarry Smith if (ngmres->monitor) { 408f109b39eSPeter Brune ierr = PetscViewerASCIIPrintf(ngmres->monitor, "picked X_M, ||F_A||_2 = %e, ||F_M||_2 = %e\n", fAnorm, fnorm);CHKERRQ(ierr); 409dfbf837cSBarry Smith } 4109f425c49SPeter Brune fnorm = fMnorm; 4119f425c49SPeter Brune nu = fnorm*fnorm; 4129f425c49SPeter Brune ierr = VecCopy(FM, F);CHKERRQ(ierr); 4139f425c49SPeter Brune ierr = VecCopy(XM, X);CHKERRQ(ierr); 4149f425c49SPeter Brune } 41519653cdaSPeter Brune } 416211b2d39SPeter Brune 41719653cdaSPeter Brune selectRestart = PETSC_FALSE; 41819653cdaSPeter Brune 41998b3e84cSPeter Brune /* difference stagnation restart */ 420*cf22de31SPeter Brune if((ngmres->epsilonB*dnorm > dminnorm) && (PetscSqrtReal(fAnorm) > ngmres->deltaB*PetscSqrtReal(fminnorm))) { 421dfbf837cSBarry Smith if (ngmres->monitor) { 422f109b39eSPeter Brune ierr = PetscViewerASCIIPrintf(ngmres->monitor, "difference restart: %e > %e\n", ngmres->epsilonB*dnorm, dminnorm);CHKERRQ(ierr); 423dfbf837cSBarry Smith } 42419653cdaSPeter Brune selectRestart = PETSC_TRUE; 42519653cdaSPeter Brune } 42698b3e84cSPeter Brune /* residual stagnation restart */ 427*cf22de31SPeter Brune if (PetscSqrtReal(fAnorm) > ngmres->gammaC*PetscSqrtReal(fminnorm)) { 428dfbf837cSBarry Smith if (ngmres->monitor) { 429*cf22de31SPeter Brune ierr = PetscViewerASCIIPrintf(ngmres->monitor, "residual restart: %e > %e\n", PetscSqrtReal(fAnorm), ngmres->gammaC*PetscSqrtReal(fminnorm));CHKERRQ(ierr); 430dfbf837cSBarry Smith } 43119653cdaSPeter Brune selectRestart = PETSC_TRUE; 43219653cdaSPeter Brune } 43319653cdaSPeter Brune 43428ed4a04SPeter Brune /* if the restart conditions persist for more than restart_it iterations, restart. */ 43519653cdaSPeter Brune if (selectRestart) { 43628ed4a04SPeter Brune restart_count++; 43728ed4a04SPeter Brune } else { 43828ed4a04SPeter Brune restart_count = 0; 43928ed4a04SPeter Brune } 44028ed4a04SPeter Brune 44128ed4a04SPeter Brune /* restart after restart conditions have persisted for a fixed number of iterations */ 44228ed4a04SPeter Brune if (restart_count >= ngmres->restart_it) { 443dfbf837cSBarry Smith if (ngmres->monitor){ 444dfbf837cSBarry Smith ierr = PetscViewerASCIIPrintf(ngmres->monitor, "Restarted at iteration %d\n", k_restart);CHKERRQ(ierr); 445dfbf837cSBarry Smith } 44628ed4a04SPeter Brune restart_count = 0; 44719653cdaSPeter Brune k_restart = 1; 44819653cdaSPeter Brune l = 1; 44998b3e84cSPeter Brune /* q_{00} = nu */ 450f109b39eSPeter Brune ngmres->fnorms[0] = fnorm; 451f109b39eSPeter Brune nu = fnorm*fnorm; 45219653cdaSPeter Brune Q(0,0) = nu; 453f109b39eSPeter Brune /* Fdot[0] = F */ 454f109b39eSPeter Brune ierr = VecCopy(X, Xdot[0]);CHKERRQ(ierr); 455f109b39eSPeter Brune ierr = VecCopy(F, Fdot[0]);CHKERRQ(ierr); 45619653cdaSPeter Brune } else { 45798b3e84cSPeter Brune /* select the current size of the subspace */ 4581e633543SBarry Smith if (l < ngmres->msize) l++; 45919653cdaSPeter Brune k_restart++; 46098b3e84cSPeter Brune /* place the current entry in the list of previous entries */ 461f109b39eSPeter Brune ierr = VecCopy(F, Fdot[ivec]);CHKERRQ(ierr); 462f109b39eSPeter Brune ierr = VecCopy(X, Xdot[ivec]);CHKERRQ(ierr); 463f109b39eSPeter Brune ngmres->fnorms[ivec] = fnorm; 464f109b39eSPeter Brune if (fminnorm > fnorm) fminnorm = fnorm; /* the minimum norm is now of r^A */ 46519653cdaSPeter Brune for (i = 0; i < l; i++) { 466f109b39eSPeter Brune ierr = VecDot(F, Fdot[i], &qentry);CHKERRQ(ierr); 46719653cdaSPeter Brune Q(i, ivec) = qentry; 46819653cdaSPeter Brune Q(ivec, i) = qentry; 46919653cdaSPeter Brune } 47019653cdaSPeter Brune } 47119653cdaSPeter Brune 4728409ca45SMatthew G Knepley ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 473087dfb9eSxuemin snes->iter = k; 474f109b39eSPeter Brune snes->norm = fnorm; 4758409ca45SMatthew G Knepley ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 4768409ca45SMatthew G Knepley SNESLogConvHistory(snes, snes->norm, snes->iter); 4778409ca45SMatthew G Knepley ierr = SNESMonitor(snes, snes->iter, snes->norm);CHKERRQ(ierr); 4788409ca45SMatthew G Knepley 479f109b39eSPeter Brune ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 480087dfb9eSxuemin if (snes->reason) PetscFunctionReturn(0); 481a312c225SMatthew G Knepley } 482a312c225SMatthew G Knepley snes->reason = SNES_DIVERGED_MAX_IT; 483a312c225SMatthew G Knepley PetscFunctionReturn(0); 484a312c225SMatthew G Knepley } 485a312c225SMatthew G Knepley 486dfbf837cSBarry Smith /*MC 487dfbf837cSBarry Smith SNESNGMRES - The Nonlinear Generalized Minimum Residual (NGMRES) method of Oosterlee and Washio. 488a312c225SMatthew G Knepley 489dfbf837cSBarry Smith Level: beginner 490dfbf837cSBarry Smith 491dfbf837cSBarry Smith "Krylov Subspace Acceleration of Nonlinear Multigrid with Application to Recirculating Flows", C. W. Oosterlee and T. Washio, 492dfbf837cSBarry Smith SIAM Journal on Scientific Computing, 21(5), 2000. 493dfbf837cSBarry Smith 494dfbf837cSBarry Smith This is also the same as the algorithm called Anderson acceleration introduced in "D. G. Anderson. Iterative procedures for nonlinear integral equations. 495dfbf837cSBarry Smith J. Assoc. Comput. Mach., 12:547–560, 1965." 496dfbf837cSBarry Smith 497dfbf837cSBarry Smith .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types) 498dfbf837cSBarry Smith M*/ 499a312c225SMatthew G Knepley 500a312c225SMatthew G Knepley EXTERN_C_BEGIN 501a312c225SMatthew G Knepley #undef __FUNCT__ 502a312c225SMatthew G Knepley #define __FUNCT__ "SNESCreate_NGMRES" 503a312c225SMatthew G Knepley PetscErrorCode SNESCreate_NGMRES(SNES snes) 504a312c225SMatthew G Knepley { 505a312c225SMatthew G Knepley SNES_NGMRES *ngmres; 506a312c225SMatthew G Knepley PetscErrorCode ierr; 507a312c225SMatthew G Knepley 508a312c225SMatthew G Knepley PetscFunctionBegin; 509a312c225SMatthew G Knepley snes->ops->destroy = SNESDestroy_NGMRES; 510a312c225SMatthew G Knepley snes->ops->setup = SNESSetUp_NGMRES; 511a312c225SMatthew G Knepley snes->ops->setfromoptions = SNESSetFromOptions_NGMRES; 512a312c225SMatthew G Knepley snes->ops->view = SNESView_NGMRES; 513a312c225SMatthew G Knepley snes->ops->solve = SNESSolve_NGMRES; 514a312c225SMatthew G Knepley snes->ops->reset = SNESReset_NGMRES; 515a312c225SMatthew G Knepley 51642f4f86dSBarry Smith snes->usespc = PETSC_TRUE; 5172c155ee1SBarry Smith snes->usesksp = PETSC_FALSE; 5182c155ee1SBarry Smith 519a312c225SMatthew G Knepley ierr = PetscNewLog(snes, SNES_NGMRES, &ngmres);CHKERRQ(ierr); 520a312c225SMatthew G Knepley snes->data = (void*) ngmres; 52119653cdaSPeter Brune ngmres->msize = 10; 52219653cdaSPeter Brune 52328ed4a04SPeter Brune ngmres->restart_it = 2; 524f109b39eSPeter Brune ngmres->gammaA = 2.0; 525f109b39eSPeter Brune ngmres->gammaC = 2.0; 526cac108bcSPeter Brune ngmres->deltaB = 0.9; 527cac108bcSPeter Brune ngmres->epsilonB = 0.1; 528f109b39eSPeter Brune ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetType_C","SNESLineSearchSetType_NGMRES",SNESLineSearchSetType_NGMRES);CHKERRQ(ierr); 529f109b39eSPeter Brune ierr = SNESLineSearchSetType(snes, SNES_LS_QUADRATIC);CHKERRQ(ierr); 530a312c225SMatthew G Knepley PetscFunctionReturn(0); 531a312c225SMatthew G Knepley } 532a312c225SMatthew G Knepley EXTERN_C_END 533