Actual source code: snesngmres.c
petsc-3.6.4 2016-04-12
1: #include <../src/snes/impls/ngmres/snesngmres.h> /*I "petscsnes.h" I*/
2: #include <petscblaslapack.h>
3: #include <petscdm.h>
5: const char *const SNESNGMRESRestartTypes[] = {"NONE","PERIODIC","DIFFERENCE","SNESNGMRESRestartType","SNES_NGMRES_RESTART_",0};
6: const char *const SNESNGMRESSelectTypes[] = {"NONE","DIFFERENCE","LINESEARCH","SNESNGMRESSelectType","SNES_NGMRES_SELECT_",0};
10: PetscErrorCode SNESReset_NGMRES(SNES snes)
11: {
12: SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data;
16: VecDestroyVecs(ngmres->msize,&ngmres->Fdot);
17: VecDestroyVecs(ngmres->msize,&ngmres->Xdot);
18: SNESLineSearchDestroy(&ngmres->additive_linesearch);
19: return(0);
20: }
24: PetscErrorCode SNESDestroy_NGMRES(SNES snes)
25: {
27: SNES_NGMRES *ngmres = (SNES_NGMRES*)snes->data;
30: SNESReset_NGMRES(snes);
31: PetscFree5(ngmres->h,ngmres->beta,ngmres->xi,ngmres->fnorms,ngmres->q);
32: PetscFree(ngmres->s);
33: PetscFree(ngmres->xnorms);
34: #if PETSC_USE_COMPLEX
35: PetscFree(ngmres->rwork);
36: #endif
37: PetscFree(ngmres->work);
38: PetscFree(snes->data);
39: return(0);
40: }
44: PetscErrorCode SNESSetUp_NGMRES(SNES snes)
45: {
46: SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data;
47: const char *optionsprefix;
48: PetscInt msize,hsize;
50: DM dm;
53: if (snes->pc && snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
54: SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"SNESNGMRES does not support left preconditioning with unpreconditioned function");
55: }
56: if (snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_PRECONDITIONED;
57: SNESSetWorkVecs(snes,5);
59: if (!snes->vec_sol) {
60: SNESGetDM(snes,&dm);
61: DMCreateGlobalVector(dm,&snes->vec_sol);
62: }
64: if (!ngmres->Xdot) {VecDuplicateVecs(snes->vec_sol,ngmres->msize,&ngmres->Xdot);}
65: if (!ngmres->Fdot) {VecDuplicateVecs(snes->vec_sol,ngmres->msize,&ngmres->Fdot);}
66: if (!ngmres->setup_called) {
67: msize = ngmres->msize; /* restart size */
68: hsize = msize * msize;
70: /* explicit least squares minimization solve */
71: PetscMalloc5(hsize,&ngmres->h, msize,&ngmres->beta, msize,&ngmres->xi, msize,&ngmres->fnorms, hsize,&ngmres->q);
72: PetscMalloc1(msize,&ngmres->xnorms);
73: ngmres->nrhs = 1;
74: ngmres->lda = msize;
75: ngmres->ldb = msize;
76: PetscMalloc1(msize,&ngmres->s);
77: PetscMemzero(ngmres->h, hsize*sizeof(PetscScalar));
78: PetscMemzero(ngmres->q, hsize*sizeof(PetscScalar));
79: PetscMemzero(ngmres->xi, msize*sizeof(PetscScalar));
80: PetscMemzero(ngmres->beta,msize*sizeof(PetscScalar));
81: ngmres->lwork = 12*msize;
82: #if PETSC_USE_COMPLEX
83: PetscMalloc1(ngmres->lwork,&ngmres->rwork);
84: #endif
85: PetscMalloc1(ngmres->lwork,&ngmres->work);
86: }
88: /* linesearch setup */
89: SNESGetOptionsPrefix(snes,&optionsprefix);
91: if (ngmres->select_type == SNES_NGMRES_SELECT_LINESEARCH) {
92: SNESLineSearchCreate(PetscObjectComm((PetscObject)snes),&ngmres->additive_linesearch);
93: SNESLineSearchSetSNES(ngmres->additive_linesearch,snes);
94: SNESLineSearchSetType(ngmres->additive_linesearch,SNESLINESEARCHL2);
95: SNESLineSearchAppendOptionsPrefix(ngmres->additive_linesearch,"additive_");
96: SNESLineSearchAppendOptionsPrefix(ngmres->additive_linesearch,optionsprefix);
97: SNESLineSearchSetFromOptions(ngmres->additive_linesearch);
98: }
100: ngmres->setup_called = PETSC_TRUE;
101: return(0);
102: }
106: PetscErrorCode SNESSetFromOptions_NGMRES(PetscOptions *PetscOptionsObject,SNES snes)
107: {
108: SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data;
110: PetscBool debug = PETSC_FALSE;
111: SNESLineSearch linesearch;
114: PetscOptionsHead(PetscOptionsObject,"SNES NGMRES options");
115: PetscOptionsEnum("-snes_ngmres_select_type","Select type","SNESNGMRESSetSelectType",SNESNGMRESSelectTypes,
116: (PetscEnum)ngmres->select_type,(PetscEnum*)&ngmres->select_type,NULL);
117: PetscOptionsEnum("-snes_ngmres_restart_type","Restart type","SNESNGMRESSetRestartType",SNESNGMRESRestartTypes,
118: (PetscEnum)ngmres->restart_type,(PetscEnum*)&ngmres->restart_type,NULL);
119: PetscOptionsBool("-snes_ngmres_candidate", "Use candidate storage", "SNES",ngmres->candidate,&ngmres->candidate,NULL);
120: PetscOptionsBool("-snes_ngmres_approxfunc","Linearly approximate the function", "SNES",ngmres->approxfunc,&ngmres->approxfunc,NULL);
121: PetscOptionsInt("-snes_ngmres_m", "Number of directions", "SNES",ngmres->msize,&ngmres->msize,NULL);
122: PetscOptionsInt("-snes_ngmres_restart", "Iterations before forced restart", "SNES",ngmres->restart_periodic,&ngmres->restart_periodic,NULL);
123: PetscOptionsInt("-snes_ngmres_restart_it", "Tolerance iterations before restart","SNES",ngmres->restart_it,&ngmres->restart_it,NULL);
124: PetscOptionsBool("-snes_ngmres_monitor", "Monitor actions of NGMRES", "SNES",ngmres->monitor ? PETSC_TRUE : PETSC_FALSE,&debug,NULL);
125: if (debug) {
126: ngmres->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));
127: }
128: PetscOptionsReal("-snes_ngmres_gammaA", "Residual selection constant", "SNES",ngmres->gammaA,&ngmres->gammaA,NULL);
129: PetscOptionsReal("-snes_ngmres_gammaC", "Residual restart constant", "SNES",ngmres->gammaC,&ngmres->gammaC,NULL);
130: PetscOptionsReal("-snes_ngmres_epsilonB", "Difference selection constant", "SNES",ngmres->epsilonB,&ngmres->epsilonB,NULL);
131: PetscOptionsReal("-snes_ngmres_deltaB", "Difference residual selection constant", "SNES",ngmres->deltaB,&ngmres->deltaB,NULL);
132: PetscOptionsBool("-snes_ngmres_single_reduction", "Aggregate reductions", "SNES",ngmres->singlereduction,&ngmres->singlereduction,NULL);
133: PetscOptionsTail();
134: if ((ngmres->gammaA > ngmres->gammaC) && (ngmres->gammaC > 2.)) ngmres->gammaC = ngmres->gammaA;
136: /* set the default type of the line search if the user hasn't already. */
137: if (!snes->linesearch) {
138: SNESGetLineSearch(snes,&linesearch);
139: SNESLineSearchSetType(linesearch,SNESLINESEARCHBASIC);
140: }
141: return(0);
142: }
146: PetscErrorCode SNESView_NGMRES(SNES snes,PetscViewer viewer)
147: {
148: SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data;
149: PetscBool iascii;
153: PetscObjectTypeCompare((PetscObject) viewer,PETSCVIEWERASCII,&iascii);
154: if (iascii) {
155: PetscViewerASCIIPrintf(viewer," Number of stored past updates: %d\n", ngmres->msize);
156: PetscViewerASCIIPrintf(viewer," Residual selection: gammaA=%1.0e, gammaC=%1.0e\n",ngmres->gammaA,ngmres->gammaC);
157: PetscViewerASCIIPrintf(viewer," Difference restart: epsilonB=%1.0e, deltaB=%1.0e\n",ngmres->epsilonB,ngmres->deltaB);
158: }
159: return(0);
160: }
164: PetscErrorCode SNESSolve_NGMRES(SNES snes)
165: {
166: SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data;
167: /* present solution, residual, and preconditioned residual */
168: Vec X,F,B,D,Y;
170: /* candidate linear combination answers */
171: Vec XA,FA,XM,FM;
173: /* coefficients and RHS to the minimization problem */
174: PetscReal fnorm,fMnorm,fAnorm;
175: PetscReal xnorm,xMnorm,xAnorm;
176: PetscReal ynorm,yMnorm,yAnorm;
177: PetscInt k,k_restart,l,ivec,restart_count = 0;
179: /* solution selection data */
180: PetscBool selectRestart;
181: PetscReal dnorm,dminnorm = 0.0;
182: PetscReal fminnorm;
184: SNESConvergedReason reason;
185: SNESLineSearchReason lssucceed;
186: PetscErrorCode ierr;
190: if (snes->xl || snes->xu || snes->ops->computevariablebounds) {
191: SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name);
192: }
194: PetscCitationsRegister(SNESCitation,&SNEScite);
195: /* variable initialization */
196: snes->reason = SNES_CONVERGED_ITERATING;
197: X = snes->vec_sol;
198: F = snes->vec_func;
199: B = snes->vec_rhs;
200: XA = snes->vec_sol_update;
201: FA = snes->work[0];
202: D = snes->work[1];
204: /* work for the line search */
205: Y = snes->work[2];
206: XM = snes->work[3];
207: FM = snes->work[4];
209: PetscObjectSAWsTakeAccess((PetscObject)snes);
210: snes->iter = 0;
211: snes->norm = 0.;
212: PetscObjectSAWsGrantAccess((PetscObject)snes);
214: /* initialization */
216: if (snes->pc && snes->pcside == PC_LEFT) {
217: SNESApplyNPC(snes,X,NULL,F);
218: SNESGetConvergedReason(snes->pc,&reason);
219: if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
220: snes->reason = SNES_DIVERGED_INNER;
221: return(0);
222: }
223: VecNorm(F,NORM_2,&fnorm);
224: } else {
225: if (!snes->vec_func_init_set) {
226: SNESComputeFunction(snes,X,F);
227: } else snes->vec_func_init_set = PETSC_FALSE;
229: VecNorm(F,NORM_2,&fnorm);
230: SNESCheckFunctionNorm(snes,fnorm);
231: }
232: fminnorm = fnorm;
234: PetscObjectSAWsTakeAccess((PetscObject)snes);
235: snes->norm = fnorm;
236: PetscObjectSAWsGrantAccess((PetscObject)snes);
237: SNESLogConvergenceHistory(snes,fnorm,0);
238: SNESMonitor(snes,0,fnorm);
239: (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);
240: if (snes->reason) return(0);
241: SNESNGMRESUpdateSubspace_Private(snes,0,0,F,fnorm,X);
243: k_restart = 1;
244: l = 1;
245: ivec = 0;
246: for (k=1; k < snes->max_its+1; k++) {
247: /* Computation of x^M */
248: if (snes->pc && snes->pcside == PC_RIGHT) {
249: VecCopy(X,XM);
250: SNESSetInitialFunction(snes->pc,F);
252: PetscLogEventBegin(SNES_NPCSolve,snes->pc,XM,B,0);
253: SNESSolve(snes->pc,B,XM);
254: PetscLogEventEnd(SNES_NPCSolve,snes->pc,XM,B,0);
256: SNESGetConvergedReason(snes->pc,&reason);
257: if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
258: snes->reason = SNES_DIVERGED_INNER;
259: return(0);
260: }
261: SNESGetNPCFunction(snes,FM,&fMnorm);
262: } else {
263: /* no preconditioner -- just take gradient descent with line search */
264: VecCopy(F,Y);
265: VecCopy(F,FM);
266: VecCopy(X,XM);
268: fMnorm = fnorm;
270: SNESLineSearchApply(snes->linesearch,XM,FM,&fMnorm,Y);
271: SNESLineSearchGetReason(snes->linesearch,&lssucceed);
272: if (lssucceed) {
273: if (++snes->numFailures >= snes->maxFailures) {
274: snes->reason = SNES_DIVERGED_LINE_SEARCH;
275: return(0);
276: }
277: }
278: }
279: SNESNGMRESFormCombinedSolution_Private(snes,ivec,l,XM,FM,fMnorm,X,XA,FA);
280: /* r = F(x) */
281: if (fminnorm > fMnorm) fminnorm = fMnorm; /* the minimum norm is now of F^M */
283: /* differences for selection and restart */
284: if (ngmres->restart_type == SNES_NGMRES_RESTART_DIFFERENCE || ngmres->select_type == SNES_NGMRES_SELECT_DIFFERENCE) {
285: SNESNGMRESNorms_Private(snes,l,X,F,XM,FM,XA,FA,D,&dnorm,&dminnorm,&xMnorm,NULL,&yMnorm,&xAnorm,&fAnorm,&yAnorm);
286: } else {
287: SNESNGMRESNorms_Private(snes,l,X,F,XM,FM,XA,FA,D,NULL,NULL,&xMnorm,NULL,&yMnorm,&xAnorm,&fAnorm,&yAnorm);
288: }
289: SNESCheckFunctionNorm(snes,fnorm);
291: /* combination (additive) or selection (multiplicative) of the N-GMRES solution */
292: SNESNGMRESSelect_Private(snes,k_restart,XM,FM,xMnorm,fMnorm,yMnorm,XA,FA,xAnorm,fAnorm,yAnorm,dnorm,fminnorm,dminnorm,X,F,Y,&xnorm,&fnorm,&ynorm);
293: selectRestart = PETSC_FALSE;
294: if (ngmres->restart_type == SNES_NGMRES_RESTART_DIFFERENCE) {
295: SNESNGMRESSelectRestart_Private(snes,l,fAnorm,dnorm,fminnorm,dminnorm,&selectRestart);
296: /* if the restart conditions persist for more than restart_it iterations, restart. */
297: if (selectRestart) restart_count++;
298: else restart_count = 0;
299: } else if (ngmres->restart_type == SNES_NGMRES_RESTART_PERIODIC) {
300: if (k_restart > ngmres->restart_periodic) {
301: if (ngmres->monitor) PetscViewerASCIIPrintf(ngmres->monitor,"periodic restart after %D iterations\n",k_restart);
302: restart_count = ngmres->restart_it;
303: }
304: }
305: ivec = k_restart % ngmres->msize; /* replace the last used part of the subspace */
306: /* restart after restart conditions have persisted for a fixed number of iterations */
307: if (restart_count >= ngmres->restart_it) {
308: if (ngmres->monitor) {
309: PetscViewerASCIIPrintf(ngmres->monitor,"Restarted at iteration %d\n",k_restart);
310: }
311: restart_count = 0;
312: k_restart = 1;
313: l = 1;
314: ivec = 0;
315: /* q_{00} = nu */
316: SNESNGMRESUpdateSubspace_Private(snes,0,0,FM,fMnorm,XM);
317: } else {
318: /* select the current size of the subspace */
319: if (l < ngmres->msize) l++;
320: k_restart++;
321: /* place the current entry in the list of previous entries */
322: if (ngmres->candidate) {
323: if (fminnorm > fMnorm) fminnorm = fMnorm;
324: SNESNGMRESUpdateSubspace_Private(snes,ivec,l,FM,fMnorm,XM);
325: } else {
326: if (fminnorm > fnorm) fminnorm = fnorm;
327: SNESNGMRESUpdateSubspace_Private(snes,ivec,l,F,fnorm,X);
328: }
329: }
331: PetscObjectSAWsTakeAccess((PetscObject)snes);
332: snes->iter = k;
333: snes->norm = fnorm;
334: PetscObjectSAWsGrantAccess((PetscObject)snes);
335: SNESLogConvergenceHistory(snes,snes->norm,snes->iter);
336: SNESMonitor(snes,snes->iter,snes->norm);
337: (*snes->ops->converged)(snes,snes->iter,0,0,fnorm,&snes->reason,snes->cnvP);
338: if (snes->reason) return(0);
339: }
340: snes->reason = SNES_DIVERGED_MAX_IT;
341: return(0);
342: }
346: /*@
347: SNESNGMRESSetRestartType - Sets the restart type for SNESNGMRES.
349: Logically Collective on SNES
351: Input Parameters:
352: + snes - the iterative context
353: - rtype - restart type
355: Options Database:
356: + -snes_ngmres_restart_type<difference,periodic,none> - set the restart type
357: - -snes_ngmres_restart[30] - sets the number of iterations before restart for periodic
359: Level: intermediate
361: SNESNGMRESRestartTypes:
362: + SNES_NGMRES_RESTART_NONE - never restart
363: . SNES_NGMRES_RESTART_DIFFERENCE - restart based upon difference criteria
364: - SNES_NGMRES_RESTART_PERIODIC - restart after a fixed number of iterations
366: Notes:
367: The default line search used is the L2 line search and it requires two additional function evaluations.
369: .keywords: SNES, SNESNGMRES, restart, type, set SNESLineSearch
370: @*/
371: PetscErrorCode SNESNGMRESSetRestartType(SNES snes,SNESNGMRESRestartType rtype)
372: {
377: PetscTryMethod(snes,"SNESNGMRESSetRestartType_C",(SNES,SNESNGMRESRestartType),(snes,rtype));
378: return(0);
379: }
383: /*@
384: SNESNGMRESSetSelectType - Sets the selection type for SNESNGMRES. This determines how the candidate solution and
385: combined solution are used to create the next iterate.
387: Logically Collective on SNES
389: Input Parameters:
390: + snes - the iterative context
391: - stype - selection type
393: Options Database:
394: . -snes_ngmres_select_type<difference,none,linesearch>
396: Level: intermediate
398: SNESNGMRESSelectTypes:
399: + SNES_NGMRES_SELECT_NONE - choose the combined solution all the time
400: . SNES_NGMRES_SELECT_DIFFERENCE - choose based upon the selection criteria
401: - SNES_NGMRES_SELECT_LINESEARCH - choose based upon line search combination
403: Notes:
404: The default line search used is the L2 line search and it requires two additional function evaluations.
406: .keywords: SNES, SNESNGMRES, selection, type, set SNESLineSearch
407: @*/
408: PetscErrorCode SNESNGMRESSetSelectType(SNES snes,SNESNGMRESSelectType stype)
409: {
414: PetscTryMethod(snes,"SNESNGMRESSetSelectType_C",(SNES,SNESNGMRESSelectType),(snes,stype));
415: return(0);
416: }
420: PetscErrorCode SNESNGMRESSetSelectType_NGMRES(SNES snes,SNESNGMRESSelectType stype)
421: {
422: SNES_NGMRES *ngmres = (SNES_NGMRES*)snes->data;
425: ngmres->select_type = stype;
426: return(0);
427: }
431: PetscErrorCode SNESNGMRESSetRestartType_NGMRES(SNES snes,SNESNGMRESRestartType rtype)
432: {
433: SNES_NGMRES *ngmres = (SNES_NGMRES*)snes->data;
436: ngmres->restart_type = rtype;
437: return(0);
438: }
440: /*MC
441: SNESNGMRES - The Nonlinear Generalized Minimum Residual method.
443: Level: beginner
445: Options Database:
446: + -snes_ngmres_select_type<difference,none,linesearch> - choose the select between candidate and combined solution
447: . -snes_ngmres_restart_type<difference,none,periodic> - choose the restart conditions
448: . -snes_ngmres_candidate - Use NGMRES variant which combines candidate solutions instead of actual solutions
449: . -snes_ngmres_m - Number of stored previous solutions and residuals
450: . -snes_ngmres_restart_it - Number of iterations the restart conditions hold before restart
451: . -snes_ngmres_gammaA - Residual tolerance for solution select between the candidate and combination
452: . -snes_ngmres_gammaC - Residual tolerance for restart
453: . -snes_ngmres_epsilonB - Difference tolerance between subsequent solutions triggering restart
454: . -snes_ngmres_deltaB - Difference tolerance between residuals triggering restart
455: . -snes_ngmres_monitor - Prints relevant information about the ngmres iteration
456: . -snes_linesearch_type <basic,l2,cp> - Line search type used for the default smoother
457: - -additive_snes_linesearch_type - linesearch type used to select between the candidate and combined solution with additive select type
459: Notes:
461: The N-GMRES method combines m previous solutions into a minimum-residual solution by solving a small linearized
462: optimization problem at each iteration.
464: References:
466: "Krylov Subspace Acceleration of Nonlinear Multigrid with Application to Recirculating Flows", C. W. Oosterlee and T. Washio,
467: SIAM Journal on Scientific Computing, 21(5), 2000.
469: .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types)
470: M*/
474: PETSC_EXTERN PetscErrorCode SNESCreate_NGMRES(SNES snes)
475: {
476: SNES_NGMRES *ngmres;
480: snes->ops->destroy = SNESDestroy_NGMRES;
481: snes->ops->setup = SNESSetUp_NGMRES;
482: snes->ops->setfromoptions = SNESSetFromOptions_NGMRES;
483: snes->ops->view = SNESView_NGMRES;
484: snes->ops->solve = SNESSolve_NGMRES;
485: snes->ops->reset = SNESReset_NGMRES;
487: snes->usespc = PETSC_TRUE;
488: snes->usesksp = PETSC_FALSE;
489: snes->pcside = PC_RIGHT;
491: PetscNewLog(snes,&ngmres);
492: snes->data = (void*) ngmres;
493: ngmres->msize = 30;
495: if (!snes->tolerancesset) {
496: snes->max_funcs = 30000;
497: snes->max_its = 10000;
498: }
500: ngmres->candidate = PETSC_FALSE;
502: ngmres->additive_linesearch = NULL;
503: ngmres->approxfunc = PETSC_FALSE;
504: ngmres->restart_it = 2;
505: ngmres->restart_periodic = 30;
506: ngmres->gammaA = 2.0;
507: ngmres->gammaC = 2.0;
508: ngmres->deltaB = 0.9;
509: ngmres->epsilonB = 0.1;
511: ngmres->restart_type = SNES_NGMRES_RESTART_DIFFERENCE;
512: ngmres->select_type = SNES_NGMRES_SELECT_DIFFERENCE;
514: PetscObjectComposeFunction((PetscObject)snes,"SNESNGMRESSetSelectType_C",SNESNGMRESSetSelectType_NGMRES);
515: PetscObjectComposeFunction((PetscObject)snes,"SNESNGMRESSetRestartType_C",SNESNGMRESSetRestartType_NGMRES);
516: return(0);
517: }