113a62661SPeter Brune #include <../src/snes/impls/ngmres/snesngmres.h> /*I "petscsnes.h" I*/ 219653cdaSPeter Brune #include <petscblaslapack.h> 3a312c225SMatthew G Knepley 46a6fc655SJed Brown const char *const SNESNGMRESRestartTypes[] = {"NONE","PERIODIC","DIFFERENCE","SNESNGMRESRestartType","SNES_NGMRES_RESTART_",0}; 56a6fc655SJed Brown const char *const SNESNGMRESSelectTypes[] = {"NONE","DIFFERENCE","LINESEARCH","SNESNGMRESSelectType","SNES_NGMRES_SELECT_",0}; 613a62661SPeter Brune 7a312c225SMatthew G Knepley #undef __FUNCT__ 8a312c225SMatthew G Knepley #define __FUNCT__ "SNESReset_NGMRES" 9a312c225SMatthew G Knepley PetscErrorCode SNESReset_NGMRES(SNES snes) 10a312c225SMatthew G Knepley { 11a312c225SMatthew G Knepley SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data; 12a312c225SMatthew G Knepley PetscErrorCode ierr; 13a312c225SMatthew G Knepley 14a312c225SMatthew G Knepley PetscFunctionBegin; 15f109b39eSPeter Brune ierr = VecDestroyVecs(ngmres->msize,&ngmres->Fdot);CHKERRQ(ierr); 16f109b39eSPeter Brune ierr = VecDestroyVecs(ngmres->msize,&ngmres->Xdot);CHKERRQ(ierr); 17f1c6b773SPeter Brune ierr = SNESLineSearchDestroy(&ngmres->additive_linesearch);CHKERRQ(ierr); 18a312c225SMatthew G Knepley PetscFunctionReturn(0); 19a312c225SMatthew G Knepley } 20a312c225SMatthew G Knepley 21a312c225SMatthew G Knepley #undef __FUNCT__ 22a312c225SMatthew G Knepley #define __FUNCT__ "SNESDestroy_NGMRES" 23a312c225SMatthew G Knepley PetscErrorCode SNESDestroy_NGMRES(SNES snes) 24a312c225SMatthew G Knepley { 25a312c225SMatthew G Knepley PetscErrorCode ierr; 2678440776SJed Brown SNES_NGMRES *ngmres = (SNES_NGMRES*)snes->data; 27a312c225SMatthew G Knepley 28a312c225SMatthew G Knepley PetscFunctionBegin; 29a312c225SMatthew G Knepley ierr = SNESReset_NGMRES(snes);CHKERRQ(ierr); 30f109b39eSPeter Brune ierr = PetscFree5(ngmres->h,ngmres->beta,ngmres->xi,ngmres->fnorms,ngmres->q);CHKERRQ(ierr); 3119653cdaSPeter Brune ierr = PetscFree(ngmres->s);CHKERRQ(ierr); 3218aa0c0cSPeter Brune ierr = PetscFree(ngmres->xnorms);CHKERRQ(ierr); 3319653cdaSPeter Brune #if PETSC_USE_COMPLEX 3422d28d08SBarry Smith ierr = PetscFree(ngmres->rwork);CHKERRQ(ierr); 3519653cdaSPeter Brune #endif 3622d28d08SBarry Smith ierr = PetscFree(ngmres->work);CHKERRQ(ierr); 3722d28d08SBarry Smith ierr = PetscFree(snes->data);CHKERRQ(ierr); 38a312c225SMatthew G Knepley PetscFunctionReturn(0); 39a312c225SMatthew G Knepley } 40a312c225SMatthew G Knepley 41a312c225SMatthew G Knepley #undef __FUNCT__ 42a312c225SMatthew G Knepley #define __FUNCT__ "SNESSetUp_NGMRES" 43a312c225SMatthew G Knepley PetscErrorCode SNESSetUp_NGMRES(SNES snes) 44a312c225SMatthew G Knepley { 45a312c225SMatthew G Knepley SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data; 46e7058c64SPeter Brune const char *optionsprefix; 4719653cdaSPeter Brune PetscInt msize,hsize; 48a312c225SMatthew G Knepley PetscErrorCode ierr; 49a312c225SMatthew G Knepley 50a312c225SMatthew G Knepley PetscFunctionBegin; 5146159c86SPeter Brune if (snes->pc && snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) { 5246159c86SPeter Brune SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"SNESNGMRES does not support left preconditioning with unpreconditioned function"); 5346159c86SPeter Brune } 546c67d002SPeter Brune if (snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_PRECONDITIONED; 55fa0ddf94SBarry Smith ierr = SNESSetWorkVecs(snes,5);CHKERRQ(ierr); 5678440776SJed Brown if (!ngmres->Xdot) {ierr = VecDuplicateVecs(snes->vec_sol,ngmres->msize,&ngmres->Xdot);CHKERRQ(ierr);} 5778440776SJed Brown if (!ngmres->Fdot) {ierr = VecDuplicateVecs(snes->vec_sol,ngmres->msize,&ngmres->Fdot);CHKERRQ(ierr);} 5878440776SJed Brown if (!ngmres->setup_called) { 59087dfb9eSxuemin msize = ngmres->msize; /* restart size */ 6019653cdaSPeter Brune hsize = msize * msize; 61087dfb9eSxuemin 6298b3e84cSPeter Brune /* explicit least squares minimization solve */ 6319653cdaSPeter Brune ierr = PetscMalloc5(hsize,PetscScalar,&ngmres->h, 6419653cdaSPeter Brune msize,PetscScalar,&ngmres->beta, 6519653cdaSPeter Brune msize,PetscScalar,&ngmres->xi, 66f109b39eSPeter Brune msize,PetscReal, &ngmres->fnorms, 6719653cdaSPeter Brune hsize,PetscScalar,&ngmres->q);CHKERRQ(ierr); 6818aa0c0cSPeter Brune if (ngmres->singlereduction) { 6918aa0c0cSPeter Brune ierr = PetscMalloc(msize*sizeof(PetscReal),&ngmres->xnorms);CHKERRQ(ierr); 7018aa0c0cSPeter Brune } 7119653cdaSPeter Brune ngmres->nrhs = 1; 7219653cdaSPeter Brune ngmres->lda = msize; 7319653cdaSPeter Brune ngmres->ldb = msize; 7419653cdaSPeter Brune ierr = PetscMalloc(msize*sizeof(PetscScalar),&ngmres->s);CHKERRQ(ierr); 7519653cdaSPeter Brune ierr = PetscMemzero(ngmres->h, hsize*sizeof(PetscScalar));CHKERRQ(ierr); 7619653cdaSPeter Brune ierr = PetscMemzero(ngmres->q, hsize*sizeof(PetscScalar));CHKERRQ(ierr); 7719653cdaSPeter Brune ierr = PetscMemzero(ngmres->xi, msize*sizeof(PetscScalar));CHKERRQ(ierr); 7819653cdaSPeter Brune ierr = PetscMemzero(ngmres->beta,msize*sizeof(PetscScalar));CHKERRQ(ierr); 7919653cdaSPeter Brune ngmres->lwork = 12*msize; 8019653cdaSPeter Brune #if PETSC_USE_COMPLEX 8122d28d08SBarry Smith ierr = PetscMalloc(sizeof(PetscReal)*ngmres->lwork,&ngmres->rwork);CHKERRQ(ierr); 8219653cdaSPeter Brune #endif 8322d28d08SBarry Smith ierr = PetscMalloc(sizeof(PetscScalar)*ngmres->lwork,&ngmres->work);CHKERRQ(ierr); 8478440776SJed Brown } 85e7058c64SPeter Brune 86e7058c64SPeter Brune /* linesearch setup */ 87e7058c64SPeter Brune ierr = SNESGetOptionsPrefix(snes,&optionsprefix);CHKERRQ(ierr); 88e7058c64SPeter Brune 8913a62661SPeter Brune if (ngmres->select_type == SNES_NGMRES_SELECT_LINESEARCH) { 90ce94432eSBarry Smith ierr = SNESLineSearchCreate(PetscObjectComm((PetscObject)snes),&ngmres->additive_linesearch);CHKERRQ(ierr); 91f1c6b773SPeter Brune ierr = SNESLineSearchSetSNES(ngmres->additive_linesearch,snes);CHKERRQ(ierr); 921a4f838cSPeter Brune ierr = SNESLineSearchSetType(ngmres->additive_linesearch,SNESLINESEARCHL2);CHKERRQ(ierr); 93f1c6b773SPeter Brune ierr = SNESLineSearchAppendOptionsPrefix(ngmres->additive_linesearch,"additive_");CHKERRQ(ierr); 94f1c6b773SPeter Brune ierr = SNESLineSearchAppendOptionsPrefix(ngmres->additive_linesearch,optionsprefix);CHKERRQ(ierr); 95f1c6b773SPeter Brune ierr = SNESLineSearchSetFromOptions(ngmres->additive_linesearch);CHKERRQ(ierr); 96e7058c64SPeter Brune } 97e7058c64SPeter Brune 9878440776SJed Brown ngmres->setup_called = PETSC_TRUE; 99a312c225SMatthew G Knepley PetscFunctionReturn(0); 100a312c225SMatthew G Knepley } 101a312c225SMatthew G Knepley 102a312c225SMatthew G Knepley #undef __FUNCT__ 103a312c225SMatthew G Knepley #define __FUNCT__ "SNESSetFromOptions_NGMRES" 104a312c225SMatthew G Knepley PetscErrorCode SNESSetFromOptions_NGMRES(SNES snes) 105a312c225SMatthew G Knepley { 106a312c225SMatthew G Knepley SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data; 107a312c225SMatthew G Knepley PetscErrorCode ierr; 108dfbf837cSBarry Smith PetscBool debug; 109f1c6b773SPeter Brune SNESLineSearch linesearch; 1100adebc6cSBarry Smith 111a312c225SMatthew G Knepley PetscFunctionBegin; 112a312c225SMatthew G Knepley ierr = PetscOptionsHead("SNES NGMRES options");CHKERRQ(ierr); 11313a62661SPeter Brune ierr = PetscOptionsEnum("-snes_ngmres_select_type","Select type","SNESNGMRESSetSelectType",SNESNGMRESSelectTypes, 1140298fd71SBarry Smith (PetscEnum)ngmres->select_type,(PetscEnum*)&ngmres->select_type,NULL);CHKERRQ(ierr); 11513a62661SPeter Brune ierr = PetscOptionsEnum("-snes_ngmres_restart_type","Restart type","SNESNGMRESSetRestartType",SNESNGMRESRestartTypes, 1160298fd71SBarry Smith (PetscEnum)ngmres->restart_type,(PetscEnum*)&ngmres->restart_type,NULL);CHKERRQ(ierr); 1170298fd71SBarry Smith ierr = PetscOptionsBool("-snes_ngmres_candidate", "Use candidate storage", "SNES",ngmres->candidate,&ngmres->candidate,NULL);CHKERRQ(ierr); 118077c4231SPeter Brune ierr = PetscOptionsBool("-snes_ngmres_approxfunc","Linearly approximate the function", "SNES",ngmres->approxfunc,&ngmres->approxfunc,NULL);CHKERRQ(ierr); 1190298fd71SBarry Smith ierr = PetscOptionsInt("-snes_ngmres_m", "Number of directions", "SNES",ngmres->msize,&ngmres->msize,NULL);CHKERRQ(ierr); 1200298fd71SBarry Smith ierr = PetscOptionsInt("-snes_ngmres_restart", "Iterations before forced restart", "SNES",ngmres->restart_periodic,&ngmres->restart_periodic,NULL);CHKERRQ(ierr); 1210298fd71SBarry Smith ierr = PetscOptionsInt("-snes_ngmres_restart_it", "Tolerance iterations before restart","SNES",ngmres->restart_it,&ngmres->restart_it,NULL);CHKERRQ(ierr); 1220298fd71SBarry Smith ierr = PetscOptionsBool("-snes_ngmres_monitor", "Monitor actions of NGMRES", "SNES",ngmres->monitor ? PETSC_TRUE : PETSC_FALSE,&debug,NULL);CHKERRQ(ierr); 123dfbf837cSBarry Smith if (debug) { 124ce94432eSBarry Smith ngmres->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));CHKERRQ(ierr); 125dfbf837cSBarry Smith } 1260298fd71SBarry Smith ierr = PetscOptionsReal("-snes_ngmres_gammaA", "Residual selection constant", "SNES",ngmres->gammaA,&ngmres->gammaA,NULL);CHKERRQ(ierr); 1270298fd71SBarry Smith ierr = PetscOptionsReal("-snes_ngmres_gammaC", "Residual restart constant", "SNES",ngmres->gammaC,&ngmres->gammaC,NULL);CHKERRQ(ierr); 1280298fd71SBarry Smith ierr = PetscOptionsReal("-snes_ngmres_epsilonB", "Difference selection constant", "SNES",ngmres->epsilonB,&ngmres->epsilonB,NULL);CHKERRQ(ierr); 1290298fd71SBarry Smith ierr = PetscOptionsReal("-snes_ngmres_deltaB", "Difference residual selection constant", "SNES",ngmres->deltaB,&ngmres->deltaB,NULL);CHKERRQ(ierr); 1300298fd71SBarry Smith ierr = PetscOptionsBool("-snes_ngmres_single_reduction", "Aggregate reductions", "SNES",ngmres->singlereduction,&ngmres->singlereduction,NULL);CHKERRQ(ierr); 131a312c225SMatthew G Knepley ierr = PetscOptionsTail();CHKERRQ(ierr); 1326a7cf640SPeter Brune if ((ngmres->gammaA > ngmres->gammaC) && (ngmres->gammaC > 2.)) ngmres->gammaC = ngmres->gammaA; 1339e764e56SPeter Brune 1349e764e56SPeter Brune /* set the default type of the line search if the user hasn't already. */ 1359e764e56SPeter Brune if (!snes->linesearch) { 1367601faf0SJed Brown ierr = SNESGetLineSearch(snes,&linesearch);CHKERRQ(ierr); 1371a4f838cSPeter Brune ierr = SNESLineSearchSetType(linesearch,SNESLINESEARCHBASIC);CHKERRQ(ierr); 1389e764e56SPeter Brune } 139a312c225SMatthew G Knepley PetscFunctionReturn(0); 140a312c225SMatthew G Knepley } 141a312c225SMatthew G Knepley 142a312c225SMatthew G Knepley #undef __FUNCT__ 143a312c225SMatthew G Knepley #define __FUNCT__ "SNESView_NGMRES" 144a312c225SMatthew G Knepley PetscErrorCode SNESView_NGMRES(SNES snes,PetscViewer viewer) 145a312c225SMatthew G Knepley { 146a312c225SMatthew G Knepley SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data; 147a312c225SMatthew G Knepley PetscBool iascii; 148a312c225SMatthew G Knepley PetscErrorCode ierr; 149a312c225SMatthew G Knepley 150a312c225SMatthew G Knepley PetscFunctionBegin; 151251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject) viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 152a312c225SMatthew G Knepley if (iascii) { 153f109b39eSPeter Brune ierr = PetscViewerASCIIPrintf(viewer," Number of stored past updates: %d\n", ngmres->msize);CHKERRQ(ierr); 154f109b39eSPeter Brune ierr = PetscViewerASCIIPrintf(viewer," Residual selection: gammaA=%1.0e, gammaC=%1.0e\n",ngmres->gammaA,ngmres->gammaC);CHKERRQ(ierr); 155f109b39eSPeter Brune ierr = PetscViewerASCIIPrintf(viewer," Difference restart: epsilonB=%1.0e, deltaB=%1.0e\n",ngmres->epsilonB,ngmres->deltaB);CHKERRQ(ierr); 156a312c225SMatthew G Knepley } 157a312c225SMatthew G Knepley PetscFunctionReturn(0); 158a312c225SMatthew G Knepley } 159a312c225SMatthew G Knepley 160a312c225SMatthew G Knepley #undef __FUNCT__ 161a312c225SMatthew G Knepley #define __FUNCT__ "SNESSolve_NGMRES" 162a312c225SMatthew G Knepley PetscErrorCode SNESSolve_NGMRES(SNES snes) 163a312c225SMatthew G Knepley { 16438774f0aSPeter Brune 165087dfb9eSxuemin SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data; 16698b3e84cSPeter Brune /* present solution, residual, and preconditioned residual */ 1679f425c49SPeter Brune Vec X,F,B,D,Y; 168f109b39eSPeter Brune 169f109b39eSPeter Brune /* candidate linear combination answers */ 170ddd40ce5SPeter Brune Vec XA,FA,XM,FM; 17119653cdaSPeter Brune 17298b3e84cSPeter Brune /* coefficients and RHS to the minimization problem */ 17318aa0c0cSPeter Brune PetscReal fnorm,fMnorm,fAnorm; 17438774f0aSPeter Brune PetscInt k,k_restart,l,ivec,restart_count = 0; 17519653cdaSPeter Brune 17698b3e84cSPeter Brune /* solution selection data */ 17738774f0aSPeter Brune PetscBool selectRestart; 17838774f0aSPeter Brune PetscReal dnorm,dminnorm = 0.0; 179d484d688SPeter Brune PetscReal fminnorm,xnorm,ynorm; 18019653cdaSPeter Brune 1811e633543SBarry Smith SNESConvergedReason reason; 18238774f0aSPeter Brune PetscBool lssucceed; 183a312c225SMatthew G Knepley PetscErrorCode ierr; 184a312c225SMatthew G Knepley 185a312c225SMatthew G Knepley PetscFunctionBegin; 18698b3e84cSPeter Brune /* variable initialization */ 187a312c225SMatthew G Knepley snes->reason = SNES_CONVERGED_ITERATING; 188f109b39eSPeter Brune X = snes->vec_sol; 189f109b39eSPeter Brune F = snes->vec_func; 190f109b39eSPeter Brune B = snes->vec_rhs; 191f109b39eSPeter Brune XA = snes->vec_sol_update; 192f109b39eSPeter Brune FA = snes->work[0]; 193f109b39eSPeter Brune D = snes->work[1]; 194f109b39eSPeter Brune 195f109b39eSPeter Brune /* work for the line search */ 196f109b39eSPeter Brune Y = snes->work[2]; 1979f425c49SPeter Brune XM = snes->work[3]; 1989f425c49SPeter Brune FM = snes->work[4]; 199a312c225SMatthew G Knepley 200*e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 201a312c225SMatthew G Knepley snes->iter = 0; 202a312c225SMatthew G Knepley snes->norm = 0.; 203*e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 20419653cdaSPeter Brune 20598b3e84cSPeter Brune /* initialization */ 20619653cdaSPeter Brune 2073a2ae377SPeter Brune if (snes->pc && snes->pcside == PC_LEFT) { 2083a2ae377SPeter Brune ierr = SNESApplyPC(snes,X,NULL,NULL,F);CHKERRQ(ierr); 2093a2ae377SPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 2103a2ae377SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 2113a2ae377SPeter Brune snes->reason = SNES_DIVERGED_INNER; 2123a2ae377SPeter Brune PetscFunctionReturn(0); 2133a2ae377SPeter Brune } 2143a2ae377SPeter Brune ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); 2153a2ae377SPeter Brune } else { 216e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 217f109b39eSPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 218a312c225SMatthew G Knepley if (snes->domainerror) { 219a312c225SMatthew G Knepley snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 220a312c225SMatthew G Knepley PetscFunctionReturn(0); 221a312c225SMatthew G Knepley } 2221aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 223e4ed7901SPeter Brune if (!snes->norm_init_set) { 2243a2ae377SPeter Brune /* convergence test */ 225f109b39eSPeter Brune ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); 226189a9710SBarry Smith if (PetscIsInfOrNanReal(fnorm)) { 227189a9710SBarry Smith snes->reason = SNES_DIVERGED_FNORM_NAN; 228189a9710SBarry Smith PetscFunctionReturn(0); 229189a9710SBarry Smith } 230e4ed7901SPeter Brune } else { 231e4ed7901SPeter Brune fnorm = snes->norm_init; 232e4ed7901SPeter Brune snes->norm_init_set = PETSC_FALSE; 233e4ed7901SPeter Brune } 2343a2ae377SPeter Brune } 235e4ed7901SPeter Brune fminnorm = fnorm; 23619653cdaSPeter Brune 23798b3e84cSPeter Brune /* q_{00} = nu */ 23838774f0aSPeter Brune SNESNGMRESUpdateSubspace_Private(snes,0,0,F,fnorm,X); 239087dfb9eSxuemin 240*e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 241f109b39eSPeter Brune snes->norm = fnorm; 242*e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 243a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr); 244f109b39eSPeter Brune ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 245f109b39eSPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 246a312c225SMatthew G Knepley if (snes->reason) PetscFunctionReturn(0); 247a312c225SMatthew G Knepley 24819653cdaSPeter Brune k_restart = 1; 24919653cdaSPeter Brune l = 1; 25009c08436SPeter Brune for (k=1; k < snes->max_its+1; k++) { 25198b3e84cSPeter Brune /* select which vector of the stored subspace will be updated */ 25298b3e84cSPeter Brune ivec = k_restart % ngmres->msize; /* replace the last used part of the subspace */ 25319653cdaSPeter Brune 25498b3e84cSPeter Brune /* Computation of x^M */ 255c40d0f55SPeter Brune if (snes->pc && snes->pcside == PC_RIGHT) { 2569f425c49SPeter Brune ierr = VecCopy(X,XM);CHKERRQ(ierr); 257e4ed7901SPeter Brune ierr = SNESSetInitialFunction(snes->pc,F);CHKERRQ(ierr); 258e4ed7901SPeter Brune ierr = SNESSetInitialFunctionNorm(snes->pc,fnorm);CHKERRQ(ierr); 25963e7833aSPeter Brune 26063e7833aSPeter Brune ierr = PetscLogEventBegin(SNES_NPCSolve,snes->pc,XM,B,0);CHKERRQ(ierr); 2619f425c49SPeter Brune ierr = SNESSolve(snes->pc,B,XM);CHKERRQ(ierr); 26263e7833aSPeter Brune ierr = PetscLogEventEnd(SNES_NPCSolve,snes->pc,XM,B,0);CHKERRQ(ierr); 26363e7833aSPeter Brune 2648cc86e31SPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 2658cc86e31SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 2668cc86e31SPeter Brune snes->reason = SNES_DIVERGED_INNER; 2678cc86e31SPeter Brune PetscFunctionReturn(0); 2688cc86e31SPeter Brune } 269ddd40ce5SPeter Brune ierr = SNESGetPCFunction(snes,FM,&fMnorm);CHKERRQ(ierr); 2708cc86e31SPeter Brune } else { 271f109b39eSPeter Brune /* no preconditioner -- just take gradient descent with line search */ 272f109b39eSPeter Brune ierr = VecCopy(F,Y);CHKERRQ(ierr); 273e7058c64SPeter Brune ierr = VecCopy(F,FM);CHKERRQ(ierr); 274e7058c64SPeter Brune ierr = VecCopy(X,XM);CHKERRQ(ierr); 2751aa26658SKarl Rupp 276e7058c64SPeter Brune fMnorm = fnorm; 2771aa26658SKarl Rupp 278f1c6b773SPeter Brune ierr = SNESLineSearchApply(snes->linesearch,XM,FM,&fMnorm,Y);CHKERRQ(ierr); 279f1c6b773SPeter Brune ierr = SNESLineSearchGetSuccess(snes->linesearch,&lssucceed);CHKERRQ(ierr); 280f109b39eSPeter Brune if (!lssucceed) { 281f109b39eSPeter Brune if (++snes->numFailures >= snes->maxFailures) { 282f109b39eSPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 283f109b39eSPeter Brune PetscFunctionReturn(0); 284f109b39eSPeter Brune } 285f109b39eSPeter Brune } 2866634f59bSPeter Brune } 287f6408107SPeter Brune ierr = SNESNGMRESFormCombinedSolution_Private(snes,l,XM,FM,fMnorm,X,XA,FA);CHKERRQ(ierr); 28898b3e84cSPeter Brune /* r = F(x) */ 2899f425c49SPeter Brune if (fminnorm > fMnorm) fminnorm = fMnorm; /* the minimum norm is now of F^M */ 29019653cdaSPeter Brune 2919f425c49SPeter Brune /* differences for selection and restart */ 29213a62661SPeter Brune if (ngmres->restart_type == SNES_NGMRES_RESTART_DIFFERENCE || ngmres->select_type == SNES_NGMRES_SELECT_DIFFERENCE) { 293fa8c639aSPeter Brune ierr = SNESNGMRESCalculateDifferences_Private(snes,l,X,F,XM,FM,XA,FA,D,&dnorm,&dminnorm,&fAnorm);CHKERRQ(ierr); 29413a62661SPeter Brune } else { 29513a62661SPeter Brune ierr = VecNorm(FA,NORM_2,&fAnorm);CHKERRQ(ierr); 29613a62661SPeter Brune } 297189a9710SBarry Smith if (PetscIsInfOrNanReal(fAnorm)) { 298189a9710SBarry Smith snes->reason = SNES_DIVERGED_FNORM_NAN; 299189a9710SBarry Smith PetscFunctionReturn(0); 300189a9710SBarry Smith } 3011aa26658SKarl Rupp 3029f425c49SPeter Brune /* combination (additive) or selection (multiplicative) of the N-GMRES solution */ 303fa8c639aSPeter Brune ierr = SNESNGMRESSelect_Private(snes,k_restart,XM,FM,fMnorm,XA,FA,fAnorm,dnorm,fminnorm,dminnorm,X,F,Y,&fnorm);CHKERRQ(ierr); 30419653cdaSPeter Brune selectRestart = PETSC_FALSE; 30513a62661SPeter Brune if (ngmres->restart_type == SNES_NGMRES_RESTART_DIFFERENCE) { 30621687c63SPeter Brune ierr = SNESNGMRESSelectRestart_Private(snes,l,fAnorm,dnorm,fminnorm,dminnorm,&selectRestart);CHKERRQ(ierr); 30728ed4a04SPeter Brune /* if the restart conditions persist for more than restart_it iterations, restart. */ 3081aa26658SKarl Rupp if (selectRestart) restart_count++; 3091aa26658SKarl Rupp else restart_count = 0; 31013a62661SPeter Brune } else if (ngmres->restart_type == SNES_NGMRES_RESTART_PERIODIC) { 31113a62661SPeter Brune if (k_restart > ngmres->restart_periodic) { 31213a62661SPeter Brune if (ngmres->monitor) ierr = PetscViewerASCIIPrintf(ngmres->monitor,"periodic restart after %D iterations\n",k_restart);CHKERRQ(ierr); 31313a62661SPeter Brune restart_count = ngmres->restart_it; 31413a62661SPeter Brune } 31513a62661SPeter Brune } 31628ed4a04SPeter Brune /* restart after restart conditions have persisted for a fixed number of iterations */ 31728ed4a04SPeter Brune if (restart_count >= ngmres->restart_it) { 318dfbf837cSBarry Smith if (ngmres->monitor) { 319dfbf837cSBarry Smith ierr = PetscViewerASCIIPrintf(ngmres->monitor,"Restarted at iteration %d\n",k_restart);CHKERRQ(ierr); 320dfbf837cSBarry Smith } 32128ed4a04SPeter Brune restart_count = 0; 32219653cdaSPeter Brune k_restart = 1; 32319653cdaSPeter Brune l = 1; 32498b3e84cSPeter Brune /* q_{00} = nu */ 32538774f0aSPeter Brune if (ngmres->candidate) { 326fa8c639aSPeter Brune ierr = SNESNGMRESUpdateSubspace_Private(snes,0,0,FM,fMnorm,XM);CHKERRQ(ierr); 327d2e16ddcSPeter Brune } else { 328fa8c639aSPeter Brune ierr = SNESNGMRESUpdateSubspace_Private(snes,0,0,F,fMnorm,X);CHKERRQ(ierr); 329d2e16ddcSPeter Brune } 33019653cdaSPeter Brune } else { 33198b3e84cSPeter Brune /* select the current size of the subspace */ 3321e633543SBarry Smith if (l < ngmres->msize) l++; 33319653cdaSPeter Brune k_restart++; 33498b3e84cSPeter Brune /* place the current entry in the list of previous entries */ 33538774f0aSPeter Brune if (ngmres->candidate) { 33638774f0aSPeter Brune if (fminnorm > fMnorm) fminnorm = fMnorm; 337fa8c639aSPeter Brune ierr = SNESNGMRESUpdateSubspace_Private(snes,ivec,l,FM,fMnorm,XM);CHKERRQ(ierr); 338d2e16ddcSPeter Brune } else { 33938774f0aSPeter Brune if (fminnorm > fnorm) fminnorm = fnorm; 340fa8c639aSPeter Brune ierr = SNESNGMRESUpdateSubspace_Private(snes,ivec,l,F,fnorm,X);CHKERRQ(ierr); 34119653cdaSPeter Brune } 342d2e16ddcSPeter Brune } 34319653cdaSPeter Brune 344*e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 345087dfb9eSxuemin snes->iter = k; 346f109b39eSPeter Brune snes->norm = fnorm; 347*e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 348a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,snes->norm,snes->iter);CHKERRQ(ierr); 3498409ca45SMatthew G Knepley ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 350d484d688SPeter Brune ierr = VecNormBegin(Y,NORM_2,&ynorm);CHKERRQ(ierr); 351d484d688SPeter Brune ierr = VecNormBegin(X,NORM_2,&xnorm);CHKERRQ(ierr); 352d484d688SPeter Brune ierr = VecNormEnd(Y,NORM_2,&ynorm);CHKERRQ(ierr); 353d484d688SPeter Brune ierr = VecNormEnd(X,NORM_2,&xnorm);CHKERRQ(ierr); 354d484d688SPeter Brune ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 355087dfb9eSxuemin if (snes->reason) PetscFunctionReturn(0); 356a312c225SMatthew G Knepley } 357a312c225SMatthew G Knepley snes->reason = SNES_DIVERGED_MAX_IT; 358a312c225SMatthew G Knepley PetscFunctionReturn(0); 359a312c225SMatthew G Knepley } 360a312c225SMatthew G Knepley 36113a62661SPeter Brune #undef __FUNCT__ 36213a62661SPeter Brune #define __FUNCT__ "SNESNGMRESSetRestartType" 36313a62661SPeter Brune /*@ 36413a62661SPeter Brune SNESNGMRESSetRestartType - Sets the restart type for SNESNGMRES. 36513a62661SPeter Brune 36613a62661SPeter Brune Logically Collective on SNES 36713a62661SPeter Brune 36813a62661SPeter Brune Input Parameters: 36913a62661SPeter Brune + snes - the iterative context 37013a62661SPeter Brune - rtype - restart type 37113a62661SPeter Brune 37213a62661SPeter Brune Options Database: 37313a62661SPeter Brune + -snes_ngmres_restart_type<difference,periodic,none> - set the restart type 3740c777b0cSPeter Brune - -snes_ngmres_restart[30] - sets the number of iterations before restart for periodic 37513a62661SPeter Brune 37613a62661SPeter Brune Level: intermediate 37713a62661SPeter Brune 37813a62661SPeter Brune SNESNGMRESRestartTypes: 37913a62661SPeter Brune + SNES_NGMRES_RESTART_NONE - never restart 38013a62661SPeter Brune . SNES_NGMRES_RESTART_DIFFERENCE - restart based upon difference criteria 38113a62661SPeter Brune - SNES_NGMRES_RESTART_PERIODIC - restart after a fixed number of iterations 38213a62661SPeter Brune 38313a62661SPeter Brune Notes: 38413a62661SPeter Brune The default line search used is the L2 line search and it requires two additional function evaluations. 38513a62661SPeter Brune 38613a62661SPeter Brune .keywords: SNES, SNESNGMRES, restart, type, set SNESLineSearch 38713a62661SPeter Brune @*/ 3880adebc6cSBarry Smith PetscErrorCode SNESNGMRESSetRestartType(SNES snes,SNESNGMRESRestartType rtype) 3890adebc6cSBarry Smith { 39013a62661SPeter Brune PetscErrorCode ierr; 3915fd66863SKarl Rupp 39213a62661SPeter Brune PetscFunctionBegin; 39313a62661SPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 39413a62661SPeter Brune ierr = PetscTryMethod(snes,"SNESNGMRESSetRestartType_C",(SNES,SNESNGMRESRestartType),(snes,rtype));CHKERRQ(ierr); 39513a62661SPeter Brune PetscFunctionReturn(0); 39613a62661SPeter Brune } 39713a62661SPeter Brune 39813a62661SPeter Brune #undef __FUNCT__ 39913a62661SPeter Brune #define __FUNCT__ "SNESNGMRESSetSelectType" 40013a62661SPeter Brune /*@ 40113a62661SPeter Brune SNESNGMRESSetSelectType - Sets the selection type for SNESNGMRES. This determines how the candidate solution and 40213a62661SPeter Brune combined solution are used to create the next iterate. 40313a62661SPeter Brune 40413a62661SPeter Brune Logically Collective on SNES 40513a62661SPeter Brune 40613a62661SPeter Brune Input Parameters: 40713a62661SPeter Brune + snes - the iterative context 40813a62661SPeter Brune - stype - selection type 40913a62661SPeter Brune 41013a62661SPeter Brune Options Database: 41113a62661SPeter Brune . -snes_ngmres_select_type<difference,none,linesearch> 41213a62661SPeter Brune 41313a62661SPeter Brune Level: intermediate 41413a62661SPeter Brune 41513a62661SPeter Brune SNESNGMRESSelectTypes: 41613a62661SPeter Brune + SNES_NGMRES_SELECT_NONE - choose the combined solution all the time 41713a62661SPeter Brune . SNES_NGMRES_SELECT_DIFFERENCE - choose based upon the selection criteria 41813a62661SPeter Brune - SNES_NGMRES_SELECT_LINESEARCH - choose based upon line search combination 41913a62661SPeter Brune 42013a62661SPeter Brune Notes: 42113a62661SPeter Brune The default line search used is the L2 line search and it requires two additional function evaluations. 42213a62661SPeter Brune 42313a62661SPeter Brune .keywords: SNES, SNESNGMRES, selection, type, set SNESLineSearch 42413a62661SPeter Brune @*/ 4250adebc6cSBarry Smith PetscErrorCode SNESNGMRESSetSelectType(SNES snes,SNESNGMRESSelectType stype) 4260adebc6cSBarry Smith { 42713a62661SPeter Brune PetscErrorCode ierr; 4285fd66863SKarl Rupp 42913a62661SPeter Brune PetscFunctionBegin; 43013a62661SPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 43113a62661SPeter Brune ierr = PetscTryMethod(snes,"SNESNGMRESSetSelectType_C",(SNES,SNESNGMRESSelectType),(snes,stype));CHKERRQ(ierr); 43213a62661SPeter Brune PetscFunctionReturn(0); 43313a62661SPeter Brune } 43413a62661SPeter Brune 43513a62661SPeter Brune #undef __FUNCT__ 43613a62661SPeter Brune #define __FUNCT__ "SNESNGMRESSetSelectType_NGMRES" 4370adebc6cSBarry Smith PetscErrorCode SNESNGMRESSetSelectType_NGMRES(SNES snes,SNESNGMRESSelectType stype) 4380adebc6cSBarry Smith { 43913a62661SPeter Brune SNES_NGMRES *ngmres = (SNES_NGMRES*)snes->data; 4405fd66863SKarl Rupp 44113a62661SPeter Brune PetscFunctionBegin; 44213a62661SPeter Brune ngmres->select_type = stype; 44313a62661SPeter Brune PetscFunctionReturn(0); 44413a62661SPeter Brune } 44513a62661SPeter Brune 44613a62661SPeter Brune #undef __FUNCT__ 44713a62661SPeter Brune #define __FUNCT__ "SNESNGMRESSetRestartType_NGMRES" 4480adebc6cSBarry Smith PetscErrorCode SNESNGMRESSetRestartType_NGMRES(SNES snes,SNESNGMRESRestartType rtype) 4490adebc6cSBarry Smith { 45013a62661SPeter Brune SNES_NGMRES *ngmres = (SNES_NGMRES*)snes->data; 4515fd66863SKarl Rupp 45213a62661SPeter Brune PetscFunctionBegin; 45313a62661SPeter Brune ngmres->restart_type = rtype; 45413a62661SPeter Brune PetscFunctionReturn(0); 45513a62661SPeter Brune } 45613a62661SPeter Brune 457dfbf837cSBarry Smith /*MC 4581867fe5bSPeter Brune SNESNGMRES - The Nonlinear Generalized Minimum Residual method. 459a312c225SMatthew G Knepley 460dfbf837cSBarry Smith Level: beginner 461dfbf837cSBarry Smith 4621867fe5bSPeter Brune Options Database: 46313a62661SPeter Brune + -snes_ngmres_select_type<difference,none,linesearch> - choose the select between candidate and combined solution 46438774f0aSPeter Brune . -snes_ngmres_restart_type<difference,none,periodic> - choose the restart conditions 46538774f0aSPeter Brune . -snes_ngmres_candidate - Use NGMRES variant which combines candidate solutions instead of actual solutions 46613a62661SPeter Brune . -snes_ngmres_m - Number of stored previous solutions and residuals 46713a62661SPeter Brune . -snes_ngmres_restart_it - Number of iterations the restart conditions hold before restart 46813a62661SPeter Brune . -snes_ngmres_gammaA - Residual tolerance for solution select between the candidate and combination 46913a62661SPeter Brune . -snes_ngmres_gammaC - Residual tolerance for restart 47013a62661SPeter Brune . -snes_ngmres_epsilonB - Difference tolerance between subsequent solutions triggering restart 47113a62661SPeter Brune . -snes_ngmres_deltaB - Difference tolerance between residuals triggering restart 47213a62661SPeter Brune . -snes_ngmres_monitor - Prints relevant information about the ngmres iteration 4735c3e6ab7SPeter Brune . -snes_linesearch_type <basic,l2,cp> - Line search type used for the default smoother 47413a62661SPeter Brune - -additive_snes_linesearch_type - linesearch type used to select between the candidate and combined solution with additive select type 4751867fe5bSPeter Brune 4761867fe5bSPeter Brune Notes: 4771867fe5bSPeter Brune 4781867fe5bSPeter Brune The N-GMRES method combines m previous solutions into a minimum-residual solution by solving a small linearized 4791867fe5bSPeter Brune optimization problem at each iteration. 4801867fe5bSPeter Brune 4811867fe5bSPeter Brune References: 4821867fe5bSPeter Brune 483dfbf837cSBarry Smith "Krylov Subspace Acceleration of Nonlinear Multigrid with Application to Recirculating Flows", C. W. Oosterlee and T. Washio, 484dfbf837cSBarry Smith SIAM Journal on Scientific Computing, 21(5), 2000. 485dfbf837cSBarry Smith 486dfbf837cSBarry Smith .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types) 487dfbf837cSBarry Smith M*/ 488a312c225SMatthew G Knepley 489a312c225SMatthew G Knepley #undef __FUNCT__ 490a312c225SMatthew G Knepley #define __FUNCT__ "SNESCreate_NGMRES" 4918cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_NGMRES(SNES snes) 492a312c225SMatthew G Knepley { 493a312c225SMatthew G Knepley SNES_NGMRES *ngmres; 494a312c225SMatthew G Knepley PetscErrorCode ierr; 495a312c225SMatthew G Knepley 496a312c225SMatthew G Knepley PetscFunctionBegin; 497a312c225SMatthew G Knepley snes->ops->destroy = SNESDestroy_NGMRES; 498a312c225SMatthew G Knepley snes->ops->setup = SNESSetUp_NGMRES; 499a312c225SMatthew G Knepley snes->ops->setfromoptions = SNESSetFromOptions_NGMRES; 500a312c225SMatthew G Knepley snes->ops->view = SNESView_NGMRES; 501a312c225SMatthew G Knepley snes->ops->solve = SNESSolve_NGMRES; 502a312c225SMatthew G Knepley snes->ops->reset = SNESReset_NGMRES; 503a312c225SMatthew G Knepley 50442f4f86dSBarry Smith snes->usespc = PETSC_TRUE; 5052c155ee1SBarry Smith snes->usesksp = PETSC_FALSE; 50646159c86SPeter Brune snes->pcside = PC_RIGHT; 5072c155ee1SBarry Smith 508a312c225SMatthew G Knepley ierr = PetscNewLog(snes,SNES_NGMRES,&ngmres);CHKERRQ(ierr); 509a312c225SMatthew G Knepley snes->data = (void*) ngmres; 510d2e16ddcSPeter Brune ngmres->msize = 30; 51119653cdaSPeter Brune 51288976e71SPeter Brune if (!snes->tolerancesset) { 5130e444f03SPeter Brune snes->max_funcs = 30000; 5140e444f03SPeter Brune snes->max_its = 10000; 51588976e71SPeter Brune } 5160e444f03SPeter Brune 51738774f0aSPeter Brune ngmres->candidate = PETSC_FALSE; 518d2e16ddcSPeter Brune 5190298fd71SBarry Smith ngmres->additive_linesearch = NULL; 520077c4231SPeter Brune ngmres->approxfunc = PETSC_FALSE; 52128ed4a04SPeter Brune ngmres->restart_it = 2; 52213a62661SPeter Brune ngmres->restart_periodic = 30; 523f109b39eSPeter Brune ngmres->gammaA = 2.0; 524f109b39eSPeter Brune ngmres->gammaC = 2.0; 525cac108bcSPeter Brune ngmres->deltaB = 0.9; 526cac108bcSPeter Brune ngmres->epsilonB = 0.1; 527e7058c64SPeter Brune 52813a62661SPeter Brune ngmres->restart_type = SNES_NGMRES_RESTART_DIFFERENCE; 52913a62661SPeter Brune ngmres->select_type = SNES_NGMRES_SELECT_DIFFERENCE; 53013a62661SPeter Brune 531bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNGMRESSetSelectType_C",SNESNGMRESSetSelectType_NGMRES);CHKERRQ(ierr); 532bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNGMRESSetRestartType_C",SNESNGMRESSetRestartType_NGMRES);CHKERRQ(ierr); 533a312c225SMatthew G Knepley PetscFunctionReturn(0); 534a312c225SMatthew G Knepley } 53599e0435eSBarry Smith 536