Actual source code: gl.c

petsc-3.7.7 2017-09-25
Report Typos and Errors
  2: #include <../src/ts/impls/implicit/gl/gl.h>                /*I   "petscts.h"   I*/
  3: #include <petscdm.h>
  4: #include <petscblaslapack.h>

  6: static const char        *TSGLErrorDirections[] = {"FORWARD","BACKWARD","TSGLErrorDirection","TSGLERROR_",0};
  7: static PetscFunctionList TSGLList;
  8: static PetscFunctionList TSGLAcceptList;
  9: static PetscBool         TSGLPackageInitialized;
 10: static PetscBool         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 PetscPowRealInt(PetscRealPart(c),p)/Factorial(p);
 31: }

 35: static PetscErrorCode TSGLGetVecs(TS ts,DM dm,Vec *Z,Vec *Ydotstage)
 36: {
 37:   TS_GL          *gl = (TS_GL*)ts->data;

 41:   if (Z) {
 42:     if (dm && dm != ts->dm) {
 43:       DMGetNamedGlobalVector(dm,"TSGL_Z",Z);
 44:     } else *Z = gl->Z;
 45:   }
 46:   if (Ydotstage) {
 47:     if (dm && dm != ts->dm) {
 48:       DMGetNamedGlobalVector(dm,"TSGL_Ydot",Ydotstage);
 49:     } else *Ydotstage = gl->Ydot[gl->stage];
 50:   }
 51:   return(0);
 52: }


 57: static PetscErrorCode TSGLRestoreVecs(TS ts,DM dm,Vec *Z,Vec *Ydotstage)
 58: {

 62:   if (Z) {
 63:     if (dm && dm != ts->dm) {
 64:       DMRestoreNamedGlobalVector(dm,"TSGL_Z",Z);
 65:     }
 66:   }
 67:   if (Ydotstage) {

 69:     if (dm && dm != ts->dm) {
 70:       DMRestoreNamedGlobalVector(dm,"TSGL_Ydot",Ydotstage);
 71:     }
 72:   }
 73:   return(0);
 74: }

 78: static PetscErrorCode DMCoarsenHook_TSGL(DM fine,DM coarse,void *ctx)
 79: {
 81:   return(0);
 82: }

 86: static PetscErrorCode DMRestrictHook_TSGL(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
 87: {
 88:   TS             ts = (TS)ctx;
 90:   Vec            Ydot,Ydot_c;

 93:   TSGLGetVecs(ts,fine,NULL,&Ydot);
 94:   TSGLGetVecs(ts,coarse,NULL,&Ydot_c);
 95:   MatRestrict(restrct,Ydot,Ydot_c);
 96:   VecPointwiseMult(Ydot_c,rscale,Ydot_c);
 97:   TSGLRestoreVecs(ts,fine,NULL,&Ydot);
 98:   TSGLRestoreVecs(ts,coarse,NULL,&Ydot_c);
 99:   return(0);
100: }

104: static PetscErrorCode DMSubDomainHook_TSGL(DM dm,DM subdm,void *ctx)
105: {
107:   return(0);
108: }

112: static PetscErrorCode DMSubDomainRestrictHook_TSGL(DM dm,VecScatter gscat, VecScatter lscat,DM subdm,void *ctx)
113: {
114:   TS             ts = (TS)ctx;
116:   Vec            Ydot,Ydot_s;

119:   TSGLGetVecs(ts,dm,NULL,&Ydot);
120:   TSGLGetVecs(ts,subdm,NULL,&Ydot_s);

122:   VecScatterBegin(gscat,Ydot,Ydot_s,INSERT_VALUES,SCATTER_FORWARD);
123:   VecScatterEnd(gscat,Ydot,Ydot_s,INSERT_VALUES,SCATTER_FORWARD);

125:   TSGLRestoreVecs(ts,dm,NULL,&Ydot);
126:   TSGLRestoreVecs(ts,subdm,NULL,&Ydot_s);
127:   return(0);
128: }

132: static PetscErrorCode TSGLSchemeCreate(PetscInt p,PetscInt q,PetscInt r,PetscInt s,const PetscScalar *c,
133:                                        const PetscScalar *a,const PetscScalar *b,const PetscScalar *u,const PetscScalar *v,TSGLScheme *inscheme)
134: {
135:   TSGLScheme     scheme;
136:   PetscInt       j;

140:   if (p < 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Scheme order must be positive");
141:   if (r < 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"At least one item must be carried between steps");
142:   if (s < 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"At least one stage is required");
144:   *inscheme = 0;
145:   PetscNew(&scheme);
146:   scheme->p = p;
147:   scheme->q = q;
148:   scheme->r = r;
149:   scheme->s = s;

151:   PetscMalloc5(s,&scheme->c,s*s,&scheme->a,r*s,&scheme->b,r*s,&scheme->u,r*r,&scheme->v);
152:   PetscMemcpy(scheme->c,c,s*sizeof(PetscScalar));
153:   for (j=0; j<s*s; j++) scheme->a[j] = (PetscAbsScalar(a[j]) < 1e-12) ? 0 : a[j];
154:   for (j=0; j<r*s; j++) scheme->b[j] = (PetscAbsScalar(b[j]) < 1e-12) ? 0 : b[j];
155:   for (j=0; j<s*r; j++) scheme->u[j] = (PetscAbsScalar(u[j]) < 1e-12) ? 0 : u[j];
156:   for (j=0; j<r*r; j++) scheme->v[j] = (PetscAbsScalar(v[j]) < 1e-12) ? 0 : v[j];

158:   PetscMalloc6(r,&scheme->alpha,r,&scheme->beta,r,&scheme->gamma,3*s,&scheme->phi,3*r,&scheme->psi,r,&scheme->stage_error);
159:   {
160:     PetscInt     i,j,k,ss=s+2;
161:     PetscBLASInt m,n,one=1,*ipiv,lwork=4*((s+3)*3+3),info,ldb;
162:     PetscReal    rcond,*sing,*workreal;
163:     PetscScalar  *ImV,*H,*bmat,*workscalar,*c=scheme->c,*a=scheme->a,*b=scheme->b,*u=scheme->u,*v=scheme->v;
164: #if !defined(PETSC_MISSING_LAPACK_GELSS)
165:     PetscBLASInt rank;
166: #endif
167:     PetscMalloc7(PetscSqr(r),&ImV,3*s,&H,3*ss,&bmat,lwork,&workscalar,5*(3+r),&workreal,r+s,&sing,r+s,&ipiv);

169:     /* column-major input */
170:     for (i=0; i<r-1; i++) {
171:       for (j=0; j<r-1; j++) ImV[i+j*r] = 1.0*(i==j) - v[(i+1)*r+j+1];
172:     }
173:     /* Build right hand side for alpha (tp - glm.B(2:end,:)*(glm.c.^(p)./factorial(p))) */
174:     for (i=1; i<r; i++) {
175:       scheme->alpha[i] = 1./Factorial(p+1-i);
176:       for (j=0; j<s; j++) scheme->alpha[i] -= b[i*s+j]*CPowF(c[j],p);
177:     }
178:     PetscBLASIntCast(r-1,&m);
179:     PetscBLASIntCast(r,&n);
180:     PetscStackCallBLAS("LAPACKgesv",LAPACKgesv_(&m,&one,ImV,&n,ipiv,scheme->alpha+1,&n,&info));
181:     if (info < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad argument to GESV");
182:     if (info > 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Bad LU factorization");

184:     /* Build right hand side for beta (tp1 - glm.B(2:end,:)*(glm.c.^(p+1)./factorial(p+1)) - e.alpha) */
185:     for (i=1; i<r; i++) {
186:       scheme->beta[i] = 1./Factorial(p+2-i) - scheme->alpha[i];
187:       for (j=0; j<s; j++) scheme->beta[i] -= b[i*s+j]*CPowF(c[j],p+1);
188:     }
189:     PetscStackCallBLAS("LAPACKgetrs",LAPACKgetrs_("No transpose",&m,&one,ImV,&n,ipiv,scheme->beta+1,&n,&info));
190:     if (info < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad argument to GETRS");
191:     if (info > 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Should not happen");

193:     /* Build stage_error vector
194:            xi = glm.c.^(p+1)/factorial(p+1) - glm.A*glm.c.^p/factorial(p) + glm.U(:,2:end)*e.alpha;
195:     */
196:     for (i=0; i<s; i++) {
197:       scheme->stage_error[i] = CPowF(c[i],p+1);
198:       for (j=0; j<s; j++) scheme->stage_error[i] -= a[i*s+j]*CPowF(c[j],p);
199:       for (j=1; j<r; j++) scheme->stage_error[i] += u[i*r+j]*scheme->alpha[j];
200:     }

202:     /* alpha[0] (epsilon in B,J,W 2007)
203:            epsilon = 1/factorial(p+1) - B(1,:)*c.^p/factorial(p) + V(1,2:end)*e.alpha;
204:     */
205:     scheme->alpha[0] = 1./Factorial(p+1);
206:     for (j=0; j<s; j++) scheme->alpha[0] -= b[0*s+j]*CPowF(c[j],p);
207:     for (j=1; j<r; j++) scheme->alpha[0] += v[0*r+j]*scheme->alpha[j];

209:     /* right hand side for gamma (glm.B(2:end,:)*e.xi - e.epsilon*eye(s-1,1)) */
210:     for (i=1; i<r; i++) {
211:       scheme->gamma[i] = (i==1 ? -1. : 0)*scheme->alpha[0];
212:       for (j=0; j<s; j++) scheme->gamma[i] += b[i*s+j]*scheme->stage_error[j];
213:     }
214:     PetscStackCallBLAS("LAPACKgetrs",LAPACKgetrs_("No transpose",&m,&one,ImV,&n,ipiv,scheme->gamma+1,&n,&info));
215:     if (info < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad argument to GETRS");
216:     if (info > 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Should not happen");

218:     /* beta[0] (rho in B,J,W 2007)
219:         e.rho = 1/factorial(p+2) - glm.B(1,:)*glm.c.^(p+1)/factorial(p+1) ...
220:             + glm.V(1,2:end)*e.beta;% - e.epsilon;
221:     % Note: The paper (B,J,W 2007) includes the last term in their definition
222:     * */
223:     scheme->beta[0] = 1./Factorial(p+2);
224:     for (j=0; j<s; j++) scheme->beta[0] -= b[0*s+j]*CPowF(c[j],p+1);
225:     for (j=1; j<r; j++) scheme->beta[0] += v[0*r+j]*scheme->beta[j];

227:     /* gamma[0] (sigma in B,J,W 2007)
228:     *   e.sigma = glm.B(1,:)*e.xi + glm.V(1,2:end)*e.gamma;
229:     * */
230:     scheme->gamma[0] = 0.0;
231:     for (j=0; j<s; j++) scheme->gamma[0] += b[0*s+j]*scheme->stage_error[j];
232:     for (j=1; j<r; j++) scheme->gamma[0] += v[0*s+j]*scheme->gamma[j];

234:     /* Assemble H
235:     *    % Determine the error estimators phi
236:        H = [[cpow(glm.c,p) + C*e.alpha] [cpow(glm.c,p+1) + C*e.beta] ...
237:                [e.xi - C*(e.gamma + 0*e.epsilon*eye(s-1,1))]]';
238:     % Paper has formula above without the 0, but that term must be left
239:     % out to satisfy the conditions they propose and to make the
240:     % example schemes work
241:     e.H = H;
242:     e.phi = (H \ [1 0 0;1 1 0;0 0 -1])';
243:     e.psi = -e.phi*C;
244:     * */
245:     for (j=0; j<s; j++) {
246:       H[0+j*3] = CPowF(c[j],p);
247:       H[1+j*3] = CPowF(c[j],p+1);
248:       H[2+j*3] = scheme->stage_error[j];
249:       for (k=1; k<r; k++) {
250:         H[0+j*3] += CPowF(c[j],k-1)*scheme->alpha[k];
251:         H[1+j*3] += CPowF(c[j],k-1)*scheme->beta[k];
252:         H[2+j*3] -= CPowF(c[j],k-1)*scheme->gamma[k];
253:       }
254:     }
255:     bmat[0+0*ss] = 1.;  bmat[0+1*ss] = 0.;  bmat[0+2*ss] = 0.;
256:     bmat[1+0*ss] = 1.;  bmat[1+1*ss] = 1.;  bmat[1+2*ss] = 0.;
257:     bmat[2+0*ss] = 0.;  bmat[2+1*ss] = 0.;  bmat[2+2*ss] = -1.;
258:     m     = 3;
259:     PetscBLASIntCast(s,&n);
260:     PetscBLASIntCast(ss,&ldb);
261:     rcond = 1e-12;
262: #if defined(PETSC_MISSING_LAPACK_GELSS)
263:     /* ESSL does not have this routine */
264:     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GELSS - Lapack routine is unavailable\nNot able to run GL time stepping.");
265: #else
266: #if defined(PETSC_USE_COMPLEX)
267:     /* ZGELSS( M, N, NRHS, A, LDA, B, LDB, S, RCOND, RANK, WORK, LWORK, RWORK, INFO) */
268:     PetscStackCallBLAS("LAPACKgelss",LAPACKgelss_(&m,&n,&m,H,&m,bmat,&ldb,sing,&rcond,&rank,workscalar,&lwork,workreal,&info));
269: #else
270:     /* DGELSS( M, N, NRHS, A, LDA, B, LDB, S, RCOND, RANK, WORK, LWORK, INFO) */
271:     PetscStackCallBLAS("LAPACKgelss",LAPACKgelss_(&m,&n,&m,H,&m,bmat,&ldb,sing,&rcond,&rank,workscalar,&lwork,&info));
272: #endif
273:     if (info < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad argument to GELSS");
274:     if (info > 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"SVD failed to converge");
275: #endif

277:     for (j=0; j<3; j++) {
278:       for (k=0; k<s; k++) scheme->phi[k+j*s] = bmat[k+j*ss];
279:     }

281:     /* the other part of the error estimator, psi in B,J,W 2007 */
282:     scheme->psi[0*r+0] = 0.;
283:     scheme->psi[1*r+0] = 0.;
284:     scheme->psi[2*r+0] = 0.;
285:     for (j=1; j<r; j++) {
286:       scheme->psi[0*r+j] = 0.;
287:       scheme->psi[1*r+j] = 0.;
288:       scheme->psi[2*r+j] = 0.;
289:       for (k=0; k<s; k++) {
290:         scheme->psi[0*r+j] -= CPowF(c[k],j-1)*scheme->phi[0*s+k];
291:         scheme->psi[1*r+j] -= CPowF(c[k],j-1)*scheme->phi[1*s+k];
292:         scheme->psi[2*r+j] -= CPowF(c[k],j-1)*scheme->phi[2*s+k];
293:       }
294:     }
295:     PetscFree7(ImV,H,bmat,workscalar,workreal,sing,ipiv);
296:   }
297:   /* Check which properties are satisfied */
298:   scheme->stiffly_accurate = PETSC_TRUE;
299:   if (scheme->c[s-1] != 1.) scheme->stiffly_accurate = PETSC_FALSE;
300:   for (j=0; j<s; j++) if (a[(s-1)*s+j] != b[j]) scheme->stiffly_accurate = PETSC_FALSE;
301:   for (j=0; j<r; j++) if (u[(s-1)*r+j] != v[j]) scheme->stiffly_accurate = PETSC_FALSE;
302:   scheme->fsal = scheme->stiffly_accurate; /* FSAL is stronger */
303:   for (j=0; j<s-1; j++) if (r>1 && b[1*s+j] != 0.) scheme->fsal = PETSC_FALSE;
304:   if (b[1*s+r-1] != 1.) scheme->fsal = PETSC_FALSE;
305:   for (j=0; j<r; j++) if (r>1 && v[1*r+j] != 0.) scheme->fsal = PETSC_FALSE;

307:   *inscheme = scheme;
308:   return(0);
309: }

313: static PetscErrorCode TSGLSchemeDestroy(TSGLScheme sc)
314: {

318:   PetscFree5(sc->c,sc->a,sc->b,sc->u,sc->v);
319:   PetscFree6(sc->alpha,sc->beta,sc->gamma,sc->phi,sc->psi,sc->stage_error);
320:   PetscFree(sc);
321:   return(0);
322: }

326: static PetscErrorCode TSGLDestroy_Default(TS_GL *gl)
327: {
329:   PetscInt       i;

332:   for (i=0; i<gl->nschemes; i++) {
333:     if (gl->schemes[i]) {TSGLSchemeDestroy(gl->schemes[i]);}
334:   }
335:   PetscFree(gl->schemes);
336:   gl->nschemes = 0;
337:   PetscMemzero(gl->type_name,sizeof(gl->type_name));
338:   return(0);
339: }

343: static PetscErrorCode TSGLViewTable_Private(PetscViewer viewer,PetscInt m,PetscInt n,const PetscScalar a[],const char name[])
344: {
346:   PetscBool      iascii;
347:   PetscInt       i,j;

350:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
351:   if (iascii) {
352:     PetscViewerASCIIPrintf(viewer,"%30s = [",name);
353:     for (i=0; i<m; i++) {
354:       if (i) {PetscViewerASCIIPrintf(viewer,"%30s   [","");}
355:       PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
356:       for (j=0; j<n; j++) {
357:         PetscViewerASCIIPrintf(viewer," %12.8g",PetscRealPart(a[i*n+j]));
358:       }
359:       PetscViewerASCIIPrintf(viewer,"]\n");
360:       PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
361:     }
362:   }
363:   return(0);
364: }


369: static PetscErrorCode TSGLSchemeView(TSGLScheme sc,PetscBool view_details,PetscViewer viewer)
370: {
372:   PetscBool      iascii;

375:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
376:   if (iascii) {
377:     PetscViewerASCIIPrintf(viewer,"GL scheme p,q,r,s = %d,%d,%d,%d\n",sc->p,sc->q,sc->r,sc->s);
378:     PetscViewerASCIIPushTab(viewer);
379:     PetscViewerASCIIPrintf(viewer,"Stiffly accurate: %s,  FSAL: %s\n",sc->stiffly_accurate ? "yes" : "no",sc->fsal ? "yes" : "no");
380:     PetscViewerASCIIPrintf(viewer,"Leading error constants: %10.3e  %10.3e  %10.3e\n",
381:                                   PetscRealPart(sc->alpha[0]),PetscRealPart(sc->beta[0]),PetscRealPart(sc->gamma[0]));
382:     TSGLViewTable_Private(viewer,1,sc->s,sc->c,"Abscissas c");
383:     if (view_details) {
384:       TSGLViewTable_Private(viewer,sc->s,sc->s,sc->a,"A");
385:       TSGLViewTable_Private(viewer,sc->r,sc->s,sc->b,"B");
386:       TSGLViewTable_Private(viewer,sc->s,sc->r,sc->u,"U");
387:       TSGLViewTable_Private(viewer,sc->r,sc->r,sc->v,"V");

389:       TSGLViewTable_Private(viewer,3,sc->s,sc->phi,"Error estimate phi");
390:       TSGLViewTable_Private(viewer,3,sc->r,sc->psi,"Error estimate psi");
391:       TSGLViewTable_Private(viewer,1,sc->r,sc->alpha,"Modify alpha");
392:       TSGLViewTable_Private(viewer,1,sc->r,sc->beta,"Modify beta");
393:       TSGLViewTable_Private(viewer,1,sc->r,sc->gamma,"Modify gamma");
394:       TSGLViewTable_Private(viewer,1,sc->s,sc->stage_error,"Stage error xi");
395:     }
396:     PetscViewerASCIIPopTab(viewer);
397:   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported",((PetscObject)viewer)->type_name);
398:   return(0);
399: }

403: static PetscErrorCode TSGLEstimateHigherMoments_Default(TSGLScheme sc,PetscReal h,Vec Ydot[],Vec Xold[],Vec hm[])
404: {
406:   PetscInt       i;

409:   if (sc->r > 64 || sc->s > 64) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Ridiculous number of stages or items passed between stages");
410:   /* build error vectors*/
411:   for (i=0; i<3; i++) {
412:     PetscScalar phih[64];
413:     PetscInt    j;
414:     for (j=0; j<sc->s; j++) phih[j] = sc->phi[i*sc->s+j]*h;
415:     VecZeroEntries(hm[i]);
416:     VecMAXPY(hm[i],sc->s,phih,Ydot);
417:     VecMAXPY(hm[i],sc->r,&sc->psi[i*sc->r],Xold);
418:   }
419:   return(0);
420: }

424: static PetscErrorCode TSGLCompleteStep_Rescale(TSGLScheme sc,PetscReal h,TSGLScheme next_sc,PetscReal next_h,Vec Ydot[],Vec Xold[],Vec X[])
425: {
427:   PetscScalar    brow[32],vrow[32];
428:   PetscInt       i,j,r,s;

431:   /* Build the new solution from (X,Ydot) */
432:   r = sc->r;
433:   s = sc->s;
434:   for (i=0; i<r; i++) {
435:     VecZeroEntries(X[i]);
436:     for (j=0; j<s; j++) brow[j] = h*sc->b[i*s+j];
437:     VecMAXPY(X[i],s,brow,Ydot);
438:     for (j=0; j<r; j++) vrow[j] = sc->v[i*r+j];
439:     VecMAXPY(X[i],r,vrow,Xold);
440:   }
441:   return(0);
442: }

446: static PetscErrorCode TSGLCompleteStep_RescaleAndModify(TSGLScheme sc,PetscReal h,TSGLScheme next_sc,PetscReal next_h,Vec Ydot[],Vec Xold[],Vec X[])
447: {
449:   PetscScalar    brow[32],vrow[32];
450:   PetscReal      ratio;
451:   PetscInt       i,j,p,r,s;

454:   /* Build the new solution from (X,Ydot) */
455:   p     = sc->p;
456:   r     = sc->r;
457:   s     = sc->s;
458:   ratio = next_h/h;
459:   for (i=0; i<r; i++) {
460:     VecZeroEntries(X[i]);
461:     for (j=0; j<s; j++) {
462:       brow[j] = h*(PetscPowRealInt(ratio,i)*sc->b[i*s+j]
463:                    + (PetscPowRealInt(ratio,i) - PetscPowRealInt(ratio,p+1))*(+ sc->alpha[i]*sc->phi[0*s+j])
464:                    + (PetscPowRealInt(ratio,i) - PetscPowRealInt(ratio,p+2))*(+ sc->beta [i]*sc->phi[1*s+j]
465:                                                                               + sc->gamma[i]*sc->phi[2*s+j]));
466:     }
467:     VecMAXPY(X[i],s,brow,Ydot);
468:     for (j=0; j<r; j++) {
469:       vrow[j] = (PetscPowRealInt(ratio,i)*sc->v[i*r+j]
470:                  + (PetscPowRealInt(ratio,i) - PetscPowRealInt(ratio,p+1))*(+ sc->alpha[i]*sc->psi[0*r+j])
471:                  + (PetscPowRealInt(ratio,i) - PetscPowRealInt(ratio,p+2))*(+ sc->beta [i]*sc->psi[1*r+j]
472:                                                                             + sc->gamma[i]*sc->psi[2*r+j]));
473:     }
474:     VecMAXPY(X[i],r,vrow,Xold);
475:   }
476:   if (r < next_sc->r) {
477:     if (r+1 != next_sc->r) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot accommodate jump in r greater than 1");
478:     VecZeroEntries(X[r]);
479:     for (j=0; j<s; j++) brow[j] = h*PetscPowRealInt(ratio,p+1)*sc->phi[0*s+j];
480:     VecMAXPY(X[r],s,brow,Ydot);
481:     for (j=0; j<r; j++) vrow[j] = PetscPowRealInt(ratio,p+1)*sc->psi[0*r+j];
482:     VecMAXPY(X[r],r,vrow,Xold);
483:   }
484:   return(0);
485: }

489: static PetscErrorCode TSGLCreate_IRKS(TS ts)
490: {
491:   TS_GL          *gl = (TS_GL*)ts->data;

495:   gl->Destroy               = TSGLDestroy_Default;
496:   gl->EstimateHigherMoments = TSGLEstimateHigherMoments_Default;
497:   gl->CompleteStep          = TSGLCompleteStep_RescaleAndModify;
498:   PetscMalloc1(10,&gl->schemes);
499:   gl->nschemes = 0;

501:   {
502:     /* p=1,q=1, r=s=2, A- and L-stable with error estimates of order 2 and 3
503:     * Listed in Butcher & Podhaisky 2006. On error estimation in general linear methods for stiff ODE.
504:     * irks(0.3,0,[.3,1],[1],1)
505:     * Note: can be made to have classical order (not stage order) 2 by replacing 0.3 with 1-sqrt(1/2)
506:     * but doing so would sacrifice the error estimator.
507:     */
508:     const PetscScalar c[2]    = {3./10., 1.};
509:     const PetscScalar a[2][2] = {{3./10., 0}, {7./10., 3./10.}};
510:     const PetscScalar b[2][2] = {{7./10., 3./10.}, {0,1}};
511:     const PetscScalar u[2][2] = {{1,0},{1,0}};
512:     const PetscScalar v[2][2] = {{1,0},{0,0}};
513:     TSGLSchemeCreate(1,1,2,2,c,*a,*b,*u,*v,&gl->schemes[gl->nschemes++]);
514:   }

516:   {
517:     /* p=q=2, r=s=3: irks(4/9,0,[1:3]/3,[0.33852],1) */
518:     /* http://www.math.auckland.ac.nz/~hpod/atlas/i2a.html */
519:     const PetscScalar c[3] = {1./3., 2./3., 1}
520:     ,a[3][3] = {{4./9.                ,0                      ,       0},
521:                 {1.03750643704090e+00 ,                  4./9.,       0},
522:                 {7.67024779410304e-01 ,  -3.81140216918943e-01,   4./9.}}
523:     ,b[3][3] = {{0.767024779410304,  -0.381140216918943,   4./9.},
524:                 {0.000000000000000,  0.000000000000000,   1.000000000000000},
525:                 {-2.075048385225385,   0.621728385225383,   1.277197204924873}}
526:     ,u[3][3] = {{1.0000000000000000,  -0.1111111111111109,  -0.0925925925925922},
527:                 {1.0000000000000000,  -0.8152842148186744,  -0.4199095530877056},
528:                 {1.0000000000000000,   0.1696709930641948,   0.0539741070314165}}
529:     ,v[3][3] = {{1.0000000000000000,  0.1696709930641948,   0.0539741070314165},
530:                 {0.000000000000000,   0.000000000000000,   0.000000000000000},
531:                 {0.000000000000000,   0.176122795075129,   0.000000000000000}};
532:     TSGLSchemeCreate(2,2,3,3,c,*a,*b,*u,*v,&gl->schemes[gl->nschemes++]);
533:   }
534:   {
535:     /* p=q=3, r=s=4: irks(9/40,0,[1:4]/4,[0.3312 1.0050],[0.49541 1;1 0]) */
536:     const PetscScalar c[4] = {0.25,0.5,0.75,1.0}
537:     ,a[4][4] = {{9./40.               ,                      0,                      0,                      0},
538:                 {2.11286958887701e-01 ,    9./40.             ,                      0,                      0},
539:                 {9.46338294287584e-01 ,  -3.42942861246094e-01,   9./40.              ,                      0},
540:                 {0.521490453970721    ,  -0.662474225622980,   0.490476425459734,   9./40.           }}
541:     ,b[4][4] = {{0.521490453970721    ,  -0.662474225622980,   0.490476425459734,   9./40.           },
542:                 {0.000000000000000    ,   0.000000000000000,   0.000000000000000,   1.000000000000000},
543:                 {-0.084677029310348   ,   1.390757514776085,  -1.568157386206001,   2.023192696767826},
544:                 {0.465383797936408    ,   1.478273530625148,  -1.930836081010182,   1.644872111193354}}
545:     ,u[4][4] = {{1.00000000000000000  ,   0.02500000000001035,  -0.02499999999999053,  -0.00442708333332865},
546:                 {1.00000000000000000  ,   0.06371304111232945,  -0.04032173972189845,  -0.01389438413189452},
547:                 {1.00000000000000000  ,  -0.07839543304147778,   0.04738685705116663,   0.02032603595928376},
548:                 {1.00000000000000000  ,   0.42550734619251651,   0.10800718022400080,  -0.01726712647760034}}
549:     ,v[4][4] = {{1.00000000000000000  ,   0.42550734619251651,   0.10800718022400080,  -0.01726712647760034},
550:                 {0.000000000000000    ,   0.000000000000000,   0.000000000000000,   0.000000000000000},
551:                 {0.000000000000000    ,  -1.761115796027561,  -0.521284157173780,   0.258249384305463},
552:                 {0.000000000000000    ,  -1.657693358744728,  -1.052227765232394,   0.521284157173780}};
553:     TSGLSchemeCreate(3,3,4,4,c,*a,*b,*u,*v,&gl->schemes[gl->nschemes++]);
554:   }
555:   {
556:     /* p=q=4, r=s=5:
557:           irks(3/11,0,[1:5]/5, [0.1715   -0.1238    0.6617],...
558:           [ -0.0812    0.4079    1.0000
559:              1.0000         0         0
560:              0.8270    1.0000         0])
561:     */
562:     const PetscScalar c[5] = {0.2,0.4,0.6,0.8,1.0}
563:     ,a[5][5] = {{2.72727272727352e-01 ,   0.00000000000000e+00,  0.00000000000000e+00 ,  0.00000000000000e+00  ,  0.00000000000000e+00},
564:                 {-1.03980153733431e-01,   2.72727272727405e-01,   0.00000000000000e+00,  0.00000000000000e+00  ,  0.00000000000000e+00},
565:                 {-1.58615400341492e+00,   7.44168951881122e-01,   2.72727272727309e-01,  0.00000000000000e+00  ,  0.00000000000000e+00},
566:                 {-8.73658042865628e-01,   5.37884671894595e-01,  -1.63298538799523e-01,   2.72727272726996e-01 ,  0.00000000000000e+00},
567:                 {2.95489397443992e-01 , -1.18481693910097e+00 , -6.68029812659953e-01 ,  1.00716687860943e+00  , 2.72727272727288e-01}}
568:     ,b[5][5] = {{2.95489397443992e-01 , -1.18481693910097e+00 , -6.68029812659953e-01 ,  1.00716687860943e+00  , 2.72727272727288e-01},
569:                 {0.00000000000000e+00 ,  1.11022302462516e-16 , -2.22044604925031e-16 ,  0.00000000000000e+00  , 1.00000000000000e+00},
570:                 {-4.05882503986005e+00,  -4.00924006567769e+00,  -1.38930610972481e+00,   4.45223930308488e+00 ,  6.32331093108427e-01},
571:                 {8.35690179937017e+00 , -2.26640927349732e+00 ,  6.86647884973826e+00 , -5.22595158025740e+00  , 4.50893068837431e+00},
572:                 {1.27656267027479e+01 ,  2.80882153840821e+00 ,  8.91173096522890e+00 , -1.07936444078906e+01  , 4.82534148988854e+00}}
573:     ,u[5][5] = {{1.00000000000000e+00 , -7.27272727273551e-02 , -3.45454545454419e-02 , -4.12121212119565e-03  ,-2.96969696964014e-04},
574:                 {1.00000000000000e+00 ,  2.31252881006154e-01 , -8.29487834416481e-03 , -9.07191207681020e-03  ,-1.70378403743473e-03},
575:                 {1.00000000000000e+00 ,  1.16925777880663e+00 ,  3.59268562942635e-02 , -4.09013451730615e-02  ,-1.02411119670164e-02},
576:                 {1.00000000000000e+00 ,  1.02634463704356e+00 ,  1.59375044913405e-01 ,  1.89673015035370e-03  ,-4.89987231897569e-03},
577:                 {1.00000000000000e+00 ,  1.27746320298021e+00 ,  2.37186008132728e-01 , -8.28694373940065e-02  ,-5.34396510196430e-02}}
578:     ,v[5][5] = {{1.00000000000000e+00 ,  1.27746320298021e+00 ,  2.37186008132728e-01 , -8.28694373940065e-02  ,-5.34396510196430e-02},
579:                 {0.00000000000000e+00 , -1.77635683940025e-15 , -1.99840144432528e-15 , -9.99200722162641e-16  ,-3.33066907387547e-16},
580:                 {0.00000000000000e+00 ,  4.37280081906924e+00 ,  5.49221645016377e-02 , -8.88913177394943e-02  , 1.12879077989154e-01},
581:                 {0.00000000000000e+00 , -1.22399504837280e+01 , -5.21287338448645e+00 , -8.03952325565291e-01  , 4.60298678047147e-01},
582:                 {0.00000000000000e+00 , -1.85178762883829e+01 , -5.21411849862624e+00 , -1.04283436528809e+00  , 7.49030161063651e-01}};
583:     TSGLSchemeCreate(4,4,5,5,c,*a,*b,*u,*v,&gl->schemes[gl->nschemes++]);
584:   }
585:   {
586:     /* p=q=5, r=s=6;
587:        irks(1/3,0,[1:6]/6,...
588:           [-0.0489    0.4228   -0.8814    0.9021],...
589:           [-0.3474   -0.6617    0.6294    0.2129
590:             0.0044   -0.4256   -0.1427   -0.8936
591:            -0.8267    0.4821    0.1371   -0.2557
592:            -0.4426   -0.3855   -0.7514    0.3014])
593:     */
594:     const PetscScalar c[6] = {1./6, 2./6, 3./6, 4./6, 5./6, 1.}
595:     ,a[6][6] = {{  3.33333333333940e-01,  0                   ,  0                   ,  0                   ,  0                   ,  0                   },
596:                 { -8.64423857333350e-02,  3.33333333332888e-01,  0                   ,  0                   ,  0                   ,  0                   },
597:                 { -2.16850174258252e+00, -2.23619072028839e+00,  3.33333333335204e-01,  0                   ,  0                   ,  0                   },
598:                 { -4.73160970138997e+00, -3.89265344629268e+00, -2.76318716520933e-01,  3.33333333335759e-01,  0                   ,  0                   },
599:                 { -6.75187540297338e+00, -7.90756533769377e+00,  7.90245051802259e-01, -4.48352364517632e-01,  3.33333333328483e-01,  0                   },
600:                 { -4.26488287921548e+00, -1.19320395589302e+01,  3.38924509887755e+00, -2.23969848002481e+00,  6.62807710124007e-01,  3.33333333335440e-01}}
601:     ,b[6][6] = {{ -4.26488287921548e+00, -1.19320395589302e+01,  3.38924509887755e+00, -2.23969848002481e+00,  6.62807710124007e-01,  3.33333333335440e-01},
602:                 { -8.88178419700125e-16,  4.44089209850063e-16, -1.54737334057131e-15, -8.88178419700125e-16,  0.00000000000000e+00,  1.00000000000001e+00},
603:                 { -2.87780425770651e+01, -1.13520448264971e+01,  2.62002318943161e+01,  2.56943874812797e+01, -3.06702268304488e+01,  6.68067773510103e+00},
604:                 {  5.47971245256474e+01,  6.80366875868284e+01, -6.50952588861999e+01, -8.28643975339097e+01,  8.17416943896414e+01, -1.17819043489036e+01},
605:                 { -2.33332114788869e+02,  6.12942539462634e+01, -4.91850135865944e+01,  1.82716844135480e+02, -1.29788173979395e+02,  3.09968095651099e+01},
606:                 { -1.72049132343751e+02,  8.60194713593999e+00,  7.98154219170200e-01,  1.50371386053218e+02, -1.18515423962066e+02,  2.50898277784663e+01}}
607:     ,u[6][6] = {{  1.00000000000000e+00, -1.66666666666870e-01, -4.16666666664335e-02, -3.85802469124815e-03, -2.25051440302250e-04, -9.64506172339142e-06},
608:                 {  1.00000000000000e+00,  8.64423857327162e-02, -4.11484912671353e-02, -1.11450903217645e-02, -1.47651050487126e-03, -1.34395070766826e-04},
609:                 {  1.00000000000000e+00,  4.57135912953434e+00,  1.06514719719137e+00,  1.33517564218007e-01,  1.11365952968659e-02,  6.12382756769504e-04},
610:                 {  1.00000000000000e+00,  9.23391519753404e+00,  2.22431212392095e+00,  2.91823807741891e-01,  2.52058456411084e-02,  1.22800542949647e-03},
611:                 {  1.00000000000000e+00,  1.48175480533865e+01,  3.73439117461835e+00,  5.14648336541804e-01,  4.76430038853402e-02,  2.56798515502156e-03},
612:                 {  1.00000000000000e+00,  1.50512347758335e+01,  4.10099701165164e+00,  5.66039141003603e-01,  3.91213893800891e-02, -2.99136269067853e-03}}
613:     ,v[6][6] = {{  1.00000000000000e+00,  1.50512347758335e+01,  4.10099701165164e+00,  5.66039141003603e-01,  3.91213893800891e-02, -2.99136269067853e-03},
614:                 {  0.00000000000000e+00, -4.88498130835069e-15, -6.43929354282591e-15, -3.55271367880050e-15, -1.22124532708767e-15, -3.12250225675825e-16},
615:                 {  0.00000000000000e+00,  1.22250171233141e+01, -1.77150760606169e+00,  3.54516769879390e-01,  6.22298845883398e-01,  2.31647447450276e-01},
616:                 {  0.00000000000000e+00, -4.48339457331040e+01, -3.57363126641880e-01,  5.18750173123425e-01,  6.55727990241799e-02,  1.63175368287079e-01},
617:                 {  0.00000000000000e+00,  1.37297394708005e+02, -1.60145272991317e+00, -5.05319555199441e+00,  1.55328940390990e-01,  9.16629423682464e-01},
618:                 {  0.00000000000000e+00,  1.05703241119022e+02, -1.16610260983038e+00, -2.99767252773859e+00, -1.13472315553890e-01,  1.09742849254729e+00}};
619:     TSGLSchemeCreate(5,5,6,6,c,*a,*b,*u,*v,&gl->schemes[gl->nschemes++]);
620:   }
621:   return(0);
622: }

626: /*@C
627:    TSGLSetType - sets the class of general linear method to use for time-stepping

629:    Collective on TS

631:    Input Parameters:
632: +  ts - the TS context
633: -  type - a method

635:    Options Database Key:
636: .  -ts_gl_type <type> - sets the method, use -help for a list of available method (e.g. irks)

638:    Notes:
639:    See "petsc/include/petscts.h" for available methods (for instance)
640: .    TSGL_IRKS - Diagonally implicit methods with inherent Runge-Kutta stability (for stiff problems)

642:    Normally, it is best to use the TSSetFromOptions() command and
643:    then set the TSGL type from the options database rather than by using
644:    this routine.  Using the options database provides the user with
645:    maximum flexibility in evaluating the many different solvers.
646:    The TSGLSetType() routine is provided for those situations where it
647:    is necessary to set the timestepping solver independently of the
648:    command line or options database.  This might be the case, for example,
649:    when the choice of solver changes during the execution of the
650:    program, and the user's application is taking responsibility for
651:    choosing the appropriate method.

653:    Level: intermediate

655: .keywords: TS, TSGL, set, type
656: @*/
657: PetscErrorCode  TSGLSetType(TS ts,TSGLType type)
658: {

664:   PetscTryMethod(ts,"TSGLSetType_C",(TS,TSGLType),(ts,type));
665:   return(0);
666: }

670: /*@C
671:    TSGLSetAcceptType - sets the acceptance test

673:    Time integrators that need to control error must have the option to reject a time step based on local error
674:    estimates.  This function allows different schemes to be set.

676:    Logically Collective on TS

678:    Input Parameters:
679: +  ts - the TS context
680: -  type - the type

682:    Options Database Key:
683: .  -ts_gl_accept_type <type> - sets the method used to determine whether to accept or reject a step

685:    Level: intermediate

687: .seealso: TS, TSGL, TSGLAcceptRegister(), TSGLAdapt, set type
688: @*/
689: PetscErrorCode  TSGLSetAcceptType(TS ts,TSGLAcceptType type)
690: {

696:   PetscTryMethod(ts,"TSGLSetAcceptType_C",(TS,TSGLAcceptType),(ts,type));
697:   return(0);
698: }

702: /*@C
703:    TSGLGetAdapt - gets the TSGLAdapt object from the TS

705:    Not Collective

707:    Input Parameter:
708: .  ts - the TS context

710:    Output Parameter:
711: .  adapt - the TSGLAdapt context

713:    Notes:
714:    This allows the user set options on the TSGLAdapt object.  Usually it is better to do this using the options
715:    database, so this function is rarely needed.

717:    Level: advanced

719: .seealso: TSGLAdapt, TSGLAdaptRegister()
720: @*/
721: PetscErrorCode  TSGLGetAdapt(TS ts,TSGLAdapt *adapt)
722: {

728:   PetscUseMethod(ts,"TSGLGetAdapt_C",(TS,TSGLAdapt*),(ts,adapt));
729:   return(0);
730: }

734: static PetscErrorCode TSGLAccept_Always(TS ts,PetscReal tleft,PetscReal h,const PetscReal enorms[],PetscBool  *accept)
735: {
737:   *accept = PETSC_TRUE;
738:   return(0);
739: }

743: static PetscErrorCode TSGLUpdateWRMS(TS ts)
744: {
745:   TS_GL          *gl = (TS_GL*)ts->data;
747:   PetscScalar    *x,*w;
748:   PetscInt       n,i;

751:   VecGetArray(gl->X[0],&x);
752:   VecGetArray(gl->W,&w);
753:   VecGetLocalSize(gl->W,&n);
754:   for (i=0; i<n; i++) w[i] = 1./(gl->wrms_atol + gl->wrms_rtol*PetscAbsScalar(x[i]));
755:   VecRestoreArray(gl->X[0],&x);
756:   VecRestoreArray(gl->W,&w);
757:   return(0);
758: }

762: static PetscErrorCode TSGLVecNormWRMS(TS ts,Vec X,PetscReal *nrm)
763: {
764:   TS_GL          *gl = (TS_GL*)ts->data;
766:   PetscScalar    *x,*w;
767:   PetscReal      sum = 0.0,gsum;
768:   PetscInt       n,N,i;

771:   VecGetArray(X,&x);
772:   VecGetArray(gl->W,&w);
773:   VecGetLocalSize(gl->W,&n);
774:   for (i=0; i<n; i++) sum += PetscAbsScalar(PetscSqr(x[i]*w[i]));
775:   VecRestoreArray(X,&x);
776:   VecRestoreArray(gl->W,&w);
777:   MPIU_Allreduce(&sum,&gsum,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)ts));
778:   VecGetSize(gl->W,&N);
779:   *nrm = PetscSqrtReal(gsum/(1.*N));
780:   return(0);
781: }

785: static PetscErrorCode TSGLSetType_GL(TS ts,TSGLType type)
786: {
787:   PetscErrorCode ierr,(*r)(TS);
788:   PetscBool      same;
789:   TS_GL          *gl = (TS_GL*)ts->data;

792:   if (gl->type_name[0]) {
793:     PetscStrcmp(gl->type_name,type,&same);
794:     if (same) return(0);
795:     (*gl->Destroy)(gl);
796:   }

798:   PetscFunctionListFind(TSGLList,type,&r);
799:   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown TSGL type \"%s\" given",type);
800:   (*r)(ts);
801:   PetscStrcpy(gl->type_name,type);
802:   return(0);
803: }

807: static PetscErrorCode TSGLSetAcceptType_GL(TS ts,TSGLAcceptType type)
808: {
809:   PetscErrorCode     ierr;
810:   TSGLAcceptFunction r;
811:   TS_GL              *gl = (TS_GL*)ts->data;

814:   PetscFunctionListFind(TSGLAcceptList,type,&r);
815:   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown TSGLAccept type \"%s\" given",type);
816:   gl->Accept = r;
817:   PetscStrncpy(gl->accept_name,type,sizeof(gl->accept_name));
818:   return(0);
819: }

823: static PetscErrorCode TSGLGetAdapt_GL(TS ts,TSGLAdapt *adapt)
824: {
826:   TS_GL          *gl = (TS_GL*)ts->data;

829:   if (!gl->adapt) {
830:     TSGLAdaptCreate(PetscObjectComm((PetscObject)ts),&gl->adapt);
831:     PetscObjectIncrementTabLevel((PetscObject)gl->adapt,(PetscObject)ts,1);
832:     PetscLogObjectParent((PetscObject)ts,(PetscObject)gl->adapt);
833:   }
834:   *adapt = gl->adapt;
835:   return(0);
836: }

840: static PetscErrorCode TSGLChooseNextScheme(TS ts,PetscReal h,const PetscReal hmnorm[],PetscInt *next_scheme,PetscReal *next_h,PetscBool  *finish)
841: {
843:   TS_GL          *gl = (TS_GL*)ts->data;
844:   PetscInt       i,n,cur_p,cur,next_sc,candidates[64],orders[64];
845:   PetscReal      errors[64],costs[64],tleft;

848:   cur   = -1;
849:   cur_p = gl->schemes[gl->current_scheme]->p;
850:   tleft = ts->max_time - (ts->ptime + ts->time_step);
851:   for (i=0,n=0; i<gl->nschemes; i++) {
852:     TSGLScheme sc = gl->schemes[i];
853:     if (sc->p < gl->min_order || gl->max_order < sc->p) continue;
854:     if (sc->p == cur_p - 1)    errors[n] = PetscAbsScalar(sc->alpha[0])*hmnorm[0];
855:     else if (sc->p == cur_p)   errors[n] = PetscAbsScalar(sc->alpha[0])*hmnorm[1];
856:     else if (sc->p == cur_p+1) errors[n] = PetscAbsScalar(sc->alpha[0])*(hmnorm[2]+hmnorm[3]);
857:     else continue;
858:     candidates[n] = i;
859:     orders[n]     = PetscMin(sc->p,sc->q); /* order of global truncation error */
860:     costs[n]      = sc->s;                 /* estimate the cost as the number of stages */
861:     if (i == gl->current_scheme) cur = n;
862:     n++;
863:   }
864:   if (cur < 0 || gl->nschemes <= cur) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Current scheme not found in scheme list");
865:   TSGLAdaptChoose(gl->adapt,n,orders,errors,costs,cur,h,tleft,&next_sc,next_h,finish);
866:   *next_scheme = candidates[next_sc];
867:   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);
868:   return(0);
869: }

873: static PetscErrorCode TSGLGetMaxSizes(TS ts,PetscInt *max_r,PetscInt *max_s)
874: {
875:   TS_GL *gl = (TS_GL*)ts->data;

878:   *max_r = gl->schemes[gl->nschemes-1]->r;
879:   *max_s = gl->schemes[gl->nschemes-1]->s;
880:   return(0);
881: }

885: static PetscErrorCode TSSolve_GL(TS ts)
886: {
887:   TS_GL               *gl = (TS_GL*)ts->data;
888:   PetscInt            i,k,its,lits,max_r,max_s;
889:   PetscBool           final_step,finish;
890:   SNESConvergedReason snesreason;
891:   PetscErrorCode      ierr;

894:   TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);

896:   TSGLGetMaxSizes(ts,&max_r,&max_s);
897:   VecCopy(ts->vec_sol,gl->X[0]);
898:   for (i=1; i<max_r; i++) {
899:     VecZeroEntries(gl->X[i]);
900:   }
901:   TSGLUpdateWRMS(ts);

903:   if (0) {
904:     /* Find consistent initial data for DAE */
905:     gl->stage_time = ts->ptime + ts->time_step;
906:     gl->scoeff = 1.;
907:     gl->stage  = 0;

909:     VecCopy(ts->vec_sol,gl->Z);
910:     VecCopy(ts->vec_sol,gl->Y);
911:     SNESSolve(ts->snes,NULL,gl->Y);
912:     SNESGetIterationNumber(ts->snes,&its);
913:     SNESGetLinearSolveIterations(ts->snes,&lits);
914:     SNESGetConvergedReason(ts->snes,&snesreason);

916:     ts->snes_its += its; ts->ksp_its += lits;
917:     if (snesreason < 0 && ts->max_snes_failures > 0 && ++ts->num_snes_failures >= ts->max_snes_failures) {
918:       ts->reason = TS_DIVERGED_NONLINEAR_SOLVE;
919:       PetscInfo2(ts,"Step=%D, nonlinear solve solve failures %D greater than current TS allowed, stopping solve\n",ts->steps,ts->num_snes_failures);
920:       return(0);
921:     }
922:   }

924:   if (gl->current_scheme < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"A starting scheme has not been provided");

926:   for (k=0,final_step=PETSC_FALSE,finish=PETSC_FALSE; k<ts->max_steps && !finish; k++) {
927:     PetscInt          j,r,s,next_scheme = 0,rejections;
928:     PetscReal         h,hmnorm[4],enorm[3],next_h;
929:     PetscBool         accept;
930:     const PetscScalar *c,*a,*u;
931:     Vec               *X,*Ydot,Y;
932:     TSGLScheme        scheme = gl->schemes[gl->current_scheme];

934:     r = scheme->r; s = scheme->s;
935:     c = scheme->c;
936:     a = scheme->a; u = scheme->u;
937:     h = ts->time_step;
938:     X = gl->X; Ydot = gl->Ydot; Y = gl->Y;

940:     if (ts->ptime > ts->max_time) break;

942:     /*
943:       We only call PreStep at the start of each STEP, not each STAGE.  This is because it is
944:       possible to fail (have to restart a step) after multiple stages.
945:     */
946:     TSPreStep(ts);

948:     rejections = 0;
949:     while (1) {
950:       for (i=0; i<s; i++) {
951:         PetscScalar shift;
952:         gl->scoeff     = 1./PetscRealPart(a[i*s+i]);
953:         shift          = gl->scoeff/ts->time_step;
954:         gl->stage      = i;
955:         gl->stage_time = ts->ptime + PetscRealPart(c[i])*h;

957:         /*
958:         * Stage equation: Y = h A Y' + U X
959:         * We assume that A is lower-triangular so that we can solve the stages (Y,Y') sequentially
960:         * Build the affine vector z_i = -[1/(h a_ii)](h sum_j a_ij y'_j + sum_j u_ij x_j)
961:         * Then y'_i = z + 1/(h a_ii) y_i
962:         */
963:         VecZeroEntries(gl->Z);
964:         for (j=0; j<r; j++) {
965:           VecAXPY(gl->Z,-shift*u[i*r+j],X[j]);
966:         }
967:         for (j=0; j<i; j++) {
968:           VecAXPY(gl->Z,-shift*h*a[i*s+j],Ydot[j]);
969:         }
970:         /* Note: Z is used within function evaluation, Ydot = Z + shift*Y */

972:         /* Compute an estimate of Y to start Newton iteration */
973:         if (gl->extrapolate) {
974:           if (i==0) {
975:             /* Linear extrapolation on the first stage */
976:             VecWAXPY(Y,c[i]*h,X[1],X[0]);
977:           } else {
978:             /* Linear extrapolation from the last stage */
979:             VecAXPY(Y,(c[i]-c[i-1])*h,Ydot[i-1]);
980:           }
981:         } else if (i==0) {        /* Directly use solution from the last step, otherwise reuse the last stage (do nothing) */
982:           VecCopy(X[0],Y);
983:         }

985:         /* Solve this stage (Ydot[i] is computed during function evaluation) */
986:         SNESSolve(ts->snes,NULL,Y);
987:         SNESGetIterationNumber(ts->snes,&its);
988:         SNESGetLinearSolveIterations(ts->snes,&lits);
989:         SNESGetConvergedReason(ts->snes,&snesreason);
990:         ts->snes_its += its; ts->ksp_its += lits;
991:         if (snesreason < 0 && ts->max_snes_failures > 0 && ++ts->num_snes_failures >= ts->max_snes_failures) {
992:           ts->reason = TS_DIVERGED_NONLINEAR_SOLVE;
993:           PetscInfo2(ts,"Step=%D, nonlinear solve solve failures %D greater than current TS allowed, stopping solve\n",ts->steps,ts->num_snes_failures);
994:           return(0);
995:         }
996:       }

998:       gl->stage_time = ts->ptime + ts->time_step;

1000:       (*gl->EstimateHigherMoments)(scheme,h,Ydot,gl->X,gl->himom);
1001:       /* hmnorm[i] = h^{p+i}x^{(p+i)} with i=0,1,2; hmnorm[3] = h^{p+2}(dx'/dx) x^{(p+1)} */
1002:       for (i=0; i<3; i++) {
1003:         TSGLVecNormWRMS(ts,gl->himom[i],&hmnorm[i+1]);
1004:       }
1005:       enorm[0] = PetscRealPart(scheme->alpha[0])*hmnorm[1];
1006:       enorm[1] = PetscRealPart(scheme->beta[0]) *hmnorm[2];
1007:       enorm[2] = PetscRealPart(scheme->gamma[0])*hmnorm[3];
1008:       (*gl->Accept)(ts,ts->max_time-gl->stage_time,h,enorm,&accept);
1009:       if (accept) goto accepted;
1010:       rejections++;
1011:       PetscInfo3(ts,"Step %D (t=%g) not accepted, rejections=%D\n",k,gl->stage_time,rejections);
1012:       if (rejections > gl->max_step_rejections) break;
1013:       /*
1014:         There are lots of reasons why a step might be rejected, including solvers not converging and other factors that
1015:         TSGLChooseNextScheme does not support.  Additionally, the error estimates may be very screwed up, so I'm not
1016:         convinced that it's safe to just compute a new error estimate using the same interface as the current adaptor
1017:         (the adaptor interface probably has to change).  Here we make an arbitrary and naive choice.  This assumes that
1018:         steps were written in Nordsieck form.  The "correct" method would be to re-complete the previous time step with
1019:         the correct "next" step size.  It is unclear to me whether the present ad-hoc method of rescaling X is stable.
1020:       */
1021:       h *= 0.5;
1022:       for (i=1; i<scheme->r; i++) {
1023:         VecScale(X[i],PetscPowRealInt(0.5,i));
1024:       }
1025:     }
1026:     SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_CONV_FAILED,"Time step %D (t=%g) not accepted after %D failures\n",k,gl->stage_time,rejections);

1028: accepted:
1029:     /* This term is not error, but it *would* be the leading term for a lower order method */
1030:     TSGLVecNormWRMS(ts,gl->X[scheme->r-1],&hmnorm[0]);
1031:     /* Correct scaling so that these are equivalent to norms of the Nordsieck vectors */

1033:     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]);
1034:     if (!final_step) {
1035:       TSGLChooseNextScheme(ts,h,hmnorm,&next_scheme,&next_h,&final_step);
1036:     } else {
1037:       /* Dummy values to complete the current step in a consistent manner */
1038:       next_scheme = gl->current_scheme;
1039:       next_h      = h;
1040:       finish      = PETSC_TRUE;
1041:     }

1043:     X        = gl->Xold;
1044:     gl->Xold = gl->X;
1045:     gl->X    = X;
1046:     (*gl->CompleteStep)(scheme,h,gl->schemes[next_scheme],next_h,Ydot,gl->Xold,gl->X);

1048:     TSGLUpdateWRMS(ts);

1050:     /* Post the solution for the user, we could avoid this copy with a small bit of cleverness */
1051:     VecCopy(gl->X[0],ts->vec_sol);
1052:     ts->ptime += h;
1053:     ts->steps++;
1054:     ts->total_steps++;

1056:     TSPostStep(ts);
1057:     TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);

1059:     gl->current_scheme = next_scheme;
1060:     ts->time_step      = next_h;
1061:   }
1062:   return(0);
1063: }

