Actual source code: gl.c
1: #define PETSCTS_DLL
3: #include gl.h
4: #include petscblaslapack.h
6: static const char *TSGLErrorDirections[] = {"FORWARD","BACKWARD","TSGLErrorDirection","TSGLERROR_",0};
7: static PetscFList TSGLList;
8: static PetscFList TSGLAcceptList;
9: static PetscTruth TSGLPackageInitialized;
10: static PetscTruth TSGLRegisterAllCalled;
12: /* This function is pure */
13: static PetscScalar Factorial(PetscInt n)
14: {
15: PetscInt i;
16: if (n < 12) { /* Can compute with 32-bit integers */
17: PetscInt f = 1;
18: for (i=2; i<=n; i++) f *= i;
19: return (PetscScalar)f;
20: } else {
21: PetscScalar f = 1.;
22: for (i=2; i<=n; i++) f *= (PetscScalar)i;
23: return f;
24: }
25: }
27: /* This function is pure */
28: static PetscScalar CPowF(PetscScalar c,PetscInt p)
29: {
30: return PetscPowScalar(c,p)/Factorial(p);
31: }
35: static PetscErrorCode TSGLSchemeCreate(PetscInt p,PetscInt q,PetscInt r,PetscInt s,const PetscScalar *c,
36: const PetscScalar *a,const PetscScalar *b,const PetscScalar *u,const PetscScalar *v,TSGLScheme *inscheme)
37: {
38: TSGLScheme scheme;
39: PetscInt j;
43: if (p < 1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Scheme order must be positive");
44: if (r < 1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"At least one item must be carried between steps");
45: if (s < 1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"At least one stage is required");
47: *inscheme = 0;
48: PetscNew(struct _TSGLScheme,&scheme);
49: scheme->p = p;
50: scheme->q = q;
51: scheme->r = r;
52: scheme->s = s;
54: PetscMalloc5(s,PetscScalar,&scheme->c,s*s,PetscScalar,&scheme->a,r*s,PetscScalar,&scheme->b,r*s,PetscScalar,&scheme->u,r*r,PetscScalar,&scheme->v);
55: PetscMemcpy(scheme->c,c,s*sizeof(PetscScalar));
56: for (j=0; j<s*s; j++) scheme->a[j] = (PetscAbsScalar(a[j]) < 1e-12) ? 0 : a[j];
57: for (j=0; j<r*s; j++) scheme->b[j] = (PetscAbsScalar(b[j]) < 1e-12) ? 0 : b[j];
58: for (j=0; j<s*r; j++) scheme->u[j] = (PetscAbsScalar(u[j]) < 1e-12) ? 0 : u[j];
59: for (j=0; j<r*r; j++) scheme->v[j] = (PetscAbsScalar(v[j]) < 1e-12) ? 0 : v[j];
61: PetscMalloc6(r,PetscScalar,&scheme->alpha,r,PetscScalar,&scheme->beta,r,PetscScalar,&scheme->gamma,3*s,PetscScalar,&scheme->phi,3*r,PetscScalar,&scheme->psi,r,PetscScalar,&scheme->stage_error);
62: {
63: PetscInt i,j,k,ss=s+2;
64: PetscBLASInt m,n,one=1,*ipiv,lwork=4*((s+3)*3+3),info,rank,ldb;
65: PetscReal rcond,*sing,*workreal;
66: PetscScalar *ImV,*H,*bmat,*workscalar,*c=scheme->c,*a=scheme->a,*b=scheme->b,*u=scheme->u,*v=scheme->v;
67: PetscMalloc7(PetscSqr(r),PetscScalar,&ImV,3*s,PetscScalar,&H,3*ss,PetscScalar,&bmat,lwork,PetscScalar,&workscalar,5*(3+r),PetscReal,&workreal,r+s,PetscReal,&sing,r+s,PetscBLASInt,&ipiv);
69: /* column-major input */
70: for (i=0; i<r-1; i++) {
71: for (j=0; j<r-1; j++) {
72: ImV[i+j*r] = 1.0*(i==j) - v[(i+1)*r+j+1];
73: }
74: }
75: /* Build right hand side for alpha (tp - glm.B(2:end,:)*(glm.c.^(p)./factorial(p))) */
76: for (i=1; i<r; i++) {
77: scheme->alpha[i] = 1./Factorial(p+1-i);
78: for (j=0; j<s; j++) scheme->alpha[i] -= b[i*s+j]*CPowF(c[j],p);
79: }
80: m = PetscBLASIntCast(r-1);
81: n = PetscBLASIntCast(r);
82: LAPACKgesv_(&m,&one,ImV,&n,ipiv,scheme->alpha+1,&n,&info);
83: if (info < 0) SETERRQ(PETSC_ERR_LIB,"Bad argument to GESV");
84: if (info > 0) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Bad LU factorization");
86: /* Build right hand side for beta (tp1 - glm.B(2:end,:)*(glm.c.^(p+1)./factorial(p+1)) - e.alpha) */
87: for (i=1; i<r; i++) {
88: scheme->beta[i] = 1./Factorial(p+2-i) - scheme->alpha[i];
89: for (j=0; j<s; j++) scheme->beta[i] -= b[i*s+j]*CPowF(c[j],p+1);
90: }
91: LAPACKgetrs_("No transpose",&m,&one,ImV,&n,ipiv,scheme->beta+1,&n,&info);
92: if (info < 0) SETERRQ(PETSC_ERR_LIB,"Bad argument to GETRS");
93: if (info > 0) SETERRQ(PETSC_ERR_LIB,"Should not happen");
95: /* Build stage_error vector
96: xi = glm.c.^(p+1)/factorial(p+1) - glm.A*glm.c.^p/factorial(p) + glm.U(:,2:end)*e.alpha;
97: */
98: for (i=0; i<s; i++) {
99: scheme->stage_error[i] = CPowF(c[i],p+1);
100: for (j=0; j<s; j++) scheme->stage_error[i] -= a[i*s+j]*CPowF(c[j],p);
101: for (j=1; j<r; j++) scheme->stage_error[i] += u[i*r+j]*scheme->alpha[j];
102: }
104: /* alpha[0] (epsilon in B,J,W 2007)
105: epsilon = 1/factorial(p+1) - B(1,:)*c.^p/factorial(p) + V(1,2:end)*e.alpha;
106: */
107: scheme->alpha[0] = 1./Factorial(p+1);
108: for (j=0; j<s; j++) scheme->alpha[0] -= b[0*s+j]*CPowF(c[j],p);
109: for (j=1; j<r; j++) scheme->alpha[0] += v[0*r+j]*scheme->alpha[j];
111: /* right hand side for gamma (glm.B(2:end,:)*e.xi - e.epsilon*eye(s-1,1)) */
112: for (i=1; i<r; i++) {
113: scheme->gamma[i] = (i==1 ? -1. : 0)*scheme->alpha[0];
114: for (j=0; j<s; j++) scheme->gamma[i] += b[i*s+j]*scheme->stage_error[j];
115: }
116: LAPACKgetrs_("No transpose",&m,&one,ImV,&n,ipiv,scheme->gamma+1,&n,&info);
117: if (info < 0) SETERRQ(PETSC_ERR_LIB,"Bad argument to GETRS");
118: if (info > 0) SETERRQ(PETSC_ERR_LIB,"Should not happen");
120: /* beta[0] (rho in B,J,W 2007)
121: e.rho = 1/factorial(p+2) - glm.B(1,:)*glm.c.^(p+1)/factorial(p+1) ...
122: + glm.V(1,2:end)*e.beta;% - e.epsilon;
123: % Note: The paper (B,J,W 2007) includes the last term in their definition
124: * */
125: scheme->beta[0] = 1./Factorial(p+2);
126: for (j=0; j<s; j++) scheme->beta[0] -= b[0*s+j]*CPowF(c[j],p+1);
127: for (j=1; j<r; j++) scheme->beta[0] += v[0*r+j]*scheme->beta[j];
129: /* gamma[0] (sigma in B,J,W 2007)
130: * e.sigma = glm.B(1,:)*e.xi + glm.V(1,2:end)*e.gamma;
131: * */
132: scheme->gamma[0] = 0.0;
133: for (j=0; j<s; j++) scheme->gamma[0] += b[0*s+j]*scheme->stage_error[j];
134: for (j=1; j<r; j++) scheme->gamma[0] += v[0*s+j]*scheme->gamma[j];
136: /* Assemble H
137: * % Determine the error estimators phi
138: H = [[cpow(glm.c,p) + C*e.alpha] [cpow(glm.c,p+1) + C*e.beta] ...
139: [e.xi - C*(e.gamma + 0*e.epsilon*eye(s-1,1))]]';
140: % Paper has formula above without the 0, but that term must be left
141: % out to satisfy the conditions they propose and to make the
142: % example schemes work
143: e.H = H;
144: e.phi = (H \ [1 0 0;1 1 0;0 0 -1])';
145: e.psi = -e.phi*C;
146: * */
147: for (j=0; j<s; j++) {
148: H[0+j*3] = CPowF(c[j],p);
149: H[1+j*3] = CPowF(c[j],p+1);
150: H[2+j*3] = scheme->stage_error[j];
151: for (k=1; k<r; k++) {
152: H[0+j*3] += CPowF(c[j],k-1)*scheme->alpha[k];
153: H[1+j*3] += CPowF(c[j],k-1)*scheme->beta[k];
154: H[2+j*3] -= CPowF(c[j],k-1)*scheme->gamma[k];
155: }
156: }
157: bmat[0+0*ss] = 1.; bmat[0+1*ss] = 0.; bmat[0+2*ss] = 0.;
158: bmat[1+0*ss] = 1.; bmat[1+1*ss] = 1.; bmat[1+2*ss] = 0.;
159: bmat[2+0*ss] = 0.; bmat[2+1*ss] = 0.; bmat[2+2*ss] = -1.;
160: m = 3;
161: n = PetscBLASIntCast(s);
162: ldb = PetscBLASIntCast(ss);
163: rcond = 1e-12;
164: #if defined(PETSC_USE_COMPLEX)
165: /* ZGELSS( M, N, NRHS, A, LDA, B, LDB, S, RCOND, RANK, WORK, LWORK, RWORK, INFO ) */
166: LAPACKgelss_(&m,&n,&m,H,&m,bmat,&ldb,sing,&rcond,&rank,workscalar,&lwork,workreal,&info);
167: #else
168: /* DGELSS( M, N, NRHS, A, LDA, B, LDB, S, RCOND, RANK, WORK, LWORK, INFO ) */
169: LAPACKgelss_(&m,&n,&m,H,&m,bmat,&ldb,sing,&rcond,&rank,workscalar,&lwork,&info);
170: #endif
171: if (info < 0) SETERRQ(PETSC_ERR_LIB,"Bad argument to GELSS");
172: if (info > 0) SETERRQ(PETSC_ERR_LIB,"SVD failed to converge");
174: for (j=0; j<3; j++) {
175: for (k=0; k<s; k++) {
176: scheme->phi[k+j*s] = bmat[k+j*ss];
177: }
178: }
180: /* the other part of the error estimator, psi in B,J,W 2007 */
181: scheme->psi[0*r+0] = 0.;
182: scheme->psi[1*r+0] = 0.;
183: scheme->psi[2*r+0] = 0.;
184: for (j=1; j<r; j++) {
185: scheme->psi[0*r+j] = 0.;
186: scheme->psi[1*r+j] = 0.;
187: scheme->psi[2*r+j] = 0.;
188: for (k=0; k<s; k++) {
189: scheme->psi[0*r+j] -= CPowF(c[k],j-1)*scheme->phi[0*s+k];
190: scheme->psi[1*r+j] -= CPowF(c[k],j-1)*scheme->phi[1*s+k];
191: scheme->psi[2*r+j] -= CPowF(c[k],j-1)*scheme->phi[2*s+k];
192: }
193: }
194: PetscFree7(ImV,H,bmat,workscalar,workreal,sing,ipiv);
195: }
196: /* Check which properties are satisfied */
197: scheme->stiffly_accurate = PETSC_TRUE;
198: if (scheme->c[s-1] != 1.) scheme->stiffly_accurate = PETSC_FALSE;
199: for (j=0; j<s; j++) if (a[(s-1)*s+j] != b[j]) scheme->stiffly_accurate = PETSC_FALSE;
200: for (j=0; j<r; j++) if (u[(s-1)*r+j] != v[j]) scheme->stiffly_accurate = PETSC_FALSE;
201: scheme->fsal = scheme->stiffly_accurate; /* FSAL is stronger */
202: for (j=0; j<s-1; j++) if (r>1 && b[1*s+j] != 0.) scheme->fsal = PETSC_FALSE;
203: if (b[1*s+r-1] != 1.) scheme->fsal = PETSC_FALSE;
204: for (j=0; j<r; j++) if (r>1 && v[1*r+j] != 0.) scheme->fsal = PETSC_FALSE;
206: *inscheme = scheme;
207: return(0);
208: }
212: static PetscErrorCode TSGLSchemeDestroy(TSGLScheme sc)
213: {
217: PetscFree5(sc->c,sc->a,sc->b,sc->u,sc->v);
218: PetscFree6(sc->alpha,sc->beta,sc->gamma,sc->phi,sc->psi,sc->stage_error);
219: PetscFree(sc);
220: return(0);
221: }
225: static PetscErrorCode TSGLDestroy_Default(TS_GL *gl)
226: {
228: PetscInt i;
231: for (i=0; i<gl->nschemes; i++) {
232: if (gl->schemes[i]) {TSGLSchemeDestroy(gl->schemes[i]);}
233: }
234: PetscFree(gl->schemes);
235: gl->schemes = 0;
236: gl->nschemes = 0;
237: PetscMemzero(gl->type_name,sizeof(gl->type_name));
238: return(0);
239: }
243: static PetscErrorCode ViewTable_Private(PetscViewer viewer,PetscInt m,PetscInt n,const PetscScalar a[],const char name[])
244: {
246: PetscTruth iascii;
247: PetscInt i,j;
250: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
251: if (iascii) {
252: PetscViewerASCIIPrintf(viewer,"%30s = [",name);
253: for (i=0; i<m; i++) {
254: if (i) {PetscViewerASCIIPrintf(viewer,"%30s [","");}
255: PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
256: for (j=0; j<n; j++) {
257: PetscViewerASCIIPrintf(viewer," %12.8g",PetscRealPart(a[i*n+j]));
258: }
259: PetscViewerASCIIPrintf(viewer,"]\n");
260: PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
261: }
262: } else {
263: SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for TS_GL",((PetscObject)viewer)->type_name);
264: }
265: return(0);
266: }
271: static PetscErrorCode TSGLSchemeView(TSGLScheme sc,PetscTruth view_details,PetscViewer viewer)
272: {
274: PetscTruth iascii;
277: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
278: if (iascii) {
279: PetscViewerASCIIPrintf(viewer,"GL scheme p,q,r,s = %d,%d,%d,%d\n",sc->p,sc->q,sc->r,sc->s);
280: PetscViewerASCIIPushTab(viewer);
281: PetscViewerASCIIPrintf(viewer,"Stiffly accurate: %s, FSAL: %s\n",sc->stiffly_accurate?"yes":"no",sc->fsal?"yes":"no");
282: PetscViewerASCIIPrintf(viewer,"Leading error constants: %10.3e %10.3e %10.3e\n",
283: PetscRealPart(sc->alpha[0]),PetscRealPart(sc->beta[0]),PetscRealPart(sc->gamma[0]));
284: ViewTable_Private(viewer,1,sc->s,sc->c,"Abscissas c");
285: if (view_details) {
286: ViewTable_Private(viewer,sc->s,sc->s,sc->a,"A");
287: ViewTable_Private(viewer,sc->r,sc->s,sc->b,"B");
288: ViewTable_Private(viewer,sc->s,sc->r,sc->u,"U");
289: ViewTable_Private(viewer,sc->r,sc->r,sc->v,"V");
291: ViewTable_Private(viewer,3,sc->s,sc->phi,"Error estimate phi");
292: ViewTable_Private(viewer,3,sc->r,sc->psi,"Error estimate psi");
293: ViewTable_Private(viewer,1,sc->r,sc->alpha,"Modify alpha");
294: ViewTable_Private(viewer,1,sc->r,sc->beta,"Modify beta");
295: ViewTable_Private(viewer,1,sc->r,sc->gamma,"Modify gamma");
296: ViewTable_Private(viewer,1,sc->s,sc->stage_error,"Stage error xi");
297: }
298: PetscViewerASCIIPopTab(viewer);
299: } else {
300: SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for TS_GL",((PetscObject)viewer)->type_name);
301: }
302: return(0);
303: }
307: static PetscErrorCode TSGLEstimateHigherMoments_Default(TSGLScheme sc,PetscReal h,Vec Ydot[],Vec Xold[],Vec hm[])
308: {
310: PetscInt i;
313: if (sc->r > 64 || sc->s > 64) SETERRQ(PETSC_ERR_PLIB,"Ridiculous number of stages or items passed between stages");
314: /* build error vectors*/
315: for (i=0; i<3; i++) {
316: PetscScalar phih[64];
317: PetscInt j;
318: for (j=0; j<sc->s; j++) phih[j] = sc->phi[i*sc->s+j]*h;
319: VecZeroEntries(hm[i]);
320: VecMAXPY(hm[i],sc->s,phih,Ydot);
321: VecMAXPY(hm[i],sc->r,&sc->psi[i*sc->r],Xold);
322: }
323: return(0);
324: }
328: static PetscErrorCode TSGLCompleteStep_Rescale(TSGLScheme sc,PetscReal h,TSGLScheme next_sc,PetscReal next_h,Vec Ydot[],Vec Xold[],Vec X[])
329: {
331: PetscScalar brow[32],vrow[32];
332: PetscInt i,j,r,s;
335: /* Build the new solution from (X,Ydot) */
336: r = sc->r;
337: s = sc->s;
338: for (i=0; i<r; i++) {
339: VecZeroEntries(X[i]);
340: for (j=0; j<s; j++) brow[j] = h*sc->b[i*s+j];
341: VecMAXPY(X[i],s,brow,Ydot);
342: for (j=0; j<r; j++) vrow[j] = sc->v[i*r+j];
343: VecMAXPY(X[i],r,vrow,Xold);
344: }
345: return(0);
346: }
350: static PetscErrorCode TSGLCompleteStep_RescaleAndModify(TSGLScheme sc,PetscReal h,TSGLScheme next_sc,PetscReal next_h,Vec Ydot[],Vec Xold[],Vec X[])
351: {
353: PetscScalar brow[32],vrow[32];
354: PetscReal ratio;
355: PetscInt i,j,p,r,s;
358: /* Build the new solution from (X,Ydot) */
359: p = sc->p;
360: r = sc->r;
361: s = sc->s;
362: ratio = next_h/h;
363: for (i=0; i<r; i++) {
364: VecZeroEntries(X[i]);
365: for (j=0; j<s; j++) {
366: brow[j] = h*(PetscPowScalar(ratio,i)*sc->b[i*s+j]
367: + (PetscPowScalar(ratio,i) - PetscPowScalar(ratio,p+1))*(+ sc->alpha[i]*sc->phi[0*s+j])
368: + (PetscPowScalar(ratio,i) - PetscPowScalar(ratio,p+2))*(+ sc->beta [i]*sc->phi[1*s+j]
369: + sc->gamma[i]*sc->phi[2*s+j]));
370: }
371: VecMAXPY(X[i],s,brow,Ydot);
372: for (j=0; j<r; j++) {
373: vrow[j] = (PetscPowScalar(ratio,i)*sc->v[i*r+j]
374: + (PetscPowScalar(ratio,i) - PetscPowScalar(ratio,p+1))*(+ sc->alpha[i]*sc->psi[0*r+j])
375: + (PetscPowScalar(ratio,i) - PetscPowScalar(ratio,p+2))*(+ sc->beta [i]*sc->psi[1*r+j]
376: + sc->gamma[i]*sc->psi[2*r+j]));
377: }
378: VecMAXPY(X[i],r,vrow,Xold);
379: }
380: if (r < next_sc->r) {
381: if (r+1 != next_sc->r) SETERRQ(PETSC_ERR_PLIB,"Cannot accommodate jump in r greater than 1");
382: VecZeroEntries(X[r]);
383: for (j=0; j<s; j++) brow[j] = h*PetscPowScalar(ratio,p+1)*sc->phi[0*s+j];
384: VecMAXPY(X[r],s,brow,Ydot);
385: for (j=0; j<r; j++) vrow[j] = PetscPowScalar(ratio,p+1)*sc->psi[0*r+j];
386: VecMAXPY(X[r],r,vrow,Xold);
387: }
388: return(0);
389: }
394: static PetscErrorCode TSGLCreate_IRKS(TS ts)
395: {
396: TS_GL *gl = (TS_GL*)ts->data;
400: gl->Destroy = TSGLDestroy_Default;
401: gl->EstimateHigherMoments = TSGLEstimateHigherMoments_Default;
402: gl->CompleteStep = TSGLCompleteStep_RescaleAndModify;
403: PetscMalloc(10*sizeof(TSGLScheme),&gl->schemes);
404: gl->nschemes = 0;
406: {
407: /* p=1,q=1, r=s=2, A- and L-stable with error estimates of order 2 and 3
408: * Listed in Butcher & Podhaisky 2006. On error estimation in general linear methods for stiff ODE.
409: * irks(0.3,0,[.3,1],[1],1)
410: * Note: can be made to have classical order (not stage order) 2 by replacing 0.3 with 1-sqrt(1/2)
411: * but doing so would sacrifice the error estimator.
412: */
413: const PetscScalar c[2] = {3./10., 1.}
414: ,a[2][2] = {{3./10., 0}, {7./10., 3./10.}}
415: ,b[2][2] = {{7./10., 3./10.}, {0,1}}
416: ,u[2][2] = {{1,0},{1,0}}
417: ,v[2][2] = {{1,0},{0,0}};
418: TSGLSchemeCreate(1,1,2,2,c,*a,*b,*u,*v,&gl->schemes[gl->nschemes++]);
419: }
421: {
422: /* p=q=2, r=s=3: irks(4/9,0,[1:3]/3,[0.33852],1) */
423: /* http://www.math.auckland.ac.nz/~hpod/atlas/i2a.html */
424: const PetscScalar c[3] = {1./3., 2./3., 1}
425: ,a[3][3] = {{4./9. ,0 , 0},
426: {1.03750643704090e+00 , 4./9., 0},
427: {7.67024779410304e-01 , -3.81140216918943e-01, 4./9.}}
428: ,b[3][3] = {{0.767024779410304, -0.381140216918943, 4./9.},
429: {0.000000000000000, 0.000000000000000, 1.000000000000000},
430: {-2.075048385225385, 0.621728385225383, 1.277197204924873}}
431: ,u[3][3] = {{1.0000000000000000, -0.1111111111111109, -0.0925925925925922},
432: {1.0000000000000000, -0.8152842148186744, -0.4199095530877056},
433: {1.0000000000000000, 0.1696709930641948, 0.0539741070314165}}
434: ,v[3][3] = {{1.0000000000000000, 0.1696709930641948, 0.0539741070314165},
435: {0.000000000000000, 0.000000000000000, 0.000000000000000},
436: {0.000000000000000, 0.176122795075129, 0.000000000000000}};
437: TSGLSchemeCreate(2,2,3,3,c,*a,*b,*u,*v,&gl->schemes[gl->nschemes++]);
438: }
439: {
440: /* p=q=3, r=s=4: irks(9/40,0,[1:4]/4,[0.3312 1.0050],[0.49541 1;1 0]) */
441: const PetscScalar c[4] = {0.25,0.5,0.75,1.0}
442: ,a[4][4] = {{9./40. , 0, 0, 0},
443: {2.11286958887701e-01 , 9./40. , 0, 0},
444: {9.46338294287584e-01 , -3.42942861246094e-01, 9./40. , 0},
445: {0.521490453970721 , -0.662474225622980, 0.490476425459734, 9./40. }}
446: ,b[4][4] = {{0.521490453970721 , -0.662474225622980, 0.490476425459734, 9./40. },
447: {0.000000000000000 , 0.000000000000000, 0.000000000000000, 1.000000000000000},
448: {-0.084677029310348 , 1.390757514776085, -1.568157386206001, 2.023192696767826},
449: {0.465383797936408 , 1.478273530625148, -1.930836081010182, 1.644872111193354}}
450: ,u[4][4] = {{1.00000000000000000 , 0.02500000000001035, -0.02499999999999053, -0.00442708333332865},
451: {1.00000000000000000 , 0.06371304111232945, -0.04032173972189845, -0.01389438413189452},
452: {1.00000000000000000 , -0.07839543304147778, 0.04738685705116663, 0.02032603595928376},
453: {1.00000000000000000 , 0.42550734619251651, 0.10800718022400080, -0.01726712647760034}}
454: ,v[4][4] = {{1.00000000000000000 , 0.42550734619251651, 0.10800718022400080, -0.01726712647760034},
455: {0.000000000000000 , 0.000000000000000, 0.000000000000000, 0.000000000000000},
456: {0.000000000000000 , -1.761115796027561, -0.521284157173780, 0.258249384305463},
457: {0.000000000000000 , -1.657693358744728, -1.052227765232394, 0.521284157173780}};
458: TSGLSchemeCreate(3,3,4,4,c,*a,*b,*u,*v,&gl->schemes[gl->nschemes++]);
459: }
460: {
461: /* p=q=4, r=s=5:
462: irks(3/11,0,[1:5]/5, [0.1715 -0.1238 0.6617],...
463: [ -0.0812 0.4079 1.0000
464: 1.0000 0 0
465: 0.8270 1.0000 0])
466: */
467: const PetscScalar c[5] = {0.2,0.4,0.6,0.8,1.0}
468: ,a[5][5] = {{2.72727272727352e-01 , 0.00000000000000e+00, 0.00000000000000e+00 , 0.00000000000000e+00 , 0.00000000000000e+00},
469: {-1.03980153733431e-01, 2.72727272727405e-01, 0.00000000000000e+00, 0.00000000000000e+00 , 0.00000000000000e+00},
470: {-1.58615400341492e+00, 7.44168951881122e-01, 2.72727272727309e-01, 0.00000000000000e+00 , 0.00000000000000e+00},
471: {-8.73658042865628e-01, 5.37884671894595e-01, -1.63298538799523e-01, 2.72727272726996e-01 , 0.00000000000000e+00},
472: {2.95489397443992e-01 , -1.18481693910097e+00 , -6.68029812659953e-01 , 1.00716687860943e+00 , 2.72727272727288e-01}}
473: ,b[5][5] = {{2.95489397443992e-01 , -1.18481693910097e+00 , -6.68029812659953e-01 , 1.00716687860943e+00 , 2.72727272727288e-01},
474: {0.00000000000000e+00 , 1.11022302462516e-16 , -2.22044604925031e-16 , 0.00000000000000e+00 , 1.00000000000000e+00},
475: {-4.05882503986005e+00, -4.00924006567769e+00, -1.38930610972481e+00, 4.45223930308488e+00 , 6.32331093108427e-01},
476: {8.35690179937017e+00 , -2.26640927349732e+00 , 6.86647884973826e+00 , -5.22595158025740e+00 , 4.50893068837431e+00},
477: {1.27656267027479e+01 , 2.80882153840821e+00 , 8.91173096522890e+00 , -1.07936444078906e+01 , 4.82534148988854e+00}}
478: ,u[5][5] = {{1.00000000000000e+00 , -7.27272727273551e-02 , -3.45454545454419e-02 , -4.12121212119565e-03 ,-2.96969696964014e-04},
479: {1.00000000000000e+00 , 2.31252881006154e-01 , -8.29487834416481e-03 , -9.07191207681020e-03 ,-1.70378403743473e-03},
480: {1.00000000000000e+00 , 1.16925777880663e+00 , 3.59268562942635e-02 , -4.09013451730615e-02 ,-1.02411119670164e-02},
481: {1.00000000000000e+00 , 1.02634463704356e+00 , 1.59375044913405e-01 , 1.89673015035370e-03 ,-4.89987231897569e-03},
482: {1.00000000000000e+00 , 1.27746320298021e+00 , 2.37186008132728e-01 , -8.28694373940065e-02 ,-5.34396510196430e-02}}
483: ,v[5][5] = {{1.00000000000000e+00 , 1.27746320298021e+00 , 2.37186008132728e-01 , -8.28694373940065e-02 ,-5.34396510196430e-02},
484: {0.00000000000000e+00 , -1.77635683940025e-15 , -1.99840144432528e-15 , -9.99200722162641e-16 ,-3.33066907387547e-16},
485: {0.00000000000000e+00 , 4.37280081906924e+00 , 5.49221645016377e-02 , -8.88913177394943e-02 , 1.12879077989154e-01},
486: {0.00000000000000e+00 , -1.22399504837280e+01 , -5.21287338448645e+00 , -8.03952325565291e-01 , 4.60298678047147e-01},
487: {0.00000000000000e+00 , -1.85178762883829e+01 , -5.21411849862624e+00 , -1.04283436528809e+00 , 7.49030161063651e-01}};
488: TSGLSchemeCreate(4,4,5,5,c,*a,*b,*u,*v,&gl->schemes[gl->nschemes++]);
489: }
490: {
491: /* p=q=5, r=s=6;
492: irks(1/3,0,[1:6]/6,...
493: [-0.0489 0.4228 -0.8814 0.9021],...
494: [-0.3474 -0.6617 0.6294 0.2129
495: 0.0044 -0.4256 -0.1427 -0.8936
496: -0.8267 0.4821 0.1371 -0.2557
497: -0.4426 -0.3855 -0.7514 0.3014])
498: */
499: const PetscScalar c[6] = {1./6, 2./6, 3./6, 4./6, 5./6, 1.}
500: ,a[6][6] = {{ 3.33333333333940e-01, 0 , 0 , 0 , 0 , 0 },
501: { -8.64423857333350e-02, 3.33333333332888e-01, 0 , 0 , 0 , 0 },
502: { -2.16850174258252e+00, -2.23619072028839e+00, 3.33333333335204e-01, 0 , 0 , 0 },
503: { -4.73160970138997e+00, -3.89265344629268e+00, -2.76318716520933e-01, 3.33333333335759e-01, 0 , 0 },
504: { -6.75187540297338e+00, -7.90756533769377e+00, 7.90245051802259e-01, -4.48352364517632e-01, 3.33333333328483e-01, 0 },
505: { -4.26488287921548e+00, -1.19320395589302e+01, 3.38924509887755e+00, -2.23969848002481e+00, 6.62807710124007e-01, 3.33333333335440e-01}}
506: ,b[6][6] = {{ -4.26488287921548e+00, -1.19320395589302e+01, 3.38924509887755e+00, -2.23969848002481e+00, 6.62807710124007e-01, 3.33333333335440e-01},
507: { -8.88178419700125e-16, 4.44089209850063e-16, -1.54737334057131e-15, -8.88178419700125e-16, 0.00000000000000e+00, 1.00000000000001e+00},
508: { -2.87780425770651e+01, -1.13520448264971e+01, 2.62002318943161e+01, 2.56943874812797e+01, -3.06702268304488e+01, 6.68067773510103e+00},
509: { 5.47971245256474e+01, 6.80366875868284e+01, -6.50952588861999e+01, -8.28643975339097e+01, 8.17416943896414e+01, -1.17819043489036e+01},
510: { -2.33332114788869e+02, 6.12942539462634e+01, -4.91850135865944e+01, 1.82716844135480e+02, -1.29788173979395e+02, 3.09968095651099e+01},
511: { -1.72049132343751e+02, 8.60194713593999e+00, 7.98154219170200e-01, 1.50371386053218e+02, -1.18515423962066e+02, 2.50898277784663e+01}}
512: ,u[6][6] = {{ 1.00000000000000e+00, -1.66666666666870e-01, -4.16666666664335e-02, -3.85802469124815e-03, -2.25051440302250e-04, -9.64506172339142e-06},
513: { 1.00000000000000e+00, 8.64423857327162e-02, -4.11484912671353e-02, -1.11450903217645e-02, -1.47651050487126e-03, -1.34395070766826e-04},
514: { 1.00000000000000e+00, 4.57135912953434e+00, 1.06514719719137e+00, 1.33517564218007e-01, 1.11365952968659e-02, 6.12382756769504e-04},
515: { 1.00000000000000e+00, 9.23391519753404e+00, 2.22431212392095e+00, 2.91823807741891e-01, 2.52058456411084e-02, 1.22800542949647e-03},
516: { 1.00000000000000e+00, 1.48175480533865e+01, 3.73439117461835e+00, 5.14648336541804e-01, 4.76430038853402e-02, 2.56798515502156e-03},
517: { 1.00000000000000e+00, 1.50512347758335e+01, 4.10099701165164e+00, 5.66039141003603e-01, 3.91213893800891e-02, -2.99136269067853e-03}}
518: ,v[6][6] = {{ 1.00000000000000e+00, 1.50512347758335e+01, 4.10099701165164e+00, 5.66039141003603e-01, 3.91213893800891e-02, -2.99136269067853e-03},
519: { 0.00000000000000e+00, -4.88498130835069e-15, -6.43929354282591e-15, -3.55271367880050e-15, -1.22124532708767e-15, -3.12250225675825e-16},
520: { 0.00000000000000e+00, 1.22250171233141e+01, -1.77150760606169e+00, 3.54516769879390e-01, 6.22298845883398e-01, 2.31647447450276e-01},
521: { 0.00000000000000e+00, -4.48339457331040e+01, -3.57363126641880e-01, 5.18750173123425e-01, 6.55727990241799e-02, 1.63175368287079e-01},
522: { 0.00000000000000e+00, 1.37297394708005e+02, -1.60145272991317e+00, -5.05319555199441e+00, 1.55328940390990e-01, 9.16629423682464e-01},
523: { 0.00000000000000e+00, 1.05703241119022e+02, -1.16610260983038e+00, -2.99767252773859e+00, -1.13472315553890e-01, 1.09742849254729e+00}};
524: TSGLSchemeCreate(5,5,6,6,c,*a,*b,*u,*v,&gl->schemes[gl->nschemes++]);
525: }
526: return(0);
527: }
533: /*@C
534: TSGLSetType - sets the class of general linear method to use for time-stepping
536: Collective on TS
538: Input Parameters:
539: + ts - the TS context
540: - type - a method
542: Options Database Key:
543: . -ts_gl_type <type> - sets the method, use -help for a list of available method (e.g. irks)
545: Notes:
546: See "petsc/include/petscts.h" for available methods (for instance)
547: . TSGL_IRKS - Diagonally implicit methods with inherent Runge-Kutta stability (for stiff problems)
549: Normally, it is best to use the TSSetFromOptions() command and
550: then set the TSGL type from the options database rather than by using
551: this routine. Using the options database provides the user with
552: maximum flexibility in evaluating the many different solvers.
553: The TSGLSetType() routine is provided for those situations where it
554: is necessary to set the timestepping solver independently of the
555: command line or options database. This might be the case, for example,
556: when the choice of solver changes during the execution of the
557: program, and the user's application is taking responsibility for
558: choosing the appropriate method.
560: Level: intermediate
562: .keywords: TS, TSGL, set, type
563: @*/
564: PetscErrorCode TSGLSetType(TS ts,const TSGLType type)
565: {
566: PetscErrorCode ierr,(*r)(TS,const TSGLType);
572: PetscObjectQueryFunction((PetscObject)ts,"TSGLSetType_C",(void(**)(void))&r);
573: if (r) {
574: (*r)(ts,type);
575: }
576: return(0);
577: }
581: /*@C
582: TSGLSetAcceptType - sets the acceptance test
584: Time integrators that need to control error must have the option to reject a time step based on local error
585: estimates. This function allows different schemes to be set.
587: Collective on TS
589: Input Parameters:
590: + ts - the TS context
591: - type - the type
593: Options Database Key:
594: . -ts_gl_accept_type <type> - sets the method used to determine whether to accept or reject a step
596: Level: intermediate
598: .seealso: TS, TSGL, TSGLAcceptRegisterDynamic(), TSGLAdapt, set type
599: @*/
600: PetscErrorCode TSGLSetAcceptType(TS ts,const TSGLAcceptType type)
601: {
602: PetscErrorCode ierr,(*r)(TS,const TSGLAcceptType);
608: PetscObjectQueryFunction((PetscObject)ts,"TSGLSetAcceptType_C",(void(**)(void))&r);
609: if (r) {
610: (*r)(ts,type);
611: }
612: return(0);
613: }
617: /*@C
618: TSGLGetAdapt - gets the TSGLAdapt object from the TS
620: Not Collective
622: Input Parameter:
623: . ts - the TS context
625: Output Parameter:
626: . adapt - the TSGLAdapt context
628: Notes:
629: This allows the user set options on the TSGLAdapt object. Usually it is better to do this using the options
630: database, so this function is rarely needed.
632: Level: advanced
634: .seealso: TSGLAdapt, TSGLAdaptRegisterDynamic()
635: @*/
636: PetscErrorCode TSGLGetAdapt(TS ts,TSGLAdapt *adapt)
637: {
638: PetscErrorCode ierr,(*r)(TS,TSGLAdapt*);
644: PetscObjectQueryFunction((PetscObject)ts,"TSGLGetAdapt_C",(void(**)(void))&r);
645: if (r) {
646: (*r)(ts,adapt);
647: }
648: return(0);
649: }
653: static PetscErrorCode TSGLAccept_Always(TS ts,PetscReal tleft,PetscReal h,const PetscReal enorms[],PetscTruth *accept)
654: {
656: *accept = PETSC_TRUE;
657: return(0);
658: }
662: static PetscErrorCode TSGLUpdateWRMS(TS ts)
663: {
664: TS_GL *gl = (TS_GL*)ts->data;
666: PetscScalar *x,*w;
667: PetscInt n,i;
670: VecGetArray(gl->X[0],&x);
671: VecGetArray(gl->W,&w);
672: VecGetLocalSize(gl->W,&n);
673: for (i=0; i<n; i++) {
674: w[i] = 1./(gl->wrms_atol + gl->wrms_rtol*PetscAbsScalar(x[i]));
675: }
676: VecRestoreArray(gl->X[0],&x);
677: VecRestoreArray(gl->W,&w);
678: return(0);
679: }
683: static PetscErrorCode TSGLVecNormWRMS(TS ts,Vec X,PetscReal *nrm)
684: {
685: TS_GL *gl = (TS_GL*)ts->data;
687: PetscScalar *x,*w,sum = 0.0;
688: PetscInt n,i;
691: VecGetArray(X,&x);
692: VecGetArray(gl->W,&w);
693: VecGetLocalSize(gl->W,&n);
694: for (i=0; i<n; i++) {
695: sum += PetscAbsScalar(PetscSqr(x[i]*w[i]));
696: }
697: VecRestoreArray(X,&x);
698: VecRestoreArray(gl->W,&w);
699: *nrm = PetscAbsScalar(PetscSqrtScalar(sum/(1.*n)));
700: return(0);
701: }
705: static PetscErrorCode TSGLSetType_GL(TS ts,const TSGLType type)
706: {
707: PetscErrorCode ierr,(*r)(TS);
708: PetscTruth same;
709: TS_GL *gl = (TS_GL*)ts->data;
712: if (gl->type_name[0]) {
713: PetscStrcmp(gl->type_name,type,&same);
714: if (same) return(0);
715: (*gl->Destroy)(gl);
716: }
718: PetscFListFind(TSGLList,((PetscObject)ts)->comm,type,(void(**)(void))&r);
719: if (!r) SETERRQ1(PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown TSGL type \"%s\" given",type);
720: (*r)(ts);
721: PetscStrcpy(gl->type_name,type);
722: return(0);
723: }
727: static PetscErrorCode TSGLSetAcceptType_GL(TS ts,const TSGLAcceptType type)
728: {
730: TSGLAcceptFunction r;
731: TS_GL *gl = (TS_GL*)ts->data;
734: PetscFListFind(TSGLAcceptList,((PetscObject)ts)->comm,type,(void(**)(void))&r);
735: if (!r) SETERRQ1(PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown TSGLAccept type \"%s\" given",type);
736: gl->Accept = r;
737: PetscStrncpy(gl->accept_name,type,sizeof(gl->accept_name));
738: return(0);
739: }
743: static PetscErrorCode TSGLGetAdapt_GL(TS ts,TSGLAdapt *adapt)
744: {
746: TS_GL *gl = (TS_GL*)ts->data;
749: if (!gl->adapt) {
750: TSGLAdaptCreate(((PetscObject)ts)->comm,&gl->adapt);
751: PetscObjectIncrementTabLevel((PetscObject)gl->adapt,(PetscObject)ts,1);
752: PetscLogObjectParent(ts,gl->adapt);
753: }
754: *adapt = gl->adapt;
755: return(0);
756: }
760: static PetscErrorCode TSGLChooseNextScheme(TS ts,PetscReal h,const PetscReal hmnorm[],PetscInt *next_scheme,PetscReal *next_h,PetscTruth *finish)
761: {
763: TS_GL *gl = (TS_GL*)ts->data;
764: PetscInt i,n,cur_p,cur,next_sc,candidates[64],orders[64];
765: PetscReal errors[64],costs[64],tleft;
768: cur = -1;
769: cur_p = gl->schemes[gl->current_scheme]->p;
770: tleft = ts->max_time - (ts->ptime + ts->time_step);
771: for (i=0,n=0; i<gl->nschemes; i++) {
772: TSGLScheme sc = gl->schemes[i];
773: if (sc->p < gl->min_order || gl->max_order < sc->p) continue;
774: if (sc->p == cur_p - 1) {
775: errors[n] = PetscAbsScalar(sc->alpha[0])*hmnorm[0];
776: } else if (sc->p == cur_p) {
777: errors[n] = PetscAbsScalar(sc->alpha[0])*hmnorm[1];
778: } else if (sc->p == cur_p+1) {
779: errors[n] = PetscAbsScalar(sc->alpha[0])*(hmnorm[2]+hmnorm[3]);
780: } else continue;
781: candidates[n] = i;
782: orders[n] = PetscMin(sc->p,sc->q); /* order of global truncation error */
783: costs[n] = sc->s; /* estimate the cost as the number of stages */
784: if (i == gl->current_scheme) cur = n;
785: n++;
786: }
787: if (cur < 0 || gl->nschemes <= cur) SETERRQ(PETSC_ERR_PLIB,"Current scheme not found in scheme list");
788: TSGLAdaptChoose(gl->adapt,n,orders,errors,costs,cur,h,tleft,&next_sc,next_h,finish);
789: *next_scheme = candidates[next_sc];
790: PetscInfo7(ts,"Adapt chose scheme %d (%d,%d,%d,%d) with step size %6.2e, finish=%d\n",*next_scheme,gl->schemes[*next_scheme]->p,gl->schemes[*next_scheme]->q,gl->schemes[*next_scheme]->r,gl->schemes[*next_scheme]->s,*next_h,*finish);
791: return(0);
792: }
796: static PetscErrorCode TSGLGetMaxSizes(TS ts,PetscInt *max_r,PetscInt *max_s)
797: {
798: TS_GL *gl = (TS_GL*)ts->data;
801: *max_r = gl->schemes[gl->nschemes-1]->r;
802: *max_s = gl->schemes[gl->nschemes-1]->s;
803: return(0);
804: }
808: static PetscErrorCode TSStep_GL(TS ts,PetscInt *steps,PetscReal *ptime)
809: {
810: TS_GL *gl = (TS_GL*)ts->data;
811: PetscInt i,k,max_steps = ts->max_steps,its,lits,max_r,max_s;
812: PetscTruth final_step,finish;
816: *steps = -ts->steps;
817: *ptime = ts->ptime;
819: TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);
821: TSGLGetMaxSizes(ts,&max_r,&max_s);
822: VecCopy(ts->vec_sol,gl->X[0]);
823: for (i=1; i<max_r; i++) {
824: VecZeroEntries(gl->X[i]);
825: }
826: TSGLUpdateWRMS(ts);
828: if (0) {
829: /* Find consistent initial data for DAE */
830: gl->stage_time = ts->ptime + ts->time_step;
831: gl->shift = 1./ts->time_step;
832: gl->stage = 0;
833: VecCopy(ts->vec_sol,gl->Z);
834: VecCopy(ts->vec_sol,gl->Y);
835: SNESSolve(ts->snes,PETSC_NULL,gl->Y);
836: SNESGetIterationNumber(ts->snes,&its);
837: SNESGetLinearSolveIterations(ts->snes,&lits);
838: ts->nonlinear_its += its; ts->linear_its += lits;
839: }
841: if (gl->current_scheme < 0) SETERRQ(PETSC_ERR_ORDER,"A starting scheme has not been provided");
843: for (k=0,final_step=PETSC_FALSE,finish=PETSC_FALSE; k<max_steps && !finish; k++) {
844: PetscInt j,r,s,next_scheme = 0,rejections;
845: PetscReal h,hmnorm[4],enorm[3],next_h;
846: PetscTruth accept;
847: const PetscScalar *c,*a,*b,*u,*v;
848: Vec *X,*Ydot,Y;
849: TSGLScheme scheme = gl->schemes[gl->current_scheme];
851: r = scheme->r; s = scheme->s;
852: c = scheme->c;
853: a = scheme->a; u = scheme->u;
854: b = scheme->b; v = scheme->v;
855: h = ts->time_step;
856: X = gl->X; Ydot = gl->Ydot; Y = gl->Y;
858: if (ts->ptime > ts->max_time) break;
860: /*
861: We only call PreStep at the start of each STEP, not each STAGE. This is because it is
862: possible to fail (have to restart a step) after multiple stages.
863: */
864: TSPreStep(ts);
866: rejections = 0;
867: while (1) {
868: for (i=0; i<s; i++) {
869: PetscScalar shift = gl->shift = 1./PetscRealPart(h*a[i*s+i]);
870: gl->stage = i;
871: gl->stage_time = ts->ptime + PetscRealPart(c[i])*h;
873: /*
874: * Stage equation: Y = h A Y' + U X
875: * We assume that A is lower-triangular so that we can solve the stages (Y,Y') sequentially
876: * Build the affine vector z_i = -[1/(h a_ii)](h sum_j a_ij y'_j + sum_j u_ij x_j)
877: * Then y'_i = z + 1/(h a_ii) y_i
878: */
879: VecZeroEntries(gl->Z);
880: for (j=0; j<r; j++) {
881: VecAXPY(gl->Z,-shift*u[i*r+j],X[j]);
882: }
883: for (j=0; j<i; j++) {
884: VecAXPY(gl->Z,-shift*h*a[i*s+j],Ydot[j]);
885: }
886: /* Note: Z is used within function evaluation, Ydot = Z + shift*Y */
888: /* Compute an estimate of Y to start Newton iteration */
889: if (gl->extrapolate) {
890: if (i==0) {
891: /* Linear extrapolation on the first stage */
892: VecWAXPY(Y,c[i]*h,X[1],X[0]);
893: } else {
894: /* Linear extrapolation from the last stage */
895: VecAXPY(Y,(c[i]-c[i-1])*h,Ydot[i-1]);
896: }
897: } else if (i==0) { /* Directly use solution from the last step, otherwise reuse the last stage (do nothing) */
898: VecCopy(X[0],Y);
899: }
901: /* Solve this stage (Ydot[i] is computed during function evaluation) */
902: SNESSolve(ts->snes,PETSC_NULL,Y);
903: SNESGetIterationNumber(ts->snes,&its);
904: SNESGetLinearSolveIterations(ts->snes,&lits);
905: ts->nonlinear_its += its; ts->linear_its += lits;
906: }
908: gl->stage_time = ts->ptime + ts->time_step;
910: (*gl->EstimateHigherMoments)(scheme,h,Ydot,gl->X,gl->himom);
911: /* hmnorm[i] = h^{p+i}x^{(p+i)} with i=0,1,2; hmnorm[3] = h^{p+2}(dx'/dx) x^{(p+1)} */
912: for (i=0; i<3; i++) {
913: TSGLVecNormWRMS(ts,gl->himom[i],&hmnorm[i+1]);
914: }
915: enorm[0] = PetscRealPart(scheme->alpha[0])*hmnorm[1];
916: enorm[1] = PetscRealPart(scheme->beta[0]) *hmnorm[2];
917: enorm[2] = PetscRealPart(scheme->gamma[0])*hmnorm[3];
918: (*gl->Accept)(ts,ts->max_time-gl->stage_time,h,enorm,&accept);
919: if (accept) goto accepted;
920: rejections++;
921: PetscInfo3(ts,"Step %D (t=%g) not accepted, rejections=%D\n",k,gl->stage_time,rejections);
922: if (rejections > gl->max_step_rejections) break;
923: /*
924: There are lots of reasons why a step might be rejected, including solvers not converging and other factors that
925: TSGLChooseNextScheme does not support. Additionally, the error estimates may be very screwed up, so I'm not
926: convinced that it's safe to just compute a new error estimate using the same interface as the current adaptor
927: (the adaptor interface probably has to change). Here we make an arbitrary and naive choice. This assumes that
928: steps were written in Nordsieck form. The "correct" method would be to re-complete the previous time step with
929: the correct "next" step size. It is unclear to me whether the present ad-hoc method of rescaling X is stable.
930: */
931: h *= 0.5;
932: for (i=1; i<scheme->r; i++) {
933: VecScale(X[i],PetscPowScalar(0.5,i));
934: }
935: }
936: SETERRQ3(PETSC_ERR_CONV_FAILED,"Time step %D (t=%g) not accepted after %D failures\n",k,gl->stage_time,rejections);
938: accepted:
939: /* This term is not error, but it *would* be the leading term for a lower order method */
940: TSGLVecNormWRMS(ts,gl->X[scheme->r-1],&hmnorm[0]);
941: /* Correct scaling so that these are equivalent to norms of the Nordsieck vectors */
943: PetscInfo4(ts,"Last moment norm %10.2e, estimated error norms %10.2e %10.2e %10.2e\n",hmnorm[0],enorm[0],enorm[1],enorm[2]);
944: if (!final_step) {
945: TSGLChooseNextScheme(ts,h,hmnorm,&next_scheme,&next_h,&final_step);
946: } else {
947: /* Dummy values to complete the current step in a consistent manner */
948: next_scheme = gl->current_scheme;
949: next_h = h;
950: finish = PETSC_TRUE;
951: }
953: X = gl->Xold;
954: gl->Xold = gl->X;
955: gl->X = X;
956: (*gl->CompleteStep)(scheme,h,gl->schemes[next_scheme],next_h,Ydot,gl->Xold,gl->X);
958: TSGLUpdateWRMS(ts);
960: /* Post the solution for the user, we could avoid this copy with a small bit of cleverness */
961: VecCopy(gl->X[0],ts->vec_sol);
962: ts->ptime += h;
963: ts->steps++;
965: TSPostStep(ts);
966: TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);
968: gl->current_scheme = next_scheme;
969: ts->time_step = next_h;
970: }
972: *steps += ts->steps;
973: *ptime = ts->ptime;
974: return(0);
975: }
977: /*------------------------------------------------------------*/
980: static PetscErrorCode TSDestroy_GL(TS ts)
981: {
982: TS_GL *gl = (TS_GL*)ts->data;
983: PetscInt max_r,max_s;
984: PetscErrorCode ierr;
987: if (gl->setupcalled) {
988: TSGLGetMaxSizes(ts,&max_r,&max_s);
989: VecDestroyVecs(gl->Xold,max_r);
990: VecDestroyVecs(gl->X,max_r);
991: VecDestroyVecs(gl->Ydot,max_s);
992: VecDestroyVecs(gl->himom,3);
993: VecDestroy(gl->W);
994: VecDestroy(gl->Y);
995: VecDestroy(gl->Z);
996: }
997: if (gl->adapt) {TSGLAdaptDestroy(gl->adapt);}
998: if (gl->Destroy) {(*gl->Destroy)(gl);}
999: PetscFree(gl);
1000: return(0);
1001: }
1003: /*
1004: This defines the nonlinear equation that is to be solved with SNES
1005: g(x) = f(t,x,z+shift*x) = 0
1006: */
1009: static PetscErrorCode TSGLFunction(SNES snes,Vec x,Vec f,void *ctx)
1010: {
1011: TS ts = (TS)ctx;
1012: TS_GL *gl = (TS_GL*)ts->data;
1013: PetscErrorCode ierr;
1016: VecWAXPY(gl->Ydot[gl->stage],gl->shift,x,gl->Z);
1017: TSComputeIFunction(ts,gl->stage_time,x,gl->Ydot[gl->stage],f);
1018: return(0);
1019: }
1023: static PetscErrorCode TSGLJacobian(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *str,void *ctx)
1024: {
1025: TS ts = (TS)ctx;
1026: TS_GL *gl = (TS_GL*)ts->data;
1027: PetscErrorCode ierr;
1030: /* gl->Xdot will have already been computed in TSGLFunction */
1031: TSComputeIJacobian(ts,gl->stage_time,x,gl->Ydot[gl->stage],gl->shift,A,B,str);
1032: return(0);
1033: }
1038: static PetscErrorCode TSSetUp_GL(TS ts)
1039: {
1040: TS_GL *gl = (TS_GL*)ts->data;
1041: Vec res;
1042: PetscInt max_r,max_s;
1043: PetscErrorCode ierr;
1046: if (ts->problem_type == TS_LINEAR) {
1047: SETERRQ(PETSC_ERR_ARG_WRONG,"Only for nonlinear problems");
1048: }
1049: gl->setupcalled = PETSC_TRUE;
1050: TSGLGetMaxSizes(ts,&max_r,&max_s);
1051: VecDuplicateVecs(ts->vec_sol,max_r,&gl->X);
1052: VecDuplicateVecs(ts->vec_sol,max_r,&gl->Xold);
1053: VecDuplicateVecs(ts->vec_sol,max_s,&gl->Ydot);
1054: VecDuplicateVecs(ts->vec_sol,3,&gl->himom);
1055: VecDuplicate(ts->vec_sol,&gl->W);
1056: VecDuplicate(ts->vec_sol,&gl->Y);
1057: VecDuplicate(ts->vec_sol,&gl->Z);
1058: VecDuplicate(ts->vec_sol,&res);
1059: SNESSetFunction(ts->snes,res,&TSGLFunction,ts);
1060: VecDestroy(res); /* Give ownership to SNES */
1061: /* This is nasty. SNESSetFromOptions() is usually called in TSSetFromOptions(). With -snes_mf_operator, it will
1062: replace A and we don't want to mess with that. With -snes_mf, A and B will be replaced as well as the function and
1063: context. Note that SNESSetFunction() normally has not been called before SNESSetFromOptions(), so when -snes_mf sets
1064: the Jacobian user context to snes->funP, it will actually be NULL. This is not a problem because both snes->funP and
1065: snes->jacP should be the TS. */
1066: {
1067: Mat A,B;
1068: PetscErrorCode (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*);
1069: void *ctx;
1070: SNESGetJacobian(ts->snes,&A,&B,&func,&ctx);
1071: SNESSetJacobian(ts->snes,A?A:ts->A,B?B:ts->B,func?func:&TSGLJacobian,ctx?ctx:ts);
1072: }
1074: /* Default acceptance tests and adaptivity */
1075: if (!gl->Accept) {TSGLSetAcceptType(ts,TSGLACCEPT_ALWAYS);}
1076: if (!gl->adapt) {TSGLGetAdapt(ts,&gl->adapt);}
1078: if (gl->current_scheme < 0) {
1079: PetscInt i;
1080: for (i=0; ; i++) {
1081: if (gl->schemes[i]->p == gl->start_order) break;
1082: if (i+1 == gl->nschemes) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"No schemes available with requested start order %d",i);
1083: }
1084: gl->current_scheme = i;
1085: }
1086: return(0);
1087: }
1088: /*------------------------------------------------------------*/
1092: static PetscErrorCode TSSetFromOptions_GL(TS ts)
1093: {
1094: TS_GL *gl = (TS_GL*)ts->data;
1095: char tname[256] = TSGL_IRKS,completef[256] = "rescale-and-modify";
1099: PetscOptionsHead("General Linear ODE solver options");
1100: {
1101: PetscTruth flg;
1102: PetscOptionsList("-ts_gl_type","Type of GL method","TSGLSetType",TSGLList,gl->type_name[0]?gl->type_name:tname,tname,sizeof(tname),&flg);
1103: if (flg || !gl->type_name[0]) {
1104: TSGLSetType(ts,tname);
1105: }
1106: PetscOptionsInt("-ts_gl_max_step_rejections","Maximum number of times to attempt a step","None",gl->max_step_rejections,&gl->max_step_rejections,PETSC_NULL);
1107: PetscOptionsInt("-ts_gl_max_order","Maximum order to try","TSGLSetMaxOrder",gl->max_order,&gl->max_order,PETSC_NULL);
1108: PetscOptionsInt("-ts_gl_min_order","Minimum order to try","TSGLSetMinOrder",gl->min_order,&gl->min_order,PETSC_NULL);
1109: PetscOptionsInt("-ts_gl_start_order","Initial order to try","TSGLSetMinOrder",gl->start_order,&gl->start_order,PETSC_NULL);
1110: PetscOptionsEnum("-ts_gl_error_direction","Which direction to look when estimating error","TSGLSetErrorDirection",TSGLErrorDirections,(PetscEnum)gl->error_direction,(PetscEnum*)&gl->error_direction,PETSC_NULL);
1111: PetscOptionsTruth("-ts_gl_extrapolate","Extrapolate stage solution from previous solution (sometimes unstable)","TSGLSetExtrapolate",gl->extrapolate,&gl->extrapolate,PETSC_NULL);
1112: PetscOptionsReal("-ts_gl_atol","Absolute tolerance","TSGLSetTolerances",gl->wrms_atol,&gl->wrms_atol,PETSC_NULL);
1113: PetscOptionsReal("-ts_gl_rtol","Relative tolerance","TSGLSetTolerances",gl->wrms_rtol,&gl->wrms_rtol,PETSC_NULL);
1114: PetscOptionsString("-ts_gl_complete","Method to use for completing the step","none",completef,completef,sizeof completef,&flg);
1115: if (flg) {
1116: PetscTruth match1,match2;
1117: PetscStrcmp(completef,"rescale",&match1);
1118: PetscStrcmp(completef,"rescale-and-modify",&match2);
1119: if (match1) gl->CompleteStep = TSGLCompleteStep_Rescale;
1120: else if (match2) gl->CompleteStep = TSGLCompleteStep_RescaleAndModify;
1121: else SETERRQ1(PETSC_ERR_ARG_UNKNOWN_TYPE,"%s",completef);
1122: }
1123: {
1124: char type[256] = TSGLACCEPT_ALWAYS;
1125: PetscOptionsList("-ts_gl_accept_type","Method to use for determining whether to accept a step","TSGLSetAcceptType",TSGLAcceptList,gl->accept_name[0]?gl->accept_name:type,type,sizeof type,&flg);
1126: if (flg || !gl->accept_name[0]) {
1127: TSGLSetAcceptType(ts,type);
1128: }
1129: }
1130: }
1131: PetscOptionsTail();
1132: {
1133: TSGLAdapt adapt;
1134: TSGLGetAdapt(ts,&adapt);
1135: TSGLAdaptSetFromOptions(adapt);
1136: }
1137: return(0);
1138: }
1142: static PetscErrorCode TSView_GL(TS ts,PetscViewer viewer)
1143: {
1144: TS_GL *gl = (TS_GL*)ts->data;
1145: PetscInt i;
1146: PetscTruth iascii,details;
1147: PetscErrorCode ierr;
1150: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
1151: if (iascii) {
1152: PetscViewerASCIIPrintf(viewer," min order %D, max order %D, current order %D\n",gl->min_order,gl->max_order,gl->schemes[gl->current_scheme]->p);
1153: PetscViewerASCIIPrintf(viewer," Error estimation: %s\n",TSGLErrorDirections[gl->error_direction]);
1154: PetscViewerASCIIPrintf(viewer," Extrapolation: %s\n",gl->extrapolate?"yes":"no");
1155: PetscViewerASCIIPrintf(viewer," Acceptance test: %s\n",gl->accept_name[0]?gl->accept_name:"(not yet set)");
1156: PetscViewerASCIIPushTab(viewer);
1157: TSGLAdaptView(gl->adapt,viewer);
1158: PetscViewerASCIIPopTab(viewer);
1159: PetscViewerASCIIPrintf(viewer," type: %s\n",gl->type_name[0]?gl->type_name:"(not yet set)");
1160: PetscViewerASCIIPrintf(viewer,"Schemes within family (%d):\n",gl->nschemes);
1161: details = PETSC_FALSE;
1162: PetscOptionsGetTruth(((PetscObject)ts)->prefix,"-ts_gl_view_detailed",&details,PETSC_NULL);
1163: PetscViewerASCIIPushTab(viewer);
1164: for (i=0; i<gl->nschemes; i++) {
1165: TSGLSchemeView(gl->schemes[i],details,viewer);
1166: }
1167: if (gl->View) {
1168: (*gl->View)(gl,viewer);
1169: }
1170: PetscViewerASCIIPopTab(viewer);
1171: } else {
1172: SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for TSGL",((PetscObject)viewer)->type_name);
1173: }
1174: return(0);
1175: }
1179: /*@C
1180: TSGLRegister - see TSGLRegisterDynamic()
1182: Level: advanced
1183: @*/
1184: PetscErrorCode TSGLRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(TS))
1185: {
1187: char fullname[PETSC_MAX_PATH_LEN];
1190: PetscFListConcat(path,name,fullname);
1191: PetscFListAdd(&TSGLList,sname,fullname,(void(*)(void))function);
1192: return(0);
1193: }
1197: /*@C
1198: TSGLAcceptRegister - see TSGLAcceptRegisterDynamic()
1200: Level: advanced
1201: @*/
1202: PetscErrorCode TSGLAcceptRegister(const char sname[],const char path[],const char name[],TSGLAcceptFunction function)
1203: {
1205: char fullname[PETSC_MAX_PATH_LEN];
1208: PetscFListConcat(path,name,fullname);
1209: PetscFListAdd(&TSGLAcceptList,sname,fullname,(void(*)(void))function);
1210: return(0);
1211: }
1215: /*@C
1216: TSGLRegisterAll - Registers all of the general linear methods in TSGL
1218: Not Collective
1220: Level: advanced
1222: .keywords: TS, TSGL, register, all
1224: .seealso: TSGLRegisterDestroy()
1225: @*/
1226: PetscErrorCode TSGLRegisterAll(const char path[])
1227: {
1231: TSGLRegisterAllCalled = PETSC_TRUE;
1233: TSGLRegisterDynamic(TSGL_IRKS,path,"TSGLCreate_IRKS",TSGLCreate_IRKS);
1234: TSGLAcceptRegisterDynamic(TSGLACCEPT_ALWAYS,path,"TSGLAccept_Always",TSGLAccept_Always);
1235: return(0);
1236: }
1240: /*@C
1241: TSGLRegisterDestroy - Frees the list of schemes that were registered by TSGLRegister()/TSGLRegisterDynamic().
1243: Not Collective
1245: Level: advanced
1247: .keywords: TSGL, register, destroy
1248: .seealso: TSGLRegister(), TSGLRegisterAll(), TSGLRegisterDynamic()
1249: @*/
1250: PetscErrorCode TSGLRegisterDestroy(void)
1251: {
1255: PetscFListDestroy(&TSGLList);
1256: TSGLRegisterAllCalled = PETSC_FALSE;
1257: return(0);
1258: }
1263: /*@C
1264: TSGLInitializePackage - This function initializes everything in the TSGL package. It is called
1265: from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_GL()
1266: when using static libraries.
1268: Input Parameter:
1269: path - The dynamic library path, or PETSC_NULL
1271: Level: developer
1273: .keywords: TS, TSGL, initialize, package
1274: .seealso: PetscInitialize()
1275: @*/
1276: PetscErrorCode TSGLInitializePackage(const char path[])
1277: {
1281: if (TSGLPackageInitialized) return(0);
1282: TSGLPackageInitialized = PETSC_TRUE;
1283: TSGLRegisterAll(PETSC_NULL);
1284: PetscRegisterFinalize(TSGLFinalizePackage);
1285: return(0);
1286: }
1290: /*@C
1291: TSGLFinalizePackage - This function destroys everything in the TSGL package. It is
1292: called from PetscFinalize().
1294: Level: developer
1296: .keywords: Petsc, destroy, package
1297: .seealso: PetscFinalize()
1298: @*/
1299: PetscErrorCode TSGLFinalizePackage(void)
1300: {
1302: TSGLPackageInitialized = PETSC_FALSE;
1303: TSGLRegisterAllCalled = PETSC_FALSE;
1304: TSGLList = PETSC_NULL;
1305: TSGLAcceptList = PETSC_NULL;
1306: return(0);
1307: }
1309: /* ------------------------------------------------------------ */
1310: /*MC
1311: TSGL - DAE solver using implicit General Linear methods
1313: These methods contain Runge-Kutta and multistep schemes as special cases. These special cases have some fundamental
1314: limitations. For example, diagonally implicit Runge-Kutta cannot have stage order greater than 1 which limits their
1315: applicability to very stiff systems. Meanwhile, multistep methods cannot be A-stable for order greater than 2 and BDF
1316: are not 0-stable for order greater than 6. GL methods can be A- and L-stable with arbitrarily high stage order and
1317: reliable error estimates for both 1 and 2 orders higher to facilitate adaptive step sizes and adaptive order schemes.
1318: All this is possible while preserving a singly diagonally implicit structure.
1320: Options database keys:
1321: + -ts_gl_type <type> - the class of general linear method (irks)
1322: . -ts_gl_rtol <tol> - relative error
1323: . -ts_gl_atol <tol> - absolute error
1324: . -ts_gl_min_order <p> - minimum order method to consider (default=1)
1325: . -ts_gl_max_order <p> - maximum order method to consider (default=3)
1326: . -ts_gl_start_order <p> - order of starting method (default=1)
1327: . -ts_gl_complete <method> - method to use for completing the step (rescale-and-modify or rescale)
1328: - -ts_adapt_type <method> - adaptive controller to use (none step both)
1330: Notes:
1331: This integrator can be applied to DAE.
1333: Diagonally implicit general linear (DIGL) methods are a generalization of diagonally implicit Runge-Kutta (DIRK).
1334: They are represented by the tableau
1336: .vb
1337: A | U
1338: -------
1339: B | V
1340: .ve
1342: combined with a vector c of abscissa. "Diagonally implicit" means that A is lower triangular.
1343: A step of the general method reads
1345: .vb
1346: [ Y ] = [A U] [ Y' ]
1347: [X^k] = [B V] [X^{k-1}]
1348: .ve
1350: where Y is the multivector of stage values, Y' is the multivector of stage derivatives, X^k is the Nordsieck vector of
1351: the solution at step k. The Nordsieck vector consists of the first r moments of the solution, given by
1353: .vb
1354: X = [x_0,x_1,...,x_{r-1}] = [x, h x', h^2 x'', ..., h^{r-1} x^{(r-1)} ]
1355: .ve
1357: If A is lower triangular, we can solve the stages (Y,Y') sequentially
1359: .vb
1360: y_i = h sum_{j=0}^{s-1} (a_ij y'_j) + sum_{j=0}^{r-1} u_ij x_j, i=0,...,{s-1}
1361: .ve
1363: and then construct the pieces to carry to the next step
1365: .vb
1366: xx_i = h sum_{j=0}^{s-1} b_ij y'_j + sum_{j=0}^{r-1} v_ij x_j, i=0,...,{r-1}
1367: .ve
1369: Note that when the equations are cast in implicit form, we are using the stage equation to define y'_i
1370: in terms of y_i and known stuff (y_j for j<i and x_j for all j).
1373: Error estimation
1375: At present, the most attractive GL methods for stiff problems are singly diagonally implicit schemes which posses
1376: Inherent Runge-Kutta Stability (IRKS). These methods have r=s, the number of items passed between steps is equal to
1377: the number of stages. The order and stage-order are one less than the number of stages. We use the error estimates
1378: in the 2007 paper which provide the following estimates
1380: .vb
1381: h^{p+1} X^{(p+1)} = phi_0^T Y' + [0 psi_0^T] Xold
1382: h^{p+2} X^{(p+2)} = phi_1^T Y' + [0 psi_1^T] Xold
1383: h^{p+2} (dx'/dx) X^{(p+1)} = phi_2^T Y' + [0 psi_2^T] Xold
1384: .ve
1386: These estimates are accurate to O(h^{p+3}).
1388: Changing the step size
1390: We use the generalized "rescale and modify" scheme, see equation (4.5) of the 2007 paper.
1392: Level: beginner
1394: References:
1395: John Butcher and Z. Jackieweicz and W. Wright, On error propagation in general linear methods for
1396: ordinary differential equations, Journal of Complexity, Vol 23 (4-6), 2007.
1398: John Butcher, Numerical methods for ordinary differential equations, second edition, Wiley, 2009.
1400: .seealso: TSCreate(), TS, TSSetType()
1402: M*/
1406: PetscErrorCode TSCreate_GL(TS ts)
1407: {
1408: TS_GL *gl;
1412: #if !defined(PETSC_USE_DYNAMIC_LIBRARIES)
1413: TSGLInitializePackage(PETSC_NULL);
1414: #endif
1416: PetscNewLog(ts,TS_GL,&gl);
1417: ts->data = (void*)gl;
1419: ts->ops->destroy = TSDestroy_GL;
1420: ts->ops->view = TSView_GL;
1421: ts->ops->setup = TSSetUp_GL;
1422: ts->ops->step = TSStep_GL;
1423: ts->ops->setfromoptions = TSSetFromOptions_GL;
1425: ts->problem_type = TS_NONLINEAR;
1426: SNESCreate(((PetscObject)ts)->comm,&ts->snes);
1427: PetscObjectIncrementTabLevel((PetscObject)ts->snes,(PetscObject)ts,1);
1429: gl->max_step_rejections = 1;
1430: gl->min_order = 1;
1431: gl->max_order = 3;
1432: gl->start_order = 1;
1433: gl->current_scheme = -1;
1434: gl->extrapolate = PETSC_FALSE;
1436: gl->wrms_atol = 1e-8;
1437: gl->wrms_rtol = 1e-5;
1439: PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSGLSetType_C", "TSGLSetType_GL", &TSGLSetType_GL);
1440: PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSGLSetAcceptType_C","TSGLSetAcceptType_GL",&TSGLSetAcceptType_GL);
1441: PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSGLGetAdapt_C", "TSGLGetAdapt_GL", &TSGLGetAdapt_GL);
1442: return(0);
1443: }