Actual source code: ms.c
petsc-3.3-p7 2013-05-11
1: #include <petsc-private/snesimpl.h> /*I "petscsnes.h" I*/
3: static const SNESMSType SNESMSDefault = SNESMSM62;
4: static PetscBool SNESMSRegisterAllCalled;
5: static PetscBool SNESMSPackageInitialized;
7: typedef struct _SNESMSTableau *SNESMSTableau;
8: struct _SNESMSTableau {
9: char *name;
10: PetscInt nstages; /* Number of stages */
11: PetscInt nregisters; /* Number of registers */
12: PetscReal stability; /* Scaled stability region */
13: PetscReal *gamma; /* Coefficients of 3S* method */
14: PetscReal *delta; /* Coefficients of 3S* method */
15: PetscReal *betasub; /* Subdiagonal of beta in Shu-Osher form */
16: };
18: typedef struct _SNESMSTableauLink *SNESMSTableauLink;
19: struct _SNESMSTableauLink {
20: struct _SNESMSTableau tab;
21: SNESMSTableauLink next;
22: };
23: static SNESMSTableauLink SNESMSTableauList;
25: typedef struct {
26: SNESMSTableau tableau; /* Tableau in low-storage form */
27: PetscReal damping; /* Damping parameter, like length of (pseudo) time step */
28: PetscBool norms; /* Compute norms, usually only for monitoring purposes */
29: } SNES_MS;
33: /*@C
34: SNESMSRegisterAll - Registers all of the multi-stage methods in SNESMS
36: Not Collective, but should be called by all processes which will need the schemes to be registered
38: Level: advanced
40: .keywords: SNES, SNESMS, register, all
42: .seealso: SNESMSRegisterDestroy()
43: @*/
44: PetscErrorCode SNESMSRegisterAll(void)
45: {
49: if (SNESMSRegisterAllCalled) return(0);
50: SNESMSRegisterAllCalled = PETSC_TRUE;
52: {
53: PetscReal gamma[3][1] = {{1.0},{0.0},{0.0}};
54: PetscReal delta[1] = {0.0};
55: PetscReal betasub[1] = {1.0};
56: SNESMSRegister(SNESMSEULER,1,3,1.0,&gamma[0][0],delta,betasub);
57: }
59: {
60: PetscReal gamma[3][6] = {
61: {0.0000000000000000E+00, -7.0304722367110606E-01, -1.9836719667506464E-01, -1.6023843981863788E+00, 9.4483822882855284E-02, -1.4204296130641869E-01},
62: {1.0000000000000000E+00, 1.1111025767083920E+00, 5.6150921583923230E-01, 7.4151723494934041E-01, 3.1714538168600587E-01, 4.6479276238548706E-01},
63: {0.0000000000000000E+00, 0.0000000000000000E+00, 0.0000000000000000E+00, 6.7968174970583317E-01, -4.1755042846051737E-03, -1.9115668129923846E-01}
64: };
65: PetscReal delta[6] = {1.0000000000000000E+00, 5.3275427433201750E-01, 6.0143544663985238E-01, 4.5874077053842177E-01, 2.7544386906104651E-01, 0.0000000000000000E+00};
66: PetscReal betasub[6] = {8.4753115429481929E-01, 7.4018896368655618E-01, 6.5963574086583309E-03, 4.6747795645517759E-01, 1.3314545813643919E-01, 5.3260800028018784E-01};
67: SNESMSRegister(SNESMSM62,6,3,1.0,&gamma[0][0],delta,betasub);
68: }
69: {
70: PetscReal gamma[3][4] = {{0,0,0,0},{0,0,0,0},{1,1,1,1}};
71: PetscReal delta[4] = {0,0,0,0};
72: PetscReal betasub[4] = {0.25, 0.5, 0.55, 1.0};
73: SNESMSRegister(SNESMSJAMESON83,4,3,1.0,&gamma[0][0],delta,betasub);
74: }
75: { /* Van Leer, Tai, and Powell (1989) 2 stage, order 1 */
76: PetscReal gamma[3][2] = {{0,0},{0,0},{1,1}};
77: PetscReal delta[2] = {0,0};
78: PetscReal betasub[2] = {0.3333,1.0};
79: SNESMSRegister(SNESMSVLTP21,2,3,1.0,&gamma[0][0],delta,betasub);
80: }
81: { /* Van Leer, Tai, and Powell (1989) 3 stage, order 1 */
82: PetscReal gamma[3][3] = {{0,0,0},{0,0,0},{1,1,1}};
83: PetscReal delta[3] = {0,0,0};
84: PetscReal betasub[3] = {0.1481,0.4000,1.0};
85: SNESMSRegister(SNESMSVLTP31,3,3,1.5,&gamma[0][0],delta,betasub);
86: }
87: { /* Van Leer, Tai, and Powell (1989) 4 stage, order 1 */
88: PetscReal gamma[3][4] = {{0,0,0,0},{0,0,0,0},{1,1,1,1}};
89: PetscReal delta[4] = {0,0,0,0};
90: PetscReal betasub[4] = {0.0833,0.2069,0.4265,1.0};
91: SNESMSRegister(SNESMSVLTP41,4,3,2.0,&gamma[0][0],delta,betasub);
92: }
93: { /* Van Leer, Tai, and Powell (1989) 5 stage, order 1 */
94: PetscReal gamma[3][5] = {{0,0,0,0,0},{0,0,0,0,0},{1,1,1,1,1}};
95: PetscReal delta[5] = {0,0,0,0,0};
96: PetscReal betasub[5] = {0.0533,0.1263,0.2375,0.4414,1.0};
97: SNESMSRegister(SNESMSVLTP51,5,3,2.5,&gamma[0][0],delta,betasub);
98: }
99: { /* Van Leer, Tai, and Powell (1989) 6 stage, order 1 */
100: PetscReal gamma[3][6] = {{0,0,0,0,0,0},{0,0,0,0,0,0},{1,1,1,1,1,1}};
101: PetscReal delta[6] = {0,0,0,0,0,0};
102: PetscReal betasub[6] = {0.0370,0.0851,0.1521,0.2562,0.4512,1.0};
103: SNESMSRegister(SNESMSVLTP61,6,3,3.0,&gamma[0][0],delta,betasub);
104: }
105: return(0);
106: }
110: /*@C
111: SNESMSRegisterDestroy - Frees the list of schemes that were registered by TSRosWRegister().
113: Not Collective
115: Level: advanced
117: .keywords: TSRosW, register, destroy
118: .seealso: TSRosWRegister(), TSRosWRegisterAll(), TSRosWRegisterDynamic()
119: @*/
120: PetscErrorCode SNESMSRegisterDestroy(void)
121: {
123: SNESMSTableauLink link;
126: while ((link = SNESMSTableauList)) {
127: SNESMSTableau t = &link->tab;
128: SNESMSTableauList = link->next;
129: PetscFree3(t->gamma,t->delta,t->betasub);
130: PetscFree(t->name);
131: PetscFree(link);
132: }
133: SNESMSRegisterAllCalled = PETSC_FALSE;
134: return(0);
135: }
139: /*@C
140: SNESMSInitializePackage - This function initializes everything in the SNESMS package. It is called
141: from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to SNESCreate_MS()
142: when using static libraries.
144: Input Parameter:
145: path - The dynamic library path, or PETSC_NULL
147: Level: developer
149: .keywords: SNES, SNESMS, initialize, package
150: .seealso: PetscInitialize()
151: @*/
152: PetscErrorCode SNESMSInitializePackage(const char path[])
153: {
157: if (SNESMSPackageInitialized) return(0);
158: SNESMSPackageInitialized = PETSC_TRUE;
159: SNESMSRegisterAll();
160: PetscRegisterFinalize(SNESMSFinalizePackage);
161: return(0);
162: }
166: /*@C
167: SNESMSFinalizePackage - This function destroys everything in the SNESMS package. It is
168: called from PetscFinalize().
170: Level: developer
172: .keywords: Petsc, destroy, package
173: .seealso: PetscFinalize()
174: @*/
175: PetscErrorCode SNESMSFinalizePackage(void)
176: {
180: SNESMSPackageInitialized = PETSC_FALSE;
181: SNESMSRegisterDestroy();
182: return(0);
183: }
187: /*@C
188: SNESMSRegister - register a multistage scheme
190: Not Collective, but the same schemes should be registered on all processes on which they will be used
192: Input Parameters:
193: + name - identifier for method
194: . nstages - number of stages
195: . nregisters - number of registers used by low-storage implementation
196: . gamma - coefficients, see Ketcheson's paper
197: . delta - coefficients, see Ketcheson's paper
198: - betasub - subdiagonal of Shu-Osher form
200: Notes:
201: The notation is described in Ketcheson (2010) Runge-Kutta methods with minimum storage implementations.
203: Level: advanced
205: .keywords: SNES, register
207: .seealso: SNESMS
208: @*/
209: PetscErrorCode SNESMSRegister(const SNESMSType name,PetscInt nstages,PetscInt nregisters,PetscReal stability,const PetscReal gamma[],const PetscReal delta[],const PetscReal betasub[])
210: {
212: SNESMSTableauLink link;
213: SNESMSTableau t;
217: if (nstages < 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Must have at least one stage");
218: if (nregisters != 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only support for methods written in 3-register form");
223: PetscMalloc(sizeof(*link),&link);
224: PetscMemzero(link,sizeof(*link));
225: t = &link->tab;
226: PetscStrallocpy(name,&t->name);
227: t->nstages = nstages;
228: t->nregisters = nregisters;
229: t->stability = stability;
230: PetscMalloc3(nstages*nregisters,PetscReal,&t->gamma,nstages,PetscReal,&t->delta,nstages,PetscReal,&t->betasub);
231: PetscMemcpy(t->gamma,gamma,nstages*nregisters*sizeof(PetscReal));
232: PetscMemcpy(t->delta,delta,nstages*sizeof(PetscReal));
233: PetscMemcpy(t->betasub,betasub,nstages*sizeof(PetscReal));
235: link->next = SNESMSTableauList;
236: SNESMSTableauList = link;
237: return(0);
238: }
242: /*
243: X - initial state, updated in-place.
244: F - residual, computed at the initial X on input
245: */
246: static PetscErrorCode SNESMSStep_3Sstar(SNES snes,Vec X,Vec F)
247: {
249: SNES_MS *ms = (SNES_MS*)snes->data;
250: SNESMSTableau t = ms->tableau;
251: const PetscReal *gamma = t->gamma,*delta = t->delta,*betasub = t->betasub;
252: Vec S1,S2,S3,Y;
253: PetscInt i,nstages = t->nstages;;
257: Y = snes->work[0];
258: S1 = X;
259: S2 = snes->work[1];
260: S3 = snes->work[2];
261: VecZeroEntries(S2);
262: VecCopy(X,S3);
263: for (i=0; i<nstages; i++) {
264: Vec Ss[4] = {S1,S2,S3,Y};
265: PetscScalar scoeff[4] = {gamma[0*nstages+i]-1.0,gamma[1*nstages+i],gamma[2*nstages+i],-betasub[i]*ms->damping};
266: VecAXPY(S2,delta[i],S1);
267: if (i>0) {
268: SNESComputeFunction(snes,S1,F);
269: if (snes->domainerror) {
270: snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
271: return(0);
272: }
273: }
274: SNES_KSPSolve(snes,snes->ksp,F,Y);
275: VecMAXPY(S1,4,scoeff,Ss);
276: }
277: return(0);
278: }
282: static PetscErrorCode SNESSolve_MS(SNES snes)
283: {
284: SNES_MS *ms = (SNES_MS*)snes->data;
285: Vec X = snes->vec_sol,F = snes->vec_func;
286: PetscReal fnorm;
287: MatStructure mstruct;
288: PetscInt i;
292: snes->reason = SNES_CONVERGED_ITERATING;
293: PetscObjectTakeAccess(snes);
294: snes->iter = 0;
295: snes->norm = 0.;
296: PetscObjectGrantAccess(snes);
297: if (!snes->vec_func_init_set) {
298: SNESComputeFunction(snes,X,F);
299: if (snes->domainerror) {
300: snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
301: return(0);
302: }
303: } else {
304: snes->vec_func_init_set = PETSC_FALSE;
305: }
306: if (snes->jacobian) { /* This method does not require a Jacobian, but it is usually preconditioned by PBJacobi */
307: SNESComputeJacobian(snes,snes->vec_sol,&snes->jacobian,&snes->jacobian_pre,&mstruct);
308: KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,mstruct);
309: }
310: if (ms->norms) {
311: if (!snes->norm_init_set) {
312: VecNorm(F,NORM_2,&fnorm); /* fnorm <- ||F|| */
313: if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
314: } else {
315: fnorm = snes->norm_init;
316: snes->norm_init_set = PETSC_FALSE;
317: }
318: /* Monitor convergence */
319: PetscObjectTakeAccess(snes);
320: snes->iter = 0;
321: snes->norm = fnorm;
322: PetscObjectGrantAccess(snes);
323: SNESLogConvHistory(snes,snes->norm,0);
324: SNESMonitor(snes,snes->iter,snes->norm);
326: /* set parameter for default relative tolerance convergence test */
327: snes->ttol = fnorm*snes->rtol;
328: /* Test for convergence */
329: (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);
330: if (snes->reason) return(0);
331: }
333: /* Call general purpose update function */
334: if (snes->ops->update) {
335: (*snes->ops->update)(snes,snes->iter);
336: }
337: for (i = 0; i < snes->max_its; i++) {
338: SNESMSStep_3Sstar(snes,X,F);
340: if (i+1 < snes->max_its || ms->norms) {
341: SNESComputeFunction(snes,X,F);
342: if (snes->domainerror) {
343: snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
344: return(0);
345: }
346: }
348: if (ms->norms) {
349: VecNorm(F,NORM_2,&fnorm); /* fnorm <- ||F|| */
350: if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
351: /* Monitor convergence */
352: PetscObjectTakeAccess(snes);
353: snes->iter = i+1;
354: snes->norm = fnorm;
355: PetscObjectGrantAccess(snes);
356: SNESLogConvHistory(snes,snes->norm,0);
357: SNESMonitor(snes,snes->iter,snes->norm);
359: /* Test for convergence */
360: (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);
361: if (snes->reason) return(0);
362: }
364: /* Call general purpose update function */
365: if (snes->ops->update) {
366: (*snes->ops->update)(snes, snes->iter);
367: }
368: }
369: if (!snes->reason) snes->reason = SNES_CONVERGED_ITS;
370: return(0);
371: }
375: static PetscErrorCode SNESSetUp_MS(SNES snes)
376: {
377: SNES_MS * ms = (SNES_MS *)snes->data;
381: if (!ms->tableau) {SNESMSSetType(snes,SNESMSDefault);}
382: SNESDefaultGetWork(snes,3);
383: SNESSetUpMatrices(snes);
384: return(0);
385: }
389: static PetscErrorCode SNESReset_MS(SNES snes)
390: {
393: return(0);
394: }
398: static PetscErrorCode SNESDestroy_MS(SNES snes)
399: {
403: PetscFree(snes->data);
404: PetscObjectComposeFunctionDynamic((PetscObject)snes,"","",PETSC_NULL);
405: return(0);
406: }
410: static PetscErrorCode SNESView_MS(SNES snes,PetscViewer viewer)
411: {
412: PetscBool iascii;
413: PetscErrorCode ierr;
414: SNES_MS *ms = (SNES_MS*)snes->data;
417: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);
418: if (iascii) {
419: SNESMSTableau tab = ms->tableau;
420: PetscViewerASCIIPrintf(viewer," multi-stage method type: %s\n",tab?tab->name:"not yet set");
421: }
422: return(0);
423: }
427: static PetscErrorCode SNESSetFromOptions_MS(SNES snes)
428: {
429: SNES_MS *ms = (SNES_MS*)snes->data;
433: PetscOptionsHead("SNES MS options");
434: {
435: SNESMSTableauLink link;
436: PetscInt count,choice;
437: PetscBool flg;
438: const char **namelist;
439: char mstype[256];
441: PetscStrncpy(mstype,SNESMSDefault,sizeof mstype);
442: for (link=SNESMSTableauList,count=0; link; link=link->next,count++) ;
443: PetscMalloc(count*sizeof(char*),&namelist);
444: for (link=SNESMSTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name;
445: PetscOptionsEList("-snes_ms_type","Multistage smoother type","SNESMSSetType",(const char*const*)namelist,count,mstype,&choice,&flg);
446: SNESMSSetType(snes,flg ? namelist[choice] : mstype);
447: PetscFree(namelist);
448: PetscOptionsReal("-snes_ms_damping","Damping for multistage method","SNESMSSetDamping",ms->damping,&ms->damping,PETSC_NULL);
449: PetscOptionsBool("-snes_ms_norms","Compute norms for monitoring","none",ms->norms,&ms->norms,PETSC_NULL);
450: }
451: PetscOptionsTail();
452: return(0);
453: }
455: EXTERN_C_BEGIN
458: PetscErrorCode SNESMSSetType_MS(SNES snes,SNESMSType mstype)
459: {
460: PetscErrorCode ierr;
461: SNES_MS *ms = (SNES_MS*)snes->data;
462: SNESMSTableauLink link;
463: PetscBool match;
466: if (ms->tableau) {
467: PetscStrcmp(ms->tableau->name,mstype,&match);
468: if (match) return(0);
469: }
470: for (link = SNESMSTableauList; link; link=link->next) {
471: PetscStrcmp(link->tab.name,mstype,&match);
472: if (match) {
473: SNESReset_MS(snes);
474: ms->tableau = &link->tab;
475: return(0);
476: }
477: }
478: SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",mstype);
479: return(0);
480: }
481: EXTERN_C_END
486: /*@C
487: SNESMSSetType - Set the type of multistage smoother
489: Logically collective
491: Input Parameter:
492: + snes - nonlinear solver context
493: - mstype - type of multistage method
495: Level: beginner
497: .seealso: SNESMSGetType(), SNESMS
498: @*/
499: PetscErrorCode SNESMSSetType(SNES snes,const SNESMSType rostype)
500: {
505: PetscTryMethod(snes,"SNESMSSetType_C",(SNES,const SNESMSType),(snes,rostype));
506: return(0);
507: }
509: /* -------------------------------------------------------------------------- */
510: /*MC
511: SNESMS - multi-stage smoothers
513: Options Database:
515: + -snes_ms_type - type of multi-stage smoother
516: - -snes_ms_damping - damping for multi-stage method
518: Notes:
519: These multistage methods are explicit Runge-Kutta methods that are often used as smoothers for
520: FAS multigrid for transport problems. In the linear case, these are equivalent to polynomial smoothers (such as Chebyshev).
522: Multi-stage smoothers should usually be preconditioned by point-block Jacobi to ensure proper scaling and to normalize the wave speeds.
524: The methods are specified in low storage form (Ketcheson 2010). New methods can be registered with SNESMSRegister().
526: References:
528: Ketcheson (2010) Runge-Kutta methods with minimum storage implementations.
530: Jameson (1983) Solution of the Euler equations for two dimensional transonic flow by a multigrid method.
532: Pierce and Giles (1997) Preconditioned multigrid methods for compressible flow calculations on stretched meshes.
534: Level: beginner
536: .seealso: SNESCreate(), SNES, SNESSetType(), SNESMS, SNESFAS, KSPCHEBYSHEV
538: M*/
539: EXTERN_C_BEGIN
542: PetscErrorCode SNESCreate_MS(SNES snes)
543: {
545: SNES_MS *ms;
548: #if !defined(PETSC_USE_DYNAMIC_LIBRARIES)
549: SNESMSInitializePackage(PETSC_NULL);
550: #endif
552: snes->ops->setup = SNESSetUp_MS;
553: snes->ops->solve = SNESSolve_MS;
554: snes->ops->destroy = SNESDestroy_MS;
555: snes->ops->setfromoptions = SNESSetFromOptions_MS;
556: snes->ops->view = SNESView_MS;
557: snes->ops->reset = SNESReset_MS;
559: snes->usespc = PETSC_FALSE;
560: snes->usesksp = PETSC_TRUE;
562: PetscNewLog(snes,SNES_MS,&ms);
563: snes->data = (void*)ms;
564: ms->damping = 0.9;
565: ms->norms = PETSC_FALSE;
567: PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESMSSetType_C","SNESMSSetType_MS",SNESMSSetType_MS);
568: return(0);
569: }
570: EXTERN_C_END