1065: /*------------------------------------------------------------*/

1069: static PetscErrorCode TSReset_GL(TS ts)
1070: {
1071:   TS_GL          *gl = (TS_GL*)ts->data;
1072:   PetscInt       max_r,max_s;

1076:   if (gl->setupcalled) {
1077:     TSGLGetMaxSizes(ts,&max_r,&max_s);
1078:     VecDestroyVecs(max_r,&gl->Xold);
1079:     VecDestroyVecs(max_r,&gl->X);
1080:     VecDestroyVecs(max_s,&gl->Ydot);
1081:     VecDestroyVecs(3,&gl->himom);
1082:     VecDestroy(&gl->W);
1083:     VecDestroy(&gl->Y);
1084:     VecDestroy(&gl->Z);
1085:   }
1086:   gl->setupcalled = PETSC_FALSE;
1087:   return(0);
1088: }

1092: static PetscErrorCode TSDestroy_GL(TS ts)
1093: {
1094:   TS_GL          *gl = (TS_GL*)ts->data;

1098:   TSReset_GL(ts);
1099:   if (gl->adapt) {TSGLAdaptDestroy(&gl->adapt);}
1100:   if (gl->Destroy) {(*gl->Destroy)(gl);}
1101:   PetscFree(ts->data);
1102:   PetscObjectComposeFunction((PetscObject)ts,"TSGLSetType_C",NULL);
1103:   PetscObjectComposeFunction((PetscObject)ts,"TSGLSetAcceptType_C",NULL);
1104:   PetscObjectComposeFunction((PetscObject)ts,"TSGLGetAdapt_C",NULL);
1105:   return(0);
1106: }

