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); 8919653cdaSPeter Brune ierr = PetscOptionsInt("-snes_ngmres_m", "Number of directions", "SNES", ngmres->msize, &ngmres->msize, PETSC_NULL);CHKERRQ(ierr); 9028ed4a04SPeter Brune ierr = PetscOptionsInt("-snes_ngmres_restart_it","Tolerance iterations before restart","SNES", ngmres->restart_it, &ngmres->restart_it, PETSC_NULL);CHKERRQ(ierr); 91dfbf837cSBarry Smith ierr = PetscOptionsBool("-snes_ngmres_monitor", "Monitor actions of NGMRES", "SNES", ngmres->monitor ? PETSC_TRUE: PETSC_FALSE, &debug, PETSC_NULL);CHKERRQ(ierr); 92dfbf837cSBarry Smith if (debug) { 93dfbf837cSBarry Smith ngmres->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr); 94dfbf837cSBarry Smith } 956a7cf640SPeter Brune ierr = PetscOptionsReal("-snes_ngmres_gammaA", "Residual selection constant", "SNES", ngmres->gammaA, &ngmres->gammaA, PETSC_NULL);CHKERRQ(ierr); 966a7cf640SPeter Brune ierr = PetscOptionsReal("-snes_ngmres_gammaC", " Residual restart constant", "SNES", ngmres->gammaC, &ngmres->gammaC, PETSC_NULL);CHKERRQ(ierr); 976a7cf640SPeter Brune ierr = PetscOptionsReal("-snes_ngmres_epsilonB", "Difference selection constant", "SNES", ngmres->epsilonB, &ngmres->epsilonB, PETSC_NULL);CHKERRQ(ierr); 986a7cf640SPeter Brune ierr = PetscOptionsReal("-snes_ngmres_deltaB", "Difference residual selection constant", "SNES", ngmres->deltaB, &ngmres->deltaB, PETSC_NULL);CHKERRQ(ierr); 99a312c225SMatthew G Knepley ierr = PetscOptionsTail();CHKERRQ(ierr); 1006a7cf640SPeter Brune if ((ngmres->gammaA > ngmres->gammaC) && (ngmres->gammaC > 2.)) ngmres->gammaC = ngmres->gammaA; 101a312c225SMatthew G Knepley PetscFunctionReturn(0); 102a312c225SMatthew G Knepley } 103a312c225SMatthew G Knepley 104a312c225SMatthew G Knepley #undef __FUNCT__ 105a312c225SMatthew G Knepley #define __FUNCT__ "SNESView_NGMRES" 106a312c225SMatthew G Knepley PetscErrorCode SNESView_NGMRES(SNES snes, PetscViewer viewer) 107a312c225SMatthew G Knepley { 108a312c225SMatthew G Knepley SNES_NGMRES *ngmres = (SNES_NGMRES *) snes->data; 109a312c225SMatthew G Knepley PetscBool iascii; 110f109b39eSPeter Brune const char *cstr; 111a312c225SMatthew G Knepley PetscErrorCode ierr; 112a312c225SMatthew G Knepley 113a312c225SMatthew G Knepley PetscFunctionBegin; 114a312c225SMatthew G Knepley ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 115a312c225SMatthew G Knepley if (iascii) { 116f109b39eSPeter Brune cstr = SNESLineSearchTypeName(snes->ls_type); 117f109b39eSPeter Brune ierr = PetscViewerASCIIPrintf(viewer, " line search variant: %s\n",cstr);CHKERRQ(ierr); 118f109b39eSPeter Brune ierr = PetscViewerASCIIPrintf(viewer, " Number of stored past updates: %d\n", ngmres->msize);CHKERRQ(ierr); 119f109b39eSPeter Brune ierr = PetscViewerASCIIPrintf(viewer, " Residual selection: gammaA=%1.0e, gammaC=%1.0e\n", ngmres->gammaA, ngmres->gammaC);CHKERRQ(ierr); 120f109b39eSPeter Brune ierr = PetscViewerASCIIPrintf(viewer, " Difference restart: epsilonB=%1.0e, deltaB=%1.0e\n", ngmres->epsilonB, ngmres->deltaB);CHKERRQ(ierr); 121a312c225SMatthew G Knepley } 122a312c225SMatthew G Knepley PetscFunctionReturn(0); 123a312c225SMatthew G Knepley } 124a312c225SMatthew G Knepley 125f109b39eSPeter Brune EXTERN_C_BEGIN 126f109b39eSPeter Brune #undef __FUNCT__ 127f109b39eSPeter Brune #define __FUNCT__ "SNESLineSearchSetType_NGMRES" 128f109b39eSPeter Brune PetscErrorCode SNESLineSearchSetType_NGMRES(SNES snes, SNESLineSearchType type) 129f109b39eSPeter Brune { 130f109b39eSPeter Brune PetscErrorCode ierr; 131f109b39eSPeter Brune PetscFunctionBegin; 132f109b39eSPeter Brune 133f109b39eSPeter Brune switch (type) { 134f109b39eSPeter Brune case SNES_LS_BASIC: 135f109b39eSPeter Brune ierr = SNESLineSearchSet(snes,SNESLineSearchNo,PETSC_NULL);CHKERRQ(ierr); 136f109b39eSPeter Brune break; 137f109b39eSPeter Brune case SNES_LS_BASIC_NONORMS: 138f109b39eSPeter Brune ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr); 139f109b39eSPeter Brune break; 140f109b39eSPeter Brune case SNES_LS_QUADRATIC: 141f109b39eSPeter Brune ierr = SNESLineSearchSet(snes,SNESLineSearchQuadraticSecant,PETSC_NULL);CHKERRQ(ierr); 142f109b39eSPeter Brune break; 143f109b39eSPeter Brune default: 144f109b39eSPeter Brune SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP,"Unknown line search type."); 145f109b39eSPeter Brune break; 146f109b39eSPeter Brune } 147f109b39eSPeter Brune snes->ls_type = type; 148f109b39eSPeter Brune PetscFunctionReturn(0); 149f109b39eSPeter Brune } 150f109b39eSPeter Brune EXTERN_C_END 151f109b39eSPeter Brune 15219653cdaSPeter Brune 153a312c225SMatthew G Knepley #undef __FUNCT__ 154a312c225SMatthew G Knepley #define __FUNCT__ "SNESSolve_NGMRES" 155211b2d39SPeter Brune 156a312c225SMatthew G Knepley PetscErrorCode SNESSolve_NGMRES(SNES snes) 157a312c225SMatthew G Knepley { 158087dfb9eSxuemin SNES_NGMRES *ngmres = (SNES_NGMRES *) snes->data; 15998b3e84cSPeter Brune /* present solution, residual, and preconditioned residual */ 160f109b39eSPeter Brune Vec X, F, B, D, G, W, Y; 161f109b39eSPeter Brune 162f109b39eSPeter Brune /* candidate linear combination answers */ 163f109b39eSPeter Brune Vec XA, FA; 16419653cdaSPeter Brune 16598b3e84cSPeter Brune /* previous iterations to construct the subspace */ 166f109b39eSPeter Brune Vec *Fdot = ngmres->Fdot; 167f109b39eSPeter Brune Vec *Xdot = ngmres->Xdot; 16819653cdaSPeter Brune 16998b3e84cSPeter Brune /* coefficients and RHS to the minimization problem */ 17019653cdaSPeter Brune PetscScalar *beta = ngmres->beta; 17119653cdaSPeter Brune PetscScalar *xi = ngmres->xi; 172f109b39eSPeter Brune PetscReal fnorm, fAnorm, gnorm, ynorm, xnorm = 0.0; 17319653cdaSPeter Brune PetscReal nu; 17419653cdaSPeter Brune PetscScalar alph_total = 0.; 17519653cdaSPeter Brune PetscScalar qentry; 17628ed4a04SPeter Brune PetscInt i, j, k, k_restart, l, ivec, restart_count = 0; 17719653cdaSPeter Brune 17898b3e84cSPeter Brune /* solution selection data */ 17919653cdaSPeter Brune PetscBool selectA, selectRestart; 180f109b39eSPeter Brune PetscReal dnorm, dminnorm, dcurnorm; 181f109b39eSPeter Brune PetscReal fminnorm; 18219653cdaSPeter Brune 1831e633543SBarry Smith SNESConvergedReason reason; 184f109b39eSPeter Brune PetscBool lssucceed; 185a312c225SMatthew G Knepley PetscErrorCode ierr; 186a312c225SMatthew G Knepley 187a312c225SMatthew G Knepley PetscFunctionBegin; 18898b3e84cSPeter Brune /* variable initialization */ 189a312c225SMatthew G Knepley snes->reason = SNES_CONVERGED_ITERATING; 190f109b39eSPeter Brune X = snes->vec_sol; 191f109b39eSPeter Brune F = snes->vec_func; 192f109b39eSPeter Brune B = snes->vec_rhs; 193f109b39eSPeter Brune XA = snes->vec_sol_update; 194f109b39eSPeter Brune FA = snes->work[0]; 195f109b39eSPeter Brune D = snes->work[1]; 196f109b39eSPeter Brune 197f109b39eSPeter Brune /* work for the line search */ 198f109b39eSPeter Brune Y = snes->work[2]; 199f109b39eSPeter Brune G = snes->work[3]; 200f109b39eSPeter Brune W = snes->work[4]; 201a312c225SMatthew G Knepley 202a312c225SMatthew G Knepley ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 203a312c225SMatthew G Knepley snes->iter = 0; 204a312c225SMatthew G Knepley snes->norm = 0.; 205a312c225SMatthew G Knepley ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 20619653cdaSPeter Brune 20798b3e84cSPeter Brune /* initialization */ 20819653cdaSPeter Brune 20998b3e84cSPeter Brune /* r = F(x) */ 210f109b39eSPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 211a312c225SMatthew G Knepley if (snes->domainerror) { 212a312c225SMatthew G Knepley snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 213a312c225SMatthew G Knepley PetscFunctionReturn(0); 214a312c225SMatthew G Knepley } 21519653cdaSPeter Brune 21698b3e84cSPeter Brune /* nu = (r, r) */ 217f109b39eSPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 218f109b39eSPeter Brune fminnorm = fnorm; 219f109b39eSPeter Brune nu = fnorm*fnorm; 220f109b39eSPeter Brune if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FP, "Infinite or not-a-number generated in function evaluation"); 22119653cdaSPeter Brune 22298b3e84cSPeter Brune /* q_{00} = nu */ 22319653cdaSPeter Brune Q(0,0) = nu; 224f109b39eSPeter Brune ngmres->fnorms[0] = fnorm; 225f109b39eSPeter Brune /* Fdot[0] = F */ 226f109b39eSPeter Brune ierr = VecCopy(X, Xdot[0]);CHKERRQ(ierr); 227f109b39eSPeter Brune ierr = VecCopy(F, Fdot[0]);CHKERRQ(ierr); 228087dfb9eSxuemin 229a312c225SMatthew G Knepley ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 230f109b39eSPeter Brune snes->norm = fnorm; 231a312c225SMatthew G Knepley ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 232f109b39eSPeter Brune SNESLogConvHistory(snes, fnorm, 0); 233f109b39eSPeter Brune ierr = SNESMonitor(snes, 0, fnorm);CHKERRQ(ierr); 234f109b39eSPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 235a312c225SMatthew G Knepley if (snes->reason) PetscFunctionReturn(0); 236a312c225SMatthew G Knepley 23719653cdaSPeter Brune k_restart = 1; 23819653cdaSPeter Brune l = 1; 23909c08436SPeter Brune for (k=1; k < snes->max_its+1; k++) { 24098b3e84cSPeter Brune /* select which vector of the stored subspace will be updated */ 24198b3e84cSPeter Brune ivec = k_restart % ngmres->msize; /* replace the last used part of the subspace */ 24219653cdaSPeter Brune 24398b3e84cSPeter Brune /* Computation of x^M */ 2448cc86e31SPeter Brune if (snes->pc) { 2458cc86e31SPeter Brune ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr); 2468cc86e31SPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 2478cc86e31SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 2488cc86e31SPeter Brune snes->reason = SNES_DIVERGED_INNER; 2498cc86e31SPeter Brune PetscFunctionReturn(0); 2508cc86e31SPeter Brune } 2518cc86e31SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 2528cc86e31SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 253d28543b3SPeter Brune } else if (snes->usegs && snes->ops->computegs) { 2548cc86e31SPeter Brune /* compute the update using the supplied Gauss-Seidel routine */ 2558cc86e31SPeter Brune ierr = SNESComputeGS(snes, B, X);CHKERRQ(ierr); 2568cc86e31SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 2578cc86e31SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 2588cc86e31SPeter Brune } else { 259f109b39eSPeter Brune /* no preconditioner -- just take gradient descent with line search */ 260f109b39eSPeter Brune ierr = VecCopy(F, Y);CHKERRQ(ierr); 261f109b39eSPeter Brune ierr = VecScale(Y, -1.0);CHKERRQ(ierr); 262f109b39eSPeter Brune ierr = (*snes->ops->linesearch)(snes,snes->lsP,X,F,Y,fnorm,xnorm,G,W,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr); 263f109b39eSPeter Brune if (!lssucceed) { 264f109b39eSPeter Brune if (++snes->numFailures >= snes->maxFailures) { 265f109b39eSPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 266f109b39eSPeter Brune PetscFunctionReturn(0); 267f109b39eSPeter Brune } 268f109b39eSPeter Brune } 269f109b39eSPeter Brune fnorm = gnorm; 270f109b39eSPeter Brune ierr = VecCopy(G, F);CHKERRQ(ierr); 271f109b39eSPeter Brune ierr = VecCopy(W, X);CHKERRQ(ierr); 2726634f59bSPeter Brune } 2731e633543SBarry Smith 27498b3e84cSPeter Brune /* r = F(x) */ 275f109b39eSPeter Brune nu = fnorm*fnorm; 276f109b39eSPeter Brune if (fminnorm > fnorm) fminnorm = fnorm; /* the minimum norm is now of F^M */ 27719653cdaSPeter Brune 27898b3e84cSPeter Brune /* construct the right hand side and xi factors */ 27919653cdaSPeter Brune for (i = 0; i < l; i++) { 280f109b39eSPeter Brune ierr = VecDot(Fdot[i], F, &xi[i]);CHKERRQ(ierr); 28119653cdaSPeter Brune beta[i] = nu - xi[i]; 28219653cdaSPeter Brune } 28319653cdaSPeter Brune 28498b3e84cSPeter Brune /* construct h */ 28519653cdaSPeter Brune for (j = 0; j < l; j++) { 28619653cdaSPeter Brune for (i = 0; i < l; i++) { 28719653cdaSPeter Brune H(i, j) = Q(i, j) - xi[i] - xi[j] + nu; 28819653cdaSPeter Brune } 28919653cdaSPeter Brune } 290f109b39eSPeter Brune 291f109b39eSPeter Brune if(l == 1) { 292f109b39eSPeter Brune /* simply set alpha[0] = beta[0] / H[0, 0] */ 293f109b39eSPeter Brune beta[0] = beta[0] / H(0, 0); 294f109b39eSPeter Brune } else { 29519653cdaSPeter Brune #ifdef PETSC_MISSING_LAPACK_GELSS 29619653cdaSPeter Brune SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "NGMRES with LS requires the LAPACK GELSS routine."); 2974a0c5b0cSMatthew G Knepley #else 29819653cdaSPeter Brune ngmres->m = PetscBLASIntCast(l); 29919653cdaSPeter Brune ngmres->n = PetscBLASIntCast(l); 30019653cdaSPeter Brune ngmres->info = PetscBLASIntCast(0); 30119653cdaSPeter Brune ngmres->rcond = -1.; 30219653cdaSPeter Brune #ifdef PETSC_USE_COMPLEX 30319653cdaSPeter Brune LAPACKgelss_(&ngmres->m, 30419653cdaSPeter Brune &ngmres->n, 30519653cdaSPeter Brune &ngmres->nrhs, 30619653cdaSPeter Brune ngmres->h, 30719653cdaSPeter Brune &ngmres->lda, 30819653cdaSPeter Brune ngmres->beta, 30919653cdaSPeter Brune &ngmres->ldb, 31019653cdaSPeter Brune ngmres->s, 31119653cdaSPeter Brune &ngmres->rcond, 31219653cdaSPeter Brune &ngmres->rank, 31319653cdaSPeter Brune ngmres->work, 31419653cdaSPeter Brune &ngmres->lwork, 31519653cdaSPeter Brune ngmres->rwork, 31619653cdaSPeter Brune &ngmres->info); 31719653cdaSPeter Brune #else 31819653cdaSPeter Brune LAPACKgelss_(&ngmres->m, 31919653cdaSPeter Brune &ngmres->n, 32019653cdaSPeter Brune &ngmres->nrhs, 32119653cdaSPeter Brune ngmres->h, 32219653cdaSPeter Brune &ngmres->lda, 32319653cdaSPeter Brune ngmres->beta, 32419653cdaSPeter Brune &ngmres->ldb, 32519653cdaSPeter Brune ngmres->s, 32619653cdaSPeter Brune &ngmres->rcond, 32719653cdaSPeter Brune &ngmres->rank, 32819653cdaSPeter Brune ngmres->work, 32919653cdaSPeter Brune &ngmres->lwork, 33019653cdaSPeter Brune &ngmres->info); 33119653cdaSPeter Brune #endif 33219653cdaSPeter Brune if (ngmres->info < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad argument to GELSS"); 33319653cdaSPeter Brune if (ngmres->info > 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"SVD failed to converge"); 334a312c225SMatthew G Knepley #endif 335f109b39eSPeter Brune } 336a312c225SMatthew G Knepley 33719653cdaSPeter Brune alph_total = 0.; 33819653cdaSPeter Brune for (i = 0; i < l; i++) { 33919653cdaSPeter Brune alph_total += beta[i]; 340c0bbabecSxuemin } 341f109b39eSPeter Brune 342f109b39eSPeter Brune ierr = VecCopy(X, XA);CHKERRQ(ierr); 343f109b39eSPeter Brune ierr = VecScale(XA, 1. - alph_total);CHKERRQ(ierr); 344087dfb9eSxuemin 34519653cdaSPeter Brune for(i=0;i<l;i++){ 346f109b39eSPeter Brune ierr= VecAXPY(XA, beta[i], Xdot[i]);CHKERRQ(ierr); 34719653cdaSPeter Brune } 348f109b39eSPeter Brune ierr = SNESComputeFunction(snes, XA, FA);CHKERRQ(ierr); 349f109b39eSPeter Brune ierr = VecNorm(FA, NORM_2, &fAnorm);CHKERRQ(ierr); 35019653cdaSPeter Brune 35119653cdaSPeter Brune selectA = PETSC_TRUE; 35298b3e84cSPeter Brune /* Conditions for choosing the accelerated answer */ 35319653cdaSPeter Brune 35498b3e84cSPeter Brune /* Criterion A -- the norm of the function isn't increased above the minimum by too much */ 355f109b39eSPeter Brune if (fAnorm >= ngmres->gammaA*fminnorm) { 35619653cdaSPeter Brune selectA = PETSC_FALSE; 357a312c225SMatthew G Knepley } 358087dfb9eSxuemin 35998b3e84cSPeter Brune /* Criterion B -- the choice of x^A isn't too close to some other choice */ 360f109b39eSPeter Brune ierr=VecCopy(XA, D);CHKERRQ(ierr); 361f109b39eSPeter Brune ierr=VecAXPY(D, -1.0, X);CHKERRQ(ierr); 362f109b39eSPeter Brune ierr=VecNorm(D, NORM_2, &dnorm);CHKERRQ(ierr); 363f109b39eSPeter Brune dminnorm = -1.0; 36419653cdaSPeter Brune for(i=0;i<l;i++) { 365f109b39eSPeter Brune ierr=VecCopy(XA, D);CHKERRQ(ierr); 366f109b39eSPeter Brune ierr=VecAXPY(D, -1.0, Xdot[i]);CHKERRQ(ierr); 367f109b39eSPeter Brune ierr=VecNorm(D, NORM_2, &dcurnorm);CHKERRQ(ierr); 368f109b39eSPeter Brune if((dcurnorm < dminnorm) || (dminnorm < 0.0)) dminnorm = dcurnorm; 36919653cdaSPeter Brune } 370*cf22de31SPeter Brune if (ngmres->epsilonB*dnorm<dminnorm || PetscSqrtReal(fnorm)<ngmres->deltaB*PetscSqrtReal(fminnorm)) { 37119653cdaSPeter Brune } else { 37219653cdaSPeter Brune selectA=PETSC_FALSE; 37319653cdaSPeter Brune } 374211b2d39SPeter Brune 375087dfb9eSxuemin 37619653cdaSPeter Brune if (selectA) { 377dfbf837cSBarry Smith if (ngmres->monitor) { 378f109b39eSPeter Brune ierr = PetscViewerASCIIPrintf(ngmres->monitor, "picked X_A, ||F_A||_2 = %e, ||F_M||_2 = %e\n", fAnorm, fnorm);CHKERRQ(ierr); 379dfbf837cSBarry Smith } 38098b3e84cSPeter Brune /* copy it over */ 381f109b39eSPeter Brune fnorm = fAnorm; 382f109b39eSPeter Brune nu = fnorm*fnorm; 383f109b39eSPeter Brune ierr = VecCopy(FA, F);CHKERRQ(ierr); 384f109b39eSPeter Brune ierr = VecCopy(XA, X);CHKERRQ(ierr); 38519653cdaSPeter Brune } else { 386dfbf837cSBarry Smith if (ngmres->monitor) { 387f109b39eSPeter Brune ierr = PetscViewerASCIIPrintf(ngmres->monitor, "picked X_M, ||F_A||_2 = %e, ||F_M||_2 = %e\n", fAnorm, fnorm);CHKERRQ(ierr); 388dfbf837cSBarry Smith } 38919653cdaSPeter Brune } 390211b2d39SPeter Brune 39119653cdaSPeter Brune selectRestart = PETSC_FALSE; 39219653cdaSPeter Brune 39398b3e84cSPeter Brune /* difference stagnation restart */ 394*cf22de31SPeter Brune if((ngmres->epsilonB*dnorm > dminnorm) && (PetscSqrtReal(fAnorm) > ngmres->deltaB*PetscSqrtReal(fminnorm))) { 395dfbf837cSBarry Smith if (ngmres->monitor) { 396f109b39eSPeter Brune ierr = PetscViewerASCIIPrintf(ngmres->monitor, "difference restart: %e > %e\n", ngmres->epsilonB*dnorm, dminnorm);CHKERRQ(ierr); 397dfbf837cSBarry Smith } 39819653cdaSPeter Brune selectRestart = PETSC_TRUE; 39919653cdaSPeter Brune } 40098b3e84cSPeter Brune /* residual stagnation restart */ 401*cf22de31SPeter Brune if (PetscSqrtReal(fAnorm) > ngmres->gammaC*PetscSqrtReal(fminnorm)) { 402dfbf837cSBarry Smith if (ngmres->monitor) { 403*cf22de31SPeter Brune ierr = PetscViewerASCIIPrintf(ngmres->monitor, "residual restart: %e > %e\n", PetscSqrtReal(fAnorm), ngmres->gammaC*PetscSqrtReal(fminnorm));CHKERRQ(ierr); 404dfbf837cSBarry Smith } 40519653cdaSPeter Brune selectRestart = PETSC_TRUE; 40619653cdaSPeter Brune } 40719653cdaSPeter Brune 40828ed4a04SPeter Brune /* if the restart conditions persist for more than restart_it iterations, restart. */ 40919653cdaSPeter Brune if (selectRestart) { 41028ed4a04SPeter Brune restart_count++; 41128ed4a04SPeter Brune } else { 41228ed4a04SPeter Brune restart_count = 0; 41328ed4a04SPeter Brune } 41428ed4a04SPeter Brune 41528ed4a04SPeter Brune /* restart after restart conditions have persisted for a fixed number of iterations */ 41628ed4a04SPeter Brune if (restart_count >= ngmres->restart_it) { 417dfbf837cSBarry Smith if (ngmres->monitor){ 418dfbf837cSBarry Smith ierr = PetscViewerASCIIPrintf(ngmres->monitor, "Restarted at iteration %d\n", k_restart);CHKERRQ(ierr); 419dfbf837cSBarry Smith } 42028ed4a04SPeter Brune restart_count = 0; 42119653cdaSPeter Brune k_restart = 1; 42219653cdaSPeter Brune l = 1; 42398b3e84cSPeter Brune /* q_{00} = nu */ 424f109b39eSPeter Brune ngmres->fnorms[0] = fnorm; 425f109b39eSPeter Brune nu = fnorm*fnorm; 42619653cdaSPeter Brune Q(0,0) = nu; 427f109b39eSPeter Brune /* Fdot[0] = F */ 428f109b39eSPeter Brune ierr = VecCopy(X, Xdot[0]);CHKERRQ(ierr); 429f109b39eSPeter Brune ierr = VecCopy(F, Fdot[0]);CHKERRQ(ierr); 43019653cdaSPeter Brune } else { 43198b3e84cSPeter Brune /* select the current size of the subspace */ 4321e633543SBarry Smith if (l < ngmres->msize) l++; 43319653cdaSPeter Brune k_restart++; 43498b3e84cSPeter Brune /* place the current entry in the list of previous entries */ 435f109b39eSPeter Brune ierr = VecCopy(F, Fdot[ivec]);CHKERRQ(ierr); 436f109b39eSPeter Brune ierr = VecCopy(X, Xdot[ivec]);CHKERRQ(ierr); 437f109b39eSPeter Brune ngmres->fnorms[ivec] = fnorm; 438f109b39eSPeter Brune if (fminnorm > fnorm) fminnorm = fnorm; /* the minimum norm is now of r^A */ 43919653cdaSPeter Brune for (i = 0; i < l; i++) { 440f109b39eSPeter Brune ierr = VecDot(F, Fdot[i], &qentry);CHKERRQ(ierr); 44119653cdaSPeter Brune Q(i, ivec) = qentry; 44219653cdaSPeter Brune Q(ivec, i) = qentry; 44319653cdaSPeter Brune } 44419653cdaSPeter Brune } 44519653cdaSPeter Brune 4468409ca45SMatthew G Knepley ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 447087dfb9eSxuemin snes->iter = k; 448f109b39eSPeter Brune snes->norm = fnorm; 4498409ca45SMatthew G Knepley ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 4508409ca45SMatthew G Knepley SNESLogConvHistory(snes, snes->norm, snes->iter); 4518409ca45SMatthew G Knepley ierr = SNESMonitor(snes, snes->iter, snes->norm);CHKERRQ(ierr); 4528409ca45SMatthew G Knepley 453f109b39eSPeter Brune ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 454087dfb9eSxuemin if (snes->reason) PetscFunctionReturn(0); 455a312c225SMatthew G Knepley } 456a312c225SMatthew G Knepley snes->reason = SNES_DIVERGED_MAX_IT; 457a312c225SMatthew G Knepley PetscFunctionReturn(0); 458a312c225SMatthew G Knepley } 459a312c225SMatthew G Knepley 460dfbf837cSBarry Smith /*MC 461dfbf837cSBarry Smith SNESNGMRES - The Nonlinear Generalized Minimum Residual (NGMRES) method of Oosterlee and Washio. 462a312c225SMatthew G Knepley 463dfbf837cSBarry Smith Level: beginner 464dfbf837cSBarry Smith 465dfbf837cSBarry Smith "Krylov Subspace Acceleration of Nonlinear Multigrid with Application to Recirculating Flows", C. W. Oosterlee and T. Washio, 466dfbf837cSBarry Smith SIAM Journal on Scientific Computing, 21(5), 2000. 467dfbf837cSBarry Smith 468dfbf837cSBarry Smith This is also the same as the algorithm called Anderson acceleration introduced in "D. G. Anderson. Iterative procedures for nonlinear integral equations. 469dfbf837cSBarry Smith J. Assoc. Comput. Mach., 12:547–560, 1965." 470dfbf837cSBarry Smith 471dfbf837cSBarry Smith .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types) 472dfbf837cSBarry Smith M*/ 473a312c225SMatthew G Knepley 474a312c225SMatthew G Knepley EXTERN_C_BEGIN 475a312c225SMatthew G Knepley #undef __FUNCT__ 476a312c225SMatthew G Knepley #define __FUNCT__ "SNESCreate_NGMRES" 477a312c225SMatthew G Knepley PetscErrorCode SNESCreate_NGMRES(SNES snes) 478a312c225SMatthew G Knepley { 479a312c225SMatthew G Knepley SNES_NGMRES *ngmres; 480a312c225SMatthew G Knepley PetscErrorCode ierr; 481a312c225SMatthew G Knepley 482a312c225SMatthew G Knepley PetscFunctionBegin; 483a312c225SMatthew G Knepley snes->ops->destroy = SNESDestroy_NGMRES; 484a312c225SMatthew G Knepley snes->ops->setup = SNESSetUp_NGMRES; 485a312c225SMatthew G Knepley snes->ops->setfromoptions = SNESSetFromOptions_NGMRES; 486a312c225SMatthew G Knepley snes->ops->view = SNESView_NGMRES; 487a312c225SMatthew G Knepley snes->ops->solve = SNESSolve_NGMRES; 488a312c225SMatthew G Knepley snes->ops->reset = SNESReset_NGMRES; 489a312c225SMatthew G Knepley 49042f4f86dSBarry Smith snes->usespc = PETSC_TRUE; 4912c155ee1SBarry Smith snes->usesksp = PETSC_FALSE; 4922c155ee1SBarry Smith 493a312c225SMatthew G Knepley ierr = PetscNewLog(snes, SNES_NGMRES, &ngmres);CHKERRQ(ierr); 494a312c225SMatthew G Knepley snes->data = (void*) ngmres; 49519653cdaSPeter Brune ngmres->msize = 10; 49619653cdaSPeter Brune 49728ed4a04SPeter Brune ngmres->restart_it = 2; 498f109b39eSPeter Brune ngmres->gammaA = 2.0; 499f109b39eSPeter Brune ngmres->gammaC = 2.0; 500cac108bcSPeter Brune ngmres->deltaB = 0.9; 501cac108bcSPeter Brune ngmres->epsilonB = 0.1; 502f109b39eSPeter Brune ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetType_C","SNESLineSearchSetType_NGMRES",SNESLineSearchSetType_NGMRES);CHKERRQ(ierr); 503f109b39eSPeter Brune ierr = SNESLineSearchSetType(snes, SNES_LS_QUADRATIC);CHKERRQ(ierr); 504a312c225SMatthew G Knepley PetscFunctionReturn(0); 505a312c225SMatthew G Knepley } 506a312c225SMatthew G Knepley EXTERN_C_END 507