Actual source code: snesngmres.c
petsc-3.5.4 2015-05-23
1: #include <../src/snes/impls/ngmres/snesngmres.h> /*I "petscsnes.h" I*/
2: #include <petscblaslapack.h>
4: const char *const SNESNGMRESRestartTypes[] = {"NONE","PERIODIC","DIFFERENCE","SNESNGMRESRestartType","SNES_NGMRES_RESTART_",0};
5: const char *const SNESNGMRESSelectTypes[] = {"NONE","DIFFERENCE","LINESEARCH","SNESNGMRESSelectType","SNES_NGMRES_SELECT_",0};
9: PetscErrorCode SNESReset_NGMRES(SNES snes)
10: {
11: SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data;
15: VecDestroyVecs(ngmres->msize,&ngmres->Fdot);
16: VecDestroyVecs(ngmres->msize,&ngmres->Xdot);
17: SNESLineSearchDestroy(&ngmres->additive_linesearch);
18: return(0);
19: }
23: PetscErrorCode SNESDestroy_NGMRES(SNES snes)
24: {
26: SNES_NGMRES *ngmres = (SNES_NGMRES*)snes->data;
29: SNESReset_NGMRES(snes);
30: PetscFree5(ngmres->h,ngmres->beta,ngmres->xi,ngmres->fnorms,ngmres->q);
31: PetscFree(ngmres->s);
32: PetscFree(ngmres->xnorms);
33: #if PETSC_USE_COMPLEX
34: PetscFree(ngmres->rwork);
35: #endif
36: PetscFree(ngmres->work);
37: PetscFree(snes->data);
38: return(0);
39: }
43: PetscErrorCode SNESSetUp_NGMRES(SNES snes)
44: {
45: SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data;
46: const char *optionsprefix;
47: PetscInt msize,hsize;
51: if (snes->pc && snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
52: SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"SNESNGMRES does not support left preconditioning with unpreconditioned function");
53: }
54: if (snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_PRECONDITIONED;
55: SNESSetWorkVecs(snes,5);
56: if (!ngmres->Xdot) {VecDuplicateVecs(snes->vec_sol,ngmres->msize,&ngmres->Xdot);}
57: if (!ngmres->Fdot) {VecDuplicateVecs(snes->vec_sol,ngmres->msize,&ngmres->Fdot);}
58: if (!ngmres->setup_called) {
59: msize = ngmres->msize; /* restart size */
60: hsize = msize * msize;
62: /* explicit least squares minimization solve */
63: PetscMalloc5(hsize,&ngmres->h, msize,&ngmres->beta, msize,&ngmres->xi, msize,&ngmres->fnorms, hsize,&ngmres->q);
64: PetscMalloc1(msize,&ngmres->xnorms);
65: ngmres->nrhs = 1;
66: ngmres->lda = msize;
67: ngmres->ldb = msize;
68: PetscMalloc1(msize,&ngmres->s);
69: PetscMemzero(ngmres->h, hsize*sizeof(PetscScalar));
70: PetscMemzero(ngmres->q, hsize*sizeof(PetscScalar));
71: PetscMemzero(ngmres->xi, msize*sizeof(PetscScalar));
72: PetscMemzero(ngmres->beta,msize*sizeof(PetscScalar));
73: ngmres->lwork = 12*msize;
74: #if PETSC_USE_COMPLEX
75: PetscMalloc(sizeof(PetscReal)*ngmres->lwork,&ngmres->rwork);
76: #endif
77: PetscMalloc(sizeof(PetscScalar)*ngmres->lwork,&ngmres->work);
78: }
80: /* linesearch setup */
81: SNESGetOptionsPrefix(snes,&optionsprefix);
83: if (ngmres->select_type == SNES_NGMRES_SELECT_LINESEARCH) {
84: SNESLineSearchCreate(PetscObjectComm((PetscObject)snes),&ngmres->additive_linesearch);
85: SNESLineSearchSetSNES(ngmres->additive_linesearch,snes);
86: SNESLineSearchSetType(ngmres->additive_linesearch,SNESLINESEARCHL2);
87: SNESLineSearchAppendOptionsPrefix(ngmres->additive_linesearch,"additive_");
88: SNESLineSearchAppendOptionsPrefix(ngmres->additive_linesearch,optionsprefix);
89: SNESLineSearchSetFromOptions(ngmres->additive_linesearch);
90: }
92: ngmres->setup_called = PETSC_TRUE;
93: return(0);
94: }
98: PetscErrorCode SNESSetFromOptions_NGMRES(SNES snes)
99: {
100: SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data;
102: PetscBool debug;
103: SNESLineSearch linesearch;
106: PetscOptionsHead("SNES NGMRES options");
107: PetscOptionsEnum("-snes_ngmres_select_type","Select type","SNESNGMRESSetSelectType",SNESNGMRESSelectTypes,
108: (PetscEnum)ngmres->select_type,(PetscEnum*)&ngmres->select_type,NULL);
109: PetscOptionsEnum("-snes_ngmres_restart_type","Restart type","SNESNGMRESSetRestartType",SNESNGMRESRestartTypes,
110: (PetscEnum)ngmres->restart_type,(PetscEnum*)&ngmres->restart_type,NULL);
111: PetscOptionsBool("-snes_ngmres_candidate", "Use candidate storage", "SNES",ngmres->candidate,&ngmres->candidate,NULL);
112: PetscOptionsBool("-snes_ngmres_approxfunc","Linearly approximate the function", "SNES",ngmres->approxfunc,&ngmres->approxfunc,NULL);
113: PetscOptionsInt("-snes_ngmres_m", "Number of directions", "SNES",ngmres->msize,&ngmres->msize,NULL);
114: PetscOptionsInt("-snes_ngmres_restart", "Iterations before forced restart", "SNES",ngmres->restart_periodic,&ngmres->restart_periodic,NULL);
115: PetscOptionsInt("-snes_ngmres_restart_it", "Tolerance iterations before restart","SNES",ngmres->restart_it,&ngmres->restart_it,NULL);
116: PetscOptionsBool("-snes_ngmres_monitor", "Monitor actions of NGMRES", "SNES",ngmres->monitor ? PETSC_TRUE : PETSC_FALSE,&debug,NULL);
117: if (debug) {
118: ngmres->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));
119: }
120: PetscOptionsReal("-snes_ngmres_gammaA", "Residual selection constant", "SNES",ngmres->gammaA,&ngmres->gammaA,NULL);
121: PetscOptionsReal("-snes_ngmres_gammaC", "Residual restart constant", "SNES",ngmres->gammaC,&ngmres->gammaC,NULL);
122: PetscOptionsReal("-snes_ngmres_epsilonB", "Difference selection constant", "SNES",ngmres->epsilonB,&ngmres->epsilonB,NULL);
123: PetscOptionsReal("-snes_ngmres_deltaB", "Difference residual selection constant", "SNES",ngmres->deltaB,&ngmres->deltaB,NULL);
124: PetscOptionsBool("-snes_ngmres_single_reduction", "Aggregate reductions", "SNES",ngmres->singlereduction,&ngmres->singlereduction,NULL);
125: PetscOptionsTail();
126: if ((ngmres->gammaA > ngmres->gammaC) && (ngmres->gammaC > 2.)) ngmres->gammaC = ngmres->gammaA;
128: /* set the default type of the line search if the user hasn't already. */
129: if (!snes->linesearch) {
130: SNESGetLineSearch(snes,&linesearch);
131: SNESLineSearchSetType(linesearch,SNESLINESEARCHBASIC);
132: }
133: return(0);
134: }
138: PetscErrorCode SNESView_NGMRES(SNES snes,PetscViewer viewer)
139: {
140: SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data;
141: PetscBool iascii;
145: PetscObjectTypeCompare((PetscObject) viewer,PETSCVIEWERASCII,&iascii);
146: if (iascii) {
147: PetscViewerASCIIPrintf(viewer," Number of stored past updates: %d\n", ngmres->msize);
148: PetscViewerASCIIPrintf(viewer," Residual selection: gammaA=%1.0e, gammaC=%1.0e\n",ngmres->gammaA,ngmres->gammaC);
149: PetscViewerASCIIPrintf(viewer," Difference restart: epsilonB=%1.0e, deltaB=%1.0e\n",ngmres->epsilonB,ngmres->deltaB);
150: }
151: return(0);
152: }
156: PetscErrorCode SNESSolve_NGMRES(SNES snes)
157: {
158: SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data;
159: /* present solution, residual, and preconditioned residual */
160: Vec X,F,B,D,Y;
162: /* candidate linear combination answers */
163: Vec XA,FA,XM,FM;
165: /* coefficients and RHS to the minimization problem */
166: PetscReal fnorm,fMnorm,fAnorm;
167: PetscReal xnorm,xMnorm,xAnorm;
168: PetscReal ynorm,yMnorm,yAnorm;
169: PetscInt k,k_restart,l,ivec,restart_count = 0;
171: /* solution selection data */
172: PetscBool selectRestart;
173: PetscReal dnorm,dminnorm = 0.0;
174: PetscReal fminnorm;
176: SNESConvergedReason reason;
177: PetscBool lssucceed;
178: PetscErrorCode ierr;
181: PetscCitationsRegister(SNESCitation,&SNEScite);
182: /* variable initialization */
183: snes->reason = SNES_CONVERGED_ITERATING;
184: X = snes->vec_sol;
185: F = snes->vec_func;
186: B = snes->vec_rhs;
187: XA = snes->vec_sol_update;
188: FA = snes->work[0];
189: D = snes->work[1];
191: /* work for the line search */
192: Y = snes->work[2];
193: XM = snes->work[3];
194: FM = snes->work[4];
196: PetscObjectSAWsTakeAccess((PetscObject)snes);
197: snes->iter = 0;
198: snes->norm = 0.;
199: PetscObjectSAWsGrantAccess((PetscObject)snes);
201: /* initialization */
203: if (snes->pc && snes->pcside == PC_LEFT) {
204: SNESApplyNPC(snes,X,NULL,F);
205: SNESGetConvergedReason(snes->pc,&reason);
206: if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
207: snes->reason = SNES_DIVERGED_INNER;
208: return(0);
209: }
210: VecNorm(F,NORM_2,&fnorm);
211: } else {
212: if (!snes->vec_func_init_set) {
213: SNESComputeFunction(snes,X,F);
214: if (snes->domainerror) {
215: snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
216: return(0);
217: }
218: } else snes->vec_func_init_set = PETSC_FALSE;
220: VecNorm(F,NORM_2,&fnorm);
221: if (PetscIsInfOrNanReal(fnorm)) {
222: snes->reason = SNES_DIVERGED_FNORM_NAN;
223: return(0);
224: }
225: }
226: fminnorm = fnorm;
228: PetscObjectSAWsTakeAccess((PetscObject)snes);
229: snes->norm = fnorm;
230: PetscObjectSAWsGrantAccess((PetscObject)snes);
231: SNESLogConvergenceHistory(snes,fnorm,0);
232: SNESMonitor(snes,0,fnorm);
233: (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);
234: if (snes->reason) return(0);
235: SNESNGMRESUpdateSubspace_Private(snes,0,0,F,fnorm,X);
237: k_restart = 1;
238: l = 1;
239: ivec = 0;
240: for (k=1; k < snes->max_its+1; k++) {
241: /* Computation of x^M */
242: if (snes->pc && snes->pcside == PC_RIGHT) {
243: VecCopy(X,XM);
244: SNESSetInitialFunction(snes->pc,F);
246: PetscLogEventBegin(SNES_NPCSolve,snes->pc,XM,B,0);
247: SNESSolve(snes->pc,B,XM);
248: PetscLogEventEnd(SNES_NPCSolve,snes->pc,XM,B,0);
250: SNESGetConvergedReason(snes->pc,&reason);
251: if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
252: snes->reason = SNES_DIVERGED_INNER;
253: return(0);
254: }
255: SNESGetNPCFunction(snes,FM,&fMnorm);
256: } else {
257: /* no preconditioner -- just take gradient descent with line search */
258: VecCopy(F,Y);
259: VecCopy(F,FM);
260: VecCopy(X,XM);
262: fMnorm = fnorm;
264: SNESLineSearchApply(snes->linesearch,XM,FM,&fMnorm,Y);
265: SNESLineSearchGetSuccess(snes->linesearch,&lssucceed);
266: if (!lssucceed) {
267: if (++snes->numFailures >= snes->maxFailures) {
268: snes->reason = SNES_DIVERGED_LINE_SEARCH;
269: return(0);
270: }
271: }
272: }
273: SNESNGMRESFormCombinedSolution_Private(snes,ivec,l,XM,FM,fMnorm,X,XA,FA);
274: /* r = F(x) */
275: if (fminnorm > fMnorm) fminnorm = fMnorm; /* the minimum norm is now of F^M */
277: /* differences for selection and restart */
278: if (ngmres->restart_type == SNES_NGMRES_RESTART_DIFFERENCE || ngmres->select_type == SNES_NGMRES_SELECT_DIFFERENCE) {
279: SNESNGMRESNorms_Private(snes,l,X,F,XM,FM,XA,FA,D,&dnorm,&dminnorm,&xMnorm,NULL,&yMnorm,&xAnorm,&fAnorm,&yAnorm);
280: } else {
281: SNESNGMRESNorms_Private(snes,l,X,F,XM,FM,XA,FA,D,NULL,NULL,&xMnorm,NULL,&yMnorm,&xAnorm,&fAnorm,&yAnorm);
282: }
283: if (PetscIsInfOrNanReal(fAnorm)) {
284: snes->reason = SNES_DIVERGED_FNORM_NAN;
285: return(0);
286: }
288: /* combination (additive) or selection (multiplicative) of the N-GMRES solution */
289: SNESNGMRESSelect_Private(snes,k_restart,XM,FM,xMnorm,fMnorm,yMnorm,XA,FA,xAnorm,fAnorm,yAnorm,dnorm,fminnorm,dminnorm,X,F,Y,&xnorm,&fnorm,&ynorm);
290: selectRestart = PETSC_FALSE;
291: if (ngmres->restart_type == SNES_NGMRES_RESTART_DIFFERENCE) {
292: SNESNGMRESSelectRestart_Private(snes,l,fAnorm,dnorm,fminnorm,dminnorm,&selectRestart);
293: /* if the restart conditions persist for more than restart_it iterations, restart. */
294: if (selectRestart) restart_count++;
295: else restart_count = 0;
296: } else if (ngmres->restart_type == SNES_NGMRES_RESTART_PERIODIC) {
297: if (k_restart > ngmres->restart_periodic) {
298: if (ngmres->monitor) PetscViewerASCIIPrintf(ngmres->monitor,"periodic restart after %D iterations\n",k_restart);
299: restart_count = ngmres->restart_it;
300: }
301: }
302: ivec = k_restart % ngmres->msize; /* replace the last used part of the subspace */
303: /* restart after restart conditions have persisted for a fixed number of iterations */
304: if (restart_count >= ngmres->restart_it) {
305: if (ngmres->monitor) {
306: PetscViewerASCIIPrintf(ngmres->monitor,"Restarted at iteration %d\n",k_restart);
307: }
308: restart_count = 0;
309: k_restart = 1;
310: l = 1;
311: ivec = 0;
312: /* q_{00} = nu */
313: SNESNGMRESUpdateSubspace_Private(snes,0,0,FM,fMnorm,XM);
314: } else {
315: /* select the current size of the subspace */
316: if (l < ngmres->msize) l++;
317: k_restart++;
318: /* place the current entry in the list of previous entries */
319: if (ngmres->candidate) {
320: if (fminnorm > fMnorm) fminnorm = fMnorm;
321: SNESNGMRESUpdateSubspace_Private(snes,ivec,l,FM,fMnorm,XM);
322: } else {
323: if (fminnorm > fnorm) fminnorm = fnorm;
324: SNESNGMRESUpdateSubspace_Private(snes,ivec,l,F,fnorm,X);
325: }
326: }
328: PetscObjectSAWsTakeAccess((PetscObject)snes);
329: snes->iter = k;
330: snes->norm = fnorm;
331: PetscObjectSAWsGrantAccess((PetscObject)snes);
332: SNESLogConvergenceHistory(snes,snes->norm,snes->iter);
333: SNESMonitor(snes,snes->iter,snes->norm);
334: (*snes->ops->converged)(snes,snes->iter,0,0,fnorm,&snes->reason,snes->cnvP);
335: if (snes->reason) return(0);
336: }
337: snes->reason = SNES_DIVERGED_MAX_IT;
338: return(0);
339: }
343: /*@
344: SNESNGMRESSetRestartType - Sets the restart type for SNESNGMRES.
346: Logically Collective on SNES
348: Input Parameters:
349: + snes - the iterative context
350: - rtype - restart type
352: Options Database:
353: + -snes_ngmres_restart_type<difference,periodic,none> - set the restart type
354: - -snes_ngmres_restart[30] - sets the number of iterations before restart for periodic
356: Level: intermediate
358: SNESNGMRESRestartTypes:
359: + SNES_NGMRES_RESTART_NONE - never restart
360: . SNES_NGMRES_RESTART_DIFFERENCE - restart based upon difference criteria
361: - SNES_NGMRES_RESTART_PERIODIC - restart after a fixed number of iterations
363: Notes:
364: The default line search used is the L2 line search and it requires two additional function evaluations.
366: .keywords: SNES, SNESNGMRES, restart, type, set SNESLineSearch
367: @*/
368: PetscErrorCode SNESNGMRESSetRestartType(SNES snes,SNESNGMRESRestartType rtype)
369: {
374: PetscTryMethod(snes,"SNESNGMRESSetRestartType_C",(SNES,SNESNGMRESRestartType),(snes,rtype));
375: return(0);
376: }
380: /*@
381: SNESNGMRESSetSelectType - Sets the selection type for SNESNGMRES. This determines how the candidate solution and
382: combined solution are used to create the next iterate.
384: Logically Collective on SNES
386: Input Parameters:
387: + snes - the iterative context
388: - stype - selection type
390: Options Database:
391: . -snes_ngmres_select_type<difference,none,linesearch>
393: Level: intermediate
395: SNESNGMRESSelectTypes:
396: + SNES_NGMRES_SELECT_NONE - choose the combined solution all the time
397: . SNES_NGMRES_SELECT_DIFFERENCE - choose based upon the selection criteria
398: - SNES_NGMRES_SELECT_LINESEARCH - choose based upon line search combination
400: Notes:
401: The default line search used is the L2 line search and it requires two additional function evaluations.
403: .keywords: SNES, SNESNGMRES, selection, type, set SNESLineSearch
404: @*/
405: PetscErrorCode SNESNGMRESSetSelectType(SNES snes,SNESNGMRESSelectType stype)
406: {
411: PetscTryMethod(snes,"SNESNGMRESSetSelectType_C",(SNES,SNESNGMRESSelectType),(snes,stype));
412: return(0);
413: }
417: PetscErrorCode SNESNGMRESSetSelectType_NGMRES(SNES snes,SNESNGMRESSelectType stype)
418: {
419: SNES_NGMRES *ngmres = (SNES_NGMRES*)snes->data;
422: ngmres->select_type = stype;
423: return(0);
424: }
428: PetscErrorCode SNESNGMRESSetRestartType_NGMRES(SNES snes,SNESNGMRESRestartType rtype)
429: {
430: SNES_NGMRES *ngmres = (SNES_NGMRES*)snes->data;
433: ngmres->restart_type = rtype;
434: return(0);
435: }
437: /*MC
438: SNESNGMRES - The Nonlinear Generalized Minimum Residual method.
440: Level: beginner
442: Options Database:
443: + -snes_ngmres_select_type<difference,none,linesearch> - choose the select between candidate and combined solution
444: . -snes_ngmres_restart_type<difference,none,periodic> - choose the restart conditions
445: . -snes_ngmres_candidate - Use NGMRES variant which combines candidate solutions instead of actual solutions
446: . -snes_ngmres_m - Number of stored previous solutions and residuals
447: . -snes_ngmres_restart_it - Number of iterations the restart conditions hold before restart
448: . -snes_ngmres_gammaA - Residual tolerance for solution select between the candidate and combination
449: . -snes_ngmres_gammaC - Residual tolerance for restart
450: . -snes_ngmres_epsilonB - Difference tolerance between subsequent solutions triggering restart
451: . -snes_ngmres_deltaB - Difference tolerance between residuals triggering restart
452: . -snes_ngmres_monitor - Prints relevant information about the ngmres iteration
453: . -snes_linesearch_type <basic,l2,cp> - Line search type used for the default smoother
454: - -additive_snes_linesearch_type - linesearch type used to select between the candidate and combined solution with additive select type
456: Notes:
458: The N-GMRES method combines m previous solutions into a minimum-residual solution by solving a small linearized
459: optimization problem at each iteration.
461: References:
463: "Krylov Subspace Acceleration of Nonlinear Multigrid with Application to Recirculating Flows", C. W. Oosterlee and T. Washio,
464: SIAM Journal on Scientific Computing, 21(5), 2000.
466: .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types)
467: M*/
471: PETSC_EXTERN PetscErrorCode SNESCreate_NGMRES(SNES snes)
472: {
473: SNES_NGMRES *ngmres;
477: snes->ops->destroy = SNESDestroy_NGMRES;
478: snes->ops->setup = SNESSetUp_NGMRES;
479: snes->ops->setfromoptions = SNESSetFromOptions_NGMRES;
480: snes->ops->view = SNESView_NGMRES;
481: snes->ops->solve = SNESSolve_NGMRES;
482: snes->ops->reset = SNESReset_NGMRES;
484: snes->usespc = PETSC_TRUE;
485: snes->usesksp = PETSC_FALSE;
486: snes->pcside = PC_RIGHT;
488: PetscNewLog(snes,&ngmres);
489: snes->data = (void*) ngmres;
490: ngmres->msize = 30;
492: if (!snes->tolerancesset) {
493: snes->max_funcs = 30000;
494: snes->max_its = 10000;
495: }
497: ngmres->candidate = PETSC_FALSE;
499: ngmres->additive_linesearch = NULL;
500: ngmres->approxfunc = PETSC_FALSE;
501: ngmres->restart_it = 2;
502: ngmres->restart_periodic = 30;
503: ngmres->gammaA = 2.0;
504: ngmres->gammaC = 2.0;
505: ngmres->deltaB = 0.9;
506: ngmres->epsilonB = 0.1;
508: ngmres->restart_type = SNES_NGMRES_RESTART_DIFFERENCE;
509: ngmres->select_type = SNES_NGMRES_SELECT_DIFFERENCE;
511: PetscObjectComposeFunction((PetscObject)snes,"SNESNGMRESSetSelectType_C",SNESNGMRESSetSelectType_NGMRES);
512: PetscObjectComposeFunction((PetscObject)snes,"SNESNGMRESSetRestartType_C",SNESNGMRESSetRestartType_NGMRES);
513: return(0);
514: }