1108: /*
1109:     This defines the nonlinear equation that is to be solved with SNES
1110:     g(x) = f(t,x,z+shift*x) = 0
1111: */
1114: static PetscErrorCode SNESTSFormFunction_GL(SNES snes,Vec x,Vec f,TS ts)
1115: {
1116:   TS_GL          *gl = (TS_GL*)ts->data;
1118:   Vec            Z,Ydot;
1119:   DM             dm,dmsave;

1122:   SNESGetDM(snes,&dm);
1123:   TSGLGetVecs(ts,dm,&Z,&Ydot);
1124:   VecWAXPY(Ydot,gl->scoeff/ts->time_step,x,Z);
1125:   dmsave = ts->dm;
1126:   ts->dm = dm;
1127:   TSComputeIFunction(ts,gl->stage_time,x,Ydot,f,PETSC_FALSE);
1128:   ts->dm = dmsave;
1129:   TSGLRestoreVecs(ts,dm,&Z,&Ydot);
1130:   return(0);
1131: }

1135: static PetscErrorCode SNESTSFormJacobian_GL(SNES snes,Vec x,Mat A,Mat B,TS ts)
1136: {
1137:   TS_GL          *gl = (TS_GL*)ts->data;
1139:   Vec            Z,Ydot;
1140:   DM             dm,dmsave;

1143:   SNESGetDM(snes,&dm);
1144:   TSGLGetVecs(ts,dm,&Z,&Ydot);
1145:   dmsave = ts->dm;
1146:   ts->dm = dm;
1147:   /* gl->Xdot will have already been computed in SNESTSFormFunction_GL */
1148:   TSComputeIJacobian(ts,gl->stage_time,x,gl->Ydot[gl->stage],gl->scoeff/ts->time_step,A,B,PETSC_FALSE);
1149:   ts->dm = dmsave;
1150:   TSGLRestoreVecs(ts,dm,&Z,&Ydot);
1151:   return(0);
1152: }


