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 63d177a5cSEmil Constantinescu static const char *TSGLLEErrorDirections[] = {"FORWARD","BACKWARD","TSGLLEErrorDirection","TSGLLEERROR_",0}; 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 #undef __FUNCT__ 343d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLEGetVecs" 353d177a5cSEmil Constantinescu static PetscErrorCode TSGLLEGetVecs(TS ts,DM dm,Vec *Z,Vec *Ydotstage) 363d177a5cSEmil Constantinescu { 373d177a5cSEmil Constantinescu TS_GLLE *gl = (TS_GLLE*)ts->data; 383d177a5cSEmil Constantinescu PetscErrorCode ierr; 393d177a5cSEmil Constantinescu 403d177a5cSEmil Constantinescu PetscFunctionBegin; 413d177a5cSEmil Constantinescu if (Z) { 423d177a5cSEmil Constantinescu if (dm && dm != ts->dm) { 433d177a5cSEmil Constantinescu ierr = DMGetNamedGlobalVector(dm,"TSGLLE_Z",Z);CHKERRQ(ierr); 443d177a5cSEmil Constantinescu } else *Z = gl->Z; 453d177a5cSEmil Constantinescu } 463d177a5cSEmil Constantinescu if (Ydotstage) { 473d177a5cSEmil Constantinescu if (dm && dm != ts->dm) { 483d177a5cSEmil Constantinescu ierr = DMGetNamedGlobalVector(dm,"TSGLLE_Ydot",Ydotstage);CHKERRQ(ierr); 493d177a5cSEmil Constantinescu } else *Ydotstage = gl->Ydot[gl->stage]; 503d177a5cSEmil Constantinescu } 513d177a5cSEmil Constantinescu PetscFunctionReturn(0); 523d177a5cSEmil Constantinescu } 533d177a5cSEmil Constantinescu 543d177a5cSEmil Constantinescu 553d177a5cSEmil Constantinescu #undef __FUNCT__ 563d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLERestoreVecs" 573d177a5cSEmil Constantinescu static PetscErrorCode TSGLLERestoreVecs(TS ts,DM dm,Vec *Z,Vec *Ydotstage) 583d177a5cSEmil Constantinescu { 593d177a5cSEmil Constantinescu PetscErrorCode ierr; 603d177a5cSEmil Constantinescu 613d177a5cSEmil Constantinescu PetscFunctionBegin; 623d177a5cSEmil Constantinescu if (Z) { 633d177a5cSEmil Constantinescu if (dm && dm != ts->dm) { 643d177a5cSEmil Constantinescu ierr = DMRestoreNamedGlobalVector(dm,"TSGLLE_Z",Z);CHKERRQ(ierr); 653d177a5cSEmil Constantinescu } 663d177a5cSEmil Constantinescu } 673d177a5cSEmil Constantinescu if (Ydotstage) { 683d177a5cSEmil Constantinescu 693d177a5cSEmil Constantinescu if (dm && dm != ts->dm) { 703d177a5cSEmil Constantinescu ierr = DMRestoreNamedGlobalVector(dm,"TSGLLE_Ydot",Ydotstage);CHKERRQ(ierr); 713d177a5cSEmil Constantinescu } 723d177a5cSEmil Constantinescu } 733d177a5cSEmil Constantinescu PetscFunctionReturn(0); 743d177a5cSEmil Constantinescu } 753d177a5cSEmil Constantinescu 763d177a5cSEmil Constantinescu #undef __FUNCT__ 773d177a5cSEmil Constantinescu #define __FUNCT__ "DMCoarsenHook_TSGLLE" 783d177a5cSEmil Constantinescu static PetscErrorCode DMCoarsenHook_TSGLLE(DM fine,DM coarse,void *ctx) 793d177a5cSEmil Constantinescu { 803d177a5cSEmil Constantinescu PetscFunctionBegin; 813d177a5cSEmil Constantinescu PetscFunctionReturn(0); 823d177a5cSEmil Constantinescu } 833d177a5cSEmil Constantinescu 843d177a5cSEmil Constantinescu #undef __FUNCT__ 853d177a5cSEmil Constantinescu #define __FUNCT__ "DMRestrictHook_TSGLLE" 863d177a5cSEmil Constantinescu static PetscErrorCode DMRestrictHook_TSGLLE(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx) 873d177a5cSEmil Constantinescu { 883d177a5cSEmil Constantinescu TS ts = (TS)ctx; 893d177a5cSEmil Constantinescu PetscErrorCode ierr; 903d177a5cSEmil Constantinescu Vec Ydot,Ydot_c; 913d177a5cSEmil Constantinescu 923d177a5cSEmil Constantinescu PetscFunctionBegin; 933d177a5cSEmil Constantinescu ierr = TSGLLEGetVecs(ts,fine,NULL,&Ydot);CHKERRQ(ierr); 943d177a5cSEmil Constantinescu ierr = TSGLLEGetVecs(ts,coarse,NULL,&Ydot_c);CHKERRQ(ierr); 953d177a5cSEmil Constantinescu ierr = MatRestrict(restrct,Ydot,Ydot_c);CHKERRQ(ierr); 963d177a5cSEmil Constantinescu ierr = VecPointwiseMult(Ydot_c,rscale,Ydot_c);CHKERRQ(ierr); 973d177a5cSEmil Constantinescu ierr = TSGLLERestoreVecs(ts,fine,NULL,&Ydot);CHKERRQ(ierr); 983d177a5cSEmil Constantinescu ierr = TSGLLERestoreVecs(ts,coarse,NULL,&Ydot_c);CHKERRQ(ierr); 993d177a5cSEmil Constantinescu PetscFunctionReturn(0); 1003d177a5cSEmil Constantinescu } 1013d177a5cSEmil Constantinescu 1023d177a5cSEmil Constantinescu #undef __FUNCT__ 1033d177a5cSEmil Constantinescu #define __FUNCT__ "DMSubDomainHook_TSGLLE" 1043d177a5cSEmil Constantinescu static PetscErrorCode DMSubDomainHook_TSGLLE(DM dm,DM subdm,void *ctx) 1053d177a5cSEmil Constantinescu { 1063d177a5cSEmil Constantinescu PetscFunctionBegin; 1073d177a5cSEmil Constantinescu PetscFunctionReturn(0); 1083d177a5cSEmil Constantinescu } 1093d177a5cSEmil Constantinescu 1103d177a5cSEmil Constantinescu #undef __FUNCT__ 1113d177a5cSEmil Constantinescu #define __FUNCT__ "DMSubDomainRestrictHook_TSGLLE" 1123d177a5cSEmil Constantinescu static PetscErrorCode DMSubDomainRestrictHook_TSGLLE(DM dm,VecScatter gscat, VecScatter lscat,DM subdm,void *ctx) 1133d177a5cSEmil Constantinescu { 1143d177a5cSEmil Constantinescu TS ts = (TS)ctx; 1153d177a5cSEmil Constantinescu PetscErrorCode ierr; 1163d177a5cSEmil Constantinescu Vec Ydot,Ydot_s; 1173d177a5cSEmil Constantinescu 1183d177a5cSEmil Constantinescu PetscFunctionBegin; 1193d177a5cSEmil Constantinescu ierr = TSGLLEGetVecs(ts,dm,NULL,&Ydot);CHKERRQ(ierr); 1203d177a5cSEmil Constantinescu ierr = TSGLLEGetVecs(ts,subdm,NULL,&Ydot_s);CHKERRQ(ierr); 1213d177a5cSEmil Constantinescu 1223d177a5cSEmil Constantinescu ierr = VecScatterBegin(gscat,Ydot,Ydot_s,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1233d177a5cSEmil Constantinescu ierr = VecScatterEnd(gscat,Ydot,Ydot_s,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1243d177a5cSEmil Constantinescu 1253d177a5cSEmil Constantinescu ierr = TSGLLERestoreVecs(ts,dm,NULL,&Ydot);CHKERRQ(ierr); 1263d177a5cSEmil Constantinescu ierr = TSGLLERestoreVecs(ts,subdm,NULL,&Ydot_s);CHKERRQ(ierr); 1273d177a5cSEmil Constantinescu PetscFunctionReturn(0); 1283d177a5cSEmil Constantinescu } 1293d177a5cSEmil Constantinescu 1303d177a5cSEmil Constantinescu #undef __FUNCT__ 1313d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLESchemeCreate" 1323d177a5cSEmil Constantinescu static PetscErrorCode TSGLLESchemeCreate(PetscInt p,PetscInt q,PetscInt r,PetscInt s,const PetscScalar *c, 1333d177a5cSEmil Constantinescu const PetscScalar *a,const PetscScalar *b,const PetscScalar *u,const PetscScalar *v,TSGLLEScheme *inscheme) 1343d177a5cSEmil Constantinescu { 1353d177a5cSEmil Constantinescu TSGLLEScheme scheme; 1363d177a5cSEmil Constantinescu PetscInt j; 1373d177a5cSEmil Constantinescu PetscErrorCode ierr; 1383d177a5cSEmil Constantinescu 1393d177a5cSEmil Constantinescu PetscFunctionBegin; 1403d177a5cSEmil Constantinescu if (p < 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Scheme order must be positive"); 1413d177a5cSEmil Constantinescu if (r < 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"At least one item must be carried between steps"); 1423d177a5cSEmil Constantinescu if (s < 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"At least one stage is required"); 1433d177a5cSEmil Constantinescu PetscValidPointer(inscheme,4); 1443d177a5cSEmil Constantinescu *inscheme = 0; 1453d177a5cSEmil Constantinescu ierr = PetscNew(&scheme);CHKERRQ(ierr); 1463d177a5cSEmil Constantinescu scheme->p = p; 1473d177a5cSEmil Constantinescu scheme->q = q; 1483d177a5cSEmil Constantinescu scheme->r = r; 1493d177a5cSEmil Constantinescu scheme->s = s; 1503d177a5cSEmil Constantinescu 1513d177a5cSEmil Constantinescu ierr = PetscMalloc5(s,&scheme->c,s*s,&scheme->a,r*s,&scheme->b,r*s,&scheme->u,r*r,&scheme->v);CHKERRQ(ierr); 1523d177a5cSEmil Constantinescu ierr = PetscMemcpy(scheme->c,c,s*sizeof(PetscScalar));CHKERRQ(ierr); 1533d177a5cSEmil Constantinescu for (j=0; j<s*s; j++) scheme->a[j] = (PetscAbsScalar(a[j]) < 1e-12) ? 0 : a[j]; 1543d177a5cSEmil Constantinescu for (j=0; j<r*s; j++) scheme->b[j] = (PetscAbsScalar(b[j]) < 1e-12) ? 0 : b[j]; 1553d177a5cSEmil Constantinescu for (j=0; j<s*r; j++) scheme->u[j] = (PetscAbsScalar(u[j]) < 1e-12) ? 0 : u[j]; 1563d177a5cSEmil Constantinescu for (j=0; j<r*r; j++) scheme->v[j] = (PetscAbsScalar(v[j]) < 1e-12) ? 0 : v[j]; 1573d177a5cSEmil Constantinescu 1583d177a5cSEmil Constantinescu ierr = PetscMalloc6(r,&scheme->alpha,r,&scheme->beta,r,&scheme->gamma,3*s,&scheme->phi,3*r,&scheme->psi,r,&scheme->stage_error);CHKERRQ(ierr); 1593d177a5cSEmil Constantinescu { 1603d177a5cSEmil Constantinescu PetscInt i,j,k,ss=s+2; 1613d177a5cSEmil Constantinescu PetscBLASInt m,n,one=1,*ipiv,lwork=4*((s+3)*3+3),info,ldb; 1623d177a5cSEmil Constantinescu PetscReal rcond,*sing,*workreal; 1633d177a5cSEmil Constantinescu PetscScalar *ImV,*H,*bmat,*workscalar,*c=scheme->c,*a=scheme->a,*b=scheme->b,*u=scheme->u,*v=scheme->v; 1643d177a5cSEmil Constantinescu #if !defined(PETSC_MISSING_LAPACK_GELSS) 1653d177a5cSEmil Constantinescu PetscBLASInt rank; 1663d177a5cSEmil Constantinescu #endif 1673d177a5cSEmil Constantinescu ierr = PetscMalloc7(PetscSqr(r),&ImV,3*s,&H,3*ss,&bmat,lwork,&workscalar,5*(3+r),&workreal,r+s,&sing,r+s,&ipiv);CHKERRQ(ierr); 1683d177a5cSEmil Constantinescu 1693d177a5cSEmil Constantinescu /* column-major input */ 1703d177a5cSEmil Constantinescu for (i=0; i<r-1; i++) { 1713d177a5cSEmil Constantinescu for (j=0; j<r-1; j++) ImV[i+j*r] = 1.0*(i==j) - v[(i+1)*r+j+1]; 1723d177a5cSEmil Constantinescu } 1733d177a5cSEmil Constantinescu /* Build right hand side for alpha (tp - glm.B(2:end,:)*(glm.c.^(p)./factorial(p))) */ 1743d177a5cSEmil Constantinescu for (i=1; i<r; i++) { 1753d177a5cSEmil Constantinescu scheme->alpha[i] = 1./Factorial(p+1-i); 1763d177a5cSEmil Constantinescu for (j=0; j<s; j++) scheme->alpha[i] -= b[i*s+j]*CPowF(c[j],p); 1773d177a5cSEmil Constantinescu } 1783d177a5cSEmil Constantinescu ierr = PetscBLASIntCast(r-1,&m);CHKERRQ(ierr); 1793d177a5cSEmil Constantinescu ierr = PetscBLASIntCast(r,&n);CHKERRQ(ierr); 1803d177a5cSEmil Constantinescu PetscStackCallBLAS("LAPACKgesv",LAPACKgesv_(&m,&one,ImV,&n,ipiv,scheme->alpha+1,&n,&info)); 1813d177a5cSEmil Constantinescu if (info < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad argument to GESV"); 1823d177a5cSEmil Constantinescu if (info > 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Bad LU factorization"); 1833d177a5cSEmil Constantinescu 1843d177a5cSEmil Constantinescu /* Build right hand side for beta (tp1 - glm.B(2:end,:)*(glm.c.^(p+1)./factorial(p+1)) - e.alpha) */ 1853d177a5cSEmil Constantinescu for (i=1; i<r; i++) { 1863d177a5cSEmil Constantinescu scheme->beta[i] = 1./Factorial(p+2-i) - scheme->alpha[i]; 1873d177a5cSEmil Constantinescu for (j=0; j<s; j++) scheme->beta[i] -= b[i*s+j]*CPowF(c[j],p+1); 1883d177a5cSEmil Constantinescu } 1893d177a5cSEmil Constantinescu PetscStackCallBLAS("LAPACKgetrs",LAPACKgetrs_("No transpose",&m,&one,ImV,&n,ipiv,scheme->beta+1,&n,&info)); 1903d177a5cSEmil Constantinescu if (info < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad argument to GETRS"); 1913d177a5cSEmil Constantinescu if (info > 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Should not happen"); 1923d177a5cSEmil Constantinescu 1933d177a5cSEmil Constantinescu /* Build stage_error vector 1943d177a5cSEmil Constantinescu xi = glm.c.^(p+1)/factorial(p+1) - glm.A*glm.c.^p/factorial(p) + glm.U(:,2:end)*e.alpha; 1953d177a5cSEmil Constantinescu */ 1963d177a5cSEmil Constantinescu for (i=0; i<s; i++) { 1973d177a5cSEmil Constantinescu scheme->stage_error[i] = CPowF(c[i],p+1); 1983d177a5cSEmil Constantinescu for (j=0; j<s; j++) scheme->stage_error[i] -= a[i*s+j]*CPowF(c[j],p); 1993d177a5cSEmil Constantinescu for (j=1; j<r; j++) scheme->stage_error[i] += u[i*r+j]*scheme->alpha[j]; 2003d177a5cSEmil Constantinescu } 2013d177a5cSEmil Constantinescu 2023d177a5cSEmil Constantinescu /* alpha[0] (epsilon in B,J,W 2007) 2033d177a5cSEmil Constantinescu epsilon = 1/factorial(p+1) - B(1,:)*c.^p/factorial(p) + V(1,2:end)*e.alpha; 2043d177a5cSEmil Constantinescu */ 2053d177a5cSEmil Constantinescu scheme->alpha[0] = 1./Factorial(p+1); 2063d177a5cSEmil Constantinescu for (j=0; j<s; j++) scheme->alpha[0] -= b[0*s+j]*CPowF(c[j],p); 2073d177a5cSEmil Constantinescu for (j=1; j<r; j++) scheme->alpha[0] += v[0*r+j]*scheme->alpha[j]; 2083d177a5cSEmil Constantinescu 2093d177a5cSEmil Constantinescu /* right hand side for gamma (glm.B(2:end,:)*e.xi - e.epsilon*eye(s-1,1)) */ 2103d177a5cSEmil Constantinescu for (i=1; i<r; i++) { 2113d177a5cSEmil Constantinescu scheme->gamma[i] = (i==1 ? -1. : 0)*scheme->alpha[0]; 2123d177a5cSEmil Constantinescu for (j=0; j<s; j++) scheme->gamma[i] += b[i*s+j]*scheme->stage_error[j]; 2133d177a5cSEmil Constantinescu } 2143d177a5cSEmil Constantinescu PetscStackCallBLAS("LAPACKgetrs",LAPACKgetrs_("No transpose",&m,&one,ImV,&n,ipiv,scheme->gamma+1,&n,&info)); 2153d177a5cSEmil Constantinescu if (info < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad argument to GETRS"); 2163d177a5cSEmil Constantinescu if (info > 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Should not happen"); 2173d177a5cSEmil Constantinescu 2183d177a5cSEmil Constantinescu /* beta[0] (rho in B,J,W 2007) 2193d177a5cSEmil Constantinescu e.rho = 1/factorial(p+2) - glm.B(1,:)*glm.c.^(p+1)/factorial(p+1) ... 2203d177a5cSEmil Constantinescu + glm.V(1,2:end)*e.beta;% - e.epsilon; 2213d177a5cSEmil Constantinescu % Note: The paper (B,J,W 2007) includes the last term in their definition 2223d177a5cSEmil Constantinescu * */ 2233d177a5cSEmil Constantinescu scheme->beta[0] = 1./Factorial(p+2); 2243d177a5cSEmil Constantinescu for (j=0; j<s; j++) scheme->beta[0] -= b[0*s+j]*CPowF(c[j],p+1); 2253d177a5cSEmil Constantinescu for (j=1; j<r; j++) scheme->beta[0] += v[0*r+j]*scheme->beta[j]; 2263d177a5cSEmil Constantinescu 2273d177a5cSEmil Constantinescu /* gamma[0] (sigma in B,J,W 2007) 2283d177a5cSEmil Constantinescu * e.sigma = glm.B(1,:)*e.xi + glm.V(1,2:end)*e.gamma; 2293d177a5cSEmil Constantinescu * */ 2303d177a5cSEmil Constantinescu scheme->gamma[0] = 0.0; 2313d177a5cSEmil Constantinescu for (j=0; j<s; j++) scheme->gamma[0] += b[0*s+j]*scheme->stage_error[j]; 2323d177a5cSEmil Constantinescu for (j=1; j<r; j++) scheme->gamma[0] += v[0*s+j]*scheme->gamma[j]; 2333d177a5cSEmil Constantinescu 2343d177a5cSEmil Constantinescu /* Assemble H 2353d177a5cSEmil Constantinescu * % Determine the error estimators phi 2363d177a5cSEmil Constantinescu H = [[cpow(glm.c,p) + C*e.alpha] [cpow(glm.c,p+1) + C*e.beta] ... 2373d177a5cSEmil Constantinescu [e.xi - C*(e.gamma + 0*e.epsilon*eye(s-1,1))]]'; 2383d177a5cSEmil Constantinescu % Paper has formula above without the 0, but that term must be left 2393d177a5cSEmil Constantinescu % out to satisfy the conditions they propose and to make the 2403d177a5cSEmil Constantinescu % example schemes work 2413d177a5cSEmil Constantinescu e.H = H; 2423d177a5cSEmil Constantinescu e.phi = (H \ [1 0 0;1 1 0;0 0 -1])'; 2433d177a5cSEmil Constantinescu e.psi = -e.phi*C; 2443d177a5cSEmil Constantinescu * */ 2453d177a5cSEmil Constantinescu for (j=0; j<s; j++) { 2463d177a5cSEmil Constantinescu H[0+j*3] = CPowF(c[j],p); 2473d177a5cSEmil Constantinescu H[1+j*3] = CPowF(c[j],p+1); 2483d177a5cSEmil Constantinescu H[2+j*3] = scheme->stage_error[j]; 2493d177a5cSEmil Constantinescu for (k=1; k<r; k++) { 2503d177a5cSEmil Constantinescu H[0+j*3] += CPowF(c[j],k-1)*scheme->alpha[k]; 2513d177a5cSEmil Constantinescu H[1+j*3] += CPowF(c[j],k-1)*scheme->beta[k]; 2523d177a5cSEmil Constantinescu H[2+j*3] -= CPowF(c[j],k-1)*scheme->gamma[k]; 2533d177a5cSEmil Constantinescu } 2543d177a5cSEmil Constantinescu } 2553d177a5cSEmil Constantinescu bmat[0+0*ss] = 1.; bmat[0+1*ss] = 0.; bmat[0+2*ss] = 0.; 2563d177a5cSEmil Constantinescu bmat[1+0*ss] = 1.; bmat[1+1*ss] = 1.; bmat[1+2*ss] = 0.; 2573d177a5cSEmil Constantinescu bmat[2+0*ss] = 0.; bmat[2+1*ss] = 0.; bmat[2+2*ss] = -1.; 2583d177a5cSEmil Constantinescu m = 3; 2593d177a5cSEmil Constantinescu ierr = PetscBLASIntCast(s,&n);CHKERRQ(ierr); 2603d177a5cSEmil Constantinescu ierr = PetscBLASIntCast(ss,&ldb);CHKERRQ(ierr); 2613d177a5cSEmil Constantinescu rcond = 1e-12; 2623d177a5cSEmil Constantinescu #if defined(PETSC_MISSING_LAPACK_GELSS) 2633d177a5cSEmil Constantinescu /* ESSL does not have this routine */ 2643d177a5cSEmil Constantinescu SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GELSS - Lapack routine is unavailable\nNot able to run GL time stepping."); 2653d177a5cSEmil Constantinescu #else 2663d177a5cSEmil Constantinescu #if defined(PETSC_USE_COMPLEX) 2673d177a5cSEmil Constantinescu /* ZGELSS( M, N, NRHS, A, LDA, B, LDB, S, RCOND, RANK, WORK, LWORK, RWORK, INFO) */ 2683d177a5cSEmil Constantinescu PetscStackCallBLAS("LAPACKgelss",LAPACKgelss_(&m,&n,&m,H,&m,bmat,&ldb,sing,&rcond,&rank,workscalar,&lwork,workreal,&info)); 2693d177a5cSEmil Constantinescu #else 2703d177a5cSEmil Constantinescu /* DGELSS( M, N, NRHS, A, LDA, B, LDB, S, RCOND, RANK, WORK, LWORK, INFO) */ 2713d177a5cSEmil Constantinescu PetscStackCallBLAS("LAPACKgelss",LAPACKgelss_(&m,&n,&m,H,&m,bmat,&ldb,sing,&rcond,&rank,workscalar,&lwork,&info)); 2723d177a5cSEmil Constantinescu #endif 2733d177a5cSEmil Constantinescu if (info < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad argument to GELSS"); 2743d177a5cSEmil Constantinescu if (info > 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"SVD failed to converge"); 2753d177a5cSEmil Constantinescu #endif 2763d177a5cSEmil Constantinescu 2773d177a5cSEmil Constantinescu for (j=0; j<3; j++) { 2783d177a5cSEmil Constantinescu for (k=0; k<s; k++) scheme->phi[k+j*s] = bmat[k+j*ss]; 2793d177a5cSEmil Constantinescu } 2803d177a5cSEmil Constantinescu 2813d177a5cSEmil Constantinescu /* the other part of the error estimator, psi in B,J,W 2007 */ 2823d177a5cSEmil Constantinescu scheme->psi[0*r+0] = 0.; 2833d177a5cSEmil Constantinescu scheme->psi[1*r+0] = 0.; 2843d177a5cSEmil Constantinescu scheme->psi[2*r+0] = 0.; 2853d177a5cSEmil Constantinescu for (j=1; j<r; j++) { 2863d177a5cSEmil Constantinescu scheme->psi[0*r+j] = 0.; 2873d177a5cSEmil Constantinescu scheme->psi[1*r+j] = 0.; 2883d177a5cSEmil Constantinescu scheme->psi[2*r+j] = 0.; 2893d177a5cSEmil Constantinescu for (k=0; k<s; k++) { 2903d177a5cSEmil Constantinescu scheme->psi[0*r+j] -= CPowF(c[k],j-1)*scheme->phi[0*s+k]; 2913d177a5cSEmil Constantinescu scheme->psi[1*r+j] -= CPowF(c[k],j-1)*scheme->phi[1*s+k]; 2923d177a5cSEmil Constantinescu scheme->psi[2*r+j] -= CPowF(c[k],j-1)*scheme->phi[2*s+k]; 2933d177a5cSEmil Constantinescu } 2943d177a5cSEmil Constantinescu } 2953d177a5cSEmil Constantinescu ierr = PetscFree7(ImV,H,bmat,workscalar,workreal,sing,ipiv);CHKERRQ(ierr); 2963d177a5cSEmil Constantinescu } 2973d177a5cSEmil Constantinescu /* Check which properties are satisfied */ 2983d177a5cSEmil Constantinescu scheme->stiffly_accurate = PETSC_TRUE; 2993d177a5cSEmil Constantinescu if (scheme->c[s-1] != 1.) scheme->stiffly_accurate = PETSC_FALSE; 3003d177a5cSEmil Constantinescu for (j=0; j<s; j++) if (a[(s-1)*s+j] != b[j]) scheme->stiffly_accurate = PETSC_FALSE; 3013d177a5cSEmil Constantinescu for (j=0; j<r; j++) if (u[(s-1)*r+j] != v[j]) scheme->stiffly_accurate = PETSC_FALSE; 3023d177a5cSEmil Constantinescu scheme->fsal = scheme->stiffly_accurate; /* FSAL is stronger */ 3033d177a5cSEmil Constantinescu for (j=0; j<s-1; j++) if (r>1 && b[1*s+j] != 0.) scheme->fsal = PETSC_FALSE; 3043d177a5cSEmil Constantinescu if (b[1*s+r-1] != 1.) scheme->fsal = PETSC_FALSE; 3053d177a5cSEmil Constantinescu for (j=0; j<r; j++) if (r>1 && v[1*r+j] != 0.) scheme->fsal = PETSC_FALSE; 3063d177a5cSEmil Constantinescu 3073d177a5cSEmil Constantinescu *inscheme = scheme; 3083d177a5cSEmil Constantinescu PetscFunctionReturn(0); 3093d177a5cSEmil Constantinescu } 3103d177a5cSEmil Constantinescu 3113d177a5cSEmil Constantinescu #undef __FUNCT__ 3123d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLESchemeDestroy" 3133d177a5cSEmil Constantinescu static PetscErrorCode TSGLLESchemeDestroy(TSGLLEScheme sc) 3143d177a5cSEmil Constantinescu { 3153d177a5cSEmil Constantinescu PetscErrorCode ierr; 3163d177a5cSEmil Constantinescu 3173d177a5cSEmil Constantinescu PetscFunctionBegin; 3183d177a5cSEmil Constantinescu ierr = PetscFree5(sc->c,sc->a,sc->b,sc->u,sc->v);CHKERRQ(ierr); 3193d177a5cSEmil Constantinescu ierr = PetscFree6(sc->alpha,sc->beta,sc->gamma,sc->phi,sc->psi,sc->stage_error);CHKERRQ(ierr); 3203d177a5cSEmil Constantinescu ierr = PetscFree(sc);CHKERRQ(ierr); 3213d177a5cSEmil Constantinescu PetscFunctionReturn(0); 3223d177a5cSEmil Constantinescu } 3233d177a5cSEmil Constantinescu 3243d177a5cSEmil Constantinescu #undef __FUNCT__ 3253d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLEDestroy_Default" 3263d177a5cSEmil Constantinescu static PetscErrorCode TSGLLEDestroy_Default(TS_GLLE *gl) 3273d177a5cSEmil Constantinescu { 3283d177a5cSEmil Constantinescu PetscErrorCode ierr; 3293d177a5cSEmil Constantinescu PetscInt i; 3303d177a5cSEmil Constantinescu 3313d177a5cSEmil Constantinescu PetscFunctionBegin; 3323d177a5cSEmil Constantinescu for (i=0; i<gl->nschemes; i++) { 3333d177a5cSEmil Constantinescu if (gl->schemes[i]) {ierr = TSGLLESchemeDestroy(gl->schemes[i]);CHKERRQ(ierr);} 3343d177a5cSEmil Constantinescu } 3353d177a5cSEmil Constantinescu ierr = PetscFree(gl->schemes);CHKERRQ(ierr); 3363d177a5cSEmil Constantinescu gl->nschemes = 0; 3373d177a5cSEmil Constantinescu ierr = PetscMemzero(gl->type_name,sizeof(gl->type_name));CHKERRQ(ierr); 3383d177a5cSEmil Constantinescu PetscFunctionReturn(0); 3393d177a5cSEmil Constantinescu } 3403d177a5cSEmil Constantinescu 3413d177a5cSEmil Constantinescu #undef __FUNCT__ 3423d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLEViewTable_Private" 3433d177a5cSEmil Constantinescu static PetscErrorCode TSGLLEViewTable_Private(PetscViewer viewer,PetscInt m,PetscInt n,const PetscScalar a[],const char name[]) 3443d177a5cSEmil Constantinescu { 3453d177a5cSEmil Constantinescu PetscErrorCode ierr; 3463d177a5cSEmil Constantinescu PetscBool iascii; 3473d177a5cSEmil Constantinescu PetscInt i,j; 3483d177a5cSEmil Constantinescu 3493d177a5cSEmil Constantinescu PetscFunctionBegin; 3503d177a5cSEmil Constantinescu ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 3513d177a5cSEmil Constantinescu if (iascii) { 3523d177a5cSEmil Constantinescu ierr = PetscViewerASCIIPrintf(viewer,"%30s = [",name);CHKERRQ(ierr); 3533d177a5cSEmil Constantinescu for (i=0; i<m; i++) { 3543d177a5cSEmil Constantinescu if (i) {ierr = PetscViewerASCIIPrintf(viewer,"%30s [","");CHKERRQ(ierr);} 3553d177a5cSEmil Constantinescu ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 3563d177a5cSEmil Constantinescu for (j=0; j<n; j++) { 3573d177a5cSEmil Constantinescu ierr = PetscViewerASCIIPrintf(viewer," %12.8g",PetscRealPart(a[i*n+j]));CHKERRQ(ierr); 3583d177a5cSEmil Constantinescu } 3593d177a5cSEmil Constantinescu ierr = PetscViewerASCIIPrintf(viewer,"]\n");CHKERRQ(ierr); 3603d177a5cSEmil Constantinescu ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 3613d177a5cSEmil Constantinescu } 3623d177a5cSEmil Constantinescu } 3633d177a5cSEmil Constantinescu PetscFunctionReturn(0); 3643d177a5cSEmil Constantinescu } 3653d177a5cSEmil Constantinescu 3663d177a5cSEmil Constantinescu 3673d177a5cSEmil Constantinescu #undef __FUNCT__ 3683d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLESchemeView" 3693d177a5cSEmil Constantinescu static PetscErrorCode TSGLLESchemeView(TSGLLEScheme sc,PetscBool view_details,PetscViewer viewer) 3703d177a5cSEmil Constantinescu { 3713d177a5cSEmil Constantinescu PetscErrorCode ierr; 3723d177a5cSEmil Constantinescu PetscBool iascii; 3733d177a5cSEmil Constantinescu 3743d177a5cSEmil Constantinescu PetscFunctionBegin; 3753d177a5cSEmil Constantinescu ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 3763d177a5cSEmil Constantinescu if (iascii) { 3773d177a5cSEmil Constantinescu ierr = PetscViewerASCIIPrintf(viewer,"GL scheme p,q,r,s = %d,%d,%d,%d\n",sc->p,sc->q,sc->r,sc->s);CHKERRQ(ierr); 3783d177a5cSEmil Constantinescu ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 3793d177a5cSEmil Constantinescu ierr = PetscViewerASCIIPrintf(viewer,"Stiffly accurate: %s, FSAL: %s\n",sc->stiffly_accurate ? "yes" : "no",sc->fsal ? "yes" : "no");CHKERRQ(ierr); 3803d177a5cSEmil Constantinescu ierr = PetscViewerASCIIPrintf(viewer,"Leading error constants: %10.3e %10.3e %10.3e\n", 3813d177a5cSEmil Constantinescu PetscRealPart(sc->alpha[0]),PetscRealPart(sc->beta[0]),PetscRealPart(sc->gamma[0]));CHKERRQ(ierr); 3823d177a5cSEmil Constantinescu ierr = TSGLLEViewTable_Private(viewer,1,sc->s,sc->c,"Abscissas c");CHKERRQ(ierr); 3833d177a5cSEmil Constantinescu if (view_details) { 3843d177a5cSEmil Constantinescu ierr = TSGLLEViewTable_Private(viewer,sc->s,sc->s,sc->a,"A");CHKERRQ(ierr); 3853d177a5cSEmil Constantinescu ierr = TSGLLEViewTable_Private(viewer,sc->r,sc->s,sc->b,"B");CHKERRQ(ierr); 3863d177a5cSEmil Constantinescu ierr = TSGLLEViewTable_Private(viewer,sc->s,sc->r,sc->u,"U");CHKERRQ(ierr); 3873d177a5cSEmil Constantinescu ierr = TSGLLEViewTable_Private(viewer,sc->r,sc->r,sc->v,"V");CHKERRQ(ierr); 3883d177a5cSEmil Constantinescu 3893d177a5cSEmil Constantinescu ierr = TSGLLEViewTable_Private(viewer,3,sc->s,sc->phi,"Error estimate phi");CHKERRQ(ierr); 3903d177a5cSEmil Constantinescu ierr = TSGLLEViewTable_Private(viewer,3,sc->r,sc->psi,"Error estimate psi");CHKERRQ(ierr); 3913d177a5cSEmil Constantinescu ierr = TSGLLEViewTable_Private(viewer,1,sc->r,sc->alpha,"Modify alpha");CHKERRQ(ierr); 3923d177a5cSEmil Constantinescu ierr = TSGLLEViewTable_Private(viewer,1,sc->r,sc->beta,"Modify beta");CHKERRQ(ierr); 3933d177a5cSEmil Constantinescu ierr = TSGLLEViewTable_Private(viewer,1,sc->r,sc->gamma,"Modify gamma");CHKERRQ(ierr); 3943d177a5cSEmil Constantinescu ierr = TSGLLEViewTable_Private(viewer,1,sc->s,sc->stage_error,"Stage error xi");CHKERRQ(ierr); 3953d177a5cSEmil Constantinescu } 3963d177a5cSEmil Constantinescu ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 3973d177a5cSEmil Constantinescu } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported",((PetscObject)viewer)->type_name); 3983d177a5cSEmil Constantinescu PetscFunctionReturn(0); 3993d177a5cSEmil Constantinescu } 4003d177a5cSEmil Constantinescu 4013d177a5cSEmil Constantinescu #undef __FUNCT__ 4023d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLEEstimateHigherMoments_Default" 4033d177a5cSEmil Constantinescu static PetscErrorCode TSGLLEEstimateHigherMoments_Default(TSGLLEScheme sc,PetscReal h,Vec Ydot[],Vec Xold[],Vec hm[]) 4043d177a5cSEmil Constantinescu { 4053d177a5cSEmil Constantinescu PetscErrorCode ierr; 4063d177a5cSEmil Constantinescu PetscInt i; 4073d177a5cSEmil Constantinescu 4083d177a5cSEmil Constantinescu PetscFunctionBegin; 4093d177a5cSEmil Constantinescu if (sc->r > 64 || sc->s > 64) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Ridiculous number of stages or items passed between stages"); 4103d177a5cSEmil Constantinescu /* build error vectors*/ 4113d177a5cSEmil Constantinescu for (i=0; i<3; i++) { 4123d177a5cSEmil Constantinescu PetscScalar phih[64]; 4133d177a5cSEmil Constantinescu PetscInt j; 4143d177a5cSEmil Constantinescu for (j=0; j<sc->s; j++) phih[j] = sc->phi[i*sc->s+j]*h; 4153d177a5cSEmil Constantinescu ierr = VecZeroEntries(hm[i]);CHKERRQ(ierr); 4163d177a5cSEmil Constantinescu ierr = VecMAXPY(hm[i],sc->s,phih,Ydot);CHKERRQ(ierr); 4173d177a5cSEmil Constantinescu ierr = VecMAXPY(hm[i],sc->r,&sc->psi[i*sc->r],Xold);CHKERRQ(ierr); 4183d177a5cSEmil Constantinescu } 4193d177a5cSEmil Constantinescu PetscFunctionReturn(0); 4203d177a5cSEmil Constantinescu } 4213d177a5cSEmil Constantinescu 4223d177a5cSEmil Constantinescu #undef __FUNCT__ 4233d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLECompleteStep_Rescale" 4243d177a5cSEmil Constantinescu static PetscErrorCode TSGLLECompleteStep_Rescale(TSGLLEScheme sc,PetscReal h,TSGLLEScheme next_sc,PetscReal next_h,Vec Ydot[],Vec Xold[],Vec X[]) 4253d177a5cSEmil Constantinescu { 4263d177a5cSEmil Constantinescu PetscErrorCode ierr; 4273d177a5cSEmil Constantinescu PetscScalar brow[32],vrow[32]; 4283d177a5cSEmil Constantinescu PetscInt i,j,r,s; 4293d177a5cSEmil Constantinescu 4303d177a5cSEmil Constantinescu PetscFunctionBegin; 4313d177a5cSEmil Constantinescu /* Build the new solution from (X,Ydot) */ 4323d177a5cSEmil Constantinescu r = sc->r; 4333d177a5cSEmil Constantinescu s = sc->s; 4343d177a5cSEmil Constantinescu for (i=0; i<r; i++) { 4353d177a5cSEmil Constantinescu ierr = VecZeroEntries(X[i]);CHKERRQ(ierr); 4363d177a5cSEmil Constantinescu for (j=0; j<s; j++) brow[j] = h*sc->b[i*s+j]; 4373d177a5cSEmil Constantinescu ierr = VecMAXPY(X[i],s,brow,Ydot);CHKERRQ(ierr); 4383d177a5cSEmil Constantinescu for (j=0; j<r; j++) vrow[j] = sc->v[i*r+j]; 4393d177a5cSEmil Constantinescu ierr = VecMAXPY(X[i],r,vrow,Xold);CHKERRQ(ierr); 4403d177a5cSEmil Constantinescu } 4413d177a5cSEmil Constantinescu PetscFunctionReturn(0); 4423d177a5cSEmil Constantinescu } 4433d177a5cSEmil Constantinescu 4443d177a5cSEmil Constantinescu #undef __FUNCT__ 4453d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLECompleteStep_RescaleAndModify" 4463d177a5cSEmil Constantinescu static PetscErrorCode TSGLLECompleteStep_RescaleAndModify(TSGLLEScheme sc,PetscReal h,TSGLLEScheme next_sc,PetscReal next_h,Vec Ydot[],Vec Xold[],Vec X[]) 4473d177a5cSEmil Constantinescu { 4483d177a5cSEmil Constantinescu PetscErrorCode ierr; 4493d177a5cSEmil Constantinescu PetscScalar brow[32],vrow[32]; 4503d177a5cSEmil Constantinescu PetscReal ratio; 4513d177a5cSEmil Constantinescu PetscInt i,j,p,r,s; 4523d177a5cSEmil Constantinescu 4533d177a5cSEmil Constantinescu PetscFunctionBegin; 4543d177a5cSEmil Constantinescu /* Build the new solution from (X,Ydot) */ 4553d177a5cSEmil Constantinescu p = sc->p; 4563d177a5cSEmil Constantinescu r = sc->r; 4573d177a5cSEmil Constantinescu s = sc->s; 4583d177a5cSEmil Constantinescu ratio = next_h/h; 4593d177a5cSEmil Constantinescu for (i=0; i<r; i++) { 4603d177a5cSEmil Constantinescu ierr = VecZeroEntries(X[i]);CHKERRQ(ierr); 4613d177a5cSEmil Constantinescu for (j=0; j<s; j++) { 4623d177a5cSEmil Constantinescu brow[j] = h*(PetscPowRealInt(ratio,i)*sc->b[i*s+j] 4633d177a5cSEmil Constantinescu + (PetscPowRealInt(ratio,i) - PetscPowRealInt(ratio,p+1))*(+ sc->alpha[i]*sc->phi[0*s+j]) 4643d177a5cSEmil Constantinescu + (PetscPowRealInt(ratio,i) - PetscPowRealInt(ratio,p+2))*(+ sc->beta [i]*sc->phi[1*s+j] 4653d177a5cSEmil Constantinescu + sc->gamma[i]*sc->phi[2*s+j])); 4663d177a5cSEmil Constantinescu } 4673d177a5cSEmil Constantinescu ierr = VecMAXPY(X[i],s,brow,Ydot);CHKERRQ(ierr); 4683d177a5cSEmil Constantinescu for (j=0; j<r; j++) { 4693d177a5cSEmil Constantinescu vrow[j] = (PetscPowRealInt(ratio,i)*sc->v[i*r+j] 4703d177a5cSEmil Constantinescu + (PetscPowRealInt(ratio,i) - PetscPowRealInt(ratio,p+1))*(+ sc->alpha[i]*sc->psi[0*r+j]) 4713d177a5cSEmil Constantinescu + (PetscPowRealInt(ratio,i) - PetscPowRealInt(ratio,p+2))*(+ sc->beta [i]*sc->psi[1*r+j] 4723d177a5cSEmil Constantinescu + sc->gamma[i]*sc->psi[2*r+j])); 4733d177a5cSEmil Constantinescu } 4743d177a5cSEmil Constantinescu ierr = VecMAXPY(X[i],r,vrow,Xold);CHKERRQ(ierr); 4753d177a5cSEmil Constantinescu } 4763d177a5cSEmil Constantinescu if (r < next_sc->r) { 4773d177a5cSEmil Constantinescu if (r+1 != next_sc->r) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot accommodate jump in r greater than 1"); 4783d177a5cSEmil Constantinescu ierr = VecZeroEntries(X[r]);CHKERRQ(ierr); 4793d177a5cSEmil Constantinescu for (j=0; j<s; j++) brow[j] = h*PetscPowRealInt(ratio,p+1)*sc->phi[0*s+j]; 4803d177a5cSEmil Constantinescu ierr = VecMAXPY(X[r],s,brow,Ydot);CHKERRQ(ierr); 4813d177a5cSEmil Constantinescu for (j=0; j<r; j++) vrow[j] = PetscPowRealInt(ratio,p+1)*sc->psi[0*r+j]; 4823d177a5cSEmil Constantinescu ierr = VecMAXPY(X[r],r,vrow,Xold);CHKERRQ(ierr); 4833d177a5cSEmil Constantinescu } 4843d177a5cSEmil Constantinescu PetscFunctionReturn(0); 4853d177a5cSEmil Constantinescu } 4863d177a5cSEmil Constantinescu 4873d177a5cSEmil Constantinescu #undef __FUNCT__ 4883d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLECreate_IRKS" 4893d177a5cSEmil Constantinescu static PetscErrorCode TSGLLECreate_IRKS(TS ts) 4903d177a5cSEmil Constantinescu { 4913d177a5cSEmil Constantinescu TS_GLLE *gl = (TS_GLLE*)ts->data; 4923d177a5cSEmil Constantinescu PetscErrorCode ierr; 4933d177a5cSEmil Constantinescu 4943d177a5cSEmil Constantinescu PetscFunctionBegin; 4953d177a5cSEmil Constantinescu gl->Destroy = TSGLLEDestroy_Default; 4963d177a5cSEmil Constantinescu gl->EstimateHigherMoments = TSGLLEEstimateHigherMoments_Default; 4973d177a5cSEmil Constantinescu gl->CompleteStep = TSGLLECompleteStep_RescaleAndModify; 4983d177a5cSEmil Constantinescu ierr = PetscMalloc1(10,&gl->schemes);CHKERRQ(ierr); 4993d177a5cSEmil Constantinescu gl->nschemes = 0; 5003d177a5cSEmil Constantinescu 5013d177a5cSEmil Constantinescu { 5023d177a5cSEmil Constantinescu /* p=1,q=1, r=s=2, A- and L-stable with error estimates of order 2 and 3 5033d177a5cSEmil Constantinescu * Listed in Butcher & Podhaisky 2006. On error estimation in general linear methods for stiff ODE. 5043d177a5cSEmil Constantinescu * irks(0.3,0,[.3,1],[1],1) 5053d177a5cSEmil Constantinescu * Note: can be made to have classical order (not stage order) 2 by replacing 0.3 with 1-sqrt(1/2) 5063d177a5cSEmil Constantinescu * but doing so would sacrifice the error estimator. 5073d177a5cSEmil Constantinescu */ 5083d177a5cSEmil Constantinescu const PetscScalar c[2] = {3./10., 1.}; 5093d177a5cSEmil Constantinescu const PetscScalar a[2][2] = {{3./10., 0}, {7./10., 3./10.}}; 5103d177a5cSEmil Constantinescu const PetscScalar b[2][2] = {{7./10., 3./10.}, {0,1}}; 5113d177a5cSEmil Constantinescu const PetscScalar u[2][2] = {{1,0},{1,0}}; 5123d177a5cSEmil Constantinescu const PetscScalar v[2][2] = {{1,0},{0,0}}; 5133d177a5cSEmil Constantinescu ierr = TSGLLESchemeCreate(1,1,2,2,c,*a,*b,*u,*v,&gl->schemes[gl->nschemes++]);CHKERRQ(ierr); 5143d177a5cSEmil Constantinescu } 5153d177a5cSEmil Constantinescu 5163d177a5cSEmil Constantinescu { 5173d177a5cSEmil Constantinescu /* p=q=2, r=s=3: irks(4/9,0,[1:3]/3,[0.33852],1) */ 5183d177a5cSEmil Constantinescu /* http://www.math.auckland.ac.nz/~hpod/atlas/i2a.html */ 5193d177a5cSEmil Constantinescu const PetscScalar c[3] = {1./3., 2./3., 1} 5203d177a5cSEmil Constantinescu ,a[3][3] = {{4./9. ,0 , 0}, 5213d177a5cSEmil Constantinescu {1.03750643704090e+00 , 4./9., 0}, 5223d177a5cSEmil Constantinescu {7.67024779410304e-01 , -3.81140216918943e-01, 4./9.}} 5233d177a5cSEmil Constantinescu ,b[3][3] = {{0.767024779410304, -0.381140216918943, 4./9.}, 5243d177a5cSEmil Constantinescu {0.000000000000000, 0.000000000000000, 1.000000000000000}, 5253d177a5cSEmil Constantinescu {-2.075048385225385, 0.621728385225383, 1.277197204924873}} 5263d177a5cSEmil Constantinescu ,u[3][3] = {{1.0000000000000000, -0.1111111111111109, -0.0925925925925922}, 5273d177a5cSEmil Constantinescu {1.0000000000000000, -0.8152842148186744, -0.4199095530877056}, 5283d177a5cSEmil Constantinescu {1.0000000000000000, 0.1696709930641948, 0.0539741070314165}} 5293d177a5cSEmil Constantinescu ,v[3][3] = {{1.0000000000000000, 0.1696709930641948, 0.0539741070314165}, 5303d177a5cSEmil Constantinescu {0.000000000000000, 0.000000000000000, 0.000000000000000}, 5313d177a5cSEmil Constantinescu {0.000000000000000, 0.176122795075129, 0.000000000000000}}; 5323d177a5cSEmil Constantinescu ierr = TSGLLESchemeCreate(2,2,3,3,c,*a,*b,*u,*v,&gl->schemes[gl->nschemes++]);CHKERRQ(ierr); 5333d177a5cSEmil Constantinescu } 5343d177a5cSEmil Constantinescu { 5353d177a5cSEmil Constantinescu /* p=q=3, r=s=4: irks(9/40,0,[1:4]/4,[0.3312 1.0050],[0.49541 1;1 0]) */ 5363d177a5cSEmil Constantinescu const PetscScalar c[4] = {0.25,0.5,0.75,1.0} 5373d177a5cSEmil Constantinescu ,a[4][4] = {{9./40. , 0, 0, 0}, 5383d177a5cSEmil Constantinescu {2.11286958887701e-01 , 9./40. , 0, 0}, 5393d177a5cSEmil Constantinescu {9.46338294287584e-01 , -3.42942861246094e-01, 9./40. , 0}, 5403d177a5cSEmil Constantinescu {0.521490453970721 , -0.662474225622980, 0.490476425459734, 9./40. }} 5413d177a5cSEmil Constantinescu ,b[4][4] = {{0.521490453970721 , -0.662474225622980, 0.490476425459734, 9./40. }, 5423d177a5cSEmil Constantinescu {0.000000000000000 , 0.000000000000000, 0.000000000000000, 1.000000000000000}, 5433d177a5cSEmil Constantinescu {-0.084677029310348 , 1.390757514776085, -1.568157386206001, 2.023192696767826}, 5443d177a5cSEmil Constantinescu {0.465383797936408 , 1.478273530625148, -1.930836081010182, 1.644872111193354}} 5453d177a5cSEmil Constantinescu ,u[4][4] = {{1.00000000000000000 , 0.02500000000001035, -0.02499999999999053, -0.00442708333332865}, 5463d177a5cSEmil Constantinescu {1.00000000000000000 , 0.06371304111232945, -0.04032173972189845, -0.01389438413189452}, 5473d177a5cSEmil Constantinescu {1.00000000000000000 , -0.07839543304147778, 0.04738685705116663, 0.02032603595928376}, 5483d177a5cSEmil Constantinescu {1.00000000000000000 , 0.42550734619251651, 0.10800718022400080, -0.01726712647760034}} 5493d177a5cSEmil Constantinescu ,v[4][4] = {{1.00000000000000000 , 0.42550734619251651, 0.10800718022400080, -0.01726712647760034}, 5503d177a5cSEmil Constantinescu {0.000000000000000 , 0.000000000000000, 0.000000000000000, 0.000000000000000}, 5513d177a5cSEmil Constantinescu {0.000000000000000 , -1.761115796027561, -0.521284157173780, 0.258249384305463}, 5523d177a5cSEmil Constantinescu {0.000000000000000 , -1.657693358744728, -1.052227765232394, 0.521284157173780}}; 5533d177a5cSEmil Constantinescu ierr = TSGLLESchemeCreate(3,3,4,4,c,*a,*b,*u,*v,&gl->schemes[gl->nschemes++]);CHKERRQ(ierr); 5543d177a5cSEmil Constantinescu } 5553d177a5cSEmil Constantinescu { 5563d177a5cSEmil Constantinescu /* p=q=4, r=s=5: 5573d177a5cSEmil Constantinescu irks(3/11,0,[1:5]/5, [0.1715 -0.1238 0.6617],... 5583d177a5cSEmil Constantinescu [ -0.0812 0.4079 1.0000 5593d177a5cSEmil Constantinescu 1.0000 0 0 5603d177a5cSEmil Constantinescu 0.8270 1.0000 0]) 5613d177a5cSEmil Constantinescu */ 5623d177a5cSEmil Constantinescu const PetscScalar c[5] = {0.2,0.4,0.6,0.8,1.0} 5633d177a5cSEmil Constantinescu ,a[5][5] = {{2.72727272727352e-01 , 0.00000000000000e+00, 0.00000000000000e+00 , 0.00000000000000e+00 , 0.00000000000000e+00}, 5643d177a5cSEmil Constantinescu {-1.03980153733431e-01, 2.72727272727405e-01, 0.00000000000000e+00, 0.00000000000000e+00 , 0.00000000000000e+00}, 5653d177a5cSEmil Constantinescu {-1.58615400341492e+00, 7.44168951881122e-01, 2.72727272727309e-01, 0.00000000000000e+00 , 0.00000000000000e+00}, 5663d177a5cSEmil Constantinescu {-8.73658042865628e-01, 5.37884671894595e-01, -1.63298538799523e-01, 2.72727272726996e-01 , 0.00000000000000e+00}, 5673d177a5cSEmil Constantinescu {2.95489397443992e-01 , -1.18481693910097e+00 , -6.68029812659953e-01 , 1.00716687860943e+00 , 2.72727272727288e-01}} 5683d177a5cSEmil Constantinescu ,b[5][5] = {{2.95489397443992e-01 , -1.18481693910097e+00 , -6.68029812659953e-01 , 1.00716687860943e+00 , 2.72727272727288e-01}, 5693d177a5cSEmil Constantinescu {0.00000000000000e+00 , 1.11022302462516e-16 , -2.22044604925031e-16 , 0.00000000000000e+00 , 1.00000000000000e+00}, 5703d177a5cSEmil Constantinescu {-4.05882503986005e+00, -4.00924006567769e+00, -1.38930610972481e+00, 4.45223930308488e+00 , 6.32331093108427e-01}, 5713d177a5cSEmil Constantinescu {8.35690179937017e+00 , -2.26640927349732e+00 , 6.86647884973826e+00 , -5.22595158025740e+00 , 4.50893068837431e+00}, 5723d177a5cSEmil Constantinescu {1.27656267027479e+01 , 2.80882153840821e+00 , 8.91173096522890e+00 , -1.07936444078906e+01 , 4.82534148988854e+00}} 5733d177a5cSEmil Constantinescu ,u[5][5] = {{1.00000000000000e+00 , -7.27272727273551e-02 , -3.45454545454419e-02 , -4.12121212119565e-03 ,-2.96969696964014e-04}, 5743d177a5cSEmil Constantinescu {1.00000000000000e+00 , 2.31252881006154e-01 , -8.29487834416481e-03 , -9.07191207681020e-03 ,-1.70378403743473e-03}, 5753d177a5cSEmil Constantinescu {1.00000000000000e+00 , 1.16925777880663e+00 , 3.59268562942635e-02 , -4.09013451730615e-02 ,-1.02411119670164e-02}, 5763d177a5cSEmil Constantinescu {1.00000000000000e+00 , 1.02634463704356e+00 , 1.59375044913405e-01 , 1.89673015035370e-03 ,-4.89987231897569e-03}, 5773d177a5cSEmil Constantinescu {1.00000000000000e+00 , 1.27746320298021e+00 , 2.37186008132728e-01 , -8.28694373940065e-02 ,-5.34396510196430e-02}} 5783d177a5cSEmil Constantinescu ,v[5][5] = {{1.00000000000000e+00 , 1.27746320298021e+00 , 2.37186008132728e-01 , -8.28694373940065e-02 ,-5.34396510196430e-02}, 5793d177a5cSEmil Constantinescu {0.00000000000000e+00 , -1.77635683940025e-15 , -1.99840144432528e-15 , -9.99200722162641e-16 ,-3.33066907387547e-16}, 5803d177a5cSEmil Constantinescu {0.00000000000000e+00 , 4.37280081906924e+00 , 5.49221645016377e-02 , -8.88913177394943e-02 , 1.12879077989154e-01}, 5813d177a5cSEmil Constantinescu {0.00000000000000e+00 , -1.22399504837280e+01 , -5.21287338448645e+00 , -8.03952325565291e-01 , 4.60298678047147e-01}, 5823d177a5cSEmil Constantinescu {0.00000000000000e+00 , -1.85178762883829e+01 , -5.21411849862624e+00 , -1.04283436528809e+00 , 7.49030161063651e-01}}; 5833d177a5cSEmil Constantinescu ierr = TSGLLESchemeCreate(4,4,5,5,c,*a,*b,*u,*v,&gl->schemes[gl->nschemes++]);CHKERRQ(ierr); 5843d177a5cSEmil Constantinescu } 5853d177a5cSEmil Constantinescu { 5863d177a5cSEmil Constantinescu /* p=q=5, r=s=6; 5873d177a5cSEmil Constantinescu irks(1/3,0,[1:6]/6,... 5883d177a5cSEmil Constantinescu [-0.0489 0.4228 -0.8814 0.9021],... 5893d177a5cSEmil Constantinescu [-0.3474 -0.6617 0.6294 0.2129 5903d177a5cSEmil Constantinescu 0.0044 -0.4256 -0.1427 -0.8936 5913d177a5cSEmil Constantinescu -0.8267 0.4821 0.1371 -0.2557 5923d177a5cSEmil Constantinescu -0.4426 -0.3855 -0.7514 0.3014]) 5933d177a5cSEmil Constantinescu */ 5943d177a5cSEmil Constantinescu const PetscScalar c[6] = {1./6, 2./6, 3./6, 4./6, 5./6, 1.} 5953d177a5cSEmil Constantinescu ,a[6][6] = {{ 3.33333333333940e-01, 0 , 0 , 0 , 0 , 0 }, 5963d177a5cSEmil Constantinescu { -8.64423857333350e-02, 3.33333333332888e-01, 0 , 0 , 0 , 0 }, 5973d177a5cSEmil Constantinescu { -2.16850174258252e+00, -2.23619072028839e+00, 3.33333333335204e-01, 0 , 0 , 0 }, 5983d177a5cSEmil Constantinescu { -4.73160970138997e+00, -3.89265344629268e+00, -2.76318716520933e-01, 3.33333333335759e-01, 0 , 0 }, 5993d177a5cSEmil Constantinescu { -6.75187540297338e+00, -7.90756533769377e+00, 7.90245051802259e-01, -4.48352364517632e-01, 3.33333333328483e-01, 0 }, 6003d177a5cSEmil Constantinescu { -4.26488287921548e+00, -1.19320395589302e+01, 3.38924509887755e+00, -2.23969848002481e+00, 6.62807710124007e-01, 3.33333333335440e-01}} 6013d177a5cSEmil Constantinescu ,b[6][6] = {{ -4.26488287921548e+00, -1.19320395589302e+01, 3.38924509887755e+00, -2.23969848002481e+00, 6.62807710124007e-01, 3.33333333335440e-01}, 6023d177a5cSEmil Constantinescu { -8.88178419700125e-16, 4.44089209850063e-16, -1.54737334057131e-15, -8.88178419700125e-16, 0.00000000000000e+00, 1.00000000000001e+00}, 6033d177a5cSEmil Constantinescu { -2.87780425770651e+01, -1.13520448264971e+01, 2.62002318943161e+01, 2.56943874812797e+01, -3.06702268304488e+01, 6.68067773510103e+00}, 6043d177a5cSEmil Constantinescu { 5.47971245256474e+01, 6.80366875868284e+01, -6.50952588861999e+01, -8.28643975339097e+01, 8.17416943896414e+01, -1.17819043489036e+01}, 6053d177a5cSEmil Constantinescu { -2.33332114788869e+02, 6.12942539462634e+01, -4.91850135865944e+01, 1.82716844135480e+02, -1.29788173979395e+02, 3.09968095651099e+01}, 6063d177a5cSEmil Constantinescu { -1.72049132343751e+02, 8.60194713593999e+00, 7.98154219170200e-01, 1.50371386053218e+02, -1.18515423962066e+02, 2.50898277784663e+01}} 6073d177a5cSEmil Constantinescu ,u[6][6] = {{ 1.00000000000000e+00, -1.66666666666870e-01, -4.16666666664335e-02, -3.85802469124815e-03, -2.25051440302250e-04, -9.64506172339142e-06}, 6083d177a5cSEmil Constantinescu { 1.00000000000000e+00, 8.64423857327162e-02, -4.11484912671353e-02, -1.11450903217645e-02, -1.47651050487126e-03, -1.34395070766826e-04}, 6093d177a5cSEmil Constantinescu { 1.00000000000000e+00, 4.57135912953434e+00, 1.06514719719137e+00, 1.33517564218007e-01, 1.11365952968659e-02, 6.12382756769504e-04}, 6103d177a5cSEmil Constantinescu { 1.00000000000000e+00, 9.23391519753404e+00, 2.22431212392095e+00, 2.91823807741891e-01, 2.52058456411084e-02, 1.22800542949647e-03}, 6113d177a5cSEmil Constantinescu { 1.00000000000000e+00, 1.48175480533865e+01, 3.73439117461835e+00, 5.14648336541804e-01, 4.76430038853402e-02, 2.56798515502156e-03}, 6123d177a5cSEmil Constantinescu { 1.00000000000000e+00, 1.50512347758335e+01, 4.10099701165164e+00, 5.66039141003603e-01, 3.91213893800891e-02, -2.99136269067853e-03}} 6133d177a5cSEmil Constantinescu ,v[6][6] = {{ 1.00000000000000e+00, 1.50512347758335e+01, 4.10099701165164e+00, 5.66039141003603e-01, 3.91213893800891e-02, -2.99136269067853e-03}, 6143d177a5cSEmil Constantinescu { 0.00000000000000e+00, -4.88498130835069e-15, -6.43929354282591e-15, -3.55271367880050e-15, -1.22124532708767e-15, -3.12250225675825e-16}, 6153d177a5cSEmil Constantinescu { 0.00000000000000e+00, 1.22250171233141e+01, -1.77150760606169e+00, 3.54516769879390e-01, 6.22298845883398e-01, 2.31647447450276e-01}, 6163d177a5cSEmil Constantinescu { 0.00000000000000e+00, -4.48339457331040e+01, -3.57363126641880e-01, 5.18750173123425e-01, 6.55727990241799e-02, 1.63175368287079e-01}, 6173d177a5cSEmil Constantinescu { 0.00000000000000e+00, 1.37297394708005e+02, -1.60145272991317e+00, -5.05319555199441e+00, 1.55328940390990e-01, 9.16629423682464e-01}, 6183d177a5cSEmil Constantinescu { 0.00000000000000e+00, 1.05703241119022e+02, -1.16610260983038e+00, -2.99767252773859e+00, -1.13472315553890e-01, 1.09742849254729e+00}}; 6193d177a5cSEmil Constantinescu ierr = TSGLLESchemeCreate(5,5,6,6,c,*a,*b,*u,*v,&gl->schemes[gl->nschemes++]);CHKERRQ(ierr); 6203d177a5cSEmil Constantinescu } 6213d177a5cSEmil Constantinescu PetscFunctionReturn(0); 6223d177a5cSEmil Constantinescu } 6233d177a5cSEmil Constantinescu 6243d177a5cSEmil Constantinescu #undef __FUNCT__ 6253d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLESetType" 6263d177a5cSEmil Constantinescu /*@C 6273d177a5cSEmil Constantinescu TSGLLESetType - sets the class of general linear method to use for time-stepping 6283d177a5cSEmil Constantinescu 6293d177a5cSEmil Constantinescu Collective on TS 6303d177a5cSEmil Constantinescu 6313d177a5cSEmil Constantinescu Input Parameters: 6323d177a5cSEmil Constantinescu + ts - the TS context 6333d177a5cSEmil Constantinescu - type - a method 6343d177a5cSEmil Constantinescu 6353d177a5cSEmil Constantinescu Options Database Key: 6363d177a5cSEmil Constantinescu . -ts_gl_type <type> - sets the method, use -help for a list of available method (e.g. irks) 6373d177a5cSEmil Constantinescu 6383d177a5cSEmil Constantinescu Notes: 6393d177a5cSEmil Constantinescu See "petsc/include/petscts.h" for available methods (for instance) 6403d177a5cSEmil Constantinescu . TSGLLE_IRKS - Diagonally implicit methods with inherent Runge-Kutta stability (for stiff problems) 6413d177a5cSEmil Constantinescu 6423d177a5cSEmil Constantinescu Normally, it is best to use the TSSetFromOptions() command and 6433d177a5cSEmil Constantinescu then set the TSGLLE type from the options database rather than by using 6443d177a5cSEmil Constantinescu this routine. Using the options database provides the user with 6453d177a5cSEmil Constantinescu maximum flexibility in evaluating the many different solvers. 6463d177a5cSEmil Constantinescu The TSGLLESetType() routine is provided for those situations where it 6473d177a5cSEmil Constantinescu is necessary to set the timestepping solver independently of the 6483d177a5cSEmil Constantinescu command line or options database. This might be the case, for example, 6493d177a5cSEmil Constantinescu when the choice of solver changes during the execution of the 6503d177a5cSEmil Constantinescu program, and the user's application is taking responsibility for 6513d177a5cSEmil Constantinescu choosing the appropriate method. 6523d177a5cSEmil Constantinescu 6533d177a5cSEmil Constantinescu Level: intermediate 6543d177a5cSEmil Constantinescu 6553d177a5cSEmil Constantinescu .keywords: TS, TSGLLE, set, type 6563d177a5cSEmil Constantinescu @*/ 6573d177a5cSEmil Constantinescu PetscErrorCode TSGLLESetType(TS ts,TSGLLEType type) 6583d177a5cSEmil Constantinescu { 6593d177a5cSEmil Constantinescu PetscErrorCode ierr; 6603d177a5cSEmil Constantinescu 6613d177a5cSEmil Constantinescu PetscFunctionBegin; 6623d177a5cSEmil Constantinescu PetscValidHeaderSpecific(ts,TS_CLASSID,1); 6633d177a5cSEmil Constantinescu PetscValidCharPointer(type,2); 6643d177a5cSEmil Constantinescu ierr = PetscTryMethod(ts,"TSGLLESetType_C",(TS,TSGLLEType),(ts,type));CHKERRQ(ierr); 6653d177a5cSEmil Constantinescu PetscFunctionReturn(0); 6663d177a5cSEmil Constantinescu } 6673d177a5cSEmil Constantinescu 6683d177a5cSEmil Constantinescu #undef __FUNCT__ 6693d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLESetAcceptType" 6703d177a5cSEmil Constantinescu /*@C 6713d177a5cSEmil Constantinescu TSGLLESetAcceptType - sets the acceptance test 6723d177a5cSEmil Constantinescu 6733d177a5cSEmil Constantinescu Time integrators that need to control error must have the option to reject a time step based on local error 6743d177a5cSEmil Constantinescu estimates. This function allows different schemes to be set. 6753d177a5cSEmil Constantinescu 6763d177a5cSEmil Constantinescu Logically Collective on TS 6773d177a5cSEmil Constantinescu 6783d177a5cSEmil Constantinescu Input Parameters: 6793d177a5cSEmil Constantinescu + ts - the TS context 6803d177a5cSEmil Constantinescu - type - the type 6813d177a5cSEmil Constantinescu 6823d177a5cSEmil Constantinescu Options Database Key: 6833d177a5cSEmil Constantinescu . -ts_gl_accept_type <type> - sets the method used to determine whether to accept or reject a step 6843d177a5cSEmil Constantinescu 6853d177a5cSEmil Constantinescu Level: intermediate 6863d177a5cSEmil Constantinescu 6873d177a5cSEmil Constantinescu .seealso: TS, TSGLLE, TSGLLEAcceptRegister(), TSGLLEAdapt, set type 6883d177a5cSEmil Constantinescu @*/ 6893d177a5cSEmil Constantinescu PetscErrorCode TSGLLESetAcceptType(TS ts,TSGLLEAcceptType type) 6903d177a5cSEmil Constantinescu { 6913d177a5cSEmil Constantinescu PetscErrorCode ierr; 6923d177a5cSEmil Constantinescu 6933d177a5cSEmil Constantinescu PetscFunctionBegin; 6943d177a5cSEmil Constantinescu PetscValidHeaderSpecific(ts,TS_CLASSID,1); 6953d177a5cSEmil Constantinescu PetscValidCharPointer(type,2); 6963d177a5cSEmil Constantinescu ierr = PetscTryMethod(ts,"TSGLLESetAcceptType_C",(TS,TSGLLEAcceptType),(ts,type));CHKERRQ(ierr); 6973d177a5cSEmil Constantinescu PetscFunctionReturn(0); 6983d177a5cSEmil Constantinescu } 6993d177a5cSEmil Constantinescu 7003d177a5cSEmil Constantinescu #undef __FUNCT__ 7013d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLEGetAdapt" 7023d177a5cSEmil Constantinescu /*@C 7033d177a5cSEmil Constantinescu TSGLLEGetAdapt - gets the TSGLLEAdapt object from the TS 7043d177a5cSEmil Constantinescu 7053d177a5cSEmil Constantinescu Not Collective 7063d177a5cSEmil Constantinescu 7073d177a5cSEmil Constantinescu Input Parameter: 7083d177a5cSEmil Constantinescu . ts - the TS context 7093d177a5cSEmil Constantinescu 7103d177a5cSEmil Constantinescu Output Parameter: 7113d177a5cSEmil Constantinescu . adapt - the TSGLLEAdapt context 7123d177a5cSEmil Constantinescu 7133d177a5cSEmil Constantinescu Notes: 7143d177a5cSEmil Constantinescu This allows the user set options on the TSGLLEAdapt object. Usually it is better to do this using the options 7153d177a5cSEmil Constantinescu database, so this function is rarely needed. 7163d177a5cSEmil Constantinescu 7173d177a5cSEmil Constantinescu Level: advanced 7183d177a5cSEmil Constantinescu 7193d177a5cSEmil Constantinescu .seealso: TSGLLEAdapt, TSGLLEAdaptRegister() 7203d177a5cSEmil Constantinescu @*/ 7213d177a5cSEmil Constantinescu PetscErrorCode TSGLLEGetAdapt(TS ts,TSGLLEAdapt *adapt) 7223d177a5cSEmil Constantinescu { 7233d177a5cSEmil Constantinescu PetscErrorCode ierr; 7243d177a5cSEmil Constantinescu 7253d177a5cSEmil Constantinescu PetscFunctionBegin; 7263d177a5cSEmil Constantinescu PetscValidHeaderSpecific(ts,TS_CLASSID,1); 7273d177a5cSEmil Constantinescu PetscValidPointer(adapt,2); 7283d177a5cSEmil Constantinescu ierr = PetscUseMethod(ts,"TSGLLEGetAdapt_C",(TS,TSGLLEAdapt*),(ts,adapt));CHKERRQ(ierr); 7293d177a5cSEmil Constantinescu PetscFunctionReturn(0); 7303d177a5cSEmil Constantinescu } 7313d177a5cSEmil Constantinescu 7323d177a5cSEmil Constantinescu #undef __FUNCT__ 7333d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLEAccept_Always" 7343d177a5cSEmil Constantinescu static PetscErrorCode TSGLLEAccept_Always(TS ts,PetscReal tleft,PetscReal h,const PetscReal enorms[],PetscBool *accept) 7353d177a5cSEmil Constantinescu { 7363d177a5cSEmil Constantinescu PetscFunctionBegin; 7373d177a5cSEmil Constantinescu *accept = PETSC_TRUE; 7383d177a5cSEmil Constantinescu PetscFunctionReturn(0); 7393d177a5cSEmil Constantinescu } 7403d177a5cSEmil Constantinescu 7413d177a5cSEmil Constantinescu #undef __FUNCT__ 7423d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLEUpdateWRMS" 7433d177a5cSEmil Constantinescu static PetscErrorCode TSGLLEUpdateWRMS(TS ts) 7443d177a5cSEmil Constantinescu { 7453d177a5cSEmil Constantinescu TS_GLLE *gl = (TS_GLLE*)ts->data; 7463d177a5cSEmil Constantinescu PetscErrorCode ierr; 7473d177a5cSEmil Constantinescu PetscScalar *x,*w; 7483d177a5cSEmil Constantinescu PetscInt n,i; 7493d177a5cSEmil Constantinescu 7503d177a5cSEmil Constantinescu PetscFunctionBegin; 7513d177a5cSEmil Constantinescu ierr = VecGetArray(gl->X[0],&x);CHKERRQ(ierr); 7523d177a5cSEmil Constantinescu ierr = VecGetArray(gl->W,&w);CHKERRQ(ierr); 7533d177a5cSEmil Constantinescu ierr = VecGetLocalSize(gl->W,&n);CHKERRQ(ierr); 7543d177a5cSEmil Constantinescu for (i=0; i<n; i++) w[i] = 1./(gl->wrms_atol + gl->wrms_rtol*PetscAbsScalar(x[i])); 7553d177a5cSEmil Constantinescu ierr = VecRestoreArray(gl->X[0],&x);CHKERRQ(ierr); 7563d177a5cSEmil Constantinescu ierr = VecRestoreArray(gl->W,&w);CHKERRQ(ierr); 7573d177a5cSEmil Constantinescu PetscFunctionReturn(0); 7583d177a5cSEmil Constantinescu } 7593d177a5cSEmil Constantinescu 7603d177a5cSEmil Constantinescu #undef __FUNCT__ 7613d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLEVecNormWRMS" 7623d177a5cSEmil Constantinescu static PetscErrorCode TSGLLEVecNormWRMS(TS ts,Vec X,PetscReal *nrm) 7633d177a5cSEmil Constantinescu { 7643d177a5cSEmil Constantinescu TS_GLLE *gl = (TS_GLLE*)ts->data; 7653d177a5cSEmil Constantinescu PetscErrorCode ierr; 7663d177a5cSEmil Constantinescu PetscScalar *x,*w; 7673d177a5cSEmil Constantinescu PetscReal sum = 0.0,gsum; 7683d177a5cSEmil Constantinescu PetscInt n,N,i; 7693d177a5cSEmil Constantinescu 7703d177a5cSEmil Constantinescu PetscFunctionBegin; 7713d177a5cSEmil Constantinescu ierr = VecGetArray(X,&x);CHKERRQ(ierr); 7723d177a5cSEmil Constantinescu ierr = VecGetArray(gl->W,&w);CHKERRQ(ierr); 7733d177a5cSEmil Constantinescu ierr = VecGetLocalSize(gl->W,&n);CHKERRQ(ierr); 7743d177a5cSEmil Constantinescu for (i=0; i<n; i++) sum += PetscAbsScalar(PetscSqr(x[i]*w[i])); 7753d177a5cSEmil Constantinescu ierr = VecRestoreArray(X,&x);CHKERRQ(ierr); 7763d177a5cSEmil Constantinescu ierr = VecRestoreArray(gl->W,&w);CHKERRQ(ierr); 7773d177a5cSEmil Constantinescu ierr = MPIU_Allreduce(&sum,&gsum,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)ts));CHKERRQ(ierr); 7783d177a5cSEmil Constantinescu ierr = VecGetSize(gl->W,&N);CHKERRQ(ierr); 7793d177a5cSEmil Constantinescu *nrm = PetscSqrtReal(gsum/(1.*N)); 7803d177a5cSEmil Constantinescu PetscFunctionReturn(0); 7813d177a5cSEmil Constantinescu } 7823d177a5cSEmil Constantinescu 7833d177a5cSEmil Constantinescu #undef __FUNCT__ 7843d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLESetType_GLLE" 7853d177a5cSEmil Constantinescu static PetscErrorCode TSGLLESetType_GLLE(TS ts,TSGLLEType type) 7863d177a5cSEmil Constantinescu { 7873d177a5cSEmil Constantinescu PetscErrorCode ierr,(*r)(TS); 7883d177a5cSEmil Constantinescu PetscBool same; 7893d177a5cSEmil Constantinescu TS_GLLE *gl = (TS_GLLE*)ts->data; 7903d177a5cSEmil Constantinescu 7913d177a5cSEmil Constantinescu PetscFunctionBegin; 7923d177a5cSEmil Constantinescu if (gl->type_name[0]) { 7933d177a5cSEmil Constantinescu ierr = PetscStrcmp(gl->type_name,type,&same);CHKERRQ(ierr); 7943d177a5cSEmil Constantinescu if (same) PetscFunctionReturn(0); 7953d177a5cSEmil Constantinescu ierr = (*gl->Destroy)(gl);CHKERRQ(ierr); 7963d177a5cSEmil Constantinescu } 7973d177a5cSEmil Constantinescu 7983d177a5cSEmil Constantinescu ierr = PetscFunctionListFind(TSGLLEList,type,&r);CHKERRQ(ierr); 7993d177a5cSEmil Constantinescu if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown TSGLLE type \"%s\" given",type); 8003d177a5cSEmil Constantinescu ierr = (*r)(ts);CHKERRQ(ierr); 8013d177a5cSEmil Constantinescu ierr = PetscStrcpy(gl->type_name,type);CHKERRQ(ierr); 8023d177a5cSEmil Constantinescu PetscFunctionReturn(0); 8033d177a5cSEmil Constantinescu } 8043d177a5cSEmil Constantinescu 8053d177a5cSEmil Constantinescu #undef __FUNCT__ 8063d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLESetAcceptType_GLLE" 8073d177a5cSEmil Constantinescu static PetscErrorCode TSGLLESetAcceptType_GLLE(TS ts,TSGLLEAcceptType type) 8083d177a5cSEmil Constantinescu { 8093d177a5cSEmil Constantinescu PetscErrorCode ierr; 8103d177a5cSEmil Constantinescu TSGLLEAcceptFunction r; 8113d177a5cSEmil Constantinescu TS_GLLE *gl = (TS_GLLE*)ts->data; 8123d177a5cSEmil Constantinescu 8133d177a5cSEmil Constantinescu PetscFunctionBegin; 8143d177a5cSEmil Constantinescu ierr = PetscFunctionListFind(TSGLLEAcceptList,type,&r);CHKERRQ(ierr); 8153d177a5cSEmil Constantinescu if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown TSGLLEAccept type \"%s\" given",type); 8163d177a5cSEmil Constantinescu gl->Accept = r; 8173d177a5cSEmil Constantinescu ierr = PetscStrncpy(gl->accept_name,type,sizeof(gl->accept_name));CHKERRQ(ierr); 8183d177a5cSEmil Constantinescu PetscFunctionReturn(0); 8193d177a5cSEmil Constantinescu } 8203d177a5cSEmil Constantinescu 8213d177a5cSEmil Constantinescu #undef __FUNCT__ 8223d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLEGetAdapt_GLLE" 8233d177a5cSEmil Constantinescu static PetscErrorCode TSGLLEGetAdapt_GLLE(TS ts,TSGLLEAdapt *adapt) 8243d177a5cSEmil Constantinescu { 8253d177a5cSEmil Constantinescu PetscErrorCode ierr; 8263d177a5cSEmil Constantinescu TS_GLLE *gl = (TS_GLLE*)ts->data; 8273d177a5cSEmil Constantinescu 8283d177a5cSEmil Constantinescu PetscFunctionBegin; 8293d177a5cSEmil Constantinescu if (!gl->adapt) { 8303d177a5cSEmil Constantinescu ierr = TSGLLEAdaptCreate(PetscObjectComm((PetscObject)ts),&gl->adapt);CHKERRQ(ierr); 8313d177a5cSEmil Constantinescu ierr = PetscObjectIncrementTabLevel((PetscObject)gl->adapt,(PetscObject)ts,1);CHKERRQ(ierr); 8323d177a5cSEmil Constantinescu ierr = PetscLogObjectParent((PetscObject)ts,(PetscObject)gl->adapt);CHKERRQ(ierr); 8333d177a5cSEmil Constantinescu } 8343d177a5cSEmil Constantinescu *adapt = gl->adapt; 8353d177a5cSEmil Constantinescu PetscFunctionReturn(0); 8363d177a5cSEmil Constantinescu } 8373d177a5cSEmil Constantinescu 8383d177a5cSEmil Constantinescu #undef __FUNCT__ 8393d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLEChooseNextScheme" 8403d177a5cSEmil Constantinescu static PetscErrorCode TSGLLEChooseNextScheme(TS ts,PetscReal h,const PetscReal hmnorm[],PetscInt *next_scheme,PetscReal *next_h,PetscBool *finish) 8413d177a5cSEmil Constantinescu { 8423d177a5cSEmil Constantinescu PetscErrorCode ierr; 8433d177a5cSEmil Constantinescu TS_GLLE *gl = (TS_GLLE*)ts->data; 8443d177a5cSEmil Constantinescu PetscInt i,n,cur_p,cur,next_sc,candidates[64],orders[64]; 8453d177a5cSEmil Constantinescu PetscReal errors[64],costs[64],tleft; 8463d177a5cSEmil Constantinescu 8473d177a5cSEmil Constantinescu PetscFunctionBegin; 8483d177a5cSEmil Constantinescu cur = -1; 8493d177a5cSEmil Constantinescu cur_p = gl->schemes[gl->current_scheme]->p; 8503d177a5cSEmil Constantinescu tleft = ts->max_time - (ts->ptime + ts->time_step); 8513d177a5cSEmil Constantinescu for (i=0,n=0; i<gl->nschemes; i++) { 8523d177a5cSEmil Constantinescu TSGLLEScheme sc = gl->schemes[i]; 8533d177a5cSEmil Constantinescu if (sc->p < gl->min_order || gl->max_order < sc->p) continue; 8543d177a5cSEmil Constantinescu if (sc->p == cur_p - 1) errors[n] = PetscAbsScalar(sc->alpha[0])*hmnorm[0]; 8553d177a5cSEmil Constantinescu else if (sc->p == cur_p) errors[n] = PetscAbsScalar(sc->alpha[0])*hmnorm[1]; 8563d177a5cSEmil Constantinescu else if (sc->p == cur_p+1) errors[n] = PetscAbsScalar(sc->alpha[0])*(hmnorm[2]+hmnorm[3]); 8573d177a5cSEmil Constantinescu else continue; 8583d177a5cSEmil Constantinescu candidates[n] = i; 8593d177a5cSEmil Constantinescu orders[n] = PetscMin(sc->p,sc->q); /* order of global truncation error */ 8603d177a5cSEmil Constantinescu costs[n] = sc->s; /* estimate the cost as the number of stages */ 8613d177a5cSEmil Constantinescu if (i == gl->current_scheme) cur = n; 8623d177a5cSEmil Constantinescu n++; 8633d177a5cSEmil Constantinescu } 8643d177a5cSEmil Constantinescu if (cur < 0 || gl->nschemes <= cur) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Current scheme not found in scheme list"); 8653d177a5cSEmil Constantinescu ierr = TSGLLEAdaptChoose(gl->adapt,n,orders,errors,costs,cur,h,tleft,&next_sc,next_h,finish);CHKERRQ(ierr); 8663d177a5cSEmil Constantinescu *next_scheme = candidates[next_sc]; 8673d177a5cSEmil Constantinescu ierr = PetscInfo7(ts,"Adapt chose scheme %d (%d,%d,%d,%d) with step size %6.2e, finish=%d\n",*next_scheme,gl->schemes[*next_scheme]->p,gl->schemes[*next_scheme]->q,gl->schemes[*next_scheme]->r,gl->schemes[*next_scheme]->s,*next_h,*finish);CHKERRQ(ierr); 8683d177a5cSEmil Constantinescu PetscFunctionReturn(0); 8693d177a5cSEmil Constantinescu } 8703d177a5cSEmil Constantinescu 8713d177a5cSEmil Constantinescu #undef __FUNCT__ 8723d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLEGetMaxSizes" 8733d177a5cSEmil Constantinescu static PetscErrorCode TSGLLEGetMaxSizes(TS ts,PetscInt *max_r,PetscInt *max_s) 8743d177a5cSEmil Constantinescu { 8753d177a5cSEmil Constantinescu TS_GLLE *gl = (TS_GLLE*)ts->data; 8763d177a5cSEmil Constantinescu 8773d177a5cSEmil Constantinescu PetscFunctionBegin; 8783d177a5cSEmil Constantinescu *max_r = gl->schemes[gl->nschemes-1]->r; 8793d177a5cSEmil Constantinescu *max_s = gl->schemes[gl->nschemes-1]->s; 8803d177a5cSEmil Constantinescu PetscFunctionReturn(0); 8813d177a5cSEmil Constantinescu } 8823d177a5cSEmil Constantinescu 8833d177a5cSEmil Constantinescu #undef __FUNCT__ 8843d177a5cSEmil Constantinescu #define __FUNCT__ "TSSolve_GLLE" 8853d177a5cSEmil Constantinescu static PetscErrorCode TSSolve_GLLE(TS ts) 8863d177a5cSEmil Constantinescu { 8873d177a5cSEmil Constantinescu TS_GLLE *gl = (TS_GLLE*)ts->data; 8883d177a5cSEmil Constantinescu PetscInt i,k,its,lits,max_r,max_s; 8893d177a5cSEmil Constantinescu PetscBool final_step,finish; 8903d177a5cSEmil Constantinescu SNESConvergedReason snesreason; 8913d177a5cSEmil Constantinescu PetscErrorCode ierr; 8923d177a5cSEmil Constantinescu 8933d177a5cSEmil Constantinescu PetscFunctionBegin; 8943d177a5cSEmil Constantinescu ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr); 8953d177a5cSEmil Constantinescu 8963d177a5cSEmil Constantinescu ierr = TSGLLEGetMaxSizes(ts,&max_r,&max_s);CHKERRQ(ierr); 8973d177a5cSEmil Constantinescu ierr = VecCopy(ts->vec_sol,gl->X[0]);CHKERRQ(ierr); 8983d177a5cSEmil Constantinescu for (i=1; i<max_r; i++) { 8993d177a5cSEmil Constantinescu ierr = VecZeroEntries(gl->X[i]);CHKERRQ(ierr); 9003d177a5cSEmil Constantinescu } 9013d177a5cSEmil Constantinescu ierr = TSGLLEUpdateWRMS(ts);CHKERRQ(ierr); 9023d177a5cSEmil Constantinescu 9033d177a5cSEmil Constantinescu if (0) { 9043d177a5cSEmil Constantinescu /* Find consistent initial data for DAE */ 9053d177a5cSEmil Constantinescu gl->stage_time = ts->ptime + ts->time_step; 9063d177a5cSEmil Constantinescu gl->scoeff = 1.; 9073d177a5cSEmil Constantinescu gl->stage = 0; 9083d177a5cSEmil Constantinescu 9093d177a5cSEmil Constantinescu ierr = VecCopy(ts->vec_sol,gl->Z);CHKERRQ(ierr); 9103d177a5cSEmil Constantinescu ierr = VecCopy(ts->vec_sol,gl->Y);CHKERRQ(ierr); 9113d177a5cSEmil Constantinescu ierr = SNESSolve(ts->snes,NULL,gl->Y);CHKERRQ(ierr); 9123d177a5cSEmil Constantinescu ierr = SNESGetIterationNumber(ts->snes,&its);CHKERRQ(ierr); 9133d177a5cSEmil Constantinescu ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr); 9143d177a5cSEmil Constantinescu ierr = SNESGetConvergedReason(ts->snes,&snesreason);CHKERRQ(ierr); 9153d177a5cSEmil Constantinescu 9163d177a5cSEmil Constantinescu ts->snes_its += its; ts->ksp_its += lits; 9173d177a5cSEmil Constantinescu if (snesreason < 0 && ts->max_snes_failures > 0 && ++ts->num_snes_failures >= ts->max_snes_failures) { 9183d177a5cSEmil Constantinescu ts->reason = TS_DIVERGED_NONLINEAR_SOLVE; 9193d177a5cSEmil Constantinescu ierr = PetscInfo2(ts,"Step=%D, nonlinear solve solve failures %D greater than current TS allowed, stopping solve\n",ts->steps,ts->num_snes_failures);CHKERRQ(ierr); 9203d177a5cSEmil Constantinescu PetscFunctionReturn(0); 9213d177a5cSEmil Constantinescu } 9223d177a5cSEmil Constantinescu } 9233d177a5cSEmil Constantinescu 9243d177a5cSEmil Constantinescu if (gl->current_scheme < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"A starting scheme has not been provided"); 9253d177a5cSEmil Constantinescu 9263d177a5cSEmil Constantinescu for (k=0,final_step=PETSC_FALSE,finish=PETSC_FALSE; k<ts->max_steps && !finish; k++) { 9273d177a5cSEmil Constantinescu PetscInt j,r,s,next_scheme = 0,rejections; 9283d177a5cSEmil Constantinescu PetscReal h,hmnorm[4],enorm[3],next_h; 9293d177a5cSEmil Constantinescu PetscBool accept; 9303d177a5cSEmil Constantinescu const PetscScalar *c,*a,*u; 9313d177a5cSEmil Constantinescu Vec *X,*Ydot,Y; 9323d177a5cSEmil Constantinescu TSGLLEScheme scheme = gl->schemes[gl->current_scheme]; 9333d177a5cSEmil Constantinescu 9343d177a5cSEmil Constantinescu r = scheme->r; s = scheme->s; 9353d177a5cSEmil Constantinescu c = scheme->c; 9363d177a5cSEmil Constantinescu a = scheme->a; u = scheme->u; 9373d177a5cSEmil Constantinescu h = ts->time_step; 9383d177a5cSEmil Constantinescu X = gl->X; Ydot = gl->Ydot; Y = gl->Y; 9393d177a5cSEmil Constantinescu 9403d177a5cSEmil Constantinescu if (ts->ptime > ts->max_time) break; 9413d177a5cSEmil Constantinescu 9423d177a5cSEmil Constantinescu /* 9433d177a5cSEmil Constantinescu We only call PreStep at the start of each STEP, not each STAGE. This is because it is 9443d177a5cSEmil Constantinescu possible to fail (have to restart a step) after multiple stages. 9453d177a5cSEmil Constantinescu */ 9463d177a5cSEmil Constantinescu ierr = TSPreStep(ts);CHKERRQ(ierr); 9473d177a5cSEmil Constantinescu 9483d177a5cSEmil Constantinescu rejections = 0; 9493d177a5cSEmil Constantinescu while (1) { 9503d177a5cSEmil Constantinescu for (i=0; i<s; i++) { 9513d177a5cSEmil Constantinescu PetscScalar shift; 9523d177a5cSEmil Constantinescu gl->scoeff = 1./PetscRealPart(a[i*s+i]); 9533d177a5cSEmil Constantinescu shift = gl->scoeff/ts->time_step; 9543d177a5cSEmil Constantinescu gl->stage = i; 9553d177a5cSEmil Constantinescu gl->stage_time = ts->ptime + PetscRealPart(c[i])*h; 9563d177a5cSEmil Constantinescu 9573d177a5cSEmil Constantinescu /* 9583d177a5cSEmil Constantinescu * Stage equation: Y = h A Y' + U X 9593d177a5cSEmil Constantinescu * We assume that A is lower-triangular so that we can solve the stages (Y,Y') sequentially 9603d177a5cSEmil Constantinescu * Build the affine vector z_i = -[1/(h a_ii)](h sum_j a_ij y'_j + sum_j u_ij x_j) 9613d177a5cSEmil Constantinescu * Then y'_i = z + 1/(h a_ii) y_i 9623d177a5cSEmil Constantinescu */ 9633d177a5cSEmil Constantinescu ierr = VecZeroEntries(gl->Z);CHKERRQ(ierr); 9643d177a5cSEmil Constantinescu for (j=0; j<r; j++) { 9653d177a5cSEmil Constantinescu ierr = VecAXPY(gl->Z,-shift*u[i*r+j],X[j]);CHKERRQ(ierr); 9663d177a5cSEmil Constantinescu } 9673d177a5cSEmil Constantinescu for (j=0; j<i; j++) { 9683d177a5cSEmil Constantinescu ierr = VecAXPY(gl->Z,-shift*h*a[i*s+j],Ydot[j]);CHKERRQ(ierr); 9693d177a5cSEmil Constantinescu } 9703d177a5cSEmil Constantinescu /* Note: Z is used within function evaluation, Ydot = Z + shift*Y */ 9713d177a5cSEmil Constantinescu 9723d177a5cSEmil Constantinescu /* Compute an estimate of Y to start Newton iteration */ 9733d177a5cSEmil Constantinescu if (gl->extrapolate) { 9743d177a5cSEmil Constantinescu if (i==0) { 9753d177a5cSEmil Constantinescu /* Linear extrapolation on the first stage */ 9763d177a5cSEmil Constantinescu ierr = VecWAXPY(Y,c[i]*h,X[1],X[0]);CHKERRQ(ierr); 9773d177a5cSEmil Constantinescu } else { 9783d177a5cSEmil Constantinescu /* Linear extrapolation from the last stage */ 9793d177a5cSEmil Constantinescu ierr = VecAXPY(Y,(c[i]-c[i-1])*h,Ydot[i-1]);CHKERRQ(ierr); 9803d177a5cSEmil Constantinescu } 9813d177a5cSEmil Constantinescu } else if (i==0) { /* Directly use solution from the last step, otherwise reuse the last stage (do nothing) */ 9823d177a5cSEmil Constantinescu ierr = VecCopy(X[0],Y);CHKERRQ(ierr); 9833d177a5cSEmil Constantinescu } 9843d177a5cSEmil Constantinescu 9853d177a5cSEmil Constantinescu /* Solve this stage (Ydot[i] is computed during function evaluation) */ 9863d177a5cSEmil Constantinescu ierr = SNESSolve(ts->snes,NULL,Y);CHKERRQ(ierr); 9873d177a5cSEmil Constantinescu ierr = SNESGetIterationNumber(ts->snes,&its);CHKERRQ(ierr); 9883d177a5cSEmil Constantinescu ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr); 9893d177a5cSEmil Constantinescu ierr = SNESGetConvergedReason(ts->snes,&snesreason);CHKERRQ(ierr); 9903d177a5cSEmil Constantinescu ts->snes_its += its; ts->ksp_its += lits; 9913d177a5cSEmil Constantinescu if (snesreason < 0 && ts->max_snes_failures > 0 && ++ts->num_snes_failures >= ts->max_snes_failures) { 9923d177a5cSEmil Constantinescu ts->reason = TS_DIVERGED_NONLINEAR_SOLVE; 9933d177a5cSEmil Constantinescu ierr = PetscInfo2(ts,"Step=%D, nonlinear solve solve failures %D greater than current TS allowed, stopping solve\n",ts->steps,ts->num_snes_failures);CHKERRQ(ierr); 9943d177a5cSEmil Constantinescu PetscFunctionReturn(0); 9953d177a5cSEmil Constantinescu } 9963d177a5cSEmil Constantinescu } 9973d177a5cSEmil Constantinescu 9983d177a5cSEmil Constantinescu gl->stage_time = ts->ptime + ts->time_step; 9993d177a5cSEmil Constantinescu 10003d177a5cSEmil Constantinescu ierr = (*gl->EstimateHigherMoments)(scheme,h,Ydot,gl->X,gl->himom);CHKERRQ(ierr); 10013d177a5cSEmil Constantinescu /* hmnorm[i] = h^{p+i}x^{(p+i)} with i=0,1,2; hmnorm[3] = h^{p+2}(dx'/dx) x^{(p+1)} */ 10023d177a5cSEmil Constantinescu for (i=0; i<3; i++) { 10033d177a5cSEmil Constantinescu ierr = TSGLLEVecNormWRMS(ts,gl->himom[i],&hmnorm[i+1]);CHKERRQ(ierr); 10043d177a5cSEmil Constantinescu } 10053d177a5cSEmil Constantinescu enorm[0] = PetscRealPart(scheme->alpha[0])*hmnorm[1]; 10063d177a5cSEmil Constantinescu enorm[1] = PetscRealPart(scheme->beta[0]) *hmnorm[2]; 10073d177a5cSEmil Constantinescu enorm[2] = PetscRealPart(scheme->gamma[0])*hmnorm[3]; 10083d177a5cSEmil Constantinescu ierr = (*gl->Accept)(ts,ts->max_time-gl->stage_time,h,enorm,&accept);CHKERRQ(ierr); 10093d177a5cSEmil Constantinescu if (accept) goto accepted; 10103d177a5cSEmil Constantinescu rejections++; 10113d177a5cSEmil Constantinescu ierr = PetscInfo3(ts,"Step %D (t=%g) not accepted, rejections=%D\n",k,gl->stage_time,rejections);CHKERRQ(ierr); 10123d177a5cSEmil Constantinescu if (rejections > gl->max_step_rejections) break; 10133d177a5cSEmil Constantinescu /* 10143d177a5cSEmil Constantinescu There are lots of reasons why a step might be rejected, including solvers not converging and other factors that 10153d177a5cSEmil Constantinescu TSGLLEChooseNextScheme does not support. Additionally, the error estimates may be very screwed up, so I'm not 10163d177a5cSEmil Constantinescu convinced that it's safe to just compute a new error estimate using the same interface as the current adaptor 10173d177a5cSEmil Constantinescu (the adaptor interface probably has to change). Here we make an arbitrary and naive choice. This assumes that 10183d177a5cSEmil Constantinescu steps were written in Nordsieck form. The "correct" method would be to re-complete the previous time step with 10193d177a5cSEmil Constantinescu the correct "next" step size. It is unclear to me whether the present ad-hoc method of rescaling X is stable. 10203d177a5cSEmil Constantinescu */ 10213d177a5cSEmil Constantinescu h *= 0.5; 10223d177a5cSEmil Constantinescu for (i=1; i<scheme->r; i++) { 10233d177a5cSEmil Constantinescu ierr = VecScale(X[i],PetscPowRealInt(0.5,i));CHKERRQ(ierr); 10243d177a5cSEmil Constantinescu } 10253d177a5cSEmil Constantinescu } 10263d177a5cSEmil Constantinescu SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_CONV_FAILED,"Time step %D (t=%g) not accepted after %D failures\n",k,gl->stage_time,rejections); 10273d177a5cSEmil Constantinescu 10283d177a5cSEmil Constantinescu accepted: 10293d177a5cSEmil Constantinescu /* This term is not error, but it *would* be the leading term for a lower order method */ 10303d177a5cSEmil Constantinescu ierr = TSGLLEVecNormWRMS(ts,gl->X[scheme->r-1],&hmnorm[0]);CHKERRQ(ierr); 10313d177a5cSEmil Constantinescu /* Correct scaling so that these are equivalent to norms of the Nordsieck vectors */ 10323d177a5cSEmil Constantinescu 10333d177a5cSEmil Constantinescu ierr = PetscInfo4(ts,"Last moment norm %10.2e, estimated error norms %10.2e %10.2e %10.2e\n",hmnorm[0],enorm[0],enorm[1],enorm[2]);CHKERRQ(ierr); 10343d177a5cSEmil Constantinescu if (!final_step) { 10353d177a5cSEmil Constantinescu ierr = TSGLLEChooseNextScheme(ts,h,hmnorm,&next_scheme,&next_h,&final_step);CHKERRQ(ierr); 10363d177a5cSEmil Constantinescu } else { 10373d177a5cSEmil Constantinescu /* Dummy values to complete the current step in a consistent manner */ 10383d177a5cSEmil Constantinescu next_scheme = gl->current_scheme; 10393d177a5cSEmil Constantinescu next_h = h; 10403d177a5cSEmil Constantinescu finish = PETSC_TRUE; 10413d177a5cSEmil Constantinescu } 10423d177a5cSEmil Constantinescu 10433d177a5cSEmil Constantinescu X = gl->Xold; 10443d177a5cSEmil Constantinescu gl->Xold = gl->X; 10453d177a5cSEmil Constantinescu gl->X = X; 10463d177a5cSEmil Constantinescu ierr = (*gl->CompleteStep)(scheme,h,gl->schemes[next_scheme],next_h,Ydot,gl->Xold,gl->X);CHKERRQ(ierr); 10473d177a5cSEmil Constantinescu 10483d177a5cSEmil Constantinescu ierr = TSGLLEUpdateWRMS(ts);CHKERRQ(ierr); 10493d177a5cSEmil Constantinescu 10503d177a5cSEmil Constantinescu /* Post the solution for the user, we could avoid this copy with a small bit of cleverness */ 10513d177a5cSEmil Constantinescu ierr = VecCopy(gl->X[0],ts->vec_sol);CHKERRQ(ierr); 10523d177a5cSEmil Constantinescu ts->ptime += h; 10533d177a5cSEmil Constantinescu ts->steps++; 10543d177a5cSEmil Constantinescu ts->total_steps++; 10553d177a5cSEmil Constantinescu 1056*d0797844SShri Abhyankar ierr = TSPostEvaluate(ts);CHKERRQ(ierr); 10573d177a5cSEmil Constantinescu ierr = TSPostStep(ts);CHKERRQ(ierr); 10583d177a5cSEmil Constantinescu ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr); 10593d177a5cSEmil Constantinescu 10603d177a5cSEmil Constantinescu gl->current_scheme = next_scheme; 10613d177a5cSEmil Constantinescu ts->time_step = next_h; 10623d177a5cSEmil Constantinescu } 10633d177a5cSEmil Constantinescu PetscFunctionReturn(0); 10643d177a5cSEmil Constantinescu } 10653d177a5cSEmil Constantinescu 10663d177a5cSEmil Constantinescu /*------------------------------------------------------------*/ 10673d177a5cSEmil Constantinescu 10683d177a5cSEmil Constantinescu #undef __FUNCT__ 10693d177a5cSEmil Constantinescu #define __FUNCT__ "TSReset_GLLE" 10703d177a5cSEmil Constantinescu static PetscErrorCode TSReset_GLLE(TS ts) 10713d177a5cSEmil Constantinescu { 10723d177a5cSEmil Constantinescu TS_GLLE *gl = (TS_GLLE*)ts->data; 10733d177a5cSEmil Constantinescu PetscInt max_r,max_s; 10743d177a5cSEmil Constantinescu PetscErrorCode ierr; 10753d177a5cSEmil Constantinescu 10763d177a5cSEmil Constantinescu PetscFunctionBegin; 10773d177a5cSEmil Constantinescu if (gl->setupcalled) { 10783d177a5cSEmil Constantinescu ierr = TSGLLEGetMaxSizes(ts,&max_r,&max_s);CHKERRQ(ierr); 10793d177a5cSEmil Constantinescu ierr = VecDestroyVecs(max_r,&gl->Xold);CHKERRQ(ierr); 10803d177a5cSEmil Constantinescu ierr = VecDestroyVecs(max_r,&gl->X);CHKERRQ(ierr); 10813d177a5cSEmil Constantinescu ierr = VecDestroyVecs(max_s,&gl->Ydot);CHKERRQ(ierr); 10823d177a5cSEmil Constantinescu ierr = VecDestroyVecs(3,&gl->himom);CHKERRQ(ierr); 10833d177a5cSEmil Constantinescu ierr = VecDestroy(&gl->W);CHKERRQ(ierr); 10843d177a5cSEmil Constantinescu ierr = VecDestroy(&gl->Y);CHKERRQ(ierr); 10853d177a5cSEmil Constantinescu ierr = VecDestroy(&gl->Z);CHKERRQ(ierr); 10863d177a5cSEmil Constantinescu } 10873d177a5cSEmil Constantinescu gl->setupcalled = PETSC_FALSE; 10883d177a5cSEmil Constantinescu PetscFunctionReturn(0); 10893d177a5cSEmil Constantinescu } 10903d177a5cSEmil Constantinescu 10913d177a5cSEmil Constantinescu #undef __FUNCT__ 10923d177a5cSEmil Constantinescu #define __FUNCT__ "TSDestroy_GLLE" 10933d177a5cSEmil Constantinescu static PetscErrorCode TSDestroy_GLLE(TS ts) 10943d177a5cSEmil Constantinescu { 10953d177a5cSEmil Constantinescu TS_GLLE *gl = (TS_GLLE*)ts->data; 10963d177a5cSEmil Constantinescu PetscErrorCode ierr; 10973d177a5cSEmil Constantinescu 10983d177a5cSEmil Constantinescu PetscFunctionBegin; 10993d177a5cSEmil Constantinescu ierr = TSReset_GLLE(ts);CHKERRQ(ierr); 11003d177a5cSEmil Constantinescu if (gl->adapt) {ierr = TSGLLEAdaptDestroy(&gl->adapt);CHKERRQ(ierr);} 11013d177a5cSEmil Constantinescu if (gl->Destroy) {ierr = (*gl->Destroy)(gl);CHKERRQ(ierr);} 11023d177a5cSEmil Constantinescu ierr = PetscFree(ts->data);CHKERRQ(ierr); 11033d177a5cSEmil Constantinescu ierr = PetscObjectComposeFunction((PetscObject)ts,"TSGLLESetType_C",NULL);CHKERRQ(ierr); 11043d177a5cSEmil Constantinescu ierr = PetscObjectComposeFunction((PetscObject)ts,"TSGLLESetAcceptType_C",NULL);CHKERRQ(ierr); 11053d177a5cSEmil Constantinescu ierr = PetscObjectComposeFunction((PetscObject)ts,"TSGLLEGetAdapt_C",NULL);CHKERRQ(ierr); 11063d177a5cSEmil Constantinescu PetscFunctionReturn(0); 11073d177a5cSEmil Constantinescu } 11083d177a5cSEmil Constantinescu 11093d177a5cSEmil Constantinescu /* 11103d177a5cSEmil Constantinescu This defines the nonlinear equation that is to be solved with SNES 11113d177a5cSEmil Constantinescu g(x) = f(t,x,z+shift*x) = 0 11123d177a5cSEmil Constantinescu */ 11133d177a5cSEmil Constantinescu #undef __FUNCT__ 11143d177a5cSEmil Constantinescu #define __FUNCT__ "SNESTSFormFunction_GLLE" 11153d177a5cSEmil Constantinescu static PetscErrorCode SNESTSFormFunction_GLLE(SNES snes,Vec x,Vec f,TS ts) 11163d177a5cSEmil Constantinescu { 11173d177a5cSEmil Constantinescu TS_GLLE *gl = (TS_GLLE*)ts->data; 11183d177a5cSEmil Constantinescu PetscErrorCode ierr; 11193d177a5cSEmil Constantinescu Vec Z,Ydot; 11203d177a5cSEmil Constantinescu DM dm,dmsave; 11213d177a5cSEmil Constantinescu 11223d177a5cSEmil Constantinescu PetscFunctionBegin; 11233d177a5cSEmil Constantinescu ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 11243d177a5cSEmil Constantinescu ierr = TSGLLEGetVecs(ts,dm,&Z,&Ydot);CHKERRQ(ierr); 11253d177a5cSEmil Constantinescu ierr = VecWAXPY(Ydot,gl->scoeff/ts->time_step,x,Z);CHKERRQ(ierr); 11263d177a5cSEmil Constantinescu dmsave = ts->dm; 11273d177a5cSEmil Constantinescu ts->dm = dm; 11283d177a5cSEmil Constantinescu ierr = TSComputeIFunction(ts,gl->stage_time,x,Ydot,f,PETSC_FALSE);CHKERRQ(ierr); 11293d177a5cSEmil Constantinescu ts->dm = dmsave; 11303d177a5cSEmil Constantinescu ierr = TSGLLERestoreVecs(ts,dm,&Z,&Ydot);CHKERRQ(ierr); 11313d177a5cSEmil Constantinescu PetscFunctionReturn(0); 11323d177a5cSEmil Constantinescu } 11333d177a5cSEmil Constantinescu 11343d177a5cSEmil Constantinescu #undef __FUNCT__ 11353d177a5cSEmil Constantinescu #define __FUNCT__ "SNESTSFormJacobian_GLLE" 11363d177a5cSEmil Constantinescu static PetscErrorCode SNESTSFormJacobian_GLLE(SNES snes,Vec x,Mat A,Mat B,TS ts) 11373d177a5cSEmil Constantinescu { 11383d177a5cSEmil Constantinescu TS_GLLE *gl = (TS_GLLE*)ts->data; 11393d177a5cSEmil Constantinescu PetscErrorCode ierr; 11403d177a5cSEmil Constantinescu Vec Z,Ydot; 11413d177a5cSEmil Constantinescu DM dm,dmsave; 11423d177a5cSEmil Constantinescu 11433d177a5cSEmil Constantinescu PetscFunctionBegin; 11443d177a5cSEmil Constantinescu ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 11453d177a5cSEmil Constantinescu ierr = TSGLLEGetVecs(ts,dm,&Z,&Ydot);CHKERRQ(ierr); 11463d177a5cSEmil Constantinescu dmsave = ts->dm; 11473d177a5cSEmil Constantinescu ts->dm = dm; 11483d177a5cSEmil Constantinescu /* gl->Xdot will have already been computed in SNESTSFormFunction_GLLE */ 11493d177a5cSEmil Constantinescu ierr = TSComputeIJacobian(ts,gl->stage_time,x,gl->Ydot[gl->stage],gl->scoeff/ts->time_step,A,B,PETSC_FALSE);CHKERRQ(ierr); 11503d177a5cSEmil Constantinescu ts->dm = dmsave; 11513d177a5cSEmil Constantinescu ierr = TSGLLERestoreVecs(ts,dm,&Z,&Ydot);CHKERRQ(ierr); 11523d177a5cSEmil Constantinescu PetscFunctionReturn(0); 11533d177a5cSEmil Constantinescu } 11543d177a5cSEmil Constantinescu 11553d177a5cSEmil Constantinescu 11563d177a5cSEmil Constantinescu #undef __FUNCT__ 11573d177a5cSEmil Constantinescu #define __FUNCT__ "TSSetUp_GLLE" 11583d177a5cSEmil Constantinescu static PetscErrorCode TSSetUp_GLLE(TS ts) 11593d177a5cSEmil Constantinescu { 11603d177a5cSEmil Constantinescu TS_GLLE *gl = (TS_GLLE*)ts->data; 11613d177a5cSEmil Constantinescu PetscInt max_r,max_s; 11623d177a5cSEmil Constantinescu PetscErrorCode ierr; 11633d177a5cSEmil Constantinescu DM dm; 11643d177a5cSEmil Constantinescu 11653d177a5cSEmil Constantinescu PetscFunctionBegin; 11663d177a5cSEmil Constantinescu gl->setupcalled = PETSC_TRUE; 11673d177a5cSEmil Constantinescu ierr = TSGLLEGetMaxSizes(ts,&max_r,&max_s);CHKERRQ(ierr); 11683d177a5cSEmil Constantinescu ierr = VecDuplicateVecs(ts->vec_sol,max_r,&gl->X);CHKERRQ(ierr); 11693d177a5cSEmil Constantinescu ierr = VecDuplicateVecs(ts->vec_sol,max_r,&gl->Xold);CHKERRQ(ierr); 11703d177a5cSEmil Constantinescu ierr = VecDuplicateVecs(ts->vec_sol,max_s,&gl->Ydot);CHKERRQ(ierr); 11713d177a5cSEmil Constantinescu ierr = VecDuplicateVecs(ts->vec_sol,3,&gl->himom);CHKERRQ(ierr); 11723d177a5cSEmil Constantinescu ierr = VecDuplicate(ts->vec_sol,&gl->W);CHKERRQ(ierr); 11733d177a5cSEmil Constantinescu ierr = VecDuplicate(ts->vec_sol,&gl->Y);CHKERRQ(ierr); 11743d177a5cSEmil Constantinescu ierr = VecDuplicate(ts->vec_sol,&gl->Z);CHKERRQ(ierr); 11753d177a5cSEmil Constantinescu 11763d177a5cSEmil Constantinescu /* Default acceptance tests and adaptivity */ 11773d177a5cSEmil Constantinescu if (!gl->Accept) {ierr = TSGLLESetAcceptType(ts,TSGLLEACCEPT_ALWAYS);CHKERRQ(ierr);} 11783d177a5cSEmil Constantinescu if (!gl->adapt) {ierr = TSGLLEGetAdapt(ts,&gl->adapt);CHKERRQ(ierr);} 11793d177a5cSEmil Constantinescu 11803d177a5cSEmil Constantinescu if (gl->current_scheme < 0) { 11813d177a5cSEmil Constantinescu PetscInt i; 11823d177a5cSEmil Constantinescu for (i=0;; i++) { 11833d177a5cSEmil Constantinescu if (gl->schemes[i]->p == gl->start_order) break; 11843d177a5cSEmil Constantinescu if (i+1 == gl->nschemes) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"No schemes available with requested start order %d",i); 11853d177a5cSEmil Constantinescu } 11863d177a5cSEmil Constantinescu gl->current_scheme = i; 11873d177a5cSEmil Constantinescu } 11883d177a5cSEmil Constantinescu ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 11893d177a5cSEmil Constantinescu if (dm) { 11903d177a5cSEmil Constantinescu ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSGLLE,DMRestrictHook_TSGLLE,ts);CHKERRQ(ierr); 11913d177a5cSEmil Constantinescu ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSGLLE,DMSubDomainRestrictHook_TSGLLE,ts);CHKERRQ(ierr); 11923d177a5cSEmil Constantinescu } 11933d177a5cSEmil Constantinescu PetscFunctionReturn(0); 11943d177a5cSEmil Constantinescu } 11953d177a5cSEmil Constantinescu /*------------------------------------------------------------*/ 11963d177a5cSEmil Constantinescu 11973d177a5cSEmil Constantinescu #undef __FUNCT__ 11983d177a5cSEmil Constantinescu #define __FUNCT__ "TSSetFromOptions_GLLE" 11993d177a5cSEmil Constantinescu static PetscErrorCode TSSetFromOptions_GLLE(PetscOptionItems *PetscOptionsObject,TS ts) 12003d177a5cSEmil Constantinescu { 12013d177a5cSEmil Constantinescu TS_GLLE *gl = (TS_GLLE*)ts->data; 12023d177a5cSEmil Constantinescu char tname[256] = TSGLLE_IRKS,completef[256] = "rescale-and-modify"; 12033d177a5cSEmil Constantinescu PetscErrorCode ierr; 12043d177a5cSEmil Constantinescu 12053d177a5cSEmil Constantinescu PetscFunctionBegin; 12063d177a5cSEmil Constantinescu ierr = PetscOptionsHead(PetscOptionsObject,"General Linear ODE solver options");CHKERRQ(ierr); 12073d177a5cSEmil Constantinescu { 12083d177a5cSEmil Constantinescu PetscBool flg; 12093d177a5cSEmil Constantinescu ierr = PetscOptionsFList("-ts_gl_type","Type of GL method","TSGLLESetType",TSGLLEList,gl->type_name[0] ? gl->type_name : tname,tname,sizeof(tname),&flg);CHKERRQ(ierr); 12103d177a5cSEmil Constantinescu if (flg || !gl->type_name[0]) { 12113d177a5cSEmil Constantinescu ierr = TSGLLESetType(ts,tname);CHKERRQ(ierr); 12123d177a5cSEmil Constantinescu } 12133d177a5cSEmil Constantinescu ierr = PetscOptionsInt("-ts_gl_max_step_rejections","Maximum number of times to attempt a step","None",gl->max_step_rejections,&gl->max_step_rejections,NULL);CHKERRQ(ierr); 12143d177a5cSEmil Constantinescu ierr = PetscOptionsInt("-ts_gl_max_order","Maximum order to try","TSGLLESetMaxOrder",gl->max_order,&gl->max_order,NULL);CHKERRQ(ierr); 12153d177a5cSEmil Constantinescu ierr = PetscOptionsInt("-ts_gl_min_order","Minimum order to try","TSGLLESetMinOrder",gl->min_order,&gl->min_order,NULL);CHKERRQ(ierr); 12163d177a5cSEmil Constantinescu ierr = PetscOptionsInt("-ts_gl_start_order","Initial order to try","TSGLLESetMinOrder",gl->start_order,&gl->start_order,NULL);CHKERRQ(ierr); 12173d177a5cSEmil Constantinescu ierr = PetscOptionsEnum("-ts_gl_error_direction","Which direction to look when estimating error","TSGLLESetErrorDirection",TSGLLEErrorDirections,(PetscEnum)gl->error_direction,(PetscEnum*)&gl->error_direction,NULL);CHKERRQ(ierr); 12183d177a5cSEmil Constantinescu ierr = PetscOptionsBool("-ts_gl_extrapolate","Extrapolate stage solution from previous solution (sometimes unstable)","TSGLLESetExtrapolate",gl->extrapolate,&gl->extrapolate,NULL);CHKERRQ(ierr); 12193d177a5cSEmil Constantinescu ierr = PetscOptionsReal("-ts_gl_atol","Absolute tolerance","TSGLLESetTolerances",gl->wrms_atol,&gl->wrms_atol,NULL);CHKERRQ(ierr); 12203d177a5cSEmil Constantinescu ierr = PetscOptionsReal("-ts_gl_rtol","Relative tolerance","TSGLLESetTolerances",gl->wrms_rtol,&gl->wrms_rtol,NULL);CHKERRQ(ierr); 12213d177a5cSEmil Constantinescu ierr = PetscOptionsString("-ts_gl_complete","Method to use for completing the step","none",completef,completef,sizeof(completef),&flg);CHKERRQ(ierr); 12223d177a5cSEmil Constantinescu if (flg) { 12233d177a5cSEmil Constantinescu PetscBool match1,match2; 12243d177a5cSEmil Constantinescu ierr = PetscStrcmp(completef,"rescale",&match1);CHKERRQ(ierr); 12253d177a5cSEmil Constantinescu ierr = PetscStrcmp(completef,"rescale-and-modify",&match2);CHKERRQ(ierr); 12263d177a5cSEmil Constantinescu if (match1) gl->CompleteStep = TSGLLECompleteStep_Rescale; 12273d177a5cSEmil Constantinescu else if (match2) gl->CompleteStep = TSGLLECompleteStep_RescaleAndModify; 12283d177a5cSEmil Constantinescu else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"%s",completef); 12293d177a5cSEmil Constantinescu } 12303d177a5cSEmil Constantinescu { 12313d177a5cSEmil Constantinescu char type[256] = TSGLLEACCEPT_ALWAYS; 12323d177a5cSEmil Constantinescu ierr = 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);CHKERRQ(ierr); 12333d177a5cSEmil Constantinescu if (flg || !gl->accept_name[0]) { 12343d177a5cSEmil Constantinescu ierr = TSGLLESetAcceptType(ts,type);CHKERRQ(ierr); 12353d177a5cSEmil Constantinescu } 12363d177a5cSEmil Constantinescu } 12373d177a5cSEmil Constantinescu { 12383d177a5cSEmil Constantinescu TSGLLEAdapt adapt; 12393d177a5cSEmil Constantinescu ierr = TSGLLEGetAdapt(ts,&adapt);CHKERRQ(ierr); 12403d177a5cSEmil Constantinescu ierr = TSGLLEAdaptSetFromOptions(PetscOptionsObject,adapt);CHKERRQ(ierr); 12413d177a5cSEmil Constantinescu } 12423d177a5cSEmil Constantinescu } 12433d177a5cSEmil Constantinescu ierr = PetscOptionsTail();CHKERRQ(ierr); 12443d177a5cSEmil Constantinescu PetscFunctionReturn(0); 12453d177a5cSEmil Constantinescu } 12463d177a5cSEmil Constantinescu 12473d177a5cSEmil Constantinescu #undef __FUNCT__ 12483d177a5cSEmil Constantinescu #define __FUNCT__ "TSView_GLLE" 12493d177a5cSEmil Constantinescu static PetscErrorCode TSView_GLLE(TS ts,PetscViewer viewer) 12503d177a5cSEmil Constantinescu { 12513d177a5cSEmil Constantinescu TS_GLLE *gl = (TS_GLLE*)ts->data; 12523d177a5cSEmil Constantinescu PetscInt i; 12533d177a5cSEmil Constantinescu PetscBool iascii,details; 12543d177a5cSEmil Constantinescu PetscErrorCode ierr; 12553d177a5cSEmil Constantinescu 12563d177a5cSEmil Constantinescu PetscFunctionBegin; 12573d177a5cSEmil Constantinescu ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 12583d177a5cSEmil Constantinescu if (iascii) { 12593d177a5cSEmil Constantinescu ierr = PetscViewerASCIIPrintf(viewer," min order %D, max order %D, current order %D\n",gl->min_order,gl->max_order,gl->schemes[gl->current_scheme]->p);CHKERRQ(ierr); 12603d177a5cSEmil Constantinescu ierr = PetscViewerASCIIPrintf(viewer," Error estimation: %s\n",TSGLLEErrorDirections[gl->error_direction]);CHKERRQ(ierr); 12613d177a5cSEmil Constantinescu ierr = PetscViewerASCIIPrintf(viewer," Extrapolation: %s\n",gl->extrapolate ? "yes" : "no");CHKERRQ(ierr); 12623d177a5cSEmil Constantinescu ierr = PetscViewerASCIIPrintf(viewer," Acceptance test: %s\n",gl->accept_name[0] ? gl->accept_name : "(not yet set)");CHKERRQ(ierr); 12633d177a5cSEmil Constantinescu ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 12643d177a5cSEmil Constantinescu ierr = TSGLLEAdaptView(gl->adapt,viewer);CHKERRQ(ierr); 12653d177a5cSEmil Constantinescu ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 12663d177a5cSEmil Constantinescu ierr = PetscViewerASCIIPrintf(viewer," type: %s\n",gl->type_name[0] ? gl->type_name : "(not yet set)");CHKERRQ(ierr); 12673d177a5cSEmil Constantinescu ierr = PetscViewerASCIIPrintf(viewer,"Schemes within family (%d):\n",gl->nschemes);CHKERRQ(ierr); 12683d177a5cSEmil Constantinescu details = PETSC_FALSE; 12693d177a5cSEmil Constantinescu ierr = PetscOptionsGetBool(((PetscObject)ts)->options,((PetscObject)ts)->prefix,"-ts_gl_view_detailed",&details,NULL);CHKERRQ(ierr); 12703d177a5cSEmil Constantinescu ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 12713d177a5cSEmil Constantinescu for (i=0; i<gl->nschemes; i++) { 12723d177a5cSEmil Constantinescu ierr = TSGLLESchemeView(gl->schemes[i],details,viewer);CHKERRQ(ierr); 12733d177a5cSEmil Constantinescu } 12743d177a5cSEmil Constantinescu if (gl->View) { 12753d177a5cSEmil Constantinescu ierr = (*gl->View)(gl,viewer);CHKERRQ(ierr); 12763d177a5cSEmil Constantinescu } 12773d177a5cSEmil Constantinescu ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 12783d177a5cSEmil Constantinescu } 12793d177a5cSEmil Constantinescu ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr); 12803d177a5cSEmil Constantinescu PetscFunctionReturn(0); 12813d177a5cSEmil Constantinescu } 12823d177a5cSEmil Constantinescu 12833d177a5cSEmil Constantinescu #undef __FUNCT__ 12843d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLERegister" 12853d177a5cSEmil Constantinescu /*@C 12863d177a5cSEmil Constantinescu TSGLLERegister - adds a TSGLLE implementation 12873d177a5cSEmil Constantinescu 12883d177a5cSEmil Constantinescu Not Collective 12893d177a5cSEmil Constantinescu 12903d177a5cSEmil Constantinescu Input Parameters: 12913d177a5cSEmil Constantinescu + name_scheme - name of user-defined general linear scheme 12923d177a5cSEmil Constantinescu - routine_create - routine to create method context 12933d177a5cSEmil Constantinescu 12943d177a5cSEmil Constantinescu Notes: 12953d177a5cSEmil Constantinescu TSGLLERegister() may be called multiple times to add several user-defined families. 12963d177a5cSEmil Constantinescu 12973d177a5cSEmil Constantinescu Sample usage: 12983d177a5cSEmil Constantinescu .vb 12993d177a5cSEmil Constantinescu TSGLLERegister("my_scheme",MySchemeCreate); 13003d177a5cSEmil Constantinescu .ve 13013d177a5cSEmil Constantinescu 13023d177a5cSEmil Constantinescu Then, your scheme can be chosen with the procedural interface via 13033d177a5cSEmil Constantinescu $ TSGLLESetType(ts,"my_scheme") 13043d177a5cSEmil Constantinescu or at runtime via the option 13053d177a5cSEmil Constantinescu $ -ts_gl_type my_scheme 13063d177a5cSEmil Constantinescu 13073d177a5cSEmil Constantinescu Level: advanced 13083d177a5cSEmil Constantinescu 13093d177a5cSEmil Constantinescu .keywords: TSGLLE, register 13103d177a5cSEmil Constantinescu 13113d177a5cSEmil Constantinescu .seealso: TSGLLERegisterAll() 13123d177a5cSEmil Constantinescu @*/ 13133d177a5cSEmil Constantinescu PetscErrorCode TSGLLERegister(const char sname[],PetscErrorCode (*function)(TS)) 13143d177a5cSEmil Constantinescu { 13153d177a5cSEmil Constantinescu PetscErrorCode ierr; 13163d177a5cSEmil Constantinescu 13173d177a5cSEmil Constantinescu PetscFunctionBegin; 13183d177a5cSEmil Constantinescu ierr = PetscFunctionListAdd(&TSGLLEList,sname,function);CHKERRQ(ierr); 13193d177a5cSEmil Constantinescu PetscFunctionReturn(0); 13203d177a5cSEmil Constantinescu } 13213d177a5cSEmil Constantinescu 13223d177a5cSEmil Constantinescu #undef __FUNCT__ 13233d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLEAcceptRegister" 13243d177a5cSEmil Constantinescu /*@C 13253d177a5cSEmil Constantinescu TSGLLEAcceptRegister - adds a TSGLLE acceptance scheme 13263d177a5cSEmil Constantinescu 13273d177a5cSEmil Constantinescu Not Collective 13283d177a5cSEmil Constantinescu 13293d177a5cSEmil Constantinescu Input Parameters: 13303d177a5cSEmil Constantinescu + name_scheme - name of user-defined acceptance scheme 13313d177a5cSEmil Constantinescu - routine_create - routine to create method context 13323d177a5cSEmil Constantinescu 13333d177a5cSEmil Constantinescu Notes: 13343d177a5cSEmil Constantinescu TSGLLEAcceptRegister() may be called multiple times to add several user-defined families. 13353d177a5cSEmil Constantinescu 13363d177a5cSEmil Constantinescu Sample usage: 13373d177a5cSEmil Constantinescu .vb 13383d177a5cSEmil Constantinescu TSGLLEAcceptRegister("my_scheme",MySchemeCreate); 13393d177a5cSEmil Constantinescu .ve 13403d177a5cSEmil Constantinescu 13413d177a5cSEmil Constantinescu Then, your scheme can be chosen with the procedural interface via 13423d177a5cSEmil Constantinescu $ TSGLLESetAcceptType(ts,"my_scheme") 13433d177a5cSEmil Constantinescu or at runtime via the option 13443d177a5cSEmil Constantinescu $ -ts_gl_accept_type my_scheme 13453d177a5cSEmil Constantinescu 13463d177a5cSEmil Constantinescu Level: advanced 13473d177a5cSEmil Constantinescu 13483d177a5cSEmil Constantinescu .keywords: TSGLLE, TSGLLEAcceptType, register 13493d177a5cSEmil Constantinescu 13503d177a5cSEmil Constantinescu .seealso: TSGLLERegisterAll() 13513d177a5cSEmil Constantinescu @*/ 13523d177a5cSEmil Constantinescu PetscErrorCode TSGLLEAcceptRegister(const char sname[],TSGLLEAcceptFunction function) 13533d177a5cSEmil Constantinescu { 13543d177a5cSEmil Constantinescu PetscErrorCode ierr; 13553d177a5cSEmil Constantinescu 13563d177a5cSEmil Constantinescu PetscFunctionBegin; 13573d177a5cSEmil Constantinescu ierr = PetscFunctionListAdd(&TSGLLEAcceptList,sname,function);CHKERRQ(ierr); 13583d177a5cSEmil Constantinescu PetscFunctionReturn(0); 13593d177a5cSEmil Constantinescu } 13603d177a5cSEmil Constantinescu 13613d177a5cSEmil Constantinescu #undef __FUNCT__ 13623d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLERegisterAll" 13633d177a5cSEmil Constantinescu /*@C 13643d177a5cSEmil Constantinescu TSGLLERegisterAll - Registers all of the general linear methods in TSGLLE 13653d177a5cSEmil Constantinescu 13663d177a5cSEmil Constantinescu Not Collective 13673d177a5cSEmil Constantinescu 13683d177a5cSEmil Constantinescu Level: advanced 13693d177a5cSEmil Constantinescu 13703d177a5cSEmil Constantinescu .keywords: TS, TSGLLE, register, all 13713d177a5cSEmil Constantinescu 13723d177a5cSEmil Constantinescu .seealso: TSGLLERegisterDestroy() 13733d177a5cSEmil Constantinescu @*/ 13743d177a5cSEmil Constantinescu PetscErrorCode TSGLLERegisterAll(void) 13753d177a5cSEmil Constantinescu { 13763d177a5cSEmil Constantinescu PetscErrorCode ierr; 13773d177a5cSEmil Constantinescu 13783d177a5cSEmil Constantinescu PetscFunctionBegin; 13793d177a5cSEmil Constantinescu if (TSGLLERegisterAllCalled) PetscFunctionReturn(0); 13803d177a5cSEmil Constantinescu TSGLLERegisterAllCalled = PETSC_TRUE; 13813d177a5cSEmil Constantinescu 13823d177a5cSEmil Constantinescu ierr = TSGLLERegister(TSGLLE_IRKS, TSGLLECreate_IRKS);CHKERRQ(ierr); 13833d177a5cSEmil Constantinescu ierr = TSGLLEAcceptRegister(TSGLLEACCEPT_ALWAYS,TSGLLEAccept_Always);CHKERRQ(ierr); 13843d177a5cSEmil Constantinescu PetscFunctionReturn(0); 13853d177a5cSEmil Constantinescu } 13863d177a5cSEmil Constantinescu 13873d177a5cSEmil Constantinescu #undef __FUNCT__ 13883d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLEInitializePackage" 13893d177a5cSEmil Constantinescu /*@C 13903d177a5cSEmil Constantinescu TSGLLEInitializePackage - This function initializes everything in the TSGLLE package. It is called 13913d177a5cSEmil Constantinescu from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_GLLE() 13923d177a5cSEmil Constantinescu when using static libraries. 13933d177a5cSEmil Constantinescu 13943d177a5cSEmil Constantinescu Level: developer 13953d177a5cSEmil Constantinescu 13963d177a5cSEmil Constantinescu .keywords: TS, TSGLLE, initialize, package 13973d177a5cSEmil Constantinescu .seealso: PetscInitialize() 13983d177a5cSEmil Constantinescu @*/ 13993d177a5cSEmil Constantinescu PetscErrorCode TSGLLEInitializePackage(void) 14003d177a5cSEmil Constantinescu { 14013d177a5cSEmil Constantinescu PetscErrorCode ierr; 14023d177a5cSEmil Constantinescu 14033d177a5cSEmil Constantinescu PetscFunctionBegin; 14043d177a5cSEmil Constantinescu if (TSGLLEPackageInitialized) PetscFunctionReturn(0); 14053d177a5cSEmil Constantinescu TSGLLEPackageInitialized = PETSC_TRUE; 14063d177a5cSEmil Constantinescu ierr = TSGLLERegisterAll();CHKERRQ(ierr); 14073d177a5cSEmil Constantinescu ierr = PetscRegisterFinalize(TSGLLEFinalizePackage);CHKERRQ(ierr); 14083d177a5cSEmil Constantinescu PetscFunctionReturn(0); 14093d177a5cSEmil Constantinescu } 14103d177a5cSEmil Constantinescu 14113d177a5cSEmil Constantinescu #undef __FUNCT__ 14123d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLEFinalizePackage" 14133d177a5cSEmil Constantinescu /*@C 14143d177a5cSEmil Constantinescu TSGLLEFinalizePackage - This function destroys everything in the TSGLLE package. It is 14153d177a5cSEmil Constantinescu called from PetscFinalize(). 14163d177a5cSEmil Constantinescu 14173d177a5cSEmil Constantinescu Level: developer 14183d177a5cSEmil Constantinescu 14193d177a5cSEmil Constantinescu .keywords: Petsc, destroy, package 14203d177a5cSEmil Constantinescu .seealso: PetscFinalize() 14213d177a5cSEmil Constantinescu @*/ 14223d177a5cSEmil Constantinescu PetscErrorCode TSGLLEFinalizePackage(void) 14233d177a5cSEmil Constantinescu { 14243d177a5cSEmil Constantinescu PetscErrorCode ierr; 14253d177a5cSEmil Constantinescu 14263d177a5cSEmil Constantinescu PetscFunctionBegin; 14273d177a5cSEmil Constantinescu ierr = PetscFunctionListDestroy(&TSGLLEList);CHKERRQ(ierr); 14283d177a5cSEmil Constantinescu ierr = PetscFunctionListDestroy(&TSGLLEAcceptList);CHKERRQ(ierr); 14293d177a5cSEmil Constantinescu TSGLLEPackageInitialized = PETSC_FALSE; 14303d177a5cSEmil Constantinescu TSGLLERegisterAllCalled = PETSC_FALSE; 14313d177a5cSEmil Constantinescu PetscFunctionReturn(0); 14323d177a5cSEmil Constantinescu } 14333d177a5cSEmil Constantinescu 14343d177a5cSEmil Constantinescu /* ------------------------------------------------------------ */ 14353d177a5cSEmil Constantinescu /*MC 14363d177a5cSEmil Constantinescu TSGLLE - DAE solver using implicit General Linear methods 14373d177a5cSEmil Constantinescu 14383d177a5cSEmil Constantinescu These methods contain Runge-Kutta and multistep schemes as special cases. These special cases have some fundamental 14393d177a5cSEmil Constantinescu limitations. For example, diagonally implicit Runge-Kutta cannot have stage order greater than 1 which limits their 14403d177a5cSEmil Constantinescu applicability to very stiff systems. Meanwhile, multistep methods cannot be A-stable for order greater than 2 and BDF 14413d177a5cSEmil Constantinescu are not 0-stable for order greater than 6. GL methods can be A- and L-stable with arbitrarily high stage order and 14423d177a5cSEmil Constantinescu reliable error estimates for both 1 and 2 orders higher to facilitate adaptive step sizes and adaptive order schemes. 14433d177a5cSEmil Constantinescu All this is possible while preserving a singly diagonally implicit structure. 14443d177a5cSEmil Constantinescu 14453d177a5cSEmil Constantinescu Options database keys: 14463d177a5cSEmil Constantinescu + -ts_gl_type <type> - the class of general linear method (irks) 14473d177a5cSEmil Constantinescu . -ts_gl_rtol <tol> - relative error 14483d177a5cSEmil Constantinescu . -ts_gl_atol <tol> - absolute error 14493d177a5cSEmil Constantinescu . -ts_gl_min_order <p> - minimum order method to consider (default=1) 14503d177a5cSEmil Constantinescu . -ts_gl_max_order <p> - maximum order method to consider (default=3) 14513d177a5cSEmil Constantinescu . -ts_gl_start_order <p> - order of starting method (default=1) 14523d177a5cSEmil Constantinescu . -ts_gl_complete <method> - method to use for completing the step (rescale-and-modify or rescale) 14533d177a5cSEmil Constantinescu - -ts_adapt_type <method> - adaptive controller to use (none step both) 14543d177a5cSEmil Constantinescu 14553d177a5cSEmil Constantinescu Notes: 14563d177a5cSEmil Constantinescu This integrator can be applied to DAE. 14573d177a5cSEmil Constantinescu 14583d177a5cSEmil Constantinescu Diagonally implicit general linear (DIGL) methods are a generalization of diagonally implicit Runge-Kutta (DIRK). 14593d177a5cSEmil Constantinescu They are represented by the tableau 14603d177a5cSEmil Constantinescu 14613d177a5cSEmil Constantinescu .vb 14623d177a5cSEmil Constantinescu A | U 14633d177a5cSEmil Constantinescu ------- 14643d177a5cSEmil Constantinescu B | V 14653d177a5cSEmil Constantinescu .ve 14663d177a5cSEmil Constantinescu 14673d177a5cSEmil Constantinescu combined with a vector c of abscissa. "Diagonally implicit" means that A is lower triangular. 14683d177a5cSEmil Constantinescu A step of the general method reads 14693d177a5cSEmil Constantinescu 14703d177a5cSEmil Constantinescu .vb 14713d177a5cSEmil Constantinescu [ Y ] = [A U] [ Y' ] 14723d177a5cSEmil Constantinescu [X^k] = [B V] [X^{k-1}] 14733d177a5cSEmil Constantinescu .ve 14743d177a5cSEmil Constantinescu 14753d177a5cSEmil Constantinescu where Y is the multivector of stage values, Y' is the multivector of stage derivatives, X^k is the Nordsieck vector of 14763d177a5cSEmil Constantinescu the solution at step k. The Nordsieck vector consists of the first r moments of the solution, given by 14773d177a5cSEmil Constantinescu 14783d177a5cSEmil Constantinescu .vb 14793d177a5cSEmil Constantinescu X = [x_0,x_1,...,x_{r-1}] = [x, h x', h^2 x'', ..., h^{r-1} x^{(r-1)} ] 14803d177a5cSEmil Constantinescu .ve 14813d177a5cSEmil Constantinescu 14823d177a5cSEmil Constantinescu If A is lower triangular, we can solve the stages (Y,Y') sequentially 14833d177a5cSEmil Constantinescu 14843d177a5cSEmil Constantinescu .vb 14853d177a5cSEmil 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} 14863d177a5cSEmil Constantinescu .ve 14873d177a5cSEmil Constantinescu 14883d177a5cSEmil Constantinescu and then construct the pieces to carry to the next step 14893d177a5cSEmil Constantinescu 14903d177a5cSEmil Constantinescu .vb 14913d177a5cSEmil 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} 14923d177a5cSEmil Constantinescu .ve 14933d177a5cSEmil Constantinescu 14943d177a5cSEmil Constantinescu Note that when the equations are cast in implicit form, we are using the stage equation to define y'_i 14953d177a5cSEmil Constantinescu in terms of y_i and known stuff (y_j for j<i and x_j for all j). 14963d177a5cSEmil Constantinescu 14973d177a5cSEmil Constantinescu 14983d177a5cSEmil Constantinescu Error estimation 14993d177a5cSEmil Constantinescu 15003d177a5cSEmil Constantinescu At present, the most attractive GL methods for stiff problems are singly diagonally implicit schemes which posses 15013d177a5cSEmil Constantinescu Inherent Runge-Kutta Stability (IRKS). These methods have r=s, the number of items passed between steps is equal to 15023d177a5cSEmil Constantinescu the number of stages. The order and stage-order are one less than the number of stages. We use the error estimates 15033d177a5cSEmil Constantinescu in the 2007 paper which provide the following estimates 15043d177a5cSEmil Constantinescu 15053d177a5cSEmil Constantinescu .vb 15063d177a5cSEmil Constantinescu h^{p+1} X^{(p+1)} = phi_0^T Y' + [0 psi_0^T] Xold 15073d177a5cSEmil Constantinescu h^{p+2} X^{(p+2)} = phi_1^T Y' + [0 psi_1^T] Xold 15083d177a5cSEmil Constantinescu h^{p+2} (dx'/dx) X^{(p+1)} = phi_2^T Y' + [0 psi_2^T] Xold 15093d177a5cSEmil Constantinescu .ve 15103d177a5cSEmil Constantinescu 15113d177a5cSEmil Constantinescu These estimates are accurate to O(h^{p+3}). 15123d177a5cSEmil Constantinescu 15133d177a5cSEmil Constantinescu Changing the step size 15143d177a5cSEmil Constantinescu 15153d177a5cSEmil Constantinescu We use the generalized "rescale and modify" scheme, see equation (4.5) of the 2007 paper. 15163d177a5cSEmil Constantinescu 15173d177a5cSEmil Constantinescu Level: beginner 15183d177a5cSEmil Constantinescu 15193d177a5cSEmil Constantinescu References: 15203d177a5cSEmil Constantinescu + 1. - John Butcher and Z. Jackieweicz and W. Wright, On error propagation in general linear methods for 15213d177a5cSEmil Constantinescu ordinary differential equations, Journal of Complexity, Vol 23, 2007. 15223d177a5cSEmil Constantinescu - 2. - John Butcher, Numerical methods for ordinary differential equations, second edition, Wiley, 2009. 15233d177a5cSEmil Constantinescu 15243d177a5cSEmil Constantinescu .seealso: TSCreate(), TS, TSSetType() 15253d177a5cSEmil Constantinescu 15263d177a5cSEmil Constantinescu M*/ 15273d177a5cSEmil Constantinescu #undef __FUNCT__ 15283d177a5cSEmil Constantinescu #define __FUNCT__ "TSCreate_GLLE" 15293d177a5cSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSCreate_GLLE(TS ts) 15303d177a5cSEmil Constantinescu { 15313d177a5cSEmil Constantinescu TS_GLLE *gl; 15323d177a5cSEmil Constantinescu PetscErrorCode ierr; 15333d177a5cSEmil Constantinescu 15343d177a5cSEmil Constantinescu PetscFunctionBegin; 15353d177a5cSEmil Constantinescu ierr = TSGLLEInitializePackage();CHKERRQ(ierr); 15363d177a5cSEmil Constantinescu 15373d177a5cSEmil Constantinescu ierr = PetscNewLog(ts,&gl);CHKERRQ(ierr); 15383d177a5cSEmil Constantinescu ts->data = (void*)gl; 15393d177a5cSEmil Constantinescu 15403d177a5cSEmil Constantinescu ts->ops->reset = TSReset_GLLE; 15413d177a5cSEmil Constantinescu ts->ops->destroy = TSDestroy_GLLE; 15423d177a5cSEmil Constantinescu ts->ops->view = TSView_GLLE; 15433d177a5cSEmil Constantinescu ts->ops->setup = TSSetUp_GLLE; 15443d177a5cSEmil Constantinescu ts->ops->solve = TSSolve_GLLE; 15453d177a5cSEmil Constantinescu ts->ops->setfromoptions = TSSetFromOptions_GLLE; 15463d177a5cSEmil Constantinescu ts->ops->snesfunction = SNESTSFormFunction_GLLE; 15473d177a5cSEmil Constantinescu ts->ops->snesjacobian = SNESTSFormJacobian_GLLE; 15483d177a5cSEmil Constantinescu 15493d177a5cSEmil Constantinescu gl->max_step_rejections = 1; 15503d177a5cSEmil Constantinescu gl->min_order = 1; 15513d177a5cSEmil Constantinescu gl->max_order = 3; 15523d177a5cSEmil Constantinescu gl->start_order = 1; 15533d177a5cSEmil Constantinescu gl->current_scheme = -1; 15543d177a5cSEmil Constantinescu gl->extrapolate = PETSC_FALSE; 15553d177a5cSEmil Constantinescu 15563d177a5cSEmil Constantinescu gl->wrms_atol = 1e-8; 15573d177a5cSEmil Constantinescu gl->wrms_rtol = 1e-5; 15583d177a5cSEmil Constantinescu 15593d177a5cSEmil Constantinescu ierr = PetscObjectComposeFunction((PetscObject)ts,"TSGLLESetType_C", &TSGLLESetType_GLLE);CHKERRQ(ierr); 15603d177a5cSEmil Constantinescu ierr = PetscObjectComposeFunction((PetscObject)ts,"TSGLLESetAcceptType_C",&TSGLLESetAcceptType_GLLE);CHKERRQ(ierr); 15613d177a5cSEmil Constantinescu ierr = PetscObjectComposeFunction((PetscObject)ts,"TSGLLEGetAdapt_C", &TSGLLEGetAdapt_GLLE);CHKERRQ(ierr); 15623d177a5cSEmil Constantinescu PetscFunctionReturn(0); 15633d177a5cSEmil Constantinescu } 1564