xref: /petsc/src/ts/impls/implicit/glle/glle.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
13d177a5cSEmil Constantinescu 
23d177a5cSEmil Constantinescu #include <../src/ts/impls/implicit/glle/glle.h>                /*I   "petscts.h"   I*/
33d177a5cSEmil Constantinescu #include <petscdm.h>
43d177a5cSEmil Constantinescu #include <petscblaslapack.h>
53d177a5cSEmil Constantinescu 
6c793f718SLisandro Dalcin static const char        *TSGLLEErrorDirections[] = {"FORWARD","BACKWARD","TSGLLEErrorDirection","TSGLLEERROR_",NULL};
73d177a5cSEmil Constantinescu static PetscFunctionList TSGLLEList;
83d177a5cSEmil Constantinescu static PetscFunctionList TSGLLEAcceptList;
93d177a5cSEmil Constantinescu static PetscBool         TSGLLEPackageInitialized;
103d177a5cSEmil Constantinescu static PetscBool         TSGLLERegisterAllCalled;
113d177a5cSEmil Constantinescu 
123d177a5cSEmil Constantinescu /* This function is pure */
133d177a5cSEmil Constantinescu static PetscScalar Factorial(PetscInt n)
143d177a5cSEmil Constantinescu {
153d177a5cSEmil Constantinescu   PetscInt i;
163d177a5cSEmil Constantinescu   if (n < 12) {                 /* Can compute with 32-bit integers */
173d177a5cSEmil Constantinescu     PetscInt f = 1;
183d177a5cSEmil Constantinescu     for (i=2; i<=n; i++) f *= i;
193d177a5cSEmil Constantinescu     return (PetscScalar)f;
203d177a5cSEmil Constantinescu   } else {
213d177a5cSEmil Constantinescu     PetscScalar f = 1.;
223d177a5cSEmil Constantinescu     for (i=2; i<=n; i++) f *= (PetscScalar)i;
233d177a5cSEmil Constantinescu     return f;
243d177a5cSEmil Constantinescu   }
253d177a5cSEmil Constantinescu }
263d177a5cSEmil Constantinescu 
273d177a5cSEmil Constantinescu /* This function is pure */
283d177a5cSEmil Constantinescu static PetscScalar CPowF(PetscScalar c,PetscInt p)
293d177a5cSEmil Constantinescu {
303d177a5cSEmil Constantinescu   return PetscPowRealInt(PetscRealPart(c),p)/Factorial(p);
313d177a5cSEmil Constantinescu }
323d177a5cSEmil Constantinescu 
333d177a5cSEmil Constantinescu static PetscErrorCode TSGLLEGetVecs(TS ts,DM dm,Vec *Z,Vec *Ydotstage)
343d177a5cSEmil Constantinescu {
353d177a5cSEmil Constantinescu   TS_GLLE        *gl = (TS_GLLE*)ts->data;
363d177a5cSEmil Constantinescu 
373d177a5cSEmil Constantinescu   PetscFunctionBegin;
383d177a5cSEmil Constantinescu   if (Z) {
393d177a5cSEmil Constantinescu     if (dm && dm != ts->dm) {
40*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMGetNamedGlobalVector(dm,"TSGLLE_Z",Z));
413d177a5cSEmil Constantinescu     } else *Z = gl->Z;
423d177a5cSEmil Constantinescu   }
433d177a5cSEmil Constantinescu   if (Ydotstage) {
443d177a5cSEmil Constantinescu     if (dm && dm != ts->dm) {
45*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMGetNamedGlobalVector(dm,"TSGLLE_Ydot",Ydotstage));
463d177a5cSEmil Constantinescu     } else *Ydotstage = gl->Ydot[gl->stage];
473d177a5cSEmil Constantinescu   }
483d177a5cSEmil Constantinescu   PetscFunctionReturn(0);
493d177a5cSEmil Constantinescu }
503d177a5cSEmil Constantinescu 
513d177a5cSEmil Constantinescu static PetscErrorCode TSGLLERestoreVecs(TS ts,DM dm,Vec *Z,Vec *Ydotstage)
523d177a5cSEmil Constantinescu {
533d177a5cSEmil Constantinescu   PetscFunctionBegin;
543d177a5cSEmil Constantinescu   if (Z) {
553d177a5cSEmil Constantinescu     if (dm && dm != ts->dm) {
56*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMRestoreNamedGlobalVector(dm,"TSGLLE_Z",Z));
573d177a5cSEmil Constantinescu     }
583d177a5cSEmil Constantinescu   }
593d177a5cSEmil Constantinescu   if (Ydotstage) {
603d177a5cSEmil Constantinescu 
613d177a5cSEmil Constantinescu     if (dm && dm != ts->dm) {
62*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMRestoreNamedGlobalVector(dm,"TSGLLE_Ydot",Ydotstage));
633d177a5cSEmil Constantinescu     }
643d177a5cSEmil Constantinescu   }
653d177a5cSEmil Constantinescu   PetscFunctionReturn(0);
663d177a5cSEmil Constantinescu }
673d177a5cSEmil Constantinescu 
683d177a5cSEmil Constantinescu static PetscErrorCode DMCoarsenHook_TSGLLE(DM fine,DM coarse,void *ctx)
693d177a5cSEmil Constantinescu {
703d177a5cSEmil Constantinescu   PetscFunctionBegin;
713d177a5cSEmil Constantinescu   PetscFunctionReturn(0);
723d177a5cSEmil Constantinescu }
733d177a5cSEmil Constantinescu 
743d177a5cSEmil Constantinescu static PetscErrorCode DMRestrictHook_TSGLLE(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
753d177a5cSEmil Constantinescu {
763d177a5cSEmil Constantinescu   TS             ts = (TS)ctx;
773d177a5cSEmil Constantinescu   Vec            Ydot,Ydot_c;
783d177a5cSEmil Constantinescu 
793d177a5cSEmil Constantinescu   PetscFunctionBegin;
80*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGLLEGetVecs(ts,fine,NULL,&Ydot));
81*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGLLEGetVecs(ts,coarse,NULL,&Ydot_c));
82*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatRestrict(restrct,Ydot,Ydot_c));
83*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecPointwiseMult(Ydot_c,rscale,Ydot_c));
84*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGLLERestoreVecs(ts,fine,NULL,&Ydot));
85*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGLLERestoreVecs(ts,coarse,NULL,&Ydot_c));
863d177a5cSEmil Constantinescu   PetscFunctionReturn(0);
873d177a5cSEmil Constantinescu }
883d177a5cSEmil Constantinescu 
893d177a5cSEmil Constantinescu static PetscErrorCode DMSubDomainHook_TSGLLE(DM dm,DM subdm,void *ctx)
903d177a5cSEmil Constantinescu {
913d177a5cSEmil Constantinescu   PetscFunctionBegin;
923d177a5cSEmil Constantinescu   PetscFunctionReturn(0);
933d177a5cSEmil Constantinescu }
943d177a5cSEmil Constantinescu 
953d177a5cSEmil Constantinescu static PetscErrorCode DMSubDomainRestrictHook_TSGLLE(DM dm,VecScatter gscat, VecScatter lscat,DM subdm,void *ctx)
963d177a5cSEmil Constantinescu {
973d177a5cSEmil Constantinescu   TS             ts = (TS)ctx;
983d177a5cSEmil Constantinescu   Vec            Ydot,Ydot_s;
993d177a5cSEmil Constantinescu 
1003d177a5cSEmil Constantinescu   PetscFunctionBegin;
101*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGLLEGetVecs(ts,dm,NULL,&Ydot));
102*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGLLEGetVecs(ts,subdm,NULL,&Ydot_s));
1033d177a5cSEmil Constantinescu 
104*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterBegin(gscat,Ydot,Ydot_s,INSERT_VALUES,SCATTER_FORWARD));
105*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterEnd(gscat,Ydot,Ydot_s,INSERT_VALUES,SCATTER_FORWARD));
1063d177a5cSEmil Constantinescu 
107*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGLLERestoreVecs(ts,dm,NULL,&Ydot));
108*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGLLERestoreVecs(ts,subdm,NULL,&Ydot_s));
1093d177a5cSEmil Constantinescu   PetscFunctionReturn(0);
1103d177a5cSEmil Constantinescu }
1113d177a5cSEmil Constantinescu 
1123d177a5cSEmil Constantinescu static PetscErrorCode TSGLLESchemeCreate(PetscInt p,PetscInt q,PetscInt r,PetscInt s,const PetscScalar *c,
1133d177a5cSEmil Constantinescu                                        const PetscScalar *a,const PetscScalar *b,const PetscScalar *u,const PetscScalar *v,TSGLLEScheme *inscheme)
1143d177a5cSEmil Constantinescu {
1153d177a5cSEmil Constantinescu   TSGLLEScheme     scheme;
1163d177a5cSEmil Constantinescu   PetscInt       j;
1173d177a5cSEmil Constantinescu 
1183d177a5cSEmil Constantinescu   PetscFunctionBegin;
1192c71b3e2SJacob Faibussowitsch   PetscCheckFalse(p < 1,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Scheme order must be positive");
1202c71b3e2SJacob Faibussowitsch   PetscCheckFalse(r < 1,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"At least one item must be carried between steps");
1212c71b3e2SJacob Faibussowitsch   PetscCheckFalse(s < 1,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"At least one stage is required");
122064a246eSJacob Faibussowitsch   PetscValidPointer(inscheme,10);
123c793f718SLisandro Dalcin   *inscheme = NULL;
124*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscNew(&scheme));
1253d177a5cSEmil Constantinescu   scheme->p = p;
1263d177a5cSEmil Constantinescu   scheme->q = q;
1273d177a5cSEmil Constantinescu   scheme->r = r;
1283d177a5cSEmil Constantinescu   scheme->s = s;
1293d177a5cSEmil Constantinescu 
130*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc5(s,&scheme->c,s*s,&scheme->a,r*s,&scheme->b,r*s,&scheme->u,r*r,&scheme->v));
131*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscArraycpy(scheme->c,c,s));
1323d177a5cSEmil Constantinescu   for (j=0; j<s*s; j++) scheme->a[j] = (PetscAbsScalar(a[j]) < 1e-12) ? 0 : a[j];
1333d177a5cSEmil Constantinescu   for (j=0; j<r*s; j++) scheme->b[j] = (PetscAbsScalar(b[j]) < 1e-12) ? 0 : b[j];
1343d177a5cSEmil Constantinescu   for (j=0; j<s*r; j++) scheme->u[j] = (PetscAbsScalar(u[j]) < 1e-12) ? 0 : u[j];
1353d177a5cSEmil Constantinescu   for (j=0; j<r*r; j++) scheme->v[j] = (PetscAbsScalar(v[j]) < 1e-12) ? 0 : v[j];
1363d177a5cSEmil Constantinescu 
137*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc6(r,&scheme->alpha,r,&scheme->beta,r,&scheme->gamma,3*s,&scheme->phi,3*r,&scheme->psi,r,&scheme->stage_error));
1383d177a5cSEmil Constantinescu   {
1393d177a5cSEmil Constantinescu     PetscInt     i,j,k,ss=s+2;
1403d177a5cSEmil Constantinescu     PetscBLASInt m,n,one=1,*ipiv,lwork=4*((s+3)*3+3),info,ldb;
1413d177a5cSEmil Constantinescu     PetscReal    rcond,*sing,*workreal;
1423d177a5cSEmil Constantinescu     PetscScalar  *ImV,*H,*bmat,*workscalar,*c=scheme->c,*a=scheme->a,*b=scheme->b,*u=scheme->u,*v=scheme->v;
1433d177a5cSEmil Constantinescu     PetscBLASInt rank;
144*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc7(PetscSqr(r),&ImV,3*s,&H,3*ss,&bmat,lwork,&workscalar,5*(3+r),&workreal,r+s,&sing,r+s,&ipiv));
1453d177a5cSEmil Constantinescu 
1463d177a5cSEmil Constantinescu     /* column-major input */
1473d177a5cSEmil Constantinescu     for (i=0; i<r-1; i++) {
1483d177a5cSEmil Constantinescu       for (j=0; j<r-1; j++) ImV[i+j*r] = 1.0*(i==j) - v[(i+1)*r+j+1];
1493d177a5cSEmil Constantinescu     }
1503d177a5cSEmil Constantinescu     /* Build right hand side for alpha (tp - glm.B(2:end,:)*(glm.c.^(p)./factorial(p))) */
1513d177a5cSEmil Constantinescu     for (i=1; i<r; i++) {
1523d177a5cSEmil Constantinescu       scheme->alpha[i] = 1./Factorial(p+1-i);
1533d177a5cSEmil Constantinescu       for (j=0; j<s; j++) scheme->alpha[i] -= b[i*s+j]*CPowF(c[j],p);
1543d177a5cSEmil Constantinescu     }
155*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscBLASIntCast(r-1,&m));
156*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscBLASIntCast(r,&n));
1573d177a5cSEmil Constantinescu     PetscStackCallBLAS("LAPACKgesv",LAPACKgesv_(&m,&one,ImV,&n,ipiv,scheme->alpha+1,&n,&info));
1582c71b3e2SJacob Faibussowitsch     PetscCheckFalse(info < 0,PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad argument to GESV");
1592c71b3e2SJacob Faibussowitsch     PetscCheckFalse(info > 0,PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Bad LU factorization");
1603d177a5cSEmil Constantinescu 
1613d177a5cSEmil Constantinescu     /* Build right hand side for beta (tp1 - glm.B(2:end,:)*(glm.c.^(p+1)./factorial(p+1)) - e.alpha) */
1623d177a5cSEmil Constantinescu     for (i=1; i<r; i++) {
1633d177a5cSEmil Constantinescu       scheme->beta[i] = 1./Factorial(p+2-i) - scheme->alpha[i];
1643d177a5cSEmil Constantinescu       for (j=0; j<s; j++) scheme->beta[i] -= b[i*s+j]*CPowF(c[j],p+1);
1653d177a5cSEmil Constantinescu     }
1663d177a5cSEmil Constantinescu     PetscStackCallBLAS("LAPACKgetrs",LAPACKgetrs_("No transpose",&m,&one,ImV,&n,ipiv,scheme->beta+1,&n,&info));
1672c71b3e2SJacob Faibussowitsch     PetscCheckFalse(info < 0,PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad argument to GETRS");
1682c71b3e2SJacob Faibussowitsch     PetscCheckFalse(info > 0,PETSC_COMM_SELF,PETSC_ERR_LIB,"Should not happen");
1693d177a5cSEmil Constantinescu 
1703d177a5cSEmil Constantinescu     /* Build stage_error vector
1713d177a5cSEmil Constantinescu            xi = glm.c.^(p+1)/factorial(p+1) - glm.A*glm.c.^p/factorial(p) + glm.U(:,2:end)*e.alpha;
1723d177a5cSEmil Constantinescu     */
1733d177a5cSEmil Constantinescu     for (i=0; i<s; i++) {
1743d177a5cSEmil Constantinescu       scheme->stage_error[i] = CPowF(c[i],p+1);
1753d177a5cSEmil Constantinescu       for (j=0; j<s; j++) scheme->stage_error[i] -= a[i*s+j]*CPowF(c[j],p);
1763d177a5cSEmil Constantinescu       for (j=1; j<r; j++) scheme->stage_error[i] += u[i*r+j]*scheme->alpha[j];
1773d177a5cSEmil Constantinescu     }
1783d177a5cSEmil Constantinescu 
1793d177a5cSEmil Constantinescu     /* alpha[0] (epsilon in B,J,W 2007)
1803d177a5cSEmil Constantinescu            epsilon = 1/factorial(p+1) - B(1,:)*c.^p/factorial(p) + V(1,2:end)*e.alpha;
1813d177a5cSEmil Constantinescu     */
1823d177a5cSEmil Constantinescu     scheme->alpha[0] = 1./Factorial(p+1);
1833d177a5cSEmil Constantinescu     for (j=0; j<s; j++) scheme->alpha[0] -= b[0*s+j]*CPowF(c[j],p);
1843d177a5cSEmil Constantinescu     for (j=1; j<r; j++) scheme->alpha[0] += v[0*r+j]*scheme->alpha[j];
1853d177a5cSEmil Constantinescu 
1863d177a5cSEmil Constantinescu     /* right hand side for gamma (glm.B(2:end,:)*e.xi - e.epsilon*eye(s-1,1)) */
1873d177a5cSEmil Constantinescu     for (i=1; i<r; i++) {
1883d177a5cSEmil Constantinescu       scheme->gamma[i] = (i==1 ? -1. : 0)*scheme->alpha[0];
1893d177a5cSEmil Constantinescu       for (j=0; j<s; j++) scheme->gamma[i] += b[i*s+j]*scheme->stage_error[j];
1903d177a5cSEmil Constantinescu     }
1913d177a5cSEmil Constantinescu     PetscStackCallBLAS("LAPACKgetrs",LAPACKgetrs_("No transpose",&m,&one,ImV,&n,ipiv,scheme->gamma+1,&n,&info));
1922c71b3e2SJacob Faibussowitsch     PetscCheckFalse(info < 0,PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad argument to GETRS");
1932c71b3e2SJacob Faibussowitsch     PetscCheckFalse(info > 0,PETSC_COMM_SELF,PETSC_ERR_LIB,"Should not happen");
1943d177a5cSEmil Constantinescu 
1953d177a5cSEmil Constantinescu     /* beta[0] (rho in B,J,W 2007)
1963d177a5cSEmil Constantinescu         e.rho = 1/factorial(p+2) - glm.B(1,:)*glm.c.^(p+1)/factorial(p+1) ...
1973d177a5cSEmil Constantinescu             + glm.V(1,2:end)*e.beta;% - e.epsilon;
1983d177a5cSEmil Constantinescu     % Note: The paper (B,J,W 2007) includes the last term in their definition
1993d177a5cSEmil Constantinescu     * */
2003d177a5cSEmil Constantinescu     scheme->beta[0] = 1./Factorial(p+2);
2013d177a5cSEmil Constantinescu     for (j=0; j<s; j++) scheme->beta[0] -= b[0*s+j]*CPowF(c[j],p+1);
2023d177a5cSEmil Constantinescu     for (j=1; j<r; j++) scheme->beta[0] += v[0*r+j]*scheme->beta[j];
2033d177a5cSEmil Constantinescu 
2043d177a5cSEmil Constantinescu     /* gamma[0] (sigma in B,J,W 2007)
2053d177a5cSEmil Constantinescu     *   e.sigma = glm.B(1,:)*e.xi + glm.V(1,2:end)*e.gamma;
2063d177a5cSEmil Constantinescu     * */
2073d177a5cSEmil Constantinescu     scheme->gamma[0] = 0.0;
2083d177a5cSEmil Constantinescu     for (j=0; j<s; j++) scheme->gamma[0] += b[0*s+j]*scheme->stage_error[j];
2093d177a5cSEmil Constantinescu     for (j=1; j<r; j++) scheme->gamma[0] += v[0*s+j]*scheme->gamma[j];
2103d177a5cSEmil Constantinescu 
2113d177a5cSEmil Constantinescu     /* Assemble H
2123d177a5cSEmil Constantinescu     *    % Determine the error estimators phi
2133d177a5cSEmil Constantinescu        H = [[cpow(glm.c,p) + C*e.alpha] [cpow(glm.c,p+1) + C*e.beta] ...
2143d177a5cSEmil Constantinescu                [e.xi - C*(e.gamma + 0*e.epsilon*eye(s-1,1))]]';
2153d177a5cSEmil Constantinescu     % Paper has formula above without the 0, but that term must be left
2163d177a5cSEmil Constantinescu     % out to satisfy the conditions they propose and to make the
2173d177a5cSEmil Constantinescu     % example schemes work
2183d177a5cSEmil Constantinescu     e.H = H;
2193d177a5cSEmil Constantinescu     e.phi = (H \ [1 0 0;1 1 0;0 0 -1])';
2203d177a5cSEmil Constantinescu     e.psi = -e.phi*C;
2213d177a5cSEmil Constantinescu     * */
2223d177a5cSEmil Constantinescu     for (j=0; j<s; j++) {
2233d177a5cSEmil Constantinescu       H[0+j*3] = CPowF(c[j],p);
2243d177a5cSEmil Constantinescu       H[1+j*3] = CPowF(c[j],p+1);
2253d177a5cSEmil Constantinescu       H[2+j*3] = scheme->stage_error[j];
2263d177a5cSEmil Constantinescu       for (k=1; k<r; k++) {
2273d177a5cSEmil Constantinescu         H[0+j*3] += CPowF(c[j],k-1)*scheme->alpha[k];
2283d177a5cSEmil Constantinescu         H[1+j*3] += CPowF(c[j],k-1)*scheme->beta[k];
2293d177a5cSEmil Constantinescu         H[2+j*3] -= CPowF(c[j],k-1)*scheme->gamma[k];
2303d177a5cSEmil Constantinescu       }
2313d177a5cSEmil Constantinescu     }
2323d177a5cSEmil Constantinescu     bmat[0+0*ss] = 1.;  bmat[0+1*ss] = 0.;  bmat[0+2*ss] = 0.;
2333d177a5cSEmil Constantinescu     bmat[1+0*ss] = 1.;  bmat[1+1*ss] = 1.;  bmat[1+2*ss] = 0.;
2343d177a5cSEmil Constantinescu     bmat[2+0*ss] = 0.;  bmat[2+1*ss] = 0.;  bmat[2+2*ss] = -1.;
2353d177a5cSEmil Constantinescu     m     = 3;
236*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscBLASIntCast(s,&n));
237*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscBLASIntCast(ss,&ldb));
2383d177a5cSEmil Constantinescu     rcond = 1e-12;
2393d177a5cSEmil Constantinescu #if defined(PETSC_USE_COMPLEX)
2403d177a5cSEmil Constantinescu     /* ZGELSS( M, N, NRHS, A, LDA, B, LDB, S, RCOND, RANK, WORK, LWORK, RWORK, INFO) */
2413d177a5cSEmil Constantinescu     PetscStackCallBLAS("LAPACKgelss",LAPACKgelss_(&m,&n,&m,H,&m,bmat,&ldb,sing,&rcond,&rank,workscalar,&lwork,workreal,&info));
2423d177a5cSEmil Constantinescu #else
2433d177a5cSEmil Constantinescu     /* DGELSS( M, N, NRHS, A, LDA, B, LDB, S, RCOND, RANK, WORK, LWORK, INFO) */
2443d177a5cSEmil Constantinescu     PetscStackCallBLAS("LAPACKgelss",LAPACKgelss_(&m,&n,&m,H,&m,bmat,&ldb,sing,&rcond,&rank,workscalar,&lwork,&info));
2453d177a5cSEmil Constantinescu #endif
2462c71b3e2SJacob Faibussowitsch     PetscCheckFalse(info < 0,PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad argument to GELSS");
2472c71b3e2SJacob Faibussowitsch     PetscCheckFalse(info > 0,PETSC_COMM_SELF,PETSC_ERR_LIB,"SVD failed to converge");
2483d177a5cSEmil Constantinescu 
2493d177a5cSEmil Constantinescu     for (j=0; j<3; j++) {
2503d177a5cSEmil Constantinescu       for (k=0; k<s; k++) scheme->phi[k+j*s] = bmat[k+j*ss];
2513d177a5cSEmil Constantinescu     }
2523d177a5cSEmil Constantinescu 
2533d177a5cSEmil Constantinescu     /* the other part of the error estimator, psi in B,J,W 2007 */
2543d177a5cSEmil Constantinescu     scheme->psi[0*r+0] = 0.;
2553d177a5cSEmil Constantinescu     scheme->psi[1*r+0] = 0.;
2563d177a5cSEmil Constantinescu     scheme->psi[2*r+0] = 0.;
2573d177a5cSEmil Constantinescu     for (j=1; j<r; j++) {
2583d177a5cSEmil Constantinescu       scheme->psi[0*r+j] = 0.;
2593d177a5cSEmil Constantinescu       scheme->psi[1*r+j] = 0.;
2603d177a5cSEmil Constantinescu       scheme->psi[2*r+j] = 0.;
2613d177a5cSEmil Constantinescu       for (k=0; k<s; k++) {
2623d177a5cSEmil Constantinescu         scheme->psi[0*r+j] -= CPowF(c[k],j-1)*scheme->phi[0*s+k];
2633d177a5cSEmil Constantinescu         scheme->psi[1*r+j] -= CPowF(c[k],j-1)*scheme->phi[1*s+k];
2643d177a5cSEmil Constantinescu         scheme->psi[2*r+j] -= CPowF(c[k],j-1)*scheme->phi[2*s+k];
2653d177a5cSEmil Constantinescu       }
2663d177a5cSEmil Constantinescu     }
267*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree7(ImV,H,bmat,workscalar,workreal,sing,ipiv));
2683d177a5cSEmil Constantinescu   }
2693d177a5cSEmil Constantinescu   /* Check which properties are satisfied */
2703d177a5cSEmil Constantinescu   scheme->stiffly_accurate = PETSC_TRUE;
2713d177a5cSEmil Constantinescu   if (scheme->c[s-1] != 1.) scheme->stiffly_accurate = PETSC_FALSE;
2723d177a5cSEmil Constantinescu   for (j=0; j<s; j++) if (a[(s-1)*s+j] != b[j]) scheme->stiffly_accurate = PETSC_FALSE;
2733d177a5cSEmil Constantinescu   for (j=0; j<r; j++) if (u[(s-1)*r+j] != v[j]) scheme->stiffly_accurate = PETSC_FALSE;
2743d177a5cSEmil Constantinescu   scheme->fsal = scheme->stiffly_accurate; /* FSAL is stronger */
2753d177a5cSEmil Constantinescu   for (j=0; j<s-1; j++) if (r>1 && b[1*s+j] != 0.) scheme->fsal = PETSC_FALSE;
2763d177a5cSEmil Constantinescu   if (b[1*s+r-1] != 1.) scheme->fsal = PETSC_FALSE;
2773d177a5cSEmil Constantinescu   for (j=0; j<r; j++) if (r>1 && v[1*r+j] != 0.) scheme->fsal = PETSC_FALSE;
2783d177a5cSEmil Constantinescu 
2793d177a5cSEmil Constantinescu   *inscheme = scheme;
2803d177a5cSEmil Constantinescu   PetscFunctionReturn(0);
2813d177a5cSEmil Constantinescu }
2823d177a5cSEmil Constantinescu 
2833d177a5cSEmil Constantinescu static PetscErrorCode TSGLLESchemeDestroy(TSGLLEScheme sc)
2843d177a5cSEmil Constantinescu {
2853d177a5cSEmil Constantinescu   PetscFunctionBegin;
286*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree5(sc->c,sc->a,sc->b,sc->u,sc->v));
287*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree6(sc->alpha,sc->beta,sc->gamma,sc->phi,sc->psi,sc->stage_error));
288*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(sc));
2893d177a5cSEmil Constantinescu   PetscFunctionReturn(0);
2903d177a5cSEmil Constantinescu }
2913d177a5cSEmil Constantinescu 
2923d177a5cSEmil Constantinescu static PetscErrorCode TSGLLEDestroy_Default(TS_GLLE *gl)
2933d177a5cSEmil Constantinescu {
2943d177a5cSEmil Constantinescu   PetscInt       i;
2953d177a5cSEmil Constantinescu 
2963d177a5cSEmil Constantinescu   PetscFunctionBegin;
2973d177a5cSEmil Constantinescu   for (i=0; i<gl->nschemes; i++) {
298*5f80ce2aSJacob Faibussowitsch     if (gl->schemes[i]) CHKERRQ(TSGLLESchemeDestroy(gl->schemes[i]));
2993d177a5cSEmil Constantinescu   }
300*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(gl->schemes));
3013d177a5cSEmil Constantinescu   gl->nschemes = 0;
302*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMemzero(gl->type_name,sizeof(gl->type_name)));
3033d177a5cSEmil Constantinescu   PetscFunctionReturn(0);
3043d177a5cSEmil Constantinescu }
3053d177a5cSEmil Constantinescu 
3063d177a5cSEmil Constantinescu static PetscErrorCode TSGLLEViewTable_Private(PetscViewer viewer,PetscInt m,PetscInt n,const PetscScalar a[],const char name[])
3073d177a5cSEmil Constantinescu {
3083d177a5cSEmil Constantinescu   PetscBool      iascii;
3093d177a5cSEmil Constantinescu   PetscInt       i,j;
3103d177a5cSEmil Constantinescu 
3113d177a5cSEmil Constantinescu   PetscFunctionBegin;
312*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii));
3133d177a5cSEmil Constantinescu   if (iascii) {
314*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIPrintf(viewer,"%30s = [",name));
3153d177a5cSEmil Constantinescu     for (i=0; i<m; i++) {
316*5f80ce2aSJacob Faibussowitsch       if (i) CHKERRQ(PetscViewerASCIIPrintf(viewer,"%30s   [",""));
317*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE));
3183d177a5cSEmil Constantinescu       for (j=0; j<n; j++) {
319*5f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscViewerASCIIPrintf(viewer," %12.8g",PetscRealPart(a[i*n+j])));
3203d177a5cSEmil Constantinescu       }
321*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscViewerASCIIPrintf(viewer,"]\n"));
322*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE));
3233d177a5cSEmil Constantinescu     }
3243d177a5cSEmil Constantinescu   }
3253d177a5cSEmil Constantinescu   PetscFunctionReturn(0);
3263d177a5cSEmil Constantinescu }
3273d177a5cSEmil Constantinescu 
3283d177a5cSEmil Constantinescu static PetscErrorCode TSGLLESchemeView(TSGLLEScheme sc,PetscBool view_details,PetscViewer viewer)
3293d177a5cSEmil Constantinescu {
3303d177a5cSEmil Constantinescu   PetscErrorCode ierr;
3313d177a5cSEmil Constantinescu   PetscBool      iascii;
3323d177a5cSEmil Constantinescu 
3333d177a5cSEmil Constantinescu   PetscFunctionBegin;
334*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii));
3353d177a5cSEmil Constantinescu   if (iascii) {
336*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIPrintf(viewer,"GL scheme p,q,r,s = %d,%d,%d,%d\n",sc->p,sc->q,sc->r,sc->s));
337*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIPushTab(viewer));
338*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIPrintf(viewer,"Stiffly accurate: %s,  FSAL: %s\n",sc->stiffly_accurate ? "yes" : "no",sc->fsal ? "yes" : "no"));
3393d177a5cSEmil Constantinescu     ierr = PetscViewerASCIIPrintf(viewer,"Leading error constants: %10.3e  %10.3e  %10.3e\n",
3403d177a5cSEmil Constantinescu                                   PetscRealPart(sc->alpha[0]),PetscRealPart(sc->beta[0]),PetscRealPart(sc->gamma[0]));CHKERRQ(ierr);
341*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TSGLLEViewTable_Private(viewer,1,sc->s,sc->c,"Abscissas c"));
3423d177a5cSEmil Constantinescu     if (view_details) {
343*5f80ce2aSJacob Faibussowitsch       CHKERRQ(TSGLLEViewTable_Private(viewer,sc->s,sc->s,sc->a,"A"));
344*5f80ce2aSJacob Faibussowitsch       CHKERRQ(TSGLLEViewTable_Private(viewer,sc->r,sc->s,sc->b,"B"));
345*5f80ce2aSJacob Faibussowitsch       CHKERRQ(TSGLLEViewTable_Private(viewer,sc->s,sc->r,sc->u,"U"));
346*5f80ce2aSJacob Faibussowitsch       CHKERRQ(TSGLLEViewTable_Private(viewer,sc->r,sc->r,sc->v,"V"));
3473d177a5cSEmil Constantinescu 
348*5f80ce2aSJacob Faibussowitsch       CHKERRQ(TSGLLEViewTable_Private(viewer,3,sc->s,sc->phi,"Error estimate phi"));
349*5f80ce2aSJacob Faibussowitsch       CHKERRQ(TSGLLEViewTable_Private(viewer,3,sc->r,sc->psi,"Error estimate psi"));
350*5f80ce2aSJacob Faibussowitsch       CHKERRQ(TSGLLEViewTable_Private(viewer,1,sc->r,sc->alpha,"Modify alpha"));
351*5f80ce2aSJacob Faibussowitsch       CHKERRQ(TSGLLEViewTable_Private(viewer,1,sc->r,sc->beta,"Modify beta"));
352*5f80ce2aSJacob Faibussowitsch       CHKERRQ(TSGLLEViewTable_Private(viewer,1,sc->r,sc->gamma,"Modify gamma"));
353*5f80ce2aSJacob Faibussowitsch       CHKERRQ(TSGLLEViewTable_Private(viewer,1,sc->s,sc->stage_error,"Stage error xi"));
3543d177a5cSEmil Constantinescu     }
355*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIPopTab(viewer));
35698921bdaSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported",((PetscObject)viewer)->type_name);
3573d177a5cSEmil Constantinescu   PetscFunctionReturn(0);
3583d177a5cSEmil Constantinescu }
3593d177a5cSEmil Constantinescu 
3603d177a5cSEmil Constantinescu static PetscErrorCode TSGLLEEstimateHigherMoments_Default(TSGLLEScheme sc,PetscReal h,Vec Ydot[],Vec Xold[],Vec hm[])
3613d177a5cSEmil Constantinescu {
3623d177a5cSEmil Constantinescu   PetscInt       i;
3633d177a5cSEmil Constantinescu 
3643d177a5cSEmil Constantinescu   PetscFunctionBegin;
3652c71b3e2SJacob Faibussowitsch   PetscCheckFalse(sc->r > 64 || sc->s > 64,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Ridiculous number of stages or items passed between stages");
3663d177a5cSEmil Constantinescu   /* build error vectors*/
3673d177a5cSEmil Constantinescu   for (i=0; i<3; i++) {
3683d177a5cSEmil Constantinescu     PetscScalar phih[64];
3693d177a5cSEmil Constantinescu     PetscInt    j;
3703d177a5cSEmil Constantinescu     for (j=0; j<sc->s; j++) phih[j] = sc->phi[i*sc->s+j]*h;
371*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecZeroEntries(hm[i]));
372*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecMAXPY(hm[i],sc->s,phih,Ydot));
373*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecMAXPY(hm[i],sc->r,&sc->psi[i*sc->r],Xold));
3743d177a5cSEmil Constantinescu   }
3753d177a5cSEmil Constantinescu   PetscFunctionReturn(0);
3763d177a5cSEmil Constantinescu }
3773d177a5cSEmil Constantinescu 
3783d177a5cSEmil Constantinescu static PetscErrorCode TSGLLECompleteStep_Rescale(TSGLLEScheme sc,PetscReal h,TSGLLEScheme next_sc,PetscReal next_h,Vec Ydot[],Vec Xold[],Vec X[])
3793d177a5cSEmil Constantinescu {
3803d177a5cSEmil Constantinescu   PetscScalar    brow[32],vrow[32];
3813d177a5cSEmil Constantinescu   PetscInt       i,j,r,s;
3823d177a5cSEmil Constantinescu 
3833d177a5cSEmil Constantinescu   PetscFunctionBegin;
3843d177a5cSEmil Constantinescu   /* Build the new solution from (X,Ydot) */
3853d177a5cSEmil Constantinescu   r = sc->r;
3863d177a5cSEmil Constantinescu   s = sc->s;
3873d177a5cSEmil Constantinescu   for (i=0; i<r; i++) {
388*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecZeroEntries(X[i]));
3893d177a5cSEmil Constantinescu     for (j=0; j<s; j++) brow[j] = h*sc->b[i*s+j];
390*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecMAXPY(X[i],s,brow,Ydot));
3913d177a5cSEmil Constantinescu     for (j=0; j<r; j++) vrow[j] = sc->v[i*r+j];
392*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecMAXPY(X[i],r,vrow,Xold));
3933d177a5cSEmil Constantinescu   }
3943d177a5cSEmil Constantinescu   PetscFunctionReturn(0);
3953d177a5cSEmil Constantinescu }
3963d177a5cSEmil Constantinescu 
3973d177a5cSEmil Constantinescu static PetscErrorCode TSGLLECompleteStep_RescaleAndModify(TSGLLEScheme sc,PetscReal h,TSGLLEScheme next_sc,PetscReal next_h,Vec Ydot[],Vec Xold[],Vec X[])
3983d177a5cSEmil Constantinescu {
3993d177a5cSEmil Constantinescu   PetscScalar    brow[32],vrow[32];
4003d177a5cSEmil Constantinescu   PetscReal      ratio;
4013d177a5cSEmil Constantinescu   PetscInt       i,j,p,r,s;
4023d177a5cSEmil Constantinescu 
4033d177a5cSEmil Constantinescu   PetscFunctionBegin;
4043d177a5cSEmil Constantinescu   /* Build the new solution from (X,Ydot) */
4053d177a5cSEmil Constantinescu   p     = sc->p;
4063d177a5cSEmil Constantinescu   r     = sc->r;
4073d177a5cSEmil Constantinescu   s     = sc->s;
4083d177a5cSEmil Constantinescu   ratio = next_h/h;
4093d177a5cSEmil Constantinescu   for (i=0; i<r; i++) {
410*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecZeroEntries(X[i]));
4113d177a5cSEmil Constantinescu     for (j=0; j<s; j++) {
4123d177a5cSEmil Constantinescu       brow[j] = h*(PetscPowRealInt(ratio,i)*sc->b[i*s+j]
4133d177a5cSEmil Constantinescu                    + (PetscPowRealInt(ratio,i) - PetscPowRealInt(ratio,p+1))*(+ sc->alpha[i]*sc->phi[0*s+j])
4143d177a5cSEmil Constantinescu                    + (PetscPowRealInt(ratio,i) - PetscPowRealInt(ratio,p+2))*(+ sc->beta [i]*sc->phi[1*s+j]
4153d177a5cSEmil Constantinescu                                                                               + sc->gamma[i]*sc->phi[2*s+j]));
4163d177a5cSEmil Constantinescu     }
417*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecMAXPY(X[i],s,brow,Ydot));
4183d177a5cSEmil Constantinescu     for (j=0; j<r; j++) {
4193d177a5cSEmil Constantinescu       vrow[j] = (PetscPowRealInt(ratio,i)*sc->v[i*r+j]
4203d177a5cSEmil Constantinescu                  + (PetscPowRealInt(ratio,i) - PetscPowRealInt(ratio,p+1))*(+ sc->alpha[i]*sc->psi[0*r+j])
4213d177a5cSEmil Constantinescu                  + (PetscPowRealInt(ratio,i) - PetscPowRealInt(ratio,p+2))*(+ sc->beta [i]*sc->psi[1*r+j]
4223d177a5cSEmil Constantinescu                                                                             + sc->gamma[i]*sc->psi[2*r+j]));
4233d177a5cSEmil Constantinescu     }
424*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecMAXPY(X[i],r,vrow,Xold));
4253d177a5cSEmil Constantinescu   }
4263d177a5cSEmil Constantinescu   if (r < next_sc->r) {
4272c71b3e2SJacob Faibussowitsch     PetscCheckFalse(r+1 != next_sc->r,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot accommodate jump in r greater than 1");
428*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecZeroEntries(X[r]));
4293d177a5cSEmil Constantinescu     for (j=0; j<s; j++) brow[j] = h*PetscPowRealInt(ratio,p+1)*sc->phi[0*s+j];
430*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecMAXPY(X[r],s,brow,Ydot));
4313d177a5cSEmil Constantinescu     for (j=0; j<r; j++) vrow[j] = PetscPowRealInt(ratio,p+1)*sc->psi[0*r+j];
432*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecMAXPY(X[r],r,vrow,Xold));
4333d177a5cSEmil Constantinescu   }
4343d177a5cSEmil Constantinescu   PetscFunctionReturn(0);
4353d177a5cSEmil Constantinescu }
4363d177a5cSEmil Constantinescu 
4373d177a5cSEmil Constantinescu static PetscErrorCode TSGLLECreate_IRKS(TS ts)
4383d177a5cSEmil Constantinescu {
4393d177a5cSEmil Constantinescu   TS_GLLE        *gl = (TS_GLLE*)ts->data;
4403d177a5cSEmil Constantinescu 
4413d177a5cSEmil Constantinescu   PetscFunctionBegin;
4423d177a5cSEmil Constantinescu   gl->Destroy               = TSGLLEDestroy_Default;
4433d177a5cSEmil Constantinescu   gl->EstimateHigherMoments = TSGLLEEstimateHigherMoments_Default;
4443d177a5cSEmil Constantinescu   gl->CompleteStep          = TSGLLECompleteStep_RescaleAndModify;
445*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(10,&gl->schemes));
4463d177a5cSEmil Constantinescu   gl->nschemes = 0;
4473d177a5cSEmil Constantinescu 
4483d177a5cSEmil Constantinescu   {
4493d177a5cSEmil Constantinescu     /* p=1,q=1, r=s=2, A- and L-stable with error estimates of order 2 and 3
4503d177a5cSEmil Constantinescu     * Listed in Butcher & Podhaisky 2006. On error estimation in general linear methods for stiff ODE.
4513d177a5cSEmil Constantinescu     * irks(0.3,0,[.3,1],[1],1)
4523d177a5cSEmil Constantinescu     * Note: can be made to have classical order (not stage order) 2 by replacing 0.3 with 1-sqrt(1/2)
4533d177a5cSEmil Constantinescu     * but doing so would sacrifice the error estimator.
4543d177a5cSEmil Constantinescu     */
4553d177a5cSEmil Constantinescu     const PetscScalar c[2]    = {3./10., 1.};
4563d177a5cSEmil Constantinescu     const PetscScalar a[2][2] = {{3./10., 0}, {7./10., 3./10.}};
4573d177a5cSEmil Constantinescu     const PetscScalar b[2][2] = {{7./10., 3./10.}, {0,1}};
4583d177a5cSEmil Constantinescu     const PetscScalar u[2][2] = {{1,0},{1,0}};
4593d177a5cSEmil Constantinescu     const PetscScalar v[2][2] = {{1,0},{0,0}};
460*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TSGLLESchemeCreate(1,1,2,2,c,*a,*b,*u,*v,&gl->schemes[gl->nschemes++]));
4613d177a5cSEmil Constantinescu   }
4623d177a5cSEmil Constantinescu 
4633d177a5cSEmil Constantinescu   {
4643d177a5cSEmil Constantinescu     /* p=q=2, r=s=3: irks(4/9,0,[1:3]/3,[0.33852],1) */
4653d177a5cSEmil Constantinescu     /* http://www.math.auckland.ac.nz/~hpod/atlas/i2a.html */
4663d177a5cSEmil Constantinescu     const PetscScalar c[3] = {1./3., 2./3., 1}
4673d177a5cSEmil Constantinescu     ,a[3][3] = {{4./9.                ,0                      ,       0},
4683d177a5cSEmil Constantinescu                 {1.03750643704090e+00 ,                  4./9.,       0},
4693d177a5cSEmil Constantinescu                 {7.67024779410304e-01 ,  -3.81140216918943e-01,   4./9.}}
4703d177a5cSEmil Constantinescu     ,b[3][3] = {{0.767024779410304,  -0.381140216918943,   4./9.},
4713d177a5cSEmil Constantinescu                 {0.000000000000000,  0.000000000000000,   1.000000000000000},
4723d177a5cSEmil Constantinescu                 {-2.075048385225385,   0.621728385225383,   1.277197204924873}}
4733d177a5cSEmil Constantinescu     ,u[3][3] = {{1.0000000000000000,  -0.1111111111111109,  -0.0925925925925922},
4743d177a5cSEmil Constantinescu                 {1.0000000000000000,  -0.8152842148186744,  -0.4199095530877056},
4753d177a5cSEmil Constantinescu                 {1.0000000000000000,   0.1696709930641948,   0.0539741070314165}}
4763d177a5cSEmil Constantinescu     ,v[3][3] = {{1.0000000000000000,  0.1696709930641948,   0.0539741070314165},
4773d177a5cSEmil Constantinescu                 {0.000000000000000,   0.000000000000000,   0.000000000000000},
4783d177a5cSEmil Constantinescu                 {0.000000000000000,   0.176122795075129,   0.000000000000000}};
479*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TSGLLESchemeCreate(2,2,3,3,c,*a,*b,*u,*v,&gl->schemes[gl->nschemes++]));
4803d177a5cSEmil Constantinescu   }
4813d177a5cSEmil Constantinescu   {
4823d177a5cSEmil Constantinescu     /* p=q=3, r=s=4: irks(9/40,0,[1:4]/4,[0.3312 1.0050],[0.49541 1;1 0]) */
4833d177a5cSEmil Constantinescu     const PetscScalar c[4] = {0.25,0.5,0.75,1.0}
4843d177a5cSEmil Constantinescu     ,a[4][4] = {{9./40.               ,                      0,                      0,                      0},
4853d177a5cSEmil Constantinescu                 {2.11286958887701e-01 ,    9./40.             ,                      0,                      0},
4863d177a5cSEmil Constantinescu                 {9.46338294287584e-01 ,  -3.42942861246094e-01,   9./40.              ,                      0},
4873d177a5cSEmil Constantinescu                 {0.521490453970721    ,  -0.662474225622980,   0.490476425459734,   9./40.           }}
4883d177a5cSEmil Constantinescu     ,b[4][4] = {{0.521490453970721    ,  -0.662474225622980,   0.490476425459734,   9./40.           },
4893d177a5cSEmil Constantinescu                 {0.000000000000000    ,   0.000000000000000,   0.000000000000000,   1.000000000000000},
4903d177a5cSEmil Constantinescu                 {-0.084677029310348   ,   1.390757514776085,  -1.568157386206001,   2.023192696767826},
4913d177a5cSEmil Constantinescu                 {0.465383797936408    ,   1.478273530625148,  -1.930836081010182,   1.644872111193354}}
4923d177a5cSEmil Constantinescu     ,u[4][4] = {{1.00000000000000000  ,   0.02500000000001035,  -0.02499999999999053,  -0.00442708333332865},
4933d177a5cSEmil Constantinescu                 {1.00000000000000000  ,   0.06371304111232945,  -0.04032173972189845,  -0.01389438413189452},
4943d177a5cSEmil Constantinescu                 {1.00000000000000000  ,  -0.07839543304147778,   0.04738685705116663,   0.02032603595928376},
4953d177a5cSEmil Constantinescu                 {1.00000000000000000  ,   0.42550734619251651,   0.10800718022400080,  -0.01726712647760034}}
4963d177a5cSEmil Constantinescu     ,v[4][4] = {{1.00000000000000000  ,   0.42550734619251651,   0.10800718022400080,  -0.01726712647760034},
4973d177a5cSEmil Constantinescu                 {0.000000000000000    ,   0.000000000000000,   0.000000000000000,   0.000000000000000},
4983d177a5cSEmil Constantinescu                 {0.000000000000000    ,  -1.761115796027561,  -0.521284157173780,   0.258249384305463},
4993d177a5cSEmil Constantinescu                 {0.000000000000000    ,  -1.657693358744728,  -1.052227765232394,   0.521284157173780}};
500*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TSGLLESchemeCreate(3,3,4,4,c,*a,*b,*u,*v,&gl->schemes[gl->nschemes++]));
5013d177a5cSEmil Constantinescu   }
5023d177a5cSEmil Constantinescu   {
5033d177a5cSEmil Constantinescu     /* p=q=4, r=s=5:
5043d177a5cSEmil Constantinescu           irks(3/11,0,[1:5]/5, [0.1715   -0.1238    0.6617],...
5053d177a5cSEmil Constantinescu           [ -0.0812    0.4079    1.0000
5063d177a5cSEmil Constantinescu              1.0000         0         0
5073d177a5cSEmil Constantinescu              0.8270    1.0000         0])
5083d177a5cSEmil Constantinescu     */
5093d177a5cSEmil Constantinescu     const PetscScalar c[5] = {0.2,0.4,0.6,0.8,1.0}
5103d177a5cSEmil Constantinescu     ,a[5][5] = {{2.72727272727352e-01 ,   0.00000000000000e+00,  0.00000000000000e+00 ,  0.00000000000000e+00  ,  0.00000000000000e+00},
5113d177a5cSEmil Constantinescu                 {-1.03980153733431e-01,   2.72727272727405e-01,   0.00000000000000e+00,  0.00000000000000e+00  ,  0.00000000000000e+00},
5123d177a5cSEmil Constantinescu                 {-1.58615400341492e+00,   7.44168951881122e-01,   2.72727272727309e-01,  0.00000000000000e+00  ,  0.00000000000000e+00},
5133d177a5cSEmil Constantinescu                 {-8.73658042865628e-01,   5.37884671894595e-01,  -1.63298538799523e-01,   2.72727272726996e-01 ,  0.00000000000000e+00},
5143d177a5cSEmil Constantinescu                 {2.95489397443992e-01 , -1.18481693910097e+00 , -6.68029812659953e-01 ,  1.00716687860943e+00  , 2.72727272727288e-01}}
5153d177a5cSEmil Constantinescu     ,b[5][5] = {{2.95489397443992e-01 , -1.18481693910097e+00 , -6.68029812659953e-01 ,  1.00716687860943e+00  , 2.72727272727288e-01},
5163d177a5cSEmil Constantinescu                 {0.00000000000000e+00 ,  1.11022302462516e-16 , -2.22044604925031e-16 ,  0.00000000000000e+00  , 1.00000000000000e+00},
5173d177a5cSEmil Constantinescu                 {-4.05882503986005e+00,  -4.00924006567769e+00,  -1.38930610972481e+00,   4.45223930308488e+00 ,  6.32331093108427e-01},
5183d177a5cSEmil Constantinescu                 {8.35690179937017e+00 , -2.26640927349732e+00 ,  6.86647884973826e+00 , -5.22595158025740e+00  , 4.50893068837431e+00},
5193d177a5cSEmil Constantinescu                 {1.27656267027479e+01 ,  2.80882153840821e+00 ,  8.91173096522890e+00 , -1.07936444078906e+01  , 4.82534148988854e+00}}
5203d177a5cSEmil Constantinescu     ,u[5][5] = {{1.00000000000000e+00 , -7.27272727273551e-02 , -3.45454545454419e-02 , -4.12121212119565e-03  ,-2.96969696964014e-04},
5213d177a5cSEmil Constantinescu                 {1.00000000000000e+00 ,  2.31252881006154e-01 , -8.29487834416481e-03 , -9.07191207681020e-03  ,-1.70378403743473e-03},
5223d177a5cSEmil Constantinescu                 {1.00000000000000e+00 ,  1.16925777880663e+00 ,  3.59268562942635e-02 , -4.09013451730615e-02  ,-1.02411119670164e-02},
5233d177a5cSEmil Constantinescu                 {1.00000000000000e+00 ,  1.02634463704356e+00 ,  1.59375044913405e-01 ,  1.89673015035370e-03  ,-4.89987231897569e-03},
5243d177a5cSEmil Constantinescu                 {1.00000000000000e+00 ,  1.27746320298021e+00 ,  2.37186008132728e-01 , -8.28694373940065e-02  ,-5.34396510196430e-02}}
5253d177a5cSEmil Constantinescu     ,v[5][5] = {{1.00000000000000e+00 ,  1.27746320298021e+00 ,  2.37186008132728e-01 , -8.28694373940065e-02  ,-5.34396510196430e-02},
5263d177a5cSEmil Constantinescu                 {0.00000000000000e+00 , -1.77635683940025e-15 , -1.99840144432528e-15 , -9.99200722162641e-16  ,-3.33066907387547e-16},
5273d177a5cSEmil Constantinescu                 {0.00000000000000e+00 ,  4.37280081906924e+00 ,  5.49221645016377e-02 , -8.88913177394943e-02  , 1.12879077989154e-01},
5283d177a5cSEmil Constantinescu                 {0.00000000000000e+00 , -1.22399504837280e+01 , -5.21287338448645e+00 , -8.03952325565291e-01  , 4.60298678047147e-01},
5293d177a5cSEmil Constantinescu                 {0.00000000000000e+00 , -1.85178762883829e+01 , -5.21411849862624e+00 , -1.04283436528809e+00  , 7.49030161063651e-01}};
530*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TSGLLESchemeCreate(4,4,5,5,c,*a,*b,*u,*v,&gl->schemes[gl->nschemes++]));
5313d177a5cSEmil Constantinescu   }
5323d177a5cSEmil Constantinescu   {
5333d177a5cSEmil Constantinescu     /* p=q=5, r=s=6;
5343d177a5cSEmil Constantinescu        irks(1/3,0,[1:6]/6,...
5353d177a5cSEmil Constantinescu           [-0.0489    0.4228   -0.8814    0.9021],...
5363d177a5cSEmil Constantinescu           [-0.3474   -0.6617    0.6294    0.2129
5373d177a5cSEmil Constantinescu             0.0044   -0.4256   -0.1427   -0.8936
5383d177a5cSEmil Constantinescu            -0.8267    0.4821    0.1371   -0.2557
5393d177a5cSEmil Constantinescu            -0.4426   -0.3855   -0.7514    0.3014])
5403d177a5cSEmil Constantinescu     */
5413d177a5cSEmil Constantinescu     const PetscScalar c[6] = {1./6, 2./6, 3./6, 4./6, 5./6, 1.}
5423d177a5cSEmil Constantinescu     ,a[6][6] = {{  3.33333333333940e-01,  0                   ,  0                   ,  0                   ,  0                   ,  0                   },
5433d177a5cSEmil Constantinescu                 { -8.64423857333350e-02,  3.33333333332888e-01,  0                   ,  0                   ,  0                   ,  0                   },
5443d177a5cSEmil Constantinescu                 { -2.16850174258252e+00, -2.23619072028839e+00,  3.33333333335204e-01,  0                   ,  0                   ,  0                   },
5453d177a5cSEmil Constantinescu                 { -4.73160970138997e+00, -3.89265344629268e+00, -2.76318716520933e-01,  3.33333333335759e-01,  0                   ,  0                   },
5463d177a5cSEmil Constantinescu                 { -6.75187540297338e+00, -7.90756533769377e+00,  7.90245051802259e-01, -4.48352364517632e-01,  3.33333333328483e-01,  0                   },
5473d177a5cSEmil Constantinescu                 { -4.26488287921548e+00, -1.19320395589302e+01,  3.38924509887755e+00, -2.23969848002481e+00,  6.62807710124007e-01,  3.33333333335440e-01}}
5483d177a5cSEmil Constantinescu     ,b[6][6] = {{ -4.26488287921548e+00, -1.19320395589302e+01,  3.38924509887755e+00, -2.23969848002481e+00,  6.62807710124007e-01,  3.33333333335440e-01},
5493d177a5cSEmil Constantinescu                 { -8.88178419700125e-16,  4.44089209850063e-16, -1.54737334057131e-15, -8.88178419700125e-16,  0.00000000000000e+00,  1.00000000000001e+00},
5503d177a5cSEmil Constantinescu                 { -2.87780425770651e+01, -1.13520448264971e+01,  2.62002318943161e+01,  2.56943874812797e+01, -3.06702268304488e+01,  6.68067773510103e+00},
5513d177a5cSEmil Constantinescu                 {  5.47971245256474e+01,  6.80366875868284e+01, -6.50952588861999e+01, -8.28643975339097e+01,  8.17416943896414e+01, -1.17819043489036e+01},
5523d177a5cSEmil Constantinescu                 { -2.33332114788869e+02,  6.12942539462634e+01, -4.91850135865944e+01,  1.82716844135480e+02, -1.29788173979395e+02,  3.09968095651099e+01},
5533d177a5cSEmil Constantinescu                 { -1.72049132343751e+02,  8.60194713593999e+00,  7.98154219170200e-01,  1.50371386053218e+02, -1.18515423962066e+02,  2.50898277784663e+01}}
5543d177a5cSEmil Constantinescu     ,u[6][6] = {{  1.00000000000000e+00, -1.66666666666870e-01, -4.16666666664335e-02, -3.85802469124815e-03, -2.25051440302250e-04, -9.64506172339142e-06},
5553d177a5cSEmil Constantinescu                 {  1.00000000000000e+00,  8.64423857327162e-02, -4.11484912671353e-02, -1.11450903217645e-02, -1.47651050487126e-03, -1.34395070766826e-04},
5563d177a5cSEmil Constantinescu                 {  1.00000000000000e+00,  4.57135912953434e+00,  1.06514719719137e+00,  1.33517564218007e-01,  1.11365952968659e-02,  6.12382756769504e-04},
5573d177a5cSEmil Constantinescu                 {  1.00000000000000e+00,  9.23391519753404e+00,  2.22431212392095e+00,  2.91823807741891e-01,  2.52058456411084e-02,  1.22800542949647e-03},
5583d177a5cSEmil Constantinescu                 {  1.00000000000000e+00,  1.48175480533865e+01,  3.73439117461835e+00,  5.14648336541804e-01,  4.76430038853402e-02,  2.56798515502156e-03},
5593d177a5cSEmil Constantinescu                 {  1.00000000000000e+00,  1.50512347758335e+01,  4.10099701165164e+00,  5.66039141003603e-01,  3.91213893800891e-02, -2.99136269067853e-03}}
5603d177a5cSEmil Constantinescu     ,v[6][6] = {{  1.00000000000000e+00,  1.50512347758335e+01,  4.10099701165164e+00,  5.66039141003603e-01,  3.91213893800891e-02, -2.99136269067853e-03},
5613d177a5cSEmil Constantinescu                 {  0.00000000000000e+00, -4.88498130835069e-15, -6.43929354282591e-15, -3.55271367880050e-15, -1.22124532708767e-15, -3.12250225675825e-16},
5623d177a5cSEmil Constantinescu                 {  0.00000000000000e+00,  1.22250171233141e+01, -1.77150760606169e+00,  3.54516769879390e-01,  6.22298845883398e-01,  2.31647447450276e-01},
5633d177a5cSEmil Constantinescu                 {  0.00000000000000e+00, -4.48339457331040e+01, -3.57363126641880e-01,  5.18750173123425e-01,  6.55727990241799e-02,  1.63175368287079e-01},
5643d177a5cSEmil Constantinescu                 {  0.00000000000000e+00,  1.37297394708005e+02, -1.60145272991317e+00, -5.05319555199441e+00,  1.55328940390990e-01,  9.16629423682464e-01},
5653d177a5cSEmil Constantinescu                 {  0.00000000000000e+00,  1.05703241119022e+02, -1.16610260983038e+00, -2.99767252773859e+00, -1.13472315553890e-01,  1.09742849254729e+00}};
566*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TSGLLESchemeCreate(5,5,6,6,c,*a,*b,*u,*v,&gl->schemes[gl->nschemes++]));
5673d177a5cSEmil Constantinescu   }
5683d177a5cSEmil Constantinescu   PetscFunctionReturn(0);
5693d177a5cSEmil Constantinescu }
5703d177a5cSEmil Constantinescu 
5713d177a5cSEmil Constantinescu /*@C
5723d177a5cSEmil Constantinescu    TSGLLESetType - sets the class of general linear method to use for time-stepping
5733d177a5cSEmil Constantinescu 
5743d177a5cSEmil Constantinescu    Collective on TS
5753d177a5cSEmil Constantinescu 
5763d177a5cSEmil Constantinescu    Input Parameters:
5773d177a5cSEmil Constantinescu +  ts - the TS context
5783d177a5cSEmil Constantinescu -  type - a method
5793d177a5cSEmil Constantinescu 
5803d177a5cSEmil Constantinescu    Options Database Key:
5813d177a5cSEmil Constantinescu .  -ts_gl_type <type> - sets the method, use -help for a list of available method (e.g. irks)
5823d177a5cSEmil Constantinescu 
5833d177a5cSEmil Constantinescu    Notes:
5843d177a5cSEmil Constantinescu    See "petsc/include/petscts.h" for available methods (for instance)
5853d177a5cSEmil Constantinescu .    TSGLLE_IRKS - Diagonally implicit methods with inherent Runge-Kutta stability (for stiff problems)
5863d177a5cSEmil Constantinescu 
5873d177a5cSEmil Constantinescu    Normally, it is best to use the TSSetFromOptions() command and
5883d177a5cSEmil Constantinescu    then set the TSGLLE type from the options database rather than by using
5893d177a5cSEmil Constantinescu    this routine.  Using the options database provides the user with
5903d177a5cSEmil Constantinescu    maximum flexibility in evaluating the many different solvers.
5913d177a5cSEmil Constantinescu    The TSGLLESetType() routine is provided for those situations where it
5923d177a5cSEmil Constantinescu    is necessary to set the timestepping solver independently of the
5933d177a5cSEmil Constantinescu    command line or options database.  This might be the case, for example,
5943d177a5cSEmil Constantinescu    when the choice of solver changes during the execution of the
5953d177a5cSEmil Constantinescu    program, and the user's application is taking responsibility for
5963d177a5cSEmil Constantinescu    choosing the appropriate method.
5973d177a5cSEmil Constantinescu 
5983d177a5cSEmil Constantinescu    Level: intermediate
5993d177a5cSEmil Constantinescu 
6003d177a5cSEmil Constantinescu @*/
6013d177a5cSEmil Constantinescu PetscErrorCode  TSGLLESetType(TS ts,TSGLLEType type)
6023d177a5cSEmil Constantinescu {
6033d177a5cSEmil Constantinescu   PetscFunctionBegin;
6043d177a5cSEmil Constantinescu   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
6053d177a5cSEmil Constantinescu   PetscValidCharPointer(type,2);
606*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscTryMethod(ts,"TSGLLESetType_C",(TS,TSGLLEType),(ts,type)));
6073d177a5cSEmil Constantinescu   PetscFunctionReturn(0);
6083d177a5cSEmil Constantinescu }
6093d177a5cSEmil Constantinescu 
6103d177a5cSEmil Constantinescu /*@C
6113d177a5cSEmil Constantinescu    TSGLLESetAcceptType - sets the acceptance test
6123d177a5cSEmil Constantinescu 
6133d177a5cSEmil Constantinescu    Time integrators that need to control error must have the option to reject a time step based on local error
6143d177a5cSEmil Constantinescu    estimates.  This function allows different schemes to be set.
6153d177a5cSEmil Constantinescu 
6163d177a5cSEmil Constantinescu    Logically Collective on TS
6173d177a5cSEmil Constantinescu 
6183d177a5cSEmil Constantinescu    Input Parameters:
6193d177a5cSEmil Constantinescu +  ts - the TS context
6203d177a5cSEmil Constantinescu -  type - the type
6213d177a5cSEmil Constantinescu 
6223d177a5cSEmil Constantinescu    Options Database Key:
6233d177a5cSEmil Constantinescu .  -ts_gl_accept_type <type> - sets the method used to determine whether to accept or reject a step
6243d177a5cSEmil Constantinescu 
6253d177a5cSEmil Constantinescu    Level: intermediate
6263d177a5cSEmil Constantinescu 
6273d177a5cSEmil Constantinescu .seealso: TS, TSGLLE, TSGLLEAcceptRegister(), TSGLLEAdapt, set type
6283d177a5cSEmil Constantinescu @*/
6293d177a5cSEmil Constantinescu PetscErrorCode  TSGLLESetAcceptType(TS ts,TSGLLEAcceptType type)
6303d177a5cSEmil Constantinescu {
6313d177a5cSEmil Constantinescu   PetscFunctionBegin;
6323d177a5cSEmil Constantinescu   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
6333d177a5cSEmil Constantinescu   PetscValidCharPointer(type,2);
634*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscTryMethod(ts,"TSGLLESetAcceptType_C",(TS,TSGLLEAcceptType),(ts,type)));
6353d177a5cSEmil Constantinescu   PetscFunctionReturn(0);
6363d177a5cSEmil Constantinescu }
6373d177a5cSEmil Constantinescu 
6383d177a5cSEmil Constantinescu /*@C
6393d177a5cSEmil Constantinescu    TSGLLEGetAdapt - gets the TSGLLEAdapt object from the TS
6403d177a5cSEmil Constantinescu 
6413d177a5cSEmil Constantinescu    Not Collective
6423d177a5cSEmil Constantinescu 
6433d177a5cSEmil Constantinescu    Input Parameter:
6443d177a5cSEmil Constantinescu .  ts - the TS context
6453d177a5cSEmil Constantinescu 
6463d177a5cSEmil Constantinescu    Output Parameter:
6473d177a5cSEmil Constantinescu .  adapt - the TSGLLEAdapt context
6483d177a5cSEmil Constantinescu 
6493d177a5cSEmil Constantinescu    Notes:
6503d177a5cSEmil Constantinescu    This allows the user set options on the TSGLLEAdapt object.  Usually it is better to do this using the options
6513d177a5cSEmil Constantinescu    database, so this function is rarely needed.
6523d177a5cSEmil Constantinescu 
6533d177a5cSEmil Constantinescu    Level: advanced
6543d177a5cSEmil Constantinescu 
6553d177a5cSEmil Constantinescu .seealso: TSGLLEAdapt, TSGLLEAdaptRegister()
6563d177a5cSEmil Constantinescu @*/
6573d177a5cSEmil Constantinescu PetscErrorCode  TSGLLEGetAdapt(TS ts,TSGLLEAdapt *adapt)
6583d177a5cSEmil Constantinescu {
6593d177a5cSEmil Constantinescu   PetscFunctionBegin;
6603d177a5cSEmil Constantinescu   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
6613d177a5cSEmil Constantinescu   PetscValidPointer(adapt,2);
662*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscUseMethod(ts,"TSGLLEGetAdapt_C",(TS,TSGLLEAdapt*),(ts,adapt)));
6633d177a5cSEmil Constantinescu   PetscFunctionReturn(0);
6643d177a5cSEmil Constantinescu }
6653d177a5cSEmil Constantinescu 
6663d177a5cSEmil Constantinescu static PetscErrorCode TSGLLEAccept_Always(TS ts,PetscReal tleft,PetscReal h,const PetscReal enorms[],PetscBool  *accept)
6673d177a5cSEmil Constantinescu {
6683d177a5cSEmil Constantinescu   PetscFunctionBegin;
6693d177a5cSEmil Constantinescu   *accept = PETSC_TRUE;
6703d177a5cSEmil Constantinescu   PetscFunctionReturn(0);
6713d177a5cSEmil Constantinescu }
6723d177a5cSEmil Constantinescu 
6733d177a5cSEmil Constantinescu static PetscErrorCode TSGLLEUpdateWRMS(TS ts)
6743d177a5cSEmil Constantinescu {
6753d177a5cSEmil Constantinescu   TS_GLLE        *gl = (TS_GLLE*)ts->data;
6763d177a5cSEmil Constantinescu   PetscScalar    *x,*w;
6773d177a5cSEmil Constantinescu   PetscInt       n,i;
6783d177a5cSEmil Constantinescu 
6793d177a5cSEmil Constantinescu   PetscFunctionBegin;
680*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(gl->X[0],&x));
681*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(gl->W,&w));
682*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetLocalSize(gl->W,&n));
6833d177a5cSEmil Constantinescu   for (i=0; i<n; i++) w[i] = 1./(gl->wrms_atol + gl->wrms_rtol*PetscAbsScalar(x[i]));
684*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(gl->X[0],&x));
685*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(gl->W,&w));
6863d177a5cSEmil Constantinescu   PetscFunctionReturn(0);
6873d177a5cSEmil Constantinescu }
6883d177a5cSEmil Constantinescu 
6893d177a5cSEmil Constantinescu static PetscErrorCode TSGLLEVecNormWRMS(TS ts,Vec X,PetscReal *nrm)
6903d177a5cSEmil Constantinescu {
6913d177a5cSEmil Constantinescu   TS_GLLE        *gl = (TS_GLLE*)ts->data;
6923d177a5cSEmil Constantinescu   PetscScalar    *x,*w;
6933d177a5cSEmil Constantinescu   PetscReal      sum = 0.0,gsum;
6943d177a5cSEmil Constantinescu   PetscInt       n,N,i;
6953d177a5cSEmil Constantinescu 
6963d177a5cSEmil Constantinescu   PetscFunctionBegin;
697*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(X,&x));
698*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(gl->W,&w));
699*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetLocalSize(gl->W,&n));
7003d177a5cSEmil Constantinescu   for (i=0; i<n; i++) sum += PetscAbsScalar(PetscSqr(x[i]*w[i]));
701*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(X,&x));
702*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(gl->W,&w));
703*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPIU_Allreduce(&sum,&gsum,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)ts)));
704*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetSize(gl->W,&N));
7053d177a5cSEmil Constantinescu   *nrm = PetscSqrtReal(gsum/(1.*N));
7063d177a5cSEmil Constantinescu   PetscFunctionReturn(0);
7073d177a5cSEmil Constantinescu }
7083d177a5cSEmil Constantinescu 
7093d177a5cSEmil Constantinescu static PetscErrorCode TSGLLESetType_GLLE(TS ts,TSGLLEType type)
7103d177a5cSEmil Constantinescu {
7113d177a5cSEmil Constantinescu   PetscBool      same;
7123d177a5cSEmil Constantinescu   TS_GLLE        *gl = (TS_GLLE*)ts->data;
713*5f80ce2aSJacob Faibussowitsch   PetscErrorCode (*r)(TS);
7143d177a5cSEmil Constantinescu 
7153d177a5cSEmil Constantinescu   PetscFunctionBegin;
7163d177a5cSEmil Constantinescu   if (gl->type_name[0]) {
717*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscStrcmp(gl->type_name,type,&same));
7183d177a5cSEmil Constantinescu     if (same) PetscFunctionReturn(0);
719*5f80ce2aSJacob Faibussowitsch     CHKERRQ((*gl->Destroy)(gl));
7203d177a5cSEmil Constantinescu   }
7213d177a5cSEmil Constantinescu 
722*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFunctionListFind(TSGLLEList,type,&r));
7233c633725SBarry Smith   PetscCheck(r,PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown TSGLLE type \"%s\" given",type);
724*5f80ce2aSJacob Faibussowitsch   CHKERRQ((*r)(ts));
725*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscStrcpy(gl->type_name,type));
7263d177a5cSEmil Constantinescu   PetscFunctionReturn(0);
7273d177a5cSEmil Constantinescu }
7283d177a5cSEmil Constantinescu 
7293d177a5cSEmil Constantinescu static PetscErrorCode TSGLLESetAcceptType_GLLE(TS ts,TSGLLEAcceptType type)
7303d177a5cSEmil Constantinescu {
7313d177a5cSEmil Constantinescu   TSGLLEAcceptFunction r;
7323d177a5cSEmil Constantinescu   TS_GLLE              *gl = (TS_GLLE*)ts->data;
7333d177a5cSEmil Constantinescu 
7343d177a5cSEmil Constantinescu   PetscFunctionBegin;
735*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFunctionListFind(TSGLLEAcceptList,type,&r));
7363c633725SBarry Smith   PetscCheck(r,PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown TSGLLEAccept type \"%s\" given",type);
7373d177a5cSEmil Constantinescu   gl->Accept = r;
738*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscStrncpy(gl->accept_name,type,sizeof(gl->accept_name)));
7393d177a5cSEmil Constantinescu   PetscFunctionReturn(0);
7403d177a5cSEmil Constantinescu }
7413d177a5cSEmil Constantinescu 
7423d177a5cSEmil Constantinescu static PetscErrorCode TSGLLEGetAdapt_GLLE(TS ts,TSGLLEAdapt *adapt)
7433d177a5cSEmil Constantinescu {
7443d177a5cSEmil Constantinescu   TS_GLLE        *gl = (TS_GLLE*)ts->data;
7453d177a5cSEmil Constantinescu 
7463d177a5cSEmil Constantinescu   PetscFunctionBegin;
7473d177a5cSEmil Constantinescu   if (!gl->adapt) {
748*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TSGLLEAdaptCreate(PetscObjectComm((PetscObject)ts),&gl->adapt));
749*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectIncrementTabLevel((PetscObject)gl->adapt,(PetscObject)ts,1));
750*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscLogObjectParent((PetscObject)ts,(PetscObject)gl->adapt));
7513d177a5cSEmil Constantinescu   }
7523d177a5cSEmil Constantinescu   *adapt = gl->adapt;
7533d177a5cSEmil Constantinescu   PetscFunctionReturn(0);
7543d177a5cSEmil Constantinescu }
7553d177a5cSEmil Constantinescu 
7563d177a5cSEmil Constantinescu static PetscErrorCode TSGLLEChooseNextScheme(TS ts,PetscReal h,const PetscReal hmnorm[],PetscInt *next_scheme,PetscReal *next_h,PetscBool  *finish)
7573d177a5cSEmil Constantinescu {
7583d177a5cSEmil Constantinescu   TS_GLLE        *gl = (TS_GLLE*)ts->data;
7593d177a5cSEmil Constantinescu   PetscInt       i,n,cur_p,cur,next_sc,candidates[64],orders[64];
7603d177a5cSEmil Constantinescu   PetscReal      errors[64],costs[64],tleft;
7613d177a5cSEmil Constantinescu 
7623d177a5cSEmil Constantinescu   PetscFunctionBegin;
7633d177a5cSEmil Constantinescu   cur   = -1;
7643d177a5cSEmil Constantinescu   cur_p = gl->schemes[gl->current_scheme]->p;
7653d177a5cSEmil Constantinescu   tleft = ts->max_time - (ts->ptime + ts->time_step);
7663d177a5cSEmil Constantinescu   for (i=0,n=0; i<gl->nschemes; i++) {
7673d177a5cSEmil Constantinescu     TSGLLEScheme sc = gl->schemes[i];
7683d177a5cSEmil Constantinescu     if (sc->p < gl->min_order || gl->max_order < sc->p) continue;
7693d177a5cSEmil Constantinescu     if (sc->p == cur_p - 1)    errors[n] = PetscAbsScalar(sc->alpha[0])*hmnorm[0];
7703d177a5cSEmil Constantinescu     else if (sc->p == cur_p)   errors[n] = PetscAbsScalar(sc->alpha[0])*hmnorm[1];
7713d177a5cSEmil Constantinescu     else if (sc->p == cur_p+1) errors[n] = PetscAbsScalar(sc->alpha[0])*(hmnorm[2]+hmnorm[3]);
7723d177a5cSEmil Constantinescu     else continue;
7733d177a5cSEmil Constantinescu     candidates[n] = i;
7743d177a5cSEmil Constantinescu     orders[n]     = PetscMin(sc->p,sc->q); /* order of global truncation error */
7753d177a5cSEmil Constantinescu     costs[n]      = sc->s;                 /* estimate the cost as the number of stages */
7763d177a5cSEmil Constantinescu     if (i == gl->current_scheme) cur = n;
7773d177a5cSEmil Constantinescu     n++;
7783d177a5cSEmil Constantinescu   }
7792c71b3e2SJacob Faibussowitsch   PetscCheckFalse(cur < 0 || gl->nschemes <= cur,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Current scheme not found in scheme list");
780*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGLLEAdaptChoose(gl->adapt,n,orders,errors,costs,cur,h,tleft,&next_sc,next_h,finish));
7813d177a5cSEmil Constantinescu   *next_scheme = candidates[next_sc];
782*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscInfo(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));
7833d177a5cSEmil Constantinescu   PetscFunctionReturn(0);
7843d177a5cSEmil Constantinescu }
7853d177a5cSEmil Constantinescu 
7863d177a5cSEmil Constantinescu static PetscErrorCode TSGLLEGetMaxSizes(TS ts,PetscInt *max_r,PetscInt *max_s)
7873d177a5cSEmil Constantinescu {
7883d177a5cSEmil Constantinescu   TS_GLLE *gl = (TS_GLLE*)ts->data;
7893d177a5cSEmil Constantinescu 
7903d177a5cSEmil Constantinescu   PetscFunctionBegin;
7913d177a5cSEmil Constantinescu   *max_r = gl->schemes[gl->nschemes-1]->r;
7923d177a5cSEmil Constantinescu   *max_s = gl->schemes[gl->nschemes-1]->s;
7933d177a5cSEmil Constantinescu   PetscFunctionReturn(0);
7943d177a5cSEmil Constantinescu }
7953d177a5cSEmil Constantinescu 
7963d177a5cSEmil Constantinescu static PetscErrorCode TSSolve_GLLE(TS ts)
7973d177a5cSEmil Constantinescu {
7983d177a5cSEmil Constantinescu   TS_GLLE             *gl = (TS_GLLE*)ts->data;
7993d177a5cSEmil Constantinescu   PetscInt            i,k,its,lits,max_r,max_s;
8003d177a5cSEmil Constantinescu   PetscBool           final_step,finish;
8013d177a5cSEmil Constantinescu   SNESConvergedReason snesreason;
8023d177a5cSEmil Constantinescu 
8033d177a5cSEmil Constantinescu   PetscFunctionBegin;
804*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol));
8053d177a5cSEmil Constantinescu 
806*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGLLEGetMaxSizes(ts,&max_r,&max_s));
807*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCopy(ts->vec_sol,gl->X[0]));
8083d177a5cSEmil Constantinescu   for (i=1; i<max_r; i++) {
809*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecZeroEntries(gl->X[i]));
8103d177a5cSEmil Constantinescu   }
811*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGLLEUpdateWRMS(ts));
8123d177a5cSEmil Constantinescu 
8133d177a5cSEmil Constantinescu   if (0) {
8143d177a5cSEmil Constantinescu     /* Find consistent initial data for DAE */
8153d177a5cSEmil Constantinescu     gl->stage_time = ts->ptime + ts->time_step;
8163d177a5cSEmil Constantinescu     gl->scoeff = 1.;
8173d177a5cSEmil Constantinescu     gl->stage  = 0;
8183d177a5cSEmil Constantinescu 
819*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecCopy(ts->vec_sol,gl->Z));
820*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecCopy(ts->vec_sol,gl->Y));
821*5f80ce2aSJacob Faibussowitsch     CHKERRQ(SNESSolve(ts->snes,NULL,gl->Y));
822*5f80ce2aSJacob Faibussowitsch     CHKERRQ(SNESGetIterationNumber(ts->snes,&its));
823*5f80ce2aSJacob Faibussowitsch     CHKERRQ(SNESGetLinearSolveIterations(ts->snes,&lits));
824*5f80ce2aSJacob Faibussowitsch     CHKERRQ(SNESGetConvergedReason(ts->snes,&snesreason));
8253d177a5cSEmil Constantinescu 
8263d177a5cSEmil Constantinescu     ts->snes_its += its; ts->ksp_its += lits;
8273d177a5cSEmil Constantinescu     if (snesreason < 0 && ts->max_snes_failures > 0 && ++ts->num_snes_failures >= ts->max_snes_failures) {
8283d177a5cSEmil Constantinescu       ts->reason = TS_DIVERGED_NONLINEAR_SOLVE;
829*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscInfo(ts,"Step=%D, nonlinear solve solve failures %D greater than current TS allowed, stopping solve\n",ts->steps,ts->num_snes_failures));
8303d177a5cSEmil Constantinescu       PetscFunctionReturn(0);
8313d177a5cSEmil Constantinescu     }
8323d177a5cSEmil Constantinescu   }
8333d177a5cSEmil Constantinescu 
8342c71b3e2SJacob Faibussowitsch   PetscCheckFalse(gl->current_scheme < 0,PETSC_COMM_SELF,PETSC_ERR_ORDER,"A starting scheme has not been provided");
8353d177a5cSEmil Constantinescu 
8363d177a5cSEmil Constantinescu   for (k=0,final_step=PETSC_FALSE,finish=PETSC_FALSE; k<ts->max_steps && !finish; k++) {
8373d177a5cSEmil Constantinescu     PetscInt          j,r,s,next_scheme = 0,rejections;
8383d177a5cSEmil Constantinescu     PetscReal         h,hmnorm[4],enorm[3],next_h;
8393d177a5cSEmil Constantinescu     PetscBool         accept;
8403d177a5cSEmil Constantinescu     const PetscScalar *c,*a,*u;
8413d177a5cSEmil Constantinescu     Vec               *X,*Ydot,Y;
8423d177a5cSEmil Constantinescu     TSGLLEScheme        scheme = gl->schemes[gl->current_scheme];
8433d177a5cSEmil Constantinescu 
8443d177a5cSEmil Constantinescu     r = scheme->r; s = scheme->s;
8453d177a5cSEmil Constantinescu     c = scheme->c;
8463d177a5cSEmil Constantinescu     a = scheme->a; u = scheme->u;
8473d177a5cSEmil Constantinescu     h = ts->time_step;
8483d177a5cSEmil Constantinescu     X = gl->X; Ydot = gl->Ydot; Y = gl->Y;
8493d177a5cSEmil Constantinescu 
8503d177a5cSEmil Constantinescu     if (ts->ptime > ts->max_time) break;
8513d177a5cSEmil Constantinescu 
8523d177a5cSEmil Constantinescu     /*
8533d177a5cSEmil Constantinescu       We only call PreStep at the start of each STEP, not each STAGE.  This is because it is
8543d177a5cSEmil Constantinescu       possible to fail (have to restart a step) after multiple stages.
8553d177a5cSEmil Constantinescu     */
856*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TSPreStep(ts));
8573d177a5cSEmil Constantinescu 
8583d177a5cSEmil Constantinescu     rejections = 0;
8593d177a5cSEmil Constantinescu     while (1) {
8603d177a5cSEmil Constantinescu       for (i=0; i<s; i++) {
8613d177a5cSEmil Constantinescu         PetscScalar shift;
8623d177a5cSEmil Constantinescu         gl->scoeff     = 1./PetscRealPart(a[i*s+i]);
8633d177a5cSEmil Constantinescu         shift          = gl->scoeff/ts->time_step;
8643d177a5cSEmil Constantinescu         gl->stage      = i;
8653d177a5cSEmil Constantinescu         gl->stage_time = ts->ptime + PetscRealPart(c[i])*h;
8663d177a5cSEmil Constantinescu 
8673d177a5cSEmil Constantinescu         /*
8683d177a5cSEmil Constantinescu         * Stage equation: Y = h A Y' + U X
8693d177a5cSEmil Constantinescu         * We assume that A is lower-triangular so that we can solve the stages (Y,Y') sequentially
8703d177a5cSEmil Constantinescu         * Build the affine vector z_i = -[1/(h a_ii)](h sum_j a_ij y'_j + sum_j u_ij x_j)
8713d177a5cSEmil Constantinescu         * Then y'_i = z + 1/(h a_ii) y_i
8723d177a5cSEmil Constantinescu         */
873*5f80ce2aSJacob Faibussowitsch         CHKERRQ(VecZeroEntries(gl->Z));
8743d177a5cSEmil Constantinescu         for (j=0; j<r; j++) {
875*5f80ce2aSJacob Faibussowitsch           CHKERRQ(VecAXPY(gl->Z,-shift*u[i*r+j],X[j]));
8763d177a5cSEmil Constantinescu         }
8773d177a5cSEmil Constantinescu         for (j=0; j<i; j++) {
878*5f80ce2aSJacob Faibussowitsch           CHKERRQ(VecAXPY(gl->Z,-shift*h*a[i*s+j],Ydot[j]));
8793d177a5cSEmil Constantinescu         }
8803d177a5cSEmil Constantinescu         /* Note: Z is used within function evaluation, Ydot = Z + shift*Y */
8813d177a5cSEmil Constantinescu 
8823d177a5cSEmil Constantinescu         /* Compute an estimate of Y to start Newton iteration */
8833d177a5cSEmil Constantinescu         if (gl->extrapolate) {
8843d177a5cSEmil Constantinescu           if (i==0) {
8853d177a5cSEmil Constantinescu             /* Linear extrapolation on the first stage */
886*5f80ce2aSJacob Faibussowitsch             CHKERRQ(VecWAXPY(Y,c[i]*h,X[1],X[0]));
8873d177a5cSEmil Constantinescu           } else {
8883d177a5cSEmil Constantinescu             /* Linear extrapolation from the last stage */
889*5f80ce2aSJacob Faibussowitsch             CHKERRQ(VecAXPY(Y,(c[i]-c[i-1])*h,Ydot[i-1]));
8903d177a5cSEmil Constantinescu           }
8913d177a5cSEmil Constantinescu         } else if (i==0) {        /* Directly use solution from the last step, otherwise reuse the last stage (do nothing) */
892*5f80ce2aSJacob Faibussowitsch           CHKERRQ(VecCopy(X[0],Y));
8933d177a5cSEmil Constantinescu         }
8943d177a5cSEmil Constantinescu 
8953d177a5cSEmil Constantinescu         /* Solve this stage (Ydot[i] is computed during function evaluation) */
896*5f80ce2aSJacob Faibussowitsch         CHKERRQ(SNESSolve(ts->snes,NULL,Y));
897*5f80ce2aSJacob Faibussowitsch         CHKERRQ(SNESGetIterationNumber(ts->snes,&its));
898*5f80ce2aSJacob Faibussowitsch         CHKERRQ(SNESGetLinearSolveIterations(ts->snes,&lits));
899*5f80ce2aSJacob Faibussowitsch         CHKERRQ(SNESGetConvergedReason(ts->snes,&snesreason));
9003d177a5cSEmil Constantinescu         ts->snes_its += its; ts->ksp_its += lits;
9013d177a5cSEmil Constantinescu         if (snesreason < 0 && ts->max_snes_failures > 0 && ++ts->num_snes_failures >= ts->max_snes_failures) {
9023d177a5cSEmil Constantinescu           ts->reason = TS_DIVERGED_NONLINEAR_SOLVE;
903*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscInfo(ts,"Step=%D, nonlinear solve solve failures %D greater than current TS allowed, stopping solve\n",ts->steps,ts->num_snes_failures));
9043d177a5cSEmil Constantinescu           PetscFunctionReturn(0);
9053d177a5cSEmil Constantinescu         }
9063d177a5cSEmil Constantinescu       }
9073d177a5cSEmil Constantinescu 
9083d177a5cSEmil Constantinescu       gl->stage_time = ts->ptime + ts->time_step;
9093d177a5cSEmil Constantinescu 
910*5f80ce2aSJacob Faibussowitsch       CHKERRQ((*gl->EstimateHigherMoments)(scheme,h,Ydot,gl->X,gl->himom));
9113d177a5cSEmil Constantinescu       /* hmnorm[i] = h^{p+i}x^{(p+i)} with i=0,1,2; hmnorm[3] = h^{p+2}(dx'/dx) x^{(p+1)} */
9123d177a5cSEmil Constantinescu       for (i=0; i<3; i++) {
913*5f80ce2aSJacob Faibussowitsch         CHKERRQ(TSGLLEVecNormWRMS(ts,gl->himom[i],&hmnorm[i+1]));
9143d177a5cSEmil Constantinescu       }
9153d177a5cSEmil Constantinescu       enorm[0] = PetscRealPart(scheme->alpha[0])*hmnorm[1];
9163d177a5cSEmil Constantinescu       enorm[1] = PetscRealPart(scheme->beta[0]) *hmnorm[2];
9173d177a5cSEmil Constantinescu       enorm[2] = PetscRealPart(scheme->gamma[0])*hmnorm[3];
918*5f80ce2aSJacob Faibussowitsch       CHKERRQ((*gl->Accept)(ts,ts->max_time-gl->stage_time,h,enorm,&accept));
9193d177a5cSEmil Constantinescu       if (accept) goto accepted;
9203d177a5cSEmil Constantinescu       rejections++;
921*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscInfo(ts,"Step %D (t=%g) not accepted, rejections=%D\n",k,gl->stage_time,rejections));
9223d177a5cSEmil Constantinescu       if (rejections > gl->max_step_rejections) break;
9233d177a5cSEmil Constantinescu       /*
9243d177a5cSEmil Constantinescu         There are lots of reasons why a step might be rejected, including solvers not converging and other factors that
9253d177a5cSEmil Constantinescu         TSGLLEChooseNextScheme does not support.  Additionally, the error estimates may be very screwed up, so I'm not
9263d177a5cSEmil Constantinescu         convinced that it's safe to just compute a new error estimate using the same interface as the current adaptor
9273d177a5cSEmil Constantinescu         (the adaptor interface probably has to change).  Here we make an arbitrary and naive choice.  This assumes that
9283d177a5cSEmil Constantinescu         steps were written in Nordsieck form.  The "correct" method would be to re-complete the previous time step with
9293d177a5cSEmil Constantinescu         the correct "next" step size.  It is unclear to me whether the present ad-hoc method of rescaling X is stable.
9303d177a5cSEmil Constantinescu       */
9313d177a5cSEmil Constantinescu       h *= 0.5;
9323d177a5cSEmil Constantinescu       for (i=1; i<scheme->r; i++) {
933*5f80ce2aSJacob Faibussowitsch         CHKERRQ(VecScale(X[i],PetscPowRealInt(0.5,i)));
9343d177a5cSEmil Constantinescu       }
9353d177a5cSEmil Constantinescu     }
93698921bdaSJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_CONV_FAILED,"Time step %D (t=%g) not accepted after %D failures",k,gl->stage_time,rejections);
9373d177a5cSEmil Constantinescu 
9383d177a5cSEmil Constantinescu accepted:
9393d177a5cSEmil Constantinescu     /* This term is not error, but it *would* be the leading term for a lower order method */
940*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TSGLLEVecNormWRMS(ts,gl->X[scheme->r-1],&hmnorm[0]));
9413d177a5cSEmil Constantinescu     /* Correct scaling so that these are equivalent to norms of the Nordsieck vectors */
9423d177a5cSEmil Constantinescu 
943*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscInfo(ts,"Last moment norm %10.2e, estimated error norms %10.2e %10.2e %10.2e\n",hmnorm[0],enorm[0],enorm[1],enorm[2]));
9443d177a5cSEmil Constantinescu     if (!final_step) {
945*5f80ce2aSJacob Faibussowitsch       CHKERRQ(TSGLLEChooseNextScheme(ts,h,hmnorm,&next_scheme,&next_h,&final_step));
9463d177a5cSEmil Constantinescu     } else {
9473d177a5cSEmil Constantinescu       /* Dummy values to complete the current step in a consistent manner */
9483d177a5cSEmil Constantinescu       next_scheme = gl->current_scheme;
9493d177a5cSEmil Constantinescu       next_h      = h;
9503d177a5cSEmil Constantinescu       finish      = PETSC_TRUE;
9513d177a5cSEmil Constantinescu     }
9523d177a5cSEmil Constantinescu 
9533d177a5cSEmil Constantinescu     X        = gl->Xold;
9543d177a5cSEmil Constantinescu     gl->Xold = gl->X;
9553d177a5cSEmil Constantinescu     gl->X    = X;
956*5f80ce2aSJacob Faibussowitsch     CHKERRQ((*gl->CompleteStep)(scheme,h,gl->schemes[next_scheme],next_h,Ydot,gl->Xold,gl->X));
9573d177a5cSEmil Constantinescu 
958*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TSGLLEUpdateWRMS(ts));
9593d177a5cSEmil Constantinescu 
9603d177a5cSEmil Constantinescu     /* Post the solution for the user, we could avoid this copy with a small bit of cleverness */
961*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecCopy(gl->X[0],ts->vec_sol));
9623d177a5cSEmil Constantinescu     ts->ptime += h;
9633d177a5cSEmil Constantinescu     ts->steps++;
9643d177a5cSEmil Constantinescu 
965*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TSPostEvaluate(ts));
966*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TSPostStep(ts));
967*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol));
9683d177a5cSEmil Constantinescu 
9693d177a5cSEmil Constantinescu     gl->current_scheme = next_scheme;
9703d177a5cSEmil Constantinescu     ts->time_step      = next_h;
9713d177a5cSEmil Constantinescu   }
9723d177a5cSEmil Constantinescu   PetscFunctionReturn(0);
9733d177a5cSEmil Constantinescu }
9743d177a5cSEmil Constantinescu 
9753d177a5cSEmil Constantinescu /*------------------------------------------------------------*/
9763d177a5cSEmil Constantinescu 
9773d177a5cSEmil Constantinescu static PetscErrorCode TSReset_GLLE(TS ts)
9783d177a5cSEmil Constantinescu {
9793d177a5cSEmil Constantinescu   TS_GLLE        *gl = (TS_GLLE*)ts->data;
9803d177a5cSEmil Constantinescu   PetscInt       max_r,max_s;
9813d177a5cSEmil Constantinescu 
9823d177a5cSEmil Constantinescu   PetscFunctionBegin;
9833d177a5cSEmil Constantinescu   if (gl->setupcalled) {
984*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TSGLLEGetMaxSizes(ts,&max_r,&max_s));
985*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroyVecs(max_r,&gl->Xold));
986*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroyVecs(max_r,&gl->X));
987*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroyVecs(max_s,&gl->Ydot));
988*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroyVecs(3,&gl->himom));
989*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&gl->W));
990*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&gl->Y));
991*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&gl->Z));
9923d177a5cSEmil Constantinescu   }
9933d177a5cSEmil Constantinescu   gl->setupcalled = PETSC_FALSE;
9943d177a5cSEmil Constantinescu   PetscFunctionReturn(0);
9953d177a5cSEmil Constantinescu }
9963d177a5cSEmil Constantinescu 
9973d177a5cSEmil Constantinescu static PetscErrorCode TSDestroy_GLLE(TS ts)
9983d177a5cSEmil Constantinescu {
9993d177a5cSEmil Constantinescu   TS_GLLE        *gl = (TS_GLLE*)ts->data;
10003d177a5cSEmil Constantinescu 
10013d177a5cSEmil Constantinescu   PetscFunctionBegin;
1002*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSReset_GLLE(ts));
1003b3a6b972SJed Brown   if (ts->dm) {
1004*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSGLLE,DMRestrictHook_TSGLLE,ts));
1005*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSGLLE,DMSubDomainRestrictHook_TSGLLE,ts));
1006b3a6b972SJed Brown   }
1007*5f80ce2aSJacob Faibussowitsch   if (gl->adapt) CHKERRQ(TSGLLEAdaptDestroy(&gl->adapt));
1008*5f80ce2aSJacob Faibussowitsch   if (gl->Destroy) CHKERRQ((*gl->Destroy)(gl));
1009*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(ts->data));
1010*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectComposeFunction((PetscObject)ts,"TSGLLESetType_C",NULL));
1011*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectComposeFunction((PetscObject)ts,"TSGLLESetAcceptType_C",NULL));
1012*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectComposeFunction((PetscObject)ts,"TSGLLEGetAdapt_C",NULL));
10133d177a5cSEmil Constantinescu   PetscFunctionReturn(0);
10143d177a5cSEmil Constantinescu }
10153d177a5cSEmil Constantinescu 
10163d177a5cSEmil Constantinescu /*
10173d177a5cSEmil Constantinescu     This defines the nonlinear equation that is to be solved with SNES
10183d177a5cSEmil Constantinescu     g(x) = f(t,x,z+shift*x) = 0
10193d177a5cSEmil Constantinescu */
10203d177a5cSEmil Constantinescu static PetscErrorCode SNESTSFormFunction_GLLE(SNES snes,Vec x,Vec f,TS ts)
10213d177a5cSEmil Constantinescu {
10223d177a5cSEmil Constantinescu   TS_GLLE        *gl = (TS_GLLE*)ts->data;
10233d177a5cSEmil Constantinescu   Vec            Z,Ydot;
10243d177a5cSEmil Constantinescu   DM             dm,dmsave;
10253d177a5cSEmil Constantinescu 
10263d177a5cSEmil Constantinescu   PetscFunctionBegin;
1027*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetDM(snes,&dm));
1028*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGLLEGetVecs(ts,dm,&Z,&Ydot));
1029*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecWAXPY(Ydot,gl->scoeff/ts->time_step,x,Z));
10303d177a5cSEmil Constantinescu   dmsave = ts->dm;
10313d177a5cSEmil Constantinescu   ts->dm = dm;
1032*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSComputeIFunction(ts,gl->stage_time,x,Ydot,f,PETSC_FALSE));
10333d177a5cSEmil Constantinescu   ts->dm = dmsave;
1034*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGLLERestoreVecs(ts,dm,&Z,&Ydot));
10353d177a5cSEmil Constantinescu   PetscFunctionReturn(0);
10363d177a5cSEmil Constantinescu }
10373d177a5cSEmil Constantinescu 
10383d177a5cSEmil Constantinescu static PetscErrorCode SNESTSFormJacobian_GLLE(SNES snes,Vec x,Mat A,Mat B,TS ts)
10393d177a5cSEmil Constantinescu {
10403d177a5cSEmil Constantinescu   TS_GLLE        *gl = (TS_GLLE*)ts->data;
10413d177a5cSEmil Constantinescu   Vec            Z,Ydot;
10423d177a5cSEmil Constantinescu   DM             dm,dmsave;
10433d177a5cSEmil Constantinescu 
10443d177a5cSEmil Constantinescu   PetscFunctionBegin;
1045*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetDM(snes,&dm));
1046*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGLLEGetVecs(ts,dm,&Z,&Ydot));
10473d177a5cSEmil Constantinescu   dmsave = ts->dm;
10483d177a5cSEmil Constantinescu   ts->dm = dm;
10493d177a5cSEmil Constantinescu   /* gl->Xdot will have already been computed in SNESTSFormFunction_GLLE */
1050*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSComputeIJacobian(ts,gl->stage_time,x,gl->Ydot[gl->stage],gl->scoeff/ts->time_step,A,B,PETSC_FALSE));
10513d177a5cSEmil Constantinescu   ts->dm = dmsave;
1052*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGLLERestoreVecs(ts,dm,&Z,&Ydot));
10533d177a5cSEmil Constantinescu   PetscFunctionReturn(0);
10543d177a5cSEmil Constantinescu }
10553d177a5cSEmil Constantinescu 
10563d177a5cSEmil Constantinescu static PetscErrorCode TSSetUp_GLLE(TS ts)
10573d177a5cSEmil Constantinescu {
10583d177a5cSEmil Constantinescu   TS_GLLE        *gl = (TS_GLLE*)ts->data;
10593d177a5cSEmil Constantinescu   PetscInt       max_r,max_s;
10603d177a5cSEmil Constantinescu   DM             dm;
10613d177a5cSEmil Constantinescu 
10623d177a5cSEmil Constantinescu   PetscFunctionBegin;
10633d177a5cSEmil Constantinescu   gl->setupcalled = PETSC_TRUE;
1064*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGLLEGetMaxSizes(ts,&max_r,&max_s));
1065*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicateVecs(ts->vec_sol,max_r,&gl->X));
1066*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicateVecs(ts->vec_sol,max_r,&gl->Xold));
1067*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicateVecs(ts->vec_sol,max_s,&gl->Ydot));
1068*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicateVecs(ts->vec_sol,3,&gl->himom));
1069*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(ts->vec_sol,&gl->W));
1070*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(ts->vec_sol,&gl->Y));
1071*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(ts->vec_sol,&gl->Z));
10723d177a5cSEmil Constantinescu 
10733d177a5cSEmil Constantinescu   /* Default acceptance tests and adaptivity */
1074*5f80ce2aSJacob Faibussowitsch   if (!gl->Accept) CHKERRQ(TSGLLESetAcceptType(ts,TSGLLEACCEPT_ALWAYS));
1075*5f80ce2aSJacob Faibussowitsch   if (!gl->adapt)  CHKERRQ(TSGLLEGetAdapt(ts,&gl->adapt));
10763d177a5cSEmil Constantinescu 
10773d177a5cSEmil Constantinescu   if (gl->current_scheme < 0) {
10783d177a5cSEmil Constantinescu     PetscInt i;
10793d177a5cSEmil Constantinescu     for (i=0;; i++) {
10803d177a5cSEmil Constantinescu       if (gl->schemes[i]->p == gl->start_order) break;
10812c71b3e2SJacob Faibussowitsch       PetscCheckFalse(i+1 == gl->nschemes,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"No schemes available with requested start order %d",i);
10823d177a5cSEmil Constantinescu     }
10833d177a5cSEmil Constantinescu     gl->current_scheme = i;
10843d177a5cSEmil Constantinescu   }
1085*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetDM(ts,&dm));
1086*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCoarsenHookAdd(dm,DMCoarsenHook_TSGLLE,DMRestrictHook_TSGLLE,ts));
1087*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSubDomainHookAdd(dm,DMSubDomainHook_TSGLLE,DMSubDomainRestrictHook_TSGLLE,ts));
10883d177a5cSEmil Constantinescu   PetscFunctionReturn(0);
10893d177a5cSEmil Constantinescu }
10903d177a5cSEmil Constantinescu /*------------------------------------------------------------*/
10913d177a5cSEmil Constantinescu 
10923d177a5cSEmil Constantinescu static PetscErrorCode TSSetFromOptions_GLLE(PetscOptionItems *PetscOptionsObject,TS ts)
10933d177a5cSEmil Constantinescu {
10943d177a5cSEmil Constantinescu   TS_GLLE        *gl        = (TS_GLLE*)ts->data;
10953d177a5cSEmil Constantinescu   char           tname[256] = TSGLLE_IRKS,completef[256] = "rescale-and-modify";
10963d177a5cSEmil Constantinescu 
10973d177a5cSEmil Constantinescu   PetscFunctionBegin;
1098*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHead(PetscOptionsObject,"General Linear ODE solver options"));
10993d177a5cSEmil Constantinescu   {
11003d177a5cSEmil Constantinescu     PetscBool flg;
1101*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsFList("-ts_gl_type","Type of GL method","TSGLLESetType",TSGLLEList,gl->type_name[0] ? gl->type_name : tname,tname,sizeof(tname),&flg));
11023d177a5cSEmil Constantinescu     if (flg || !gl->type_name[0]) {
1103*5f80ce2aSJacob Faibussowitsch       CHKERRQ(TSGLLESetType(ts,tname));
11043d177a5cSEmil Constantinescu     }
1105*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsInt("-ts_gl_max_step_rejections","Maximum number of times to attempt a step","None",gl->max_step_rejections,&gl->max_step_rejections,NULL));
1106*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsInt("-ts_gl_max_order","Maximum order to try","TSGLLESetMaxOrder",gl->max_order,&gl->max_order,NULL));
1107*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsInt("-ts_gl_min_order","Minimum order to try","TSGLLESetMinOrder",gl->min_order,&gl->min_order,NULL));
1108*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsInt("-ts_gl_start_order","Initial order to try","TSGLLESetMinOrder",gl->start_order,&gl->start_order,NULL));
1109*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsEnum("-ts_gl_error_direction","Which direction to look when estimating error","TSGLLESetErrorDirection",TSGLLEErrorDirections,(PetscEnum)gl->error_direction,(PetscEnum*)&gl->error_direction,NULL));
1110*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsBool("-ts_gl_extrapolate","Extrapolate stage solution from previous solution (sometimes unstable)","TSGLLESetExtrapolate",gl->extrapolate,&gl->extrapolate,NULL));
1111*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsReal("-ts_gl_atol","Absolute tolerance","TSGLLESetTolerances",gl->wrms_atol,&gl->wrms_atol,NULL));
1112*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsReal("-ts_gl_rtol","Relative tolerance","TSGLLESetTolerances",gl->wrms_rtol,&gl->wrms_rtol,NULL));
1113*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsString("-ts_gl_complete","Method to use for completing the step","none",completef,completef,sizeof(completef),&flg));
11143d177a5cSEmil Constantinescu     if (flg) {
11153d177a5cSEmil Constantinescu       PetscBool match1,match2;
1116*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscStrcmp(completef,"rescale",&match1));
1117*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscStrcmp(completef,"rescale-and-modify",&match2));
11183d177a5cSEmil Constantinescu       if (match1)      gl->CompleteStep = TSGLLECompleteStep_Rescale;
11193d177a5cSEmil Constantinescu       else if (match2) gl->CompleteStep = TSGLLECompleteStep_RescaleAndModify;
112098921bdaSJacob Faibussowitsch       else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"%s",completef);
11213d177a5cSEmil Constantinescu     }
11223d177a5cSEmil Constantinescu     {
11233d177a5cSEmil Constantinescu       char type[256] = TSGLLEACCEPT_ALWAYS;
1124*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscOptionsFList("-ts_gl_accept_type","Method to use for determining whether to accept a step","TSGLLESetAcceptType",TSGLLEAcceptList,gl->accept_name[0] ? gl->accept_name : type,type,sizeof(type),&flg));
11253d177a5cSEmil Constantinescu       if (flg || !gl->accept_name[0]) {
1126*5f80ce2aSJacob Faibussowitsch         CHKERRQ(TSGLLESetAcceptType(ts,type));
11273d177a5cSEmil Constantinescu       }
11283d177a5cSEmil Constantinescu     }
11293d177a5cSEmil Constantinescu     {
11303d177a5cSEmil Constantinescu       TSGLLEAdapt adapt;
1131*5f80ce2aSJacob Faibussowitsch       CHKERRQ(TSGLLEGetAdapt(ts,&adapt));
1132*5f80ce2aSJacob Faibussowitsch       CHKERRQ(TSGLLEAdaptSetFromOptions(PetscOptionsObject,adapt));
11333d177a5cSEmil Constantinescu     }
11343d177a5cSEmil Constantinescu   }
1135*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsTail());
11363d177a5cSEmil Constantinescu   PetscFunctionReturn(0);
11373d177a5cSEmil Constantinescu }
11383d177a5cSEmil Constantinescu 
11393d177a5cSEmil Constantinescu static PetscErrorCode TSView_GLLE(TS ts,PetscViewer viewer)
11403d177a5cSEmil Constantinescu {
11413d177a5cSEmil Constantinescu   TS_GLLE        *gl = (TS_GLLE*)ts->data;
11423d177a5cSEmil Constantinescu   PetscInt       i;
11433d177a5cSEmil Constantinescu   PetscBool      iascii,details;
11443d177a5cSEmil Constantinescu 
11453d177a5cSEmil Constantinescu   PetscFunctionBegin;
1146*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii));
11473d177a5cSEmil Constantinescu   if (iascii) {
1148*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIPrintf(viewer,"  min order %D, max order %D, current order %D\n",gl->min_order,gl->max_order,gl->schemes[gl->current_scheme]->p));
1149*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIPrintf(viewer,"  Error estimation: %s\n",TSGLLEErrorDirections[gl->error_direction]));
1150*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIPrintf(viewer,"  Extrapolation: %s\n",gl->extrapolate ? "yes" : "no"));
1151*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIPrintf(viewer,"  Acceptance test: %s\n",gl->accept_name[0] ? gl->accept_name : "(not yet set)"));
1152*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIPushTab(viewer));
1153*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TSGLLEAdaptView(gl->adapt,viewer));
1154*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIPopTab(viewer));
1155*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIPrintf(viewer,"  type: %s\n",gl->type_name[0] ? gl->type_name : "(not yet set)"));
1156*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIPrintf(viewer,"Schemes within family (%d):\n",gl->nschemes));
11573d177a5cSEmil Constantinescu     details = PETSC_FALSE;
1158*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsGetBool(((PetscObject)ts)->options,((PetscObject)ts)->prefix,"-ts_gl_view_detailed",&details,NULL));
1159*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIPushTab(viewer));
11603d177a5cSEmil Constantinescu     for (i=0; i<gl->nschemes; i++) {
1161*5f80ce2aSJacob Faibussowitsch       CHKERRQ(TSGLLESchemeView(gl->schemes[i],details,viewer));
11623d177a5cSEmil Constantinescu     }
11633d177a5cSEmil Constantinescu     if (gl->View) {
1164*5f80ce2aSJacob Faibussowitsch       CHKERRQ((*gl->View)(gl,viewer));
11653d177a5cSEmil Constantinescu     }
1166*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIPopTab(viewer));
11673d177a5cSEmil Constantinescu   }
11683d177a5cSEmil Constantinescu   PetscFunctionReturn(0);
11693d177a5cSEmil Constantinescu }
11703d177a5cSEmil Constantinescu 
11713d177a5cSEmil Constantinescu /*@C
11723d177a5cSEmil Constantinescu    TSGLLERegister -  adds a TSGLLE implementation
11733d177a5cSEmil Constantinescu 
11743d177a5cSEmil Constantinescu    Not Collective
11753d177a5cSEmil Constantinescu 
11763d177a5cSEmil Constantinescu    Input Parameters:
11773d177a5cSEmil Constantinescu +  name_scheme - name of user-defined general linear scheme
11783d177a5cSEmil Constantinescu -  routine_create - routine to create method context
11793d177a5cSEmil Constantinescu 
11803d177a5cSEmil Constantinescu    Notes:
11813d177a5cSEmil Constantinescu    TSGLLERegister() may be called multiple times to add several user-defined families.
11823d177a5cSEmil Constantinescu 
11833d177a5cSEmil Constantinescu    Sample usage:
11843d177a5cSEmil Constantinescu .vb
11853d177a5cSEmil Constantinescu    TSGLLERegister("my_scheme",MySchemeCreate);
11863d177a5cSEmil Constantinescu .ve
11873d177a5cSEmil Constantinescu 
11883d177a5cSEmil Constantinescu    Then, your scheme can be chosen with the procedural interface via
11893d177a5cSEmil Constantinescu $     TSGLLESetType(ts,"my_scheme")
11903d177a5cSEmil Constantinescu    or at runtime via the option
11913d177a5cSEmil Constantinescu $     -ts_gl_type my_scheme
11923d177a5cSEmil Constantinescu 
11933d177a5cSEmil Constantinescu    Level: advanced
11943d177a5cSEmil Constantinescu 
11953d177a5cSEmil Constantinescu .seealso: TSGLLERegisterAll()
11963d177a5cSEmil Constantinescu @*/
11973d177a5cSEmil Constantinescu PetscErrorCode  TSGLLERegister(const char sname[],PetscErrorCode (*function)(TS))
11983d177a5cSEmil Constantinescu {
11993d177a5cSEmil Constantinescu   PetscFunctionBegin;
1200*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGLLEInitializePackage());
1201*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFunctionListAdd(&TSGLLEList,sname,function));
12023d177a5cSEmil Constantinescu   PetscFunctionReturn(0);
12033d177a5cSEmil Constantinescu }
12043d177a5cSEmil Constantinescu 
12053d177a5cSEmil Constantinescu /*@C
12063d177a5cSEmil Constantinescu    TSGLLEAcceptRegister -  adds a TSGLLE acceptance scheme
12073d177a5cSEmil Constantinescu 
12083d177a5cSEmil Constantinescu    Not Collective
12093d177a5cSEmil Constantinescu 
12103d177a5cSEmil Constantinescu    Input Parameters:
12113d177a5cSEmil Constantinescu +  name_scheme - name of user-defined acceptance scheme
12123d177a5cSEmil Constantinescu -  routine_create - routine to create method context
12133d177a5cSEmil Constantinescu 
12143d177a5cSEmil Constantinescu    Notes:
12153d177a5cSEmil Constantinescu    TSGLLEAcceptRegister() may be called multiple times to add several user-defined families.
12163d177a5cSEmil Constantinescu 
12173d177a5cSEmil Constantinescu    Sample usage:
12183d177a5cSEmil Constantinescu .vb
12193d177a5cSEmil Constantinescu    TSGLLEAcceptRegister("my_scheme",MySchemeCreate);
12203d177a5cSEmil Constantinescu .ve
12213d177a5cSEmil Constantinescu 
12223d177a5cSEmil Constantinescu    Then, your scheme can be chosen with the procedural interface via
12233d177a5cSEmil Constantinescu $     TSGLLESetAcceptType(ts,"my_scheme")
12243d177a5cSEmil Constantinescu    or at runtime via the option
12253d177a5cSEmil Constantinescu $     -ts_gl_accept_type my_scheme
12263d177a5cSEmil Constantinescu 
12273d177a5cSEmil Constantinescu    Level: advanced
12283d177a5cSEmil Constantinescu 
12293d177a5cSEmil Constantinescu .seealso: TSGLLERegisterAll()
12303d177a5cSEmil Constantinescu @*/
12313d177a5cSEmil Constantinescu PetscErrorCode  TSGLLEAcceptRegister(const char sname[],TSGLLEAcceptFunction function)
12323d177a5cSEmil Constantinescu {
12333d177a5cSEmil Constantinescu   PetscFunctionBegin;
1234*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFunctionListAdd(&TSGLLEAcceptList,sname,function));
12353d177a5cSEmil Constantinescu   PetscFunctionReturn(0);
12363d177a5cSEmil Constantinescu }
12373d177a5cSEmil Constantinescu 
12383d177a5cSEmil Constantinescu /*@C
12393d177a5cSEmil Constantinescu   TSGLLERegisterAll - Registers all of the general linear methods in TSGLLE
12403d177a5cSEmil Constantinescu 
12413d177a5cSEmil Constantinescu   Not Collective
12423d177a5cSEmil Constantinescu 
12433d177a5cSEmil Constantinescu   Level: advanced
12443d177a5cSEmil Constantinescu 
12453d177a5cSEmil Constantinescu .seealso:  TSGLLERegisterDestroy()
12463d177a5cSEmil Constantinescu @*/
12473d177a5cSEmil Constantinescu PetscErrorCode  TSGLLERegisterAll(void)
12483d177a5cSEmil Constantinescu {
12493d177a5cSEmil Constantinescu   PetscFunctionBegin;
12503d177a5cSEmil Constantinescu   if (TSGLLERegisterAllCalled) PetscFunctionReturn(0);
12513d177a5cSEmil Constantinescu   TSGLLERegisterAllCalled = PETSC_TRUE;
12523d177a5cSEmil Constantinescu 
1253*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGLLERegister(TSGLLE_IRKS,              TSGLLECreate_IRKS));
1254*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGLLEAcceptRegister(TSGLLEACCEPT_ALWAYS,TSGLLEAccept_Always));
12553d177a5cSEmil Constantinescu   PetscFunctionReturn(0);
12563d177a5cSEmil Constantinescu }
12573d177a5cSEmil Constantinescu 
12583d177a5cSEmil Constantinescu /*@C
12593d177a5cSEmil Constantinescu   TSGLLEInitializePackage - This function initializes everything in the TSGLLE package. It is called
12608a690491SBarry Smith   from TSInitializePackage().
12613d177a5cSEmil Constantinescu 
12623d177a5cSEmil Constantinescu   Level: developer
12633d177a5cSEmil Constantinescu 
12643d177a5cSEmil Constantinescu .seealso: PetscInitialize()
12653d177a5cSEmil Constantinescu @*/
12663d177a5cSEmil Constantinescu PetscErrorCode  TSGLLEInitializePackage(void)
12673d177a5cSEmil Constantinescu {
12683d177a5cSEmil Constantinescu   PetscFunctionBegin;
12693d177a5cSEmil Constantinescu   if (TSGLLEPackageInitialized) PetscFunctionReturn(0);
12703d177a5cSEmil Constantinescu   TSGLLEPackageInitialized = PETSC_TRUE;
1271*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGLLERegisterAll());
1272*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRegisterFinalize(TSGLLEFinalizePackage));
12733d177a5cSEmil Constantinescu   PetscFunctionReturn(0);
12743d177a5cSEmil Constantinescu }
12753d177a5cSEmil Constantinescu 
12763d177a5cSEmil Constantinescu /*@C
12773d177a5cSEmil Constantinescu   TSGLLEFinalizePackage - This function destroys everything in the TSGLLE package. It is
12783d177a5cSEmil Constantinescu   called from PetscFinalize().
12793d177a5cSEmil Constantinescu 
12803d177a5cSEmil Constantinescu   Level: developer
12813d177a5cSEmil Constantinescu 
12823d177a5cSEmil Constantinescu .seealso: PetscFinalize()
12833d177a5cSEmil Constantinescu @*/
12843d177a5cSEmil Constantinescu PetscErrorCode  TSGLLEFinalizePackage(void)
12853d177a5cSEmil Constantinescu {
12863d177a5cSEmil Constantinescu   PetscFunctionBegin;
1287*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFunctionListDestroy(&TSGLLEList));
1288*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFunctionListDestroy(&TSGLLEAcceptList));
12893d177a5cSEmil Constantinescu   TSGLLEPackageInitialized = PETSC_FALSE;
12903d177a5cSEmil Constantinescu   TSGLLERegisterAllCalled  = PETSC_FALSE;
12913d177a5cSEmil Constantinescu   PetscFunctionReturn(0);
12923d177a5cSEmil Constantinescu }
12933d177a5cSEmil Constantinescu 
12943d177a5cSEmil Constantinescu /* ------------------------------------------------------------ */
12953d177a5cSEmil Constantinescu /*MC
12963d177a5cSEmil Constantinescu       TSGLLE - DAE solver using implicit General Linear methods
12973d177a5cSEmil Constantinescu 
12983d177a5cSEmil Constantinescu   These methods contain Runge-Kutta and multistep schemes as special cases.  These special cases have some fundamental
12993d177a5cSEmil Constantinescu   limitations.  For example, diagonally implicit Runge-Kutta cannot have stage order greater than 1 which limits their
13003d177a5cSEmil Constantinescu   applicability to very stiff systems.  Meanwhile, multistep methods cannot be A-stable for order greater than 2 and BDF
13013d177a5cSEmil Constantinescu   are not 0-stable for order greater than 6.  GL methods can be A- and L-stable with arbitrarily high stage order and
13023d177a5cSEmil Constantinescu   reliable error estimates for both 1 and 2 orders higher to facilitate adaptive step sizes and adaptive order schemes.
13033d177a5cSEmil Constantinescu   All this is possible while preserving a singly diagonally implicit structure.
13043d177a5cSEmil Constantinescu 
13053d177a5cSEmil Constantinescu   Options database keys:
13063d177a5cSEmil Constantinescu +  -ts_gl_type <type> - the class of general linear method (irks)
13073d177a5cSEmil Constantinescu .  -ts_gl_rtol <tol>  - relative error
13083d177a5cSEmil Constantinescu .  -ts_gl_atol <tol>  - absolute error
13093d177a5cSEmil Constantinescu .  -ts_gl_min_order <p> - minimum order method to consider (default=1)
13103d177a5cSEmil Constantinescu .  -ts_gl_max_order <p> - maximum order method to consider (default=3)
13113d177a5cSEmil Constantinescu .  -ts_gl_start_order <p> - order of starting method (default=1)
13123d177a5cSEmil Constantinescu .  -ts_gl_complete <method> - method to use for completing the step (rescale-and-modify or rescale)
13133d177a5cSEmil Constantinescu -  -ts_adapt_type <method> - adaptive controller to use (none step both)
13143d177a5cSEmil Constantinescu 
13153d177a5cSEmil Constantinescu   Notes:
13163d177a5cSEmil Constantinescu   This integrator can be applied to DAE.
13173d177a5cSEmil Constantinescu 
13183d177a5cSEmil Constantinescu   Diagonally implicit general linear (DIGL) methods are a generalization of diagonally implicit Runge-Kutta (DIRK).
13193d177a5cSEmil Constantinescu   They are represented by the tableau
13203d177a5cSEmil Constantinescu 
13213d177a5cSEmil Constantinescu .vb
13223d177a5cSEmil Constantinescu   A  |  U
13233d177a5cSEmil Constantinescu   -------
13243d177a5cSEmil Constantinescu   B  |  V
13253d177a5cSEmil Constantinescu .ve
13263d177a5cSEmil Constantinescu 
13273d177a5cSEmil Constantinescu   combined with a vector c of abscissa.  "Diagonally implicit" means that A is lower triangular.
13283d177a5cSEmil Constantinescu   A step of the general method reads
13293d177a5cSEmil Constantinescu 
13303d177a5cSEmil Constantinescu .vb
13313d177a5cSEmil Constantinescu   [ Y ] = [A  U] [  Y'   ]
13323d177a5cSEmil Constantinescu   [X^k] = [B  V] [X^{k-1}]
13333d177a5cSEmil Constantinescu .ve
13343d177a5cSEmil Constantinescu 
13353d177a5cSEmil Constantinescu   where Y is the multivector of stage values, Y' is the multivector of stage derivatives, X^k is the Nordsieck vector of
13363d177a5cSEmil Constantinescu   the solution at step k.  The Nordsieck vector consists of the first r moments of the solution, given by
13373d177a5cSEmil Constantinescu 
13383d177a5cSEmil Constantinescu .vb
13393d177a5cSEmil Constantinescu   X = [x_0,x_1,...,x_{r-1}] = [x, h x', h^2 x'', ..., h^{r-1} x^{(r-1)} ]
13403d177a5cSEmil Constantinescu .ve
13413d177a5cSEmil Constantinescu 
13423d177a5cSEmil Constantinescu   If A is lower triangular, we can solve the stages (Y,Y') sequentially
13433d177a5cSEmil Constantinescu 
13443d177a5cSEmil Constantinescu .vb
13453d177a5cSEmil Constantinescu   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}
13463d177a5cSEmil Constantinescu .ve
13473d177a5cSEmil Constantinescu 
13483d177a5cSEmil Constantinescu   and then construct the pieces to carry to the next step
13493d177a5cSEmil Constantinescu 
13503d177a5cSEmil Constantinescu .vb
13513d177a5cSEmil Constantinescu   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}
13523d177a5cSEmil Constantinescu .ve
13533d177a5cSEmil Constantinescu 
13543d177a5cSEmil Constantinescu   Note that when the equations are cast in implicit form, we are using the stage equation to define y'_i
13553d177a5cSEmil Constantinescu   in terms of y_i and known stuff (y_j for j<i and x_j for all j).
13563d177a5cSEmil Constantinescu 
13573d177a5cSEmil Constantinescu   Error estimation
13583d177a5cSEmil Constantinescu 
13593d177a5cSEmil Constantinescu   At present, the most attractive GL methods for stiff problems are singly diagonally implicit schemes which posses
13603d177a5cSEmil Constantinescu   Inherent Runge-Kutta Stability (IRKS).  These methods have r=s, the number of items passed between steps is equal to
13613d177a5cSEmil Constantinescu   the number of stages.  The order and stage-order are one less than the number of stages.  We use the error estimates
13623d177a5cSEmil Constantinescu   in the 2007 paper which provide the following estimates
13633d177a5cSEmil Constantinescu 
13643d177a5cSEmil Constantinescu .vb
13653d177a5cSEmil Constantinescu   h^{p+1} X^{(p+1)}          = phi_0^T Y' + [0 psi_0^T] Xold
13663d177a5cSEmil Constantinescu   h^{p+2} X^{(p+2)}          = phi_1^T Y' + [0 psi_1^T] Xold
13673d177a5cSEmil Constantinescu   h^{p+2} (dx'/dx) X^{(p+1)} = phi_2^T Y' + [0 psi_2^T] Xold
13683d177a5cSEmil Constantinescu .ve
13693d177a5cSEmil Constantinescu 
13703d177a5cSEmil Constantinescu   These estimates are accurate to O(h^{p+3}).
13713d177a5cSEmil Constantinescu 
13723d177a5cSEmil Constantinescu   Changing the step size
13733d177a5cSEmil Constantinescu 
13743d177a5cSEmil Constantinescu   We use the generalized "rescale and modify" scheme, see equation (4.5) of the 2007 paper.
13753d177a5cSEmil Constantinescu 
13763d177a5cSEmil Constantinescu   Level: beginner
13773d177a5cSEmil Constantinescu 
13783d177a5cSEmil Constantinescu   References:
1379606c0280SSatish Balay + * - John Butcher and Z. Jackieweicz and W. Wright, On error propagation in general linear methods for
13803d177a5cSEmil Constantinescu   ordinary differential equations, Journal of Complexity, Vol 23, 2007.
1381606c0280SSatish Balay - * - John Butcher, Numerical methods for ordinary differential equations, second edition, Wiley, 2009.
13823d177a5cSEmil Constantinescu 
13833d177a5cSEmil Constantinescu .seealso:  TSCreate(), TS, TSSetType()
13843d177a5cSEmil Constantinescu 
13853d177a5cSEmil Constantinescu M*/
13863d177a5cSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSCreate_GLLE(TS ts)
13873d177a5cSEmil Constantinescu {
13883d177a5cSEmil Constantinescu   TS_GLLE        *gl;
13893d177a5cSEmil Constantinescu 
13903d177a5cSEmil Constantinescu   PetscFunctionBegin;
1391*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGLLEInitializePackage());
13923d177a5cSEmil Constantinescu 
1393*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscNewLog(ts,&gl));
13943d177a5cSEmil Constantinescu   ts->data = (void*)gl;
13953d177a5cSEmil Constantinescu 
13963d177a5cSEmil Constantinescu   ts->ops->reset          = TSReset_GLLE;
13973d177a5cSEmil Constantinescu   ts->ops->destroy        = TSDestroy_GLLE;
13983d177a5cSEmil Constantinescu   ts->ops->view           = TSView_GLLE;
13993d177a5cSEmil Constantinescu   ts->ops->setup          = TSSetUp_GLLE;
14003d177a5cSEmil Constantinescu   ts->ops->solve          = TSSolve_GLLE;
14013d177a5cSEmil Constantinescu   ts->ops->setfromoptions = TSSetFromOptions_GLLE;
14023d177a5cSEmil Constantinescu   ts->ops->snesfunction   = SNESTSFormFunction_GLLE;
14033d177a5cSEmil Constantinescu   ts->ops->snesjacobian   = SNESTSFormJacobian_GLLE;
14043d177a5cSEmil Constantinescu 
1405efd4aadfSBarry Smith   ts->usessnes = PETSC_TRUE;
1406efd4aadfSBarry Smith 
14073d177a5cSEmil Constantinescu   gl->max_step_rejections = 1;
14083d177a5cSEmil Constantinescu   gl->min_order           = 1;
14093d177a5cSEmil Constantinescu   gl->max_order           = 3;
14103d177a5cSEmil Constantinescu   gl->start_order         = 1;
14113d177a5cSEmil Constantinescu   gl->current_scheme      = -1;
14123d177a5cSEmil Constantinescu   gl->extrapolate         = PETSC_FALSE;
14133d177a5cSEmil Constantinescu 
14143d177a5cSEmil Constantinescu   gl->wrms_atol = 1e-8;
14153d177a5cSEmil Constantinescu   gl->wrms_rtol = 1e-5;
14163d177a5cSEmil Constantinescu 
1417*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectComposeFunction((PetscObject)ts,"TSGLLESetType_C",      &TSGLLESetType_GLLE));
1418*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectComposeFunction((PetscObject)ts,"TSGLLESetAcceptType_C",&TSGLLESetAcceptType_GLLE));
1419*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectComposeFunction((PetscObject)ts,"TSGLLEGetAdapt_C",     &TSGLLEGetAdapt_GLLE));
14203d177a5cSEmil Constantinescu   PetscFunctionReturn(0);
14213d177a5cSEmil Constantinescu }
1422