1157: static PetscErrorCode TSSetUp_GL(TS ts)
1158: {
1159:   TS_GL          *gl = (TS_GL*)ts->data;
1160:   PetscInt       max_r,max_s;
1162:   DM             dm;

1165:   gl->setupcalled = PETSC_TRUE;
1166:   TSGLGetMaxSizes(ts,&max_r,&max_s);
1167:   VecDuplicateVecs(ts->vec_sol,max_r,&gl->X);
1168:   VecDuplicateVecs(ts->vec_sol,max_r,&gl->Xold);
1169:   VecDuplicateVecs(ts->vec_sol,max_s,&gl->Ydot);
1170:   VecDuplicateVecs(ts->vec_sol,3,&gl->himom);
1171:   VecDuplicate(ts->vec_sol,&gl->W);
1172:   VecDuplicate(ts->vec_sol,&gl->Y);
1173:   VecDuplicate(ts->vec_sol,&gl->Z);

1175:   /* Default acceptance tests and adaptivity */
1176:   if (!gl->Accept) {TSGLSetAcceptType(ts,TSGLACCEPT_ALWAYS);}
1177:   if (!gl->adapt)  {TSGLGetAdapt(ts,&gl->adapt);}

1179:   if (gl->current_scheme < 0) {
1180:     PetscInt i;
1181:     for (i=0;; i++) {
1182:       if (gl->schemes[i]->p == gl->start_order) break;
1183:       if (i+1 == gl->nschemes) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"No schemes available with requested start order %d",i);
1184:     }
1185:     gl->current_scheme = i;
1186:   }
1187:   TSGetDM(ts,&dm);
1188:   if (dm) {
1189:     DMCoarsenHookAdd(dm,DMCoarsenHook_TSGL,DMRestrictHook_TSGL,ts);
1190:     DMSubDomainHookAdd(dm,DMSubDomainHook_TSGL,DMSubDomainRestrictHook_TSGL,ts);
1191:   }
1192:   return(0);
1193: }
1194: /*------------------------------------------------------------*/

1198: static PetscErrorCode TSSetFromOptions_GL(PetscOptionItems *PetscOptionsObject,TS ts)
1199: {
1200:   TS_GL          *gl        = (TS_GL*)ts->data;
1201:   char           tname[256] = TSGL_IRKS,completef[256] = "rescale-and-modify";

1205:   PetscOptionsHead(PetscOptionsObject,"General Linear ODE solver options");
1206:   {
1207:     PetscBool flg;
1208:     PetscOptionsFList("-ts_gl_type","Type of GL method","TSGLSetType",TSGLList,gl->type_name[0] ? gl->type_name : tname,tname,sizeof(tname),&flg);
1209:     if (flg || !gl->type_name[0]) {
1210:       TSGLSetType(ts,tname);
1211:     }
1212:     PetscOptionsInt("-ts_gl_max_step_rejections","Maximum number of times to attempt a step","None",gl->max_step_rejections,&gl->max_step_rejections,NULL);
1213:     PetscOptionsInt("-ts_gl_max_order","Maximum order to try","TSGLSetMaxOrder",gl->max_order,&gl->max_order,NULL);
1214:     PetscOptionsInt("-ts_gl_min_order","Minimum order to try","TSGLSetMinOrder",gl->min_order,&gl->min_order,NULL);
1215:     PetscOptionsInt("-ts_gl_start_order","Initial order to try","TSGLSetMinOrder",gl->start_order,&gl->start_order,NULL);
1216:     PetscOptionsEnum("-ts_gl_error_direction","Which direction to look when estimating error","TSGLSetErrorDirection",TSGLErrorDirections,(PetscEnum)gl->error_direction,(PetscEnum*)&gl->error_direction,NULL);
1217:     PetscOptionsBool("-ts_gl_extrapolate","Extrapolate stage solution from previous solution (sometimes unstable)","TSGLSetExtrapolate",gl->extrapolate,&gl->extrapolate,NULL);
1218:     PetscOptionsReal("-ts_gl_atol","Absolute tolerance","TSGLSetTolerances",gl->wrms_atol,&gl->wrms_atol,NULL);
1219:     PetscOptionsReal("-ts_gl_rtol","Relative tolerance","TSGLSetTolerances",gl->wrms_rtol,&gl->wrms_rtol,NULL);
1220:     PetscOptionsString("-ts_gl_complete","Method to use for completing the step","none",completef,completef,sizeof(completef),&flg);
1221:     if (flg) {
1222:       PetscBool match1,match2;
1223:       PetscStrcmp(completef,"rescale",&match1);
1224:       PetscStrcmp(completef,"rescale-and-modify",&match2);
1225:       if (match1)      gl->CompleteStep = TSGLCompleteStep_Rescale;
1226:       else if (match2) gl->CompleteStep = TSGLCompleteStep_RescaleAndModify;
1227:       else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"%s",completef);
1228:     }
1229:     {
1230:       char type[256] = TSGLACCEPT_ALWAYS;
1231:       PetscOptionsFList("-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);
1232:       if (flg || !gl->accept_name[0]) {
1233:         TSGLSetAcceptType(ts,type);
1234:       }
1235:     }
1236:     {
1237:       TSGLAdapt adapt;
1238:       TSGLGetAdapt(ts,&adapt);
1239:       TSGLAdaptSetFromOptions(PetscOptionsObject,adapt);
1240:     }
1241:   }
1242:   PetscOptionsTail();
1243:   return(0);
1244: }

1248: static PetscErrorCode TSView_GL(TS ts,PetscViewer viewer)
1249: {
1250:   TS_GL          *gl = (TS_GL*)ts->data;
1251:   PetscInt       i;
1252:   PetscBool      iascii,details;

1256:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
1257:   if (iascii) {
1258:     PetscViewerASCIIPrintf(viewer,"  min order %D, max order %D, current order %D\n",gl->min_order,gl->max_order,gl->schemes[gl->current_scheme]->p);
1259:     PetscViewerASCIIPrintf(viewer,"  Error estimation: %s\n",TSGLErrorDirections[gl->error_direction]);
1260:     PetscViewerASCIIPrintf(viewer,"  Extrapolation: %s\n",gl->extrapolate ? "yes" : "no");
1261:     PetscViewerASCIIPrintf(viewer,"  Acceptance test: %s\n",gl->accept_name[0] ? gl->accept_name : "(not yet set)");
1262:     PetscViewerASCIIPushTab(viewer);
1263:     TSGLAdaptView(gl->adapt,viewer);
1264:     PetscViewerASCIIPopTab(viewer);
1265:     PetscViewerASCIIPrintf(viewer,"  type: %s\n",gl->type_name[0] ? gl->type_name : "(not yet set)");
1266:     PetscViewerASCIIPrintf(viewer,"Schemes within family (%d):\n",gl->nschemes);
1267:     details = PETSC_FALSE;
1268:     PetscOptionsGetBool(((PetscObject)ts)->options,((PetscObject)ts)->prefix,"-ts_gl_view_detailed",&details,NULL);
1269:     PetscViewerASCIIPushTab(viewer);
1270:     for (i=0; i<gl->nschemes; i++) {
1271:       TSGLSchemeView(gl->schemes[i],details,viewer);
1272:     }
1273:     if (gl->View) {
1274:       (*gl->View)(gl,viewer);
1275:     }
1276:     PetscViewerASCIIPopTab(viewer);
1277:   }
1278:   SNESView(ts->snes,viewer);
1279:   return(0);
1280: }

1284: /*@C
1285:    TSGLRegister -  adds a TSGL implementation

1287:    Not Collective

1289:    Input Parameters:
1290: +  name_scheme - name of user-defined general linear scheme
1291: -  routine_create - routine to create method context

1293:    Notes:
1294:    TSGLRegister() may be called multiple times to add several user-defined families.

1296:    Sample usage:
1297: .vb
1298:    TSGLRegister("my_scheme",MySchemeCreate);
1299: .ve

1301:    Then, your scheme can be chosen with the procedural interface via
1302: $     TSGLSetType(ts,"my_scheme")
1303:    or at runtime via the option
1304: $     -ts_gl_type my_scheme

1306:    Level: advanced

1308: .keywords: TSGL, register

1310: .seealso: TSGLRegisterAll()
1311: @*/
1312: PetscErrorCode  TSGLRegister(const char sname[],PetscErrorCode (*function)(TS))
1313: {

1317:   PetscFunctionListAdd(&TSGLList,sname,function);
1318:   return(0);
1319: }

1323: /*@C
1324:    TSGLAcceptRegister -  adds a TSGL acceptance scheme

1326:    Not Collective

1328:    Input Parameters:
1329: +  name_scheme - name of user-defined acceptance scheme
1330: -  routine_create - routine to create method context

1332:    Notes:
1333:    TSGLAcceptRegister() may be called multiple times to add several user-defined families.

1335:    Sample usage:
1336: .vb
1337:    TSGLAcceptRegister("my_scheme",MySchemeCreate);
1338: .ve

1340:    Then, your scheme can be chosen with the procedural interface via
1341: $     TSGLSetAcceptType(ts,"my_scheme")
1342:    or at runtime via the option
1343: $     -ts_gl_accept_type my_scheme

1345:    Level: advanced

1347: .keywords: TSGL, TSGLAcceptType, register

1349: .seealso: TSGLRegisterAll()
1350: @*/
1351: PetscErrorCode  TSGLAcceptRegister(const char sname[],TSGLAcceptFunction function)
1352: {

1356:   PetscFunctionListAdd(&TSGLAcceptList,sname,function);
1357:   return(0);
1358: }

1362: /*@C
1363:   TSGLRegisterAll - Registers all of the general linear methods in TSGL

1365:   Not Collective

1367:   Level: advanced

1369: .keywords: TS, TSGL, register, all

1371: .seealso:  TSGLRegisterDestroy()
1372: @*/
1373: PetscErrorCode  TSGLRegisterAll(void)
1374: {

1378:   if (TSGLRegisterAllCalled) return(0);
1379:   TSGLRegisterAllCalled = PETSC_TRUE;

1381:   TSGLRegister(TSGL_IRKS,              TSGLCreate_IRKS);
1382:   TSGLAcceptRegister(TSGLACCEPT_ALWAYS,TSGLAccept_Always);
1383:   return(0);
1384: }

1388: /*@C
1389:   TSGLInitializePackage - This function initializes everything in the TSGL package. It is called
1390:   from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_GL()
1391:   when using static libraries.

1393:   Level: developer

1395: .keywords: TS, TSGL, initialize, package
1396: .seealso: PetscInitialize()
1397: @*/
1398: PetscErrorCode  TSGLInitializePackage(void)
1399: {

1403:   if (TSGLPackageInitialized) return(0);
1404:   TSGLPackageInitialized = PETSC_TRUE;
1405:   TSGLRegisterAll();
1406:   PetscRegisterFinalize(TSGLFinalizePackage);
1407:   return(0);
1408: }

1412: /*@C
1413:   TSGLFinalizePackage - This function destroys everything in the TSGL package. It is
1414:   called from PetscFinalize().

1416:   Level: developer

1418: .keywords: Petsc, destroy, package
1419: .seealso: PetscFinalize()
1420: @*/
1421: PetscErrorCode  TSGLFinalizePackage(void)
1422: {

1426:   PetscFunctionListDestroy(&TSGLList);
1427:   PetscFunctionListDestroy(&TSGLAcceptList);
1428:   TSGLPackageInitialized = PETSC_FALSE;
1429:   TSGLRegisterAllCalled  = PETSC_FALSE;
1430:   return(0);
1431: }

1433: /* ------------------------------------------------------------ */
1434: /*MC
1435:       TSGL - DAE solver using implicit General Linear methods

1437:   These methods contain Runge-Kutta and multistep schemes as special cases.  These special cases have some fundamental
1438:   limitations.  For example, diagonally implicit Runge-Kutta cannot have stage order greater than 1 which limits their
1439:   applicability to very stiff systems.  Meanwhile, multistep methods cannot be A-stable for order greater than 2 and BDF
1440:   are not 0-stable for order greater than 6.  GL methods can be A- and L-stable with arbitrarily high stage order and
1441:   reliable error estimates for both 1 and 2 orders higher to facilitate adaptive step sizes and adaptive order schemes.
1442:   All this is possible while preserving a singly diagonally implicit structure.

1444:   Options database keys:
1445: +  -ts_gl_type <type> - the class of general linear method (irks)
1446: .  -ts_gl_rtol <tol>  - relative error
1447: .  -ts_gl_atol <tol>  - absolute error
1448: .  -ts_gl_min_order <p> - minimum order method to consider (default=1)
1449: .  -ts_gl_max_order <p> - maximum order method to consider (default=3)
1450: .  -ts_gl_start_order <p> - order of starting method (default=1)
1451: .  -ts_gl_complete <method> - method to use for completing the step (rescale-and-modify or rescale)
1452: -  -ts_adapt_type <method> - adaptive controller to use (none step both)

1454:   Notes:
1455:   This integrator can be applied to DAE.

1457:   Diagonally implicit general linear (DIGL) methods are a generalization of diagonally implicit Runge-Kutta (DIRK).
1458:   They are represented by the tableau

1460: .vb
1461:   A  |  U
1462:   -------
1463:   B  |  V
1464: .ve

1466:   combined with a vector c of abscissa.  "Diagonally implicit" means that A is lower triangular.
1467:   A step of the general method reads

1469: .vb
1470:   [ Y ] = [A  U] [  Y'   ]
1471:   [X^k] = [B  V] [X^{k-1}]
1472: .ve

1474:   where Y is the multivector of stage values, Y' is the multivector of stage derivatives, X^k is the Nordsieck vector of
1475:   the solution at step k.  The Nordsieck vector consists of the first r moments of the solution, given by

1477: .vb
1478:   X = [x_0,x_1,...,x_{r-1}] = [x, h x', h^2 x'', ..., h^{r-1} x^{(r-1)} ]
1479: .ve

1481:   If A is lower triangular, we can solve the stages (Y,Y') sequentially

1483: .vb
1484:   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}
1485: .ve

1487:   and then construct the pieces to carry to the next step

1489: .vb
1490:   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}
1491: .ve

1493:   Note that when the equations are cast in implicit form, we are using the stage equation to define y'_i
1494:   in terms of y_i and known stuff (y_j for j<i and x_j for all j).


1497:   Error estimation

1499:   At present, the most attractive GL methods for stiff problems are singly diagonally implicit schemes which posses
1500:   Inherent Runge-Kutta Stability (IRKS).  These methods have r=s, the number of items passed between steps is equal to
1501:   the number of stages.  The order and stage-order are one less than the number of stages.  We use the error estimates
1502:   in the 2007 paper which provide the following estimates

1504: .vb
1505:   h^{p+1} X^{(p+1)}          = phi_0^T Y' + [0 psi_0^T] Xold
1506:   h^{p+2} X^{(p+2)}          = phi_1^T Y' + [0 psi_1^T] Xold
1507:   h^{p+2} (dx'/dx) X^{(p+1)} = phi_2^T Y' + [0 psi_2^T] Xold
1508: .ve

1510:   These estimates are accurate to O(h^{p+3}).

1512:   Changing the step size

1514:   We use the generalized "rescale and modify" scheme, see equation (4.5) of the 2007 paper.

1516:   Level: beginner

1518:   References:
1519: +  1. - John Butcher and Z. Jackieweicz and W. Wright, On error propagation in general linear methods for
1520:   ordinary differential equations, Journal of Complexity, Vol 23, 2007.
1521: -  2. - John Butcher, Numerical methods for ordinary differential equations, second edition, Wiley, 2009.

1523: .seealso:  TSCreate(), TS, TSSetType()

1525: M*/
1528: PETSC_EXTERN PetscErrorCode TSCreate_GL(TS ts)
1529: {
1530:   TS_GL          *gl;

1534:   TSGLInitializePackage();

1536:   PetscNewLog(ts,&gl);
1537:   ts->data = (void*)gl;

1539:   ts->ops->reset          = TSReset_GL;
1540:   ts->ops->destroy        = TSDestroy_GL;
1541:   ts->ops->view           = TSView_GL;
1542:   ts->ops->setup          = TSSetUp_GL;
1543:   ts->ops->solve          = TSSolve_GL;
1544:   ts->ops->setfromoptions = TSSetFromOptions_GL;
1545:   ts->ops->snesfunction   = SNESTSFormFunction_GL;
1546:   ts->ops->snesjacobian   = SNESTSFormJacobian_GL;

1548:   gl->max_step_rejections = 1;
1549:   gl->min_order           = 1;
1550:   gl->max_order           = 3;
1551:   gl->start_order         = 1;
1552:   gl->current_scheme      = -1;
1553:   gl->extrapolate         = PETSC_FALSE;

1555:   gl->wrms_atol = 1e-8;
1556:   gl->wrms_rtol = 1e-5;

1558:   PetscObjectComposeFunction((PetscObject)ts,"TSGLSetType_C",      &TSGLSetType_GL);
1559:   PetscObjectComposeFunction((PetscObject)ts,"TSGLSetAcceptType_C",&TSGLSetAcceptType_GL);
1560:   PetscObjectComposeFunction((PetscObject)ts,"TSGLGetAdapt_C",     &TSGLGetAdapt_GL);
1561:   return(0);
1562: }