13d177a5cSEmil Constantinescu 23d177a5cSEmil Constantinescu #include <../src/ts/impls/implicit/glle/glle.h> /*I "petscts.h" I*/ 33d177a5cSEmil Constantinescu #include <petscdm.h> 43d177a5cSEmil Constantinescu #include <petscblaslapack.h> 53d177a5cSEmil Constantinescu 6c793f718SLisandro Dalcin static const char *TSGLLEErrorDirections[] = {"FORWARD","BACKWARD","TSGLLEErrorDirection","TSGLLEERROR_",NULL}; 73d177a5cSEmil Constantinescu static PetscFunctionList TSGLLEList; 83d177a5cSEmil Constantinescu static PetscFunctionList TSGLLEAcceptList; 93d177a5cSEmil Constantinescu static PetscBool TSGLLEPackageInitialized; 103d177a5cSEmil Constantinescu static PetscBool TSGLLERegisterAllCalled; 113d177a5cSEmil Constantinescu 123d177a5cSEmil Constantinescu /* This function is pure */ 133d177a5cSEmil Constantinescu static PetscScalar Factorial(PetscInt n) 143d177a5cSEmil Constantinescu { 153d177a5cSEmil Constantinescu PetscInt i; 163d177a5cSEmil Constantinescu if (n < 12) { /* Can compute with 32-bit integers */ 173d177a5cSEmil Constantinescu PetscInt f = 1; 183d177a5cSEmil Constantinescu for (i=2; i<=n; i++) f *= i; 193d177a5cSEmil Constantinescu return (PetscScalar)f; 203d177a5cSEmil Constantinescu } else { 213d177a5cSEmil Constantinescu PetscScalar f = 1.; 223d177a5cSEmil Constantinescu for (i=2; i<=n; i++) f *= (PetscScalar)i; 233d177a5cSEmil Constantinescu return f; 243d177a5cSEmil Constantinescu } 253d177a5cSEmil Constantinescu } 263d177a5cSEmil Constantinescu 273d177a5cSEmil Constantinescu /* This function is pure */ 283d177a5cSEmil Constantinescu static PetscScalar CPowF(PetscScalar c,PetscInt p) 293d177a5cSEmil Constantinescu { 303d177a5cSEmil Constantinescu return PetscPowRealInt(PetscRealPart(c),p)/Factorial(p); 313d177a5cSEmil Constantinescu } 323d177a5cSEmil Constantinescu 333d177a5cSEmil Constantinescu static PetscErrorCode TSGLLEGetVecs(TS ts,DM dm,Vec *Z,Vec *Ydotstage) 343d177a5cSEmil Constantinescu { 353d177a5cSEmil Constantinescu TS_GLLE *gl = (TS_GLLE*)ts->data; 363d177a5cSEmil Constantinescu 373d177a5cSEmil Constantinescu PetscFunctionBegin; 383d177a5cSEmil Constantinescu if (Z) { 393d177a5cSEmil Constantinescu if (dm && dm != ts->dm) { 409566063dSJacob Faibussowitsch PetscCall(DMGetNamedGlobalVector(dm,"TSGLLE_Z",Z)); 413d177a5cSEmil Constantinescu } else *Z = gl->Z; 423d177a5cSEmil Constantinescu } 433d177a5cSEmil Constantinescu if (Ydotstage) { 443d177a5cSEmil Constantinescu if (dm && dm != ts->dm) { 459566063dSJacob Faibussowitsch PetscCall(DMGetNamedGlobalVector(dm,"TSGLLE_Ydot",Ydotstage)); 463d177a5cSEmil Constantinescu } else *Ydotstage = gl->Ydot[gl->stage]; 473d177a5cSEmil Constantinescu } 483d177a5cSEmil Constantinescu PetscFunctionReturn(0); 493d177a5cSEmil Constantinescu } 503d177a5cSEmil Constantinescu 513d177a5cSEmil Constantinescu static PetscErrorCode TSGLLERestoreVecs(TS ts,DM dm,Vec *Z,Vec *Ydotstage) 523d177a5cSEmil Constantinescu { 533d177a5cSEmil Constantinescu PetscFunctionBegin; 543d177a5cSEmil Constantinescu if (Z) { 553d177a5cSEmil Constantinescu if (dm && dm != ts->dm) { 569566063dSJacob Faibussowitsch PetscCall(DMRestoreNamedGlobalVector(dm,"TSGLLE_Z",Z)); 573d177a5cSEmil Constantinescu } 583d177a5cSEmil Constantinescu } 593d177a5cSEmil Constantinescu if (Ydotstage) { 603d177a5cSEmil Constantinescu 613d177a5cSEmil Constantinescu if (dm && dm != ts->dm) { 629566063dSJacob Faibussowitsch PetscCall(DMRestoreNamedGlobalVector(dm,"TSGLLE_Ydot",Ydotstage)); 633d177a5cSEmil Constantinescu } 643d177a5cSEmil Constantinescu } 653d177a5cSEmil Constantinescu PetscFunctionReturn(0); 663d177a5cSEmil Constantinescu } 673d177a5cSEmil Constantinescu 683d177a5cSEmil Constantinescu static PetscErrorCode DMCoarsenHook_TSGLLE(DM fine,DM coarse,void *ctx) 693d177a5cSEmil Constantinescu { 703d177a5cSEmil Constantinescu PetscFunctionBegin; 713d177a5cSEmil Constantinescu PetscFunctionReturn(0); 723d177a5cSEmil Constantinescu } 733d177a5cSEmil Constantinescu 743d177a5cSEmil Constantinescu static PetscErrorCode DMRestrictHook_TSGLLE(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx) 753d177a5cSEmil Constantinescu { 763d177a5cSEmil Constantinescu TS ts = (TS)ctx; 773d177a5cSEmil Constantinescu Vec Ydot,Ydot_c; 783d177a5cSEmil Constantinescu 793d177a5cSEmil Constantinescu PetscFunctionBegin; 809566063dSJacob Faibussowitsch PetscCall(TSGLLEGetVecs(ts,fine,NULL,&Ydot)); 819566063dSJacob Faibussowitsch PetscCall(TSGLLEGetVecs(ts,coarse,NULL,&Ydot_c)); 829566063dSJacob Faibussowitsch PetscCall(MatRestrict(restrct,Ydot,Ydot_c)); 839566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(Ydot_c,rscale,Ydot_c)); 849566063dSJacob Faibussowitsch PetscCall(TSGLLERestoreVecs(ts,fine,NULL,&Ydot)); 859566063dSJacob Faibussowitsch PetscCall(TSGLLERestoreVecs(ts,coarse,NULL,&Ydot_c)); 863d177a5cSEmil Constantinescu PetscFunctionReturn(0); 873d177a5cSEmil Constantinescu } 883d177a5cSEmil Constantinescu 893d177a5cSEmil Constantinescu static PetscErrorCode DMSubDomainHook_TSGLLE(DM dm,DM subdm,void *ctx) 903d177a5cSEmil Constantinescu { 913d177a5cSEmil Constantinescu PetscFunctionBegin; 923d177a5cSEmil Constantinescu PetscFunctionReturn(0); 933d177a5cSEmil Constantinescu } 943d177a5cSEmil Constantinescu 953d177a5cSEmil Constantinescu static PetscErrorCode DMSubDomainRestrictHook_TSGLLE(DM dm,VecScatter gscat, VecScatter lscat,DM subdm,void *ctx) 963d177a5cSEmil Constantinescu { 973d177a5cSEmil Constantinescu TS ts = (TS)ctx; 983d177a5cSEmil Constantinescu Vec Ydot,Ydot_s; 993d177a5cSEmil Constantinescu 1003d177a5cSEmil Constantinescu PetscFunctionBegin; 1019566063dSJacob Faibussowitsch PetscCall(TSGLLEGetVecs(ts,dm,NULL,&Ydot)); 1029566063dSJacob Faibussowitsch PetscCall(TSGLLEGetVecs(ts,subdm,NULL,&Ydot_s)); 1033d177a5cSEmil Constantinescu 1049566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(gscat,Ydot,Ydot_s,INSERT_VALUES,SCATTER_FORWARD)); 1059566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(gscat,Ydot,Ydot_s,INSERT_VALUES,SCATTER_FORWARD)); 1063d177a5cSEmil Constantinescu 1079566063dSJacob Faibussowitsch PetscCall(TSGLLERestoreVecs(ts,dm,NULL,&Ydot)); 1089566063dSJacob Faibussowitsch PetscCall(TSGLLERestoreVecs(ts,subdm,NULL,&Ydot_s)); 1093d177a5cSEmil Constantinescu PetscFunctionReturn(0); 1103d177a5cSEmil Constantinescu } 1113d177a5cSEmil Constantinescu 1123d177a5cSEmil Constantinescu static PetscErrorCode TSGLLESchemeCreate(PetscInt p,PetscInt q,PetscInt r,PetscInt s,const PetscScalar *c, 1133d177a5cSEmil Constantinescu const PetscScalar *a,const PetscScalar *b,const PetscScalar *u,const PetscScalar *v,TSGLLEScheme *inscheme) 1143d177a5cSEmil Constantinescu { 1153d177a5cSEmil Constantinescu TSGLLEScheme scheme; 1163d177a5cSEmil Constantinescu PetscInt j; 1173d177a5cSEmil Constantinescu 1183d177a5cSEmil Constantinescu PetscFunctionBegin; 11908401ef6SPierre Jolivet PetscCheck(p >= 1,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Scheme order must be positive"); 12008401ef6SPierre Jolivet PetscCheck(r >= 1,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"At least one item must be carried between steps"); 12108401ef6SPierre Jolivet PetscCheck(s >= 1,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"At least one stage is required"); 122064a246eSJacob Faibussowitsch PetscValidPointer(inscheme,10); 123c793f718SLisandro Dalcin *inscheme = NULL; 1249566063dSJacob Faibussowitsch PetscCall(PetscNew(&scheme)); 1253d177a5cSEmil Constantinescu scheme->p = p; 1263d177a5cSEmil Constantinescu scheme->q = q; 1273d177a5cSEmil Constantinescu scheme->r = r; 1283d177a5cSEmil Constantinescu scheme->s = s; 1293d177a5cSEmil Constantinescu 1309566063dSJacob Faibussowitsch PetscCall(PetscMalloc5(s,&scheme->c,s*s,&scheme->a,r*s,&scheme->b,r*s,&scheme->u,r*r,&scheme->v)); 1319566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(scheme->c,c,s)); 1323d177a5cSEmil Constantinescu for (j=0; j<s*s; j++) scheme->a[j] = (PetscAbsScalar(a[j]) < 1e-12) ? 0 : a[j]; 1333d177a5cSEmil Constantinescu for (j=0; j<r*s; j++) scheme->b[j] = (PetscAbsScalar(b[j]) < 1e-12) ? 0 : b[j]; 1343d177a5cSEmil Constantinescu for (j=0; j<s*r; j++) scheme->u[j] = (PetscAbsScalar(u[j]) < 1e-12) ? 0 : u[j]; 1353d177a5cSEmil Constantinescu for (j=0; j<r*r; j++) scheme->v[j] = (PetscAbsScalar(v[j]) < 1e-12) ? 0 : v[j]; 1363d177a5cSEmil Constantinescu 1379566063dSJacob Faibussowitsch PetscCall(PetscMalloc6(r,&scheme->alpha,r,&scheme->beta,r,&scheme->gamma,3*s,&scheme->phi,3*r,&scheme->psi,r,&scheme->stage_error)); 1383d177a5cSEmil Constantinescu { 1393d177a5cSEmil Constantinescu PetscInt i,j,k,ss=s+2; 1403d177a5cSEmil Constantinescu PetscBLASInt m,n,one=1,*ipiv,lwork=4*((s+3)*3+3),info,ldb; 1413d177a5cSEmil Constantinescu PetscReal rcond,*sing,*workreal; 1423d177a5cSEmil Constantinescu PetscScalar *ImV,*H,*bmat,*workscalar,*c=scheme->c,*a=scheme->a,*b=scheme->b,*u=scheme->u,*v=scheme->v; 1433d177a5cSEmil Constantinescu PetscBLASInt rank; 1449566063dSJacob Faibussowitsch PetscCall(PetscMalloc7(PetscSqr(r),&ImV,3*s,&H,3*ss,&bmat,lwork,&workscalar,5*(3+r),&workreal,r+s,&sing,r+s,&ipiv)); 1453d177a5cSEmil Constantinescu 1463d177a5cSEmil Constantinescu /* column-major input */ 1473d177a5cSEmil Constantinescu for (i=0; i<r-1; i++) { 1483d177a5cSEmil Constantinescu for (j=0; j<r-1; j++) ImV[i+j*r] = 1.0*(i==j) - v[(i+1)*r+j+1]; 1493d177a5cSEmil Constantinescu } 1503d177a5cSEmil Constantinescu /* Build right hand side for alpha (tp - glm.B(2:end,:)*(glm.c.^(p)./factorial(p))) */ 1513d177a5cSEmil Constantinescu for (i=1; i<r; i++) { 1523d177a5cSEmil Constantinescu scheme->alpha[i] = 1./Factorial(p+1-i); 1533d177a5cSEmil Constantinescu for (j=0; j<s; j++) scheme->alpha[i] -= b[i*s+j]*CPowF(c[j],p); 1543d177a5cSEmil Constantinescu } 1559566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(r-1,&m)); 1569566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(r,&n)); 157792fecdfSBarry Smith PetscCallBLAS("LAPACKgesv",LAPACKgesv_(&m,&one,ImV,&n,ipiv,scheme->alpha+1,&n,&info)); 15808401ef6SPierre Jolivet PetscCheck(info >= 0,PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad argument to GESV"); 15908401ef6SPierre Jolivet PetscCheck(info <= 0,PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Bad LU factorization"); 1603d177a5cSEmil Constantinescu 1613d177a5cSEmil Constantinescu /* Build right hand side for beta (tp1 - glm.B(2:end,:)*(glm.c.^(p+1)./factorial(p+1)) - e.alpha) */ 1623d177a5cSEmil Constantinescu for (i=1; i<r; i++) { 1633d177a5cSEmil Constantinescu scheme->beta[i] = 1./Factorial(p+2-i) - scheme->alpha[i]; 1643d177a5cSEmil Constantinescu for (j=0; j<s; j++) scheme->beta[i] -= b[i*s+j]*CPowF(c[j],p+1); 1653d177a5cSEmil Constantinescu } 166792fecdfSBarry Smith PetscCallBLAS("LAPACKgetrs",LAPACKgetrs_("No transpose",&m,&one,ImV,&n,ipiv,scheme->beta+1,&n,&info)); 16708401ef6SPierre Jolivet PetscCheck(info >= 0,PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad argument to GETRS"); 16808401ef6SPierre Jolivet PetscCheck(info <= 0,PETSC_COMM_SELF,PETSC_ERR_LIB,"Should not happen"); 1693d177a5cSEmil Constantinescu 1703d177a5cSEmil Constantinescu /* Build stage_error vector 1713d177a5cSEmil Constantinescu xi = glm.c.^(p+1)/factorial(p+1) - glm.A*glm.c.^p/factorial(p) + glm.U(:,2:end)*e.alpha; 1723d177a5cSEmil Constantinescu */ 1733d177a5cSEmil Constantinescu for (i=0; i<s; i++) { 1743d177a5cSEmil Constantinescu scheme->stage_error[i] = CPowF(c[i],p+1); 1753d177a5cSEmil Constantinescu for (j=0; j<s; j++) scheme->stage_error[i] -= a[i*s+j]*CPowF(c[j],p); 1763d177a5cSEmil Constantinescu for (j=1; j<r; j++) scheme->stage_error[i] += u[i*r+j]*scheme->alpha[j]; 1773d177a5cSEmil Constantinescu } 1783d177a5cSEmil Constantinescu 1793d177a5cSEmil Constantinescu /* alpha[0] (epsilon in B,J,W 2007) 1803d177a5cSEmil Constantinescu epsilon = 1/factorial(p+1) - B(1,:)*c.^p/factorial(p) + V(1,2:end)*e.alpha; 1813d177a5cSEmil Constantinescu */ 1823d177a5cSEmil Constantinescu scheme->alpha[0] = 1./Factorial(p+1); 1833d177a5cSEmil Constantinescu for (j=0; j<s; j++) scheme->alpha[0] -= b[0*s+j]*CPowF(c[j],p); 1843d177a5cSEmil Constantinescu for (j=1; j<r; j++) scheme->alpha[0] += v[0*r+j]*scheme->alpha[j]; 1853d177a5cSEmil Constantinescu 1863d177a5cSEmil Constantinescu /* right hand side for gamma (glm.B(2:end,:)*e.xi - e.epsilon*eye(s-1,1)) */ 1873d177a5cSEmil Constantinescu for (i=1; i<r; i++) { 1883d177a5cSEmil Constantinescu scheme->gamma[i] = (i==1 ? -1. : 0)*scheme->alpha[0]; 1893d177a5cSEmil Constantinescu for (j=0; j<s; j++) scheme->gamma[i] += b[i*s+j]*scheme->stage_error[j]; 1903d177a5cSEmil Constantinescu } 191792fecdfSBarry Smith PetscCallBLAS("LAPACKgetrs",LAPACKgetrs_("No transpose",&m,&one,ImV,&n,ipiv,scheme->gamma+1,&n,&info)); 19208401ef6SPierre Jolivet PetscCheck(info >= 0,PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad argument to GETRS"); 19308401ef6SPierre Jolivet PetscCheck(info <= 0,PETSC_COMM_SELF,PETSC_ERR_LIB,"Should not happen"); 1943d177a5cSEmil Constantinescu 1953d177a5cSEmil Constantinescu /* beta[0] (rho in B,J,W 2007) 1963d177a5cSEmil Constantinescu e.rho = 1/factorial(p+2) - glm.B(1,:)*glm.c.^(p+1)/factorial(p+1) ... 1973d177a5cSEmil Constantinescu + glm.V(1,2:end)*e.beta;% - e.epsilon; 1983d177a5cSEmil Constantinescu % Note: The paper (B,J,W 2007) includes the last term in their definition 1993d177a5cSEmil Constantinescu * */ 2003d177a5cSEmil Constantinescu scheme->beta[0] = 1./Factorial(p+2); 2013d177a5cSEmil Constantinescu for (j=0; j<s; j++) scheme->beta[0] -= b[0*s+j]*CPowF(c[j],p+1); 2023d177a5cSEmil Constantinescu for (j=1; j<r; j++) scheme->beta[0] += v[0*r+j]*scheme->beta[j]; 2033d177a5cSEmil Constantinescu 2043d177a5cSEmil Constantinescu /* gamma[0] (sigma in B,J,W 2007) 2053d177a5cSEmil Constantinescu * e.sigma = glm.B(1,:)*e.xi + glm.V(1,2:end)*e.gamma; 2063d177a5cSEmil Constantinescu * */ 2073d177a5cSEmil Constantinescu scheme->gamma[0] = 0.0; 2083d177a5cSEmil Constantinescu for (j=0; j<s; j++) scheme->gamma[0] += b[0*s+j]*scheme->stage_error[j]; 2093d177a5cSEmil Constantinescu for (j=1; j<r; j++) scheme->gamma[0] += v[0*s+j]*scheme->gamma[j]; 2103d177a5cSEmil Constantinescu 2113d177a5cSEmil Constantinescu /* Assemble H 21263a3b9bcSJacob Faibussowitsch * % " PetscInt_FMT "etermine the error estimators phi 2133d177a5cSEmil Constantinescu H = [[cpow(glm.c,p) + C*e.alpha] [cpow(glm.c,p+1) + C*e.beta] ... 2143d177a5cSEmil Constantinescu [e.xi - C*(e.gamma + 0*e.epsilon*eye(s-1,1))]]'; 2153d177a5cSEmil Constantinescu % Paper has formula above without the 0, but that term must be left 2163d177a5cSEmil Constantinescu % out to satisfy the conditions they propose and to make the 2173d177a5cSEmil Constantinescu % example schemes work 2183d177a5cSEmil Constantinescu e.H = H; 2193d177a5cSEmil Constantinescu e.phi = (H \ [1 0 0;1 1 0;0 0 -1])'; 2203d177a5cSEmil Constantinescu e.psi = -e.phi*C; 2213d177a5cSEmil Constantinescu * */ 2223d177a5cSEmil Constantinescu for (j=0; j<s; j++) { 2233d177a5cSEmil Constantinescu H[0+j*3] = CPowF(c[j],p); 2243d177a5cSEmil Constantinescu H[1+j*3] = CPowF(c[j],p+1); 2253d177a5cSEmil Constantinescu H[2+j*3] = scheme->stage_error[j]; 2263d177a5cSEmil Constantinescu for (k=1; k<r; k++) { 2273d177a5cSEmil Constantinescu H[0+j*3] += CPowF(c[j],k-1)*scheme->alpha[k]; 2283d177a5cSEmil Constantinescu H[1+j*3] += CPowF(c[j],k-1)*scheme->beta[k]; 2293d177a5cSEmil Constantinescu H[2+j*3] -= CPowF(c[j],k-1)*scheme->gamma[k]; 2303d177a5cSEmil Constantinescu } 2313d177a5cSEmil Constantinescu } 2323d177a5cSEmil Constantinescu bmat[0+0*ss] = 1.; bmat[0+1*ss] = 0.; bmat[0+2*ss] = 0.; 2333d177a5cSEmil Constantinescu bmat[1+0*ss] = 1.; bmat[1+1*ss] = 1.; bmat[1+2*ss] = 0.; 2343d177a5cSEmil Constantinescu bmat[2+0*ss] = 0.; bmat[2+1*ss] = 0.; bmat[2+2*ss] = -1.; 2353d177a5cSEmil Constantinescu m = 3; 2369566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(s,&n)); 2379566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(ss,&ldb)); 2383d177a5cSEmil Constantinescu rcond = 1e-12; 2393d177a5cSEmil Constantinescu #if defined(PETSC_USE_COMPLEX) 2403d177a5cSEmil Constantinescu /* ZGELSS( M, N, NRHS, A, LDA, B, LDB, S, RCOND, RANK, WORK, LWORK, RWORK, INFO) */ 241792fecdfSBarry Smith PetscCallBLAS("LAPACKgelss",LAPACKgelss_(&m,&n,&m,H,&m,bmat,&ldb,sing,&rcond,&rank,workscalar,&lwork,workreal,&info)); 2423d177a5cSEmil Constantinescu #else 2433d177a5cSEmil Constantinescu /* DGELSS( M, N, NRHS, A, LDA, B, LDB, S, RCOND, RANK, WORK, LWORK, INFO) */ 244792fecdfSBarry Smith PetscCallBLAS("LAPACKgelss",LAPACKgelss_(&m,&n,&m,H,&m,bmat,&ldb,sing,&rcond,&rank,workscalar,&lwork,&info)); 2453d177a5cSEmil Constantinescu #endif 24608401ef6SPierre Jolivet PetscCheck(info >= 0,PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad argument to GELSS"); 24708401ef6SPierre Jolivet PetscCheck(info <= 0,PETSC_COMM_SELF,PETSC_ERR_LIB,"SVD failed to converge"); 2483d177a5cSEmil Constantinescu 2493d177a5cSEmil Constantinescu for (j=0; j<3; j++) { 2503d177a5cSEmil Constantinescu for (k=0; k<s; k++) scheme->phi[k+j*s] = bmat[k+j*ss]; 2513d177a5cSEmil Constantinescu } 2523d177a5cSEmil Constantinescu 2533d177a5cSEmil Constantinescu /* the other part of the error estimator, psi in B,J,W 2007 */ 2543d177a5cSEmil Constantinescu scheme->psi[0*r+0] = 0.; 2553d177a5cSEmil Constantinescu scheme->psi[1*r+0] = 0.; 2563d177a5cSEmil Constantinescu scheme->psi[2*r+0] = 0.; 2573d177a5cSEmil Constantinescu for (j=1; j<r; j++) { 2583d177a5cSEmil Constantinescu scheme->psi[0*r+j] = 0.; 2593d177a5cSEmil Constantinescu scheme->psi[1*r+j] = 0.; 2603d177a5cSEmil Constantinescu scheme->psi[2*r+j] = 0.; 2613d177a5cSEmil Constantinescu for (k=0; k<s; k++) { 2623d177a5cSEmil Constantinescu scheme->psi[0*r+j] -= CPowF(c[k],j-1)*scheme->phi[0*s+k]; 2633d177a5cSEmil Constantinescu scheme->psi[1*r+j] -= CPowF(c[k],j-1)*scheme->phi[1*s+k]; 2643d177a5cSEmil Constantinescu scheme->psi[2*r+j] -= CPowF(c[k],j-1)*scheme->phi[2*s+k]; 2653d177a5cSEmil Constantinescu } 2663d177a5cSEmil Constantinescu } 2679566063dSJacob Faibussowitsch PetscCall(PetscFree7(ImV,H,bmat,workscalar,workreal,sing,ipiv)); 2683d177a5cSEmil Constantinescu } 2693d177a5cSEmil Constantinescu /* Check which properties are satisfied */ 2703d177a5cSEmil Constantinescu scheme->stiffly_accurate = PETSC_TRUE; 2713d177a5cSEmil Constantinescu if (scheme->c[s-1] != 1.) scheme->stiffly_accurate = PETSC_FALSE; 2723d177a5cSEmil Constantinescu for (j=0; j<s; j++) if (a[(s-1)*s+j] != b[j]) scheme->stiffly_accurate = PETSC_FALSE; 2733d177a5cSEmil Constantinescu for (j=0; j<r; j++) if (u[(s-1)*r+j] != v[j]) scheme->stiffly_accurate = PETSC_FALSE; 2743d177a5cSEmil Constantinescu scheme->fsal = scheme->stiffly_accurate; /* FSAL is stronger */ 2753d177a5cSEmil Constantinescu for (j=0; j<s-1; j++) if (r>1 && b[1*s+j] != 0.) scheme->fsal = PETSC_FALSE; 2763d177a5cSEmil Constantinescu if (b[1*s+r-1] != 1.) scheme->fsal = PETSC_FALSE; 2773d177a5cSEmil Constantinescu for (j=0; j<r; j++) if (r>1 && v[1*r+j] != 0.) scheme->fsal = PETSC_FALSE; 2783d177a5cSEmil Constantinescu 2793d177a5cSEmil Constantinescu *inscheme = scheme; 2803d177a5cSEmil Constantinescu PetscFunctionReturn(0); 2813d177a5cSEmil Constantinescu } 2823d177a5cSEmil Constantinescu 2833d177a5cSEmil Constantinescu static PetscErrorCode TSGLLESchemeDestroy(TSGLLEScheme sc) 2843d177a5cSEmil Constantinescu { 2853d177a5cSEmil Constantinescu PetscFunctionBegin; 2869566063dSJacob Faibussowitsch PetscCall(PetscFree5(sc->c,sc->a,sc->b,sc->u,sc->v)); 2879566063dSJacob Faibussowitsch PetscCall(PetscFree6(sc->alpha,sc->beta,sc->gamma,sc->phi,sc->psi,sc->stage_error)); 2889566063dSJacob Faibussowitsch PetscCall(PetscFree(sc)); 2893d177a5cSEmil Constantinescu PetscFunctionReturn(0); 2903d177a5cSEmil Constantinescu } 2913d177a5cSEmil Constantinescu 2923d177a5cSEmil Constantinescu static PetscErrorCode TSGLLEDestroy_Default(TS_GLLE *gl) 2933d177a5cSEmil Constantinescu { 2943d177a5cSEmil Constantinescu PetscInt i; 2953d177a5cSEmil Constantinescu 2963d177a5cSEmil Constantinescu PetscFunctionBegin; 2973d177a5cSEmil Constantinescu for (i=0; i<gl->nschemes; i++) { 2989566063dSJacob Faibussowitsch if (gl->schemes[i]) PetscCall(TSGLLESchemeDestroy(gl->schemes[i])); 2993d177a5cSEmil Constantinescu } 3009566063dSJacob Faibussowitsch PetscCall(PetscFree(gl->schemes)); 3013d177a5cSEmil Constantinescu gl->nschemes = 0; 3029566063dSJacob Faibussowitsch PetscCall(PetscMemzero(gl->type_name,sizeof(gl->type_name))); 3033d177a5cSEmil Constantinescu PetscFunctionReturn(0); 3043d177a5cSEmil Constantinescu } 3053d177a5cSEmil Constantinescu 3063d177a5cSEmil Constantinescu static PetscErrorCode TSGLLEViewTable_Private(PetscViewer viewer,PetscInt m,PetscInt n,const PetscScalar a[],const char name[]) 3073d177a5cSEmil Constantinescu { 3083d177a5cSEmil Constantinescu PetscBool iascii; 3093d177a5cSEmil Constantinescu PetscInt i,j; 3103d177a5cSEmil Constantinescu 3113d177a5cSEmil Constantinescu PetscFunctionBegin; 3129566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii)); 3133d177a5cSEmil Constantinescu if (iascii) { 3149566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"%30s = [",name)); 3153d177a5cSEmil Constantinescu for (i=0; i<m; i++) { 3169566063dSJacob Faibussowitsch if (i) PetscCall(PetscViewerASCIIPrintf(viewer,"%30s [","")); 3179566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE)); 3183d177a5cSEmil Constantinescu for (j=0; j<n; j++) { 31963a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," %12.8g",(double)PetscRealPart(a[i*n+j]))); 3203d177a5cSEmil Constantinescu } 3219566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"]\n")); 3229566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE)); 3233d177a5cSEmil Constantinescu } 3243d177a5cSEmil Constantinescu } 3253d177a5cSEmil Constantinescu PetscFunctionReturn(0); 3263d177a5cSEmil Constantinescu } 3273d177a5cSEmil Constantinescu 3283d177a5cSEmil Constantinescu static PetscErrorCode TSGLLESchemeView(TSGLLEScheme sc,PetscBool view_details,PetscViewer viewer) 3293d177a5cSEmil Constantinescu { 3303d177a5cSEmil Constantinescu PetscBool iascii; 3313d177a5cSEmil Constantinescu 3323d177a5cSEmil Constantinescu PetscFunctionBegin; 3339566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii)); 3343d177a5cSEmil Constantinescu if (iascii) { 33563a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"GL scheme p,q,r,s = %" PetscInt_FMT ",%" PetscInt_FMT ",%" PetscInt_FMT ",%" PetscInt_FMT "\n",sc->p,sc->q,sc->r,sc->s)); 3369566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 3379566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Stiffly accurate: %s, FSAL: %s\n",sc->stiffly_accurate ? "yes" : "no",sc->fsal ? "yes" : "no")); 3389566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Leading error constants: %10.3e %10.3e %10.3e\n", 33963a3b9bcSJacob Faibussowitsch (double)PetscRealPart(sc->alpha[0]),(double)PetscRealPart(sc->beta[0]),(double)PetscRealPart(sc->gamma[0]))); 3409566063dSJacob Faibussowitsch PetscCall(TSGLLEViewTable_Private(viewer,1,sc->s,sc->c,"Abscissas c")); 3413d177a5cSEmil Constantinescu if (view_details) { 3429566063dSJacob Faibussowitsch PetscCall(TSGLLEViewTable_Private(viewer,sc->s,sc->s,sc->a,"A")); 3439566063dSJacob Faibussowitsch PetscCall(TSGLLEViewTable_Private(viewer,sc->r,sc->s,sc->b,"B")); 3449566063dSJacob Faibussowitsch PetscCall(TSGLLEViewTable_Private(viewer,sc->s,sc->r,sc->u,"U")); 3459566063dSJacob Faibussowitsch PetscCall(TSGLLEViewTable_Private(viewer,sc->r,sc->r,sc->v,"V")); 3463d177a5cSEmil Constantinescu 3479566063dSJacob Faibussowitsch PetscCall(TSGLLEViewTable_Private(viewer,3,sc->s,sc->phi,"Error estimate phi")); 3489566063dSJacob Faibussowitsch PetscCall(TSGLLEViewTable_Private(viewer,3,sc->r,sc->psi,"Error estimate psi")); 3499566063dSJacob Faibussowitsch PetscCall(TSGLLEViewTable_Private(viewer,1,sc->r,sc->alpha,"Modify alpha")); 3509566063dSJacob Faibussowitsch PetscCall(TSGLLEViewTable_Private(viewer,1,sc->r,sc->beta,"Modify beta")); 3519566063dSJacob Faibussowitsch PetscCall(TSGLLEViewTable_Private(viewer,1,sc->r,sc->gamma,"Modify gamma")); 3529566063dSJacob Faibussowitsch PetscCall(TSGLLEViewTable_Private(viewer,1,sc->s,sc->stage_error,"Stage error xi")); 3533d177a5cSEmil Constantinescu } 3549566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 35598921bdaSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported",((PetscObject)viewer)->type_name); 3563d177a5cSEmil Constantinescu PetscFunctionReturn(0); 3573d177a5cSEmil Constantinescu } 3583d177a5cSEmil Constantinescu 3593d177a5cSEmil Constantinescu static PetscErrorCode TSGLLEEstimateHigherMoments_Default(TSGLLEScheme sc,PetscReal h,Vec Ydot[],Vec Xold[],Vec hm[]) 3603d177a5cSEmil Constantinescu { 3613d177a5cSEmil Constantinescu PetscInt i; 3623d177a5cSEmil Constantinescu 3633d177a5cSEmil Constantinescu PetscFunctionBegin; 364cad9d221SBarry Smith PetscCheck(sc->r <= 64 && sc->s <= 64,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Ridiculous number of stages or items passed between stages"); 3653d177a5cSEmil Constantinescu /* build error vectors*/ 3663d177a5cSEmil Constantinescu for (i=0; i<3; i++) { 3673d177a5cSEmil Constantinescu PetscScalar phih[64]; 3683d177a5cSEmil Constantinescu PetscInt j; 3693d177a5cSEmil Constantinescu for (j=0; j<sc->s; j++) phih[j] = sc->phi[i*sc->s+j]*h; 3709566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(hm[i])); 3719566063dSJacob Faibussowitsch PetscCall(VecMAXPY(hm[i],sc->s,phih,Ydot)); 3729566063dSJacob Faibussowitsch PetscCall(VecMAXPY(hm[i],sc->r,&sc->psi[i*sc->r],Xold)); 3733d177a5cSEmil Constantinescu } 3743d177a5cSEmil Constantinescu PetscFunctionReturn(0); 3753d177a5cSEmil Constantinescu } 3763d177a5cSEmil Constantinescu 3773d177a5cSEmil Constantinescu static PetscErrorCode TSGLLECompleteStep_Rescale(TSGLLEScheme sc,PetscReal h,TSGLLEScheme next_sc,PetscReal next_h,Vec Ydot[],Vec Xold[],Vec X[]) 3783d177a5cSEmil Constantinescu { 3793d177a5cSEmil Constantinescu PetscScalar brow[32],vrow[32]; 3803d177a5cSEmil Constantinescu PetscInt i,j,r,s; 3813d177a5cSEmil Constantinescu 3823d177a5cSEmil Constantinescu PetscFunctionBegin; 3833d177a5cSEmil Constantinescu /* Build the new solution from (X,Ydot) */ 3843d177a5cSEmil Constantinescu r = sc->r; 3853d177a5cSEmil Constantinescu s = sc->s; 3863d177a5cSEmil Constantinescu for (i=0; i<r; i++) { 3879566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(X[i])); 3883d177a5cSEmil Constantinescu for (j=0; j<s; j++) brow[j] = h*sc->b[i*s+j]; 3899566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X[i],s,brow,Ydot)); 3903d177a5cSEmil Constantinescu for (j=0; j<r; j++) vrow[j] = sc->v[i*r+j]; 3919566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X[i],r,vrow,Xold)); 3923d177a5cSEmil Constantinescu } 3933d177a5cSEmil Constantinescu PetscFunctionReturn(0); 3943d177a5cSEmil Constantinescu } 3953d177a5cSEmil Constantinescu 3963d177a5cSEmil Constantinescu static PetscErrorCode TSGLLECompleteStep_RescaleAndModify(TSGLLEScheme sc,PetscReal h,TSGLLEScheme next_sc,PetscReal next_h,Vec Ydot[],Vec Xold[],Vec X[]) 3973d177a5cSEmil Constantinescu { 3983d177a5cSEmil Constantinescu PetscScalar brow[32],vrow[32]; 3993d177a5cSEmil Constantinescu PetscReal ratio; 4003d177a5cSEmil Constantinescu PetscInt i,j,p,r,s; 4013d177a5cSEmil Constantinescu 4023d177a5cSEmil Constantinescu PetscFunctionBegin; 4033d177a5cSEmil Constantinescu /* Build the new solution from (X,Ydot) */ 4043d177a5cSEmil Constantinescu p = sc->p; 4053d177a5cSEmil Constantinescu r = sc->r; 4063d177a5cSEmil Constantinescu s = sc->s; 4073d177a5cSEmil Constantinescu ratio = next_h/h; 4083d177a5cSEmil Constantinescu for (i=0; i<r; i++) { 4099566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(X[i])); 4103d177a5cSEmil Constantinescu for (j=0; j<s; j++) { 4113d177a5cSEmil Constantinescu brow[j] = h*(PetscPowRealInt(ratio,i)*sc->b[i*s+j] 4123d177a5cSEmil Constantinescu + (PetscPowRealInt(ratio,i) - PetscPowRealInt(ratio,p+1))*(+ sc->alpha[i]*sc->phi[0*s+j]) 4133d177a5cSEmil Constantinescu + (PetscPowRealInt(ratio,i) - PetscPowRealInt(ratio,p+2))*(+ sc->beta [i]*sc->phi[1*s+j] 4143d177a5cSEmil Constantinescu + sc->gamma[i]*sc->phi[2*s+j])); 4153d177a5cSEmil Constantinescu } 4169566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X[i],s,brow,Ydot)); 4173d177a5cSEmil Constantinescu for (j=0; j<r; j++) { 4183d177a5cSEmil Constantinescu vrow[j] = (PetscPowRealInt(ratio,i)*sc->v[i*r+j] 4193d177a5cSEmil Constantinescu + (PetscPowRealInt(ratio,i) - PetscPowRealInt(ratio,p+1))*(+ sc->alpha[i]*sc->psi[0*r+j]) 4203d177a5cSEmil Constantinescu + (PetscPowRealInt(ratio,i) - PetscPowRealInt(ratio,p+2))*(+ sc->beta [i]*sc->psi[1*r+j] 4213d177a5cSEmil Constantinescu + sc->gamma[i]*sc->psi[2*r+j])); 4223d177a5cSEmil Constantinescu } 4239566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X[i],r,vrow,Xold)); 4243d177a5cSEmil Constantinescu } 4253d177a5cSEmil Constantinescu if (r < next_sc->r) { 42608401ef6SPierre Jolivet PetscCheck(r+1 == next_sc->r,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot accommodate jump in r greater than 1"); 4279566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(X[r])); 4283d177a5cSEmil Constantinescu for (j=0; j<s; j++) brow[j] = h*PetscPowRealInt(ratio,p+1)*sc->phi[0*s+j]; 4299566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X[r],s,brow,Ydot)); 4303d177a5cSEmil Constantinescu for (j=0; j<r; j++) vrow[j] = PetscPowRealInt(ratio,p+1)*sc->psi[0*r+j]; 4319566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X[r],r,vrow,Xold)); 4323d177a5cSEmil Constantinescu } 4333d177a5cSEmil Constantinescu PetscFunctionReturn(0); 4343d177a5cSEmil Constantinescu } 4353d177a5cSEmil Constantinescu 4363d177a5cSEmil Constantinescu static PetscErrorCode TSGLLECreate_IRKS(TS ts) 4373d177a5cSEmil Constantinescu { 4383d177a5cSEmil Constantinescu TS_GLLE *gl = (TS_GLLE*)ts->data; 4393d177a5cSEmil Constantinescu 4403d177a5cSEmil Constantinescu PetscFunctionBegin; 4413d177a5cSEmil Constantinescu gl->Destroy = TSGLLEDestroy_Default; 4423d177a5cSEmil Constantinescu gl->EstimateHigherMoments = TSGLLEEstimateHigherMoments_Default; 4433d177a5cSEmil Constantinescu gl->CompleteStep = TSGLLECompleteStep_RescaleAndModify; 4449566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(10,&gl->schemes)); 4453d177a5cSEmil Constantinescu gl->nschemes = 0; 4463d177a5cSEmil Constantinescu 4473d177a5cSEmil Constantinescu { 4483d177a5cSEmil Constantinescu /* p=1,q=1, r=s=2, A- and L-stable with error estimates of order 2 and 3 4493d177a5cSEmil Constantinescu * Listed in Butcher & Podhaisky 2006. On error estimation in general linear methods for stiff ODE. 4503d177a5cSEmil Constantinescu * irks(0.3,0,[.3,1],[1],1) 4513d177a5cSEmil Constantinescu * Note: can be made to have classical order (not stage order) 2 by replacing 0.3 with 1-sqrt(1/2) 4523d177a5cSEmil Constantinescu * but doing so would sacrifice the error estimator. 4533d177a5cSEmil Constantinescu */ 4543d177a5cSEmil Constantinescu const PetscScalar c[2] = {3./10., 1.}; 4553d177a5cSEmil Constantinescu const PetscScalar a[2][2] = {{3./10., 0}, {7./10., 3./10.}}; 4563d177a5cSEmil Constantinescu const PetscScalar b[2][2] = {{7./10., 3./10.}, {0,1}}; 4573d177a5cSEmil Constantinescu const PetscScalar u[2][2] = {{1,0},{1,0}}; 4583d177a5cSEmil Constantinescu const PetscScalar v[2][2] = {{1,0},{0,0}}; 4599566063dSJacob Faibussowitsch PetscCall(TSGLLESchemeCreate(1,1,2,2,c,*a,*b,*u,*v,&gl->schemes[gl->nschemes++])); 4603d177a5cSEmil Constantinescu } 4613d177a5cSEmil Constantinescu 4623d177a5cSEmil Constantinescu { 4633d177a5cSEmil Constantinescu /* p=q=2, r=s=3: irks(4/9,0,[1:3]/3,[0.33852],1) */ 4643d177a5cSEmil Constantinescu /* http://www.math.auckland.ac.nz/~hpod/atlas/i2a.html */ 4653d177a5cSEmil Constantinescu const PetscScalar c[3] = {1./3., 2./3., 1} 4663d177a5cSEmil Constantinescu ,a[3][3] = {{4./9. ,0 , 0}, 4673d177a5cSEmil Constantinescu {1.03750643704090e+00 , 4./9., 0}, 4683d177a5cSEmil Constantinescu {7.67024779410304e-01 , -3.81140216918943e-01, 4./9.}} 4693d177a5cSEmil Constantinescu ,b[3][3] = {{0.767024779410304, -0.381140216918943, 4./9.}, 4703d177a5cSEmil Constantinescu {0.000000000000000, 0.000000000000000, 1.000000000000000}, 4713d177a5cSEmil Constantinescu {-2.075048385225385, 0.621728385225383, 1.277197204924873}} 4723d177a5cSEmil Constantinescu ,u[3][3] = {{1.0000000000000000, -0.1111111111111109, -0.0925925925925922}, 4733d177a5cSEmil Constantinescu {1.0000000000000000, -0.8152842148186744, -0.4199095530877056}, 4743d177a5cSEmil Constantinescu {1.0000000000000000, 0.1696709930641948, 0.0539741070314165}} 4753d177a5cSEmil Constantinescu ,v[3][3] = {{1.0000000000000000, 0.1696709930641948, 0.0539741070314165}, 4763d177a5cSEmil Constantinescu {0.000000000000000, 0.000000000000000, 0.000000000000000}, 4773d177a5cSEmil Constantinescu {0.000000000000000, 0.176122795075129, 0.000000000000000}}; 4789566063dSJacob Faibussowitsch PetscCall(TSGLLESchemeCreate(2,2,3,3,c,*a,*b,*u,*v,&gl->schemes[gl->nschemes++])); 4793d177a5cSEmil Constantinescu } 4803d177a5cSEmil Constantinescu { 4813d177a5cSEmil Constantinescu /* p=q=3, r=s=4: irks(9/40,0,[1:4]/4,[0.3312 1.0050],[0.49541 1;1 0]) */ 4823d177a5cSEmil Constantinescu const PetscScalar c[4] = {0.25,0.5,0.75,1.0} 4833d177a5cSEmil Constantinescu ,a[4][4] = {{9./40. , 0, 0, 0}, 4843d177a5cSEmil Constantinescu {2.11286958887701e-01 , 9./40. , 0, 0}, 4853d177a5cSEmil Constantinescu {9.46338294287584e-01 , -3.42942861246094e-01, 9./40. , 0}, 4863d177a5cSEmil Constantinescu {0.521490453970721 , -0.662474225622980, 0.490476425459734, 9./40. }} 4873d177a5cSEmil Constantinescu ,b[4][4] = {{0.521490453970721 , -0.662474225622980, 0.490476425459734, 9./40. }, 4883d177a5cSEmil Constantinescu {0.000000000000000 , 0.000000000000000, 0.000000000000000, 1.000000000000000}, 4893d177a5cSEmil Constantinescu {-0.084677029310348 , 1.390757514776085, -1.568157386206001, 2.023192696767826}, 4903d177a5cSEmil Constantinescu {0.465383797936408 , 1.478273530625148, -1.930836081010182, 1.644872111193354}} 4913d177a5cSEmil Constantinescu ,u[4][4] = {{1.00000000000000000 , 0.02500000000001035, -0.02499999999999053, -0.00442708333332865}, 4923d177a5cSEmil Constantinescu {1.00000000000000000 , 0.06371304111232945, -0.04032173972189845, -0.01389438413189452}, 4933d177a5cSEmil Constantinescu {1.00000000000000000 , -0.07839543304147778, 0.04738685705116663, 0.02032603595928376}, 4943d177a5cSEmil Constantinescu {1.00000000000000000 , 0.42550734619251651, 0.10800718022400080, -0.01726712647760034}} 4953d177a5cSEmil Constantinescu ,v[4][4] = {{1.00000000000000000 , 0.42550734619251651, 0.10800718022400080, -0.01726712647760034}, 4963d177a5cSEmil Constantinescu {0.000000000000000 , 0.000000000000000, 0.000000000000000, 0.000000000000000}, 4973d177a5cSEmil Constantinescu {0.000000000000000 , -1.761115796027561, -0.521284157173780, 0.258249384305463}, 4983d177a5cSEmil Constantinescu {0.000000000000000 , -1.657693358744728, -1.052227765232394, 0.521284157173780}}; 4999566063dSJacob Faibussowitsch PetscCall(TSGLLESchemeCreate(3,3,4,4,c,*a,*b,*u,*v,&gl->schemes[gl->nschemes++])); 5003d177a5cSEmil Constantinescu } 5013d177a5cSEmil Constantinescu { 5023d177a5cSEmil Constantinescu /* p=q=4, r=s=5: 5033d177a5cSEmil Constantinescu irks(3/11,0,[1:5]/5, [0.1715 -0.1238 0.6617],... 5043d177a5cSEmil Constantinescu [ -0.0812 0.4079 1.0000 5053d177a5cSEmil Constantinescu 1.0000 0 0 5063d177a5cSEmil Constantinescu 0.8270 1.0000 0]) 5073d177a5cSEmil Constantinescu */ 5083d177a5cSEmil Constantinescu const PetscScalar c[5] = {0.2,0.4,0.6,0.8,1.0} 5093d177a5cSEmil Constantinescu ,a[5][5] = {{2.72727272727352e-01 , 0.00000000000000e+00, 0.00000000000000e+00 , 0.00000000000000e+00 , 0.00000000000000e+00}, 5103d177a5cSEmil Constantinescu {-1.03980153733431e-01, 2.72727272727405e-01, 0.00000000000000e+00, 0.00000000000000e+00 , 0.00000000000000e+00}, 5113d177a5cSEmil Constantinescu {-1.58615400341492e+00, 7.44168951881122e-01, 2.72727272727309e-01, 0.00000000000000e+00 , 0.00000000000000e+00}, 5123d177a5cSEmil Constantinescu {-8.73658042865628e-01, 5.37884671894595e-01, -1.63298538799523e-01, 2.72727272726996e-01 , 0.00000000000000e+00}, 5133d177a5cSEmil Constantinescu {2.95489397443992e-01 , -1.18481693910097e+00 , -6.68029812659953e-01 , 1.00716687860943e+00 , 2.72727272727288e-01}} 5143d177a5cSEmil Constantinescu ,b[5][5] = {{2.95489397443992e-01 , -1.18481693910097e+00 , -6.68029812659953e-01 , 1.00716687860943e+00 , 2.72727272727288e-01}, 5153d177a5cSEmil Constantinescu {0.00000000000000e+00 , 1.11022302462516e-16 , -2.22044604925031e-16 , 0.00000000000000e+00 , 1.00000000000000e+00}, 5163d177a5cSEmil Constantinescu {-4.05882503986005e+00, -4.00924006567769e+00, -1.38930610972481e+00, 4.45223930308488e+00 , 6.32331093108427e-01}, 5173d177a5cSEmil Constantinescu {8.35690179937017e+00 , -2.26640927349732e+00 , 6.86647884973826e+00 , -5.22595158025740e+00 , 4.50893068837431e+00}, 5183d177a5cSEmil Constantinescu {1.27656267027479e+01 , 2.80882153840821e+00 , 8.91173096522890e+00 , -1.07936444078906e+01 , 4.82534148988854e+00}} 5193d177a5cSEmil Constantinescu ,u[5][5] = {{1.00000000000000e+00 , -7.27272727273551e-02 , -3.45454545454419e-02 , -4.12121212119565e-03 ,-2.96969696964014e-04}, 5203d177a5cSEmil Constantinescu {1.00000000000000e+00 , 2.31252881006154e-01 , -8.29487834416481e-03 , -9.07191207681020e-03 ,-1.70378403743473e-03}, 5213d177a5cSEmil Constantinescu {1.00000000000000e+00 , 1.16925777880663e+00 , 3.59268562942635e-02 , -4.09013451730615e-02 ,-1.02411119670164e-02}, 5223d177a5cSEmil Constantinescu {1.00000000000000e+00 , 1.02634463704356e+00 , 1.59375044913405e-01 , 1.89673015035370e-03 ,-4.89987231897569e-03}, 5233d177a5cSEmil Constantinescu {1.00000000000000e+00 , 1.27746320298021e+00 , 2.37186008132728e-01 , -8.28694373940065e-02 ,-5.34396510196430e-02}} 5243d177a5cSEmil Constantinescu ,v[5][5] = {{1.00000000000000e+00 , 1.27746320298021e+00 , 2.37186008132728e-01 , -8.28694373940065e-02 ,-5.34396510196430e-02}, 5253d177a5cSEmil Constantinescu {0.00000000000000e+00 , -1.77635683940025e-15 , -1.99840144432528e-15 , -9.99200722162641e-16 ,-3.33066907387547e-16}, 5263d177a5cSEmil Constantinescu {0.00000000000000e+00 , 4.37280081906924e+00 , 5.49221645016377e-02 , -8.88913177394943e-02 , 1.12879077989154e-01}, 5273d177a5cSEmil Constantinescu {0.00000000000000e+00 , -1.22399504837280e+01 , -5.21287338448645e+00 , -8.03952325565291e-01 , 4.60298678047147e-01}, 5283d177a5cSEmil Constantinescu {0.00000000000000e+00 , -1.85178762883829e+01 , -5.21411849862624e+00 , -1.04283436528809e+00 , 7.49030161063651e-01}}; 5299566063dSJacob Faibussowitsch PetscCall(TSGLLESchemeCreate(4,4,5,5,c,*a,*b,*u,*v,&gl->schemes[gl->nschemes++])); 5303d177a5cSEmil Constantinescu } 5313d177a5cSEmil Constantinescu { 5323d177a5cSEmil Constantinescu /* p=q=5, r=s=6; 5333d177a5cSEmil Constantinescu irks(1/3,0,[1:6]/6,... 5343d177a5cSEmil Constantinescu [-0.0489 0.4228 -0.8814 0.9021],... 5353d177a5cSEmil Constantinescu [-0.3474 -0.6617 0.6294 0.2129 5363d177a5cSEmil Constantinescu 0.0044 -0.4256 -0.1427 -0.8936 5373d177a5cSEmil Constantinescu -0.8267 0.4821 0.1371 -0.2557 5383d177a5cSEmil Constantinescu -0.4426 -0.3855 -0.7514 0.3014]) 5393d177a5cSEmil Constantinescu */ 5403d177a5cSEmil Constantinescu const PetscScalar c[6] = {1./6, 2./6, 3./6, 4./6, 5./6, 1.} 5413d177a5cSEmil Constantinescu ,a[6][6] = {{ 3.33333333333940e-01, 0 , 0 , 0 , 0 , 0 }, 5423d177a5cSEmil Constantinescu { -8.64423857333350e-02, 3.33333333332888e-01, 0 , 0 , 0 , 0 }, 5433d177a5cSEmil Constantinescu { -2.16850174258252e+00, -2.23619072028839e+00, 3.33333333335204e-01, 0 , 0 , 0 }, 5443d177a5cSEmil Constantinescu { -4.73160970138997e+00, -3.89265344629268e+00, -2.76318716520933e-01, 3.33333333335759e-01, 0 , 0 }, 5453d177a5cSEmil Constantinescu { -6.75187540297338e+00, -7.90756533769377e+00, 7.90245051802259e-01, -4.48352364517632e-01, 3.33333333328483e-01, 0 }, 5463d177a5cSEmil Constantinescu { -4.26488287921548e+00, -1.19320395589302e+01, 3.38924509887755e+00, -2.23969848002481e+00, 6.62807710124007e-01, 3.33333333335440e-01}} 5473d177a5cSEmil Constantinescu ,b[6][6] = {{ -4.26488287921548e+00, -1.19320395589302e+01, 3.38924509887755e+00, -2.23969848002481e+00, 6.62807710124007e-01, 3.33333333335440e-01}, 5483d177a5cSEmil Constantinescu { -8.88178419700125e-16, 4.44089209850063e-16, -1.54737334057131e-15, -8.88178419700125e-16, 0.00000000000000e+00, 1.00000000000001e+00}, 5493d177a5cSEmil Constantinescu { -2.87780425770651e+01, -1.13520448264971e+01, 2.62002318943161e+01, 2.56943874812797e+01, -3.06702268304488e+01, 6.68067773510103e+00}, 5503d177a5cSEmil Constantinescu { 5.47971245256474e+01, 6.80366875868284e+01, -6.50952588861999e+01, -8.28643975339097e+01, 8.17416943896414e+01, -1.17819043489036e+01}, 5513d177a5cSEmil Constantinescu { -2.33332114788869e+02, 6.12942539462634e+01, -4.91850135865944e+01, 1.82716844135480e+02, -1.29788173979395e+02, 3.09968095651099e+01}, 5523d177a5cSEmil Constantinescu { -1.72049132343751e+02, 8.60194713593999e+00, 7.98154219170200e-01, 1.50371386053218e+02, -1.18515423962066e+02, 2.50898277784663e+01}} 5533d177a5cSEmil Constantinescu ,u[6][6] = {{ 1.00000000000000e+00, -1.66666666666870e-01, -4.16666666664335e-02, -3.85802469124815e-03, -2.25051440302250e-04, -9.64506172339142e-06}, 5543d177a5cSEmil Constantinescu { 1.00000000000000e+00, 8.64423857327162e-02, -4.11484912671353e-02, -1.11450903217645e-02, -1.47651050487126e-03, -1.34395070766826e-04}, 5553d177a5cSEmil Constantinescu { 1.00000000000000e+00, 4.57135912953434e+00, 1.06514719719137e+00, 1.33517564218007e-01, 1.11365952968659e-02, 6.12382756769504e-04}, 5563d177a5cSEmil Constantinescu { 1.00000000000000e+00, 9.23391519753404e+00, 2.22431212392095e+00, 2.91823807741891e-01, 2.52058456411084e-02, 1.22800542949647e-03}, 5573d177a5cSEmil Constantinescu { 1.00000000000000e+00, 1.48175480533865e+01, 3.73439117461835e+00, 5.14648336541804e-01, 4.76430038853402e-02, 2.56798515502156e-03}, 5583d177a5cSEmil Constantinescu { 1.00000000000000e+00, 1.50512347758335e+01, 4.10099701165164e+00, 5.66039141003603e-01, 3.91213893800891e-02, -2.99136269067853e-03}} 5593d177a5cSEmil Constantinescu ,v[6][6] = {{ 1.00000000000000e+00, 1.50512347758335e+01, 4.10099701165164e+00, 5.66039141003603e-01, 3.91213893800891e-02, -2.99136269067853e-03}, 5603d177a5cSEmil Constantinescu { 0.00000000000000e+00, -4.88498130835069e-15, -6.43929354282591e-15, -3.55271367880050e-15, -1.22124532708767e-15, -3.12250225675825e-16}, 5613d177a5cSEmil Constantinescu { 0.00000000000000e+00, 1.22250171233141e+01, -1.77150760606169e+00, 3.54516769879390e-01, 6.22298845883398e-01, 2.31647447450276e-01}, 5623d177a5cSEmil Constantinescu { 0.00000000000000e+00, -4.48339457331040e+01, -3.57363126641880e-01, 5.18750173123425e-01, 6.55727990241799e-02, 1.63175368287079e-01}, 5633d177a5cSEmil Constantinescu { 0.00000000000000e+00, 1.37297394708005e+02, -1.60145272991317e+00, -5.05319555199441e+00, 1.55328940390990e-01, 9.16629423682464e-01}, 5643d177a5cSEmil Constantinescu { 0.00000000000000e+00, 1.05703241119022e+02, -1.16610260983038e+00, -2.99767252773859e+00, -1.13472315553890e-01, 1.09742849254729e+00}}; 5659566063dSJacob Faibussowitsch PetscCall(TSGLLESchemeCreate(5,5,6,6,c,*a,*b,*u,*v,&gl->schemes[gl->nschemes++])); 5663d177a5cSEmil Constantinescu } 5673d177a5cSEmil Constantinescu PetscFunctionReturn(0); 5683d177a5cSEmil Constantinescu } 5693d177a5cSEmil Constantinescu 5703d177a5cSEmil Constantinescu /*@C 5713d177a5cSEmil Constantinescu TSGLLESetType - sets the class of general linear method to use for time-stepping 5723d177a5cSEmil Constantinescu 5733d177a5cSEmil Constantinescu Collective on TS 5743d177a5cSEmil Constantinescu 5753d177a5cSEmil Constantinescu Input Parameters: 5763d177a5cSEmil Constantinescu + ts - the TS context 5773d177a5cSEmil Constantinescu - type - a method 5783d177a5cSEmil Constantinescu 5793d177a5cSEmil Constantinescu Options Database Key: 5803d177a5cSEmil Constantinescu . -ts_gl_type <type> - sets the method, use -help for a list of available method (e.g. irks) 5813d177a5cSEmil Constantinescu 5823d177a5cSEmil Constantinescu Notes: 5833d177a5cSEmil Constantinescu See "petsc/include/petscts.h" for available methods (for instance) 5843d177a5cSEmil Constantinescu . TSGLLE_IRKS - Diagonally implicit methods with inherent Runge-Kutta stability (for stiff problems) 5853d177a5cSEmil Constantinescu 5863d177a5cSEmil Constantinescu Normally, it is best to use the TSSetFromOptions() command and 5873d177a5cSEmil Constantinescu then set the TSGLLE type from the options database rather than by using 5883d177a5cSEmil Constantinescu this routine. Using the options database provides the user with 5893d177a5cSEmil Constantinescu maximum flexibility in evaluating the many different solvers. 5903d177a5cSEmil Constantinescu The TSGLLESetType() routine is provided for those situations where it 5913d177a5cSEmil Constantinescu is necessary to set the timestepping solver independently of the 5923d177a5cSEmil Constantinescu command line or options database. This might be the case, for example, 5933d177a5cSEmil Constantinescu when the choice of solver changes during the execution of the 5943d177a5cSEmil Constantinescu program, and the user's application is taking responsibility for 5953d177a5cSEmil Constantinescu choosing the appropriate method. 5963d177a5cSEmil Constantinescu 5973d177a5cSEmil Constantinescu Level: intermediate 5983d177a5cSEmil Constantinescu 5993d177a5cSEmil Constantinescu @*/ 6003d177a5cSEmil Constantinescu PetscErrorCode TSGLLESetType(TS ts,TSGLLEType type) 6013d177a5cSEmil Constantinescu { 6023d177a5cSEmil Constantinescu PetscFunctionBegin; 6033d177a5cSEmil Constantinescu PetscValidHeaderSpecific(ts,TS_CLASSID,1); 6043d177a5cSEmil Constantinescu PetscValidCharPointer(type,2); 605cac4c232SBarry Smith PetscTryMethod(ts,"TSGLLESetType_C",(TS,TSGLLEType),(ts,type)); 6063d177a5cSEmil Constantinescu PetscFunctionReturn(0); 6073d177a5cSEmil Constantinescu } 6083d177a5cSEmil Constantinescu 6093d177a5cSEmil Constantinescu /*@C 6103d177a5cSEmil Constantinescu TSGLLESetAcceptType - sets the acceptance test 6113d177a5cSEmil Constantinescu 6123d177a5cSEmil Constantinescu Time integrators that need to control error must have the option to reject a time step based on local error 6133d177a5cSEmil Constantinescu estimates. This function allows different schemes to be set. 6143d177a5cSEmil Constantinescu 6153d177a5cSEmil Constantinescu Logically Collective on TS 6163d177a5cSEmil Constantinescu 6173d177a5cSEmil Constantinescu Input Parameters: 6183d177a5cSEmil Constantinescu + ts - the TS context 6193d177a5cSEmil Constantinescu - type - the type 6203d177a5cSEmil Constantinescu 6213d177a5cSEmil Constantinescu Options Database Key: 6223d177a5cSEmil Constantinescu . -ts_gl_accept_type <type> - sets the method used to determine whether to accept or reject a step 6233d177a5cSEmil Constantinescu 6243d177a5cSEmil Constantinescu Level: intermediate 6253d177a5cSEmil Constantinescu 626db781477SPatrick Sanan .seealso: `TS`, `TSGLLE`, `TSGLLEAcceptRegister()`, `TSGLLEAdapt`, `set` `type` 6273d177a5cSEmil Constantinescu @*/ 6283d177a5cSEmil Constantinescu PetscErrorCode TSGLLESetAcceptType(TS ts,TSGLLEAcceptType type) 6293d177a5cSEmil Constantinescu { 6303d177a5cSEmil Constantinescu PetscFunctionBegin; 6313d177a5cSEmil Constantinescu PetscValidHeaderSpecific(ts,TS_CLASSID,1); 6323d177a5cSEmil Constantinescu PetscValidCharPointer(type,2); 633cac4c232SBarry Smith PetscTryMethod(ts,"TSGLLESetAcceptType_C",(TS,TSGLLEAcceptType),(ts,type)); 6343d177a5cSEmil Constantinescu PetscFunctionReturn(0); 6353d177a5cSEmil Constantinescu } 6363d177a5cSEmil Constantinescu 6373d177a5cSEmil Constantinescu /*@C 6383d177a5cSEmil Constantinescu TSGLLEGetAdapt - gets the TSGLLEAdapt object from the TS 6393d177a5cSEmil Constantinescu 6403d177a5cSEmil Constantinescu Not Collective 6413d177a5cSEmil Constantinescu 6423d177a5cSEmil Constantinescu Input Parameter: 6433d177a5cSEmil Constantinescu . ts - the TS context 6443d177a5cSEmil Constantinescu 6453d177a5cSEmil Constantinescu Output Parameter: 6463d177a5cSEmil Constantinescu . adapt - the TSGLLEAdapt context 6473d177a5cSEmil Constantinescu 6483d177a5cSEmil Constantinescu Notes: 6493d177a5cSEmil Constantinescu This allows the user set options on the TSGLLEAdapt object. Usually it is better to do this using the options 6503d177a5cSEmil Constantinescu database, so this function is rarely needed. 6513d177a5cSEmil Constantinescu 6523d177a5cSEmil Constantinescu Level: advanced 6533d177a5cSEmil Constantinescu 654db781477SPatrick Sanan .seealso: `TSGLLEAdapt`, `TSGLLEAdaptRegister()` 6553d177a5cSEmil Constantinescu @*/ 6563d177a5cSEmil Constantinescu PetscErrorCode TSGLLEGetAdapt(TS ts,TSGLLEAdapt *adapt) 6573d177a5cSEmil Constantinescu { 6583d177a5cSEmil Constantinescu PetscFunctionBegin; 6593d177a5cSEmil Constantinescu PetscValidHeaderSpecific(ts,TS_CLASSID,1); 6603d177a5cSEmil Constantinescu PetscValidPointer(adapt,2); 661cac4c232SBarry Smith PetscUseMethod(ts,"TSGLLEGetAdapt_C",(TS,TSGLLEAdapt*),(ts,adapt)); 6623d177a5cSEmil Constantinescu PetscFunctionReturn(0); 6633d177a5cSEmil Constantinescu } 6643d177a5cSEmil Constantinescu 6653d177a5cSEmil Constantinescu static PetscErrorCode TSGLLEAccept_Always(TS ts,PetscReal tleft,PetscReal h,const PetscReal enorms[],PetscBool *accept) 6663d177a5cSEmil Constantinescu { 6673d177a5cSEmil Constantinescu PetscFunctionBegin; 6683d177a5cSEmil Constantinescu *accept = PETSC_TRUE; 6693d177a5cSEmil Constantinescu PetscFunctionReturn(0); 6703d177a5cSEmil Constantinescu } 6713d177a5cSEmil Constantinescu 6723d177a5cSEmil Constantinescu static PetscErrorCode TSGLLEUpdateWRMS(TS ts) 6733d177a5cSEmil Constantinescu { 6743d177a5cSEmil Constantinescu TS_GLLE *gl = (TS_GLLE*)ts->data; 6753d177a5cSEmil Constantinescu PetscScalar *x,*w; 6763d177a5cSEmil Constantinescu PetscInt n,i; 6773d177a5cSEmil Constantinescu 6783d177a5cSEmil Constantinescu PetscFunctionBegin; 6799566063dSJacob Faibussowitsch PetscCall(VecGetArray(gl->X[0],&x)); 6809566063dSJacob Faibussowitsch PetscCall(VecGetArray(gl->W,&w)); 6819566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(gl->W,&n)); 6823d177a5cSEmil Constantinescu for (i=0; i<n; i++) w[i] = 1./(gl->wrms_atol + gl->wrms_rtol*PetscAbsScalar(x[i])); 6839566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(gl->X[0],&x)); 6849566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(gl->W,&w)); 6853d177a5cSEmil Constantinescu PetscFunctionReturn(0); 6863d177a5cSEmil Constantinescu } 6873d177a5cSEmil Constantinescu 6883d177a5cSEmil Constantinescu static PetscErrorCode TSGLLEVecNormWRMS(TS ts,Vec X,PetscReal *nrm) 6893d177a5cSEmil Constantinescu { 6903d177a5cSEmil Constantinescu TS_GLLE *gl = (TS_GLLE*)ts->data; 6913d177a5cSEmil Constantinescu PetscScalar *x,*w; 6923d177a5cSEmil Constantinescu PetscReal sum = 0.0,gsum; 6933d177a5cSEmil Constantinescu PetscInt n,N,i; 6943d177a5cSEmil Constantinescu 6953d177a5cSEmil Constantinescu PetscFunctionBegin; 6969566063dSJacob Faibussowitsch PetscCall(VecGetArray(X,&x)); 6979566063dSJacob Faibussowitsch PetscCall(VecGetArray(gl->W,&w)); 6989566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(gl->W,&n)); 6993d177a5cSEmil Constantinescu for (i=0; i<n; i++) sum += PetscAbsScalar(PetscSqr(x[i]*w[i])); 7009566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(X,&x)); 7019566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(gl->W,&w)); 7021c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(&sum,&gsum,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)ts))); 7039566063dSJacob Faibussowitsch PetscCall(VecGetSize(gl->W,&N)); 7043d177a5cSEmil Constantinescu *nrm = PetscSqrtReal(gsum/(1.*N)); 7053d177a5cSEmil Constantinescu PetscFunctionReturn(0); 7063d177a5cSEmil Constantinescu } 7073d177a5cSEmil Constantinescu 7083d177a5cSEmil Constantinescu static PetscErrorCode TSGLLESetType_GLLE(TS ts,TSGLLEType type) 7093d177a5cSEmil Constantinescu { 7103d177a5cSEmil Constantinescu PetscBool same; 7113d177a5cSEmil Constantinescu TS_GLLE *gl = (TS_GLLE*)ts->data; 7125f80ce2aSJacob Faibussowitsch PetscErrorCode (*r)(TS); 7133d177a5cSEmil Constantinescu 7143d177a5cSEmil Constantinescu PetscFunctionBegin; 7153d177a5cSEmil Constantinescu if (gl->type_name[0]) { 7169566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(gl->type_name,type,&same)); 7173d177a5cSEmil Constantinescu if (same) PetscFunctionReturn(0); 7189566063dSJacob Faibussowitsch PetscCall((*gl->Destroy)(gl)); 7193d177a5cSEmil Constantinescu } 7203d177a5cSEmil Constantinescu 7219566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(TSGLLEList,type,&r)); 7223c633725SBarry Smith PetscCheck(r,PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown TSGLLE type \"%s\" given",type); 7239566063dSJacob Faibussowitsch PetscCall((*r)(ts)); 7249566063dSJacob Faibussowitsch PetscCall(PetscStrcpy(gl->type_name,type)); 7253d177a5cSEmil Constantinescu PetscFunctionReturn(0); 7263d177a5cSEmil Constantinescu } 7273d177a5cSEmil Constantinescu 7283d177a5cSEmil Constantinescu static PetscErrorCode TSGLLESetAcceptType_GLLE(TS ts,TSGLLEAcceptType type) 7293d177a5cSEmil Constantinescu { 7303d177a5cSEmil Constantinescu TSGLLEAcceptFunction r; 7313d177a5cSEmil Constantinescu TS_GLLE *gl = (TS_GLLE*)ts->data; 7323d177a5cSEmil Constantinescu 7333d177a5cSEmil Constantinescu PetscFunctionBegin; 7349566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(TSGLLEAcceptList,type,&r)); 7353c633725SBarry Smith PetscCheck(r,PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown TSGLLEAccept type \"%s\" given",type); 7363d177a5cSEmil Constantinescu gl->Accept = r; 7379566063dSJacob Faibussowitsch PetscCall(PetscStrncpy(gl->accept_name,type,sizeof(gl->accept_name))); 7383d177a5cSEmil Constantinescu PetscFunctionReturn(0); 7393d177a5cSEmil Constantinescu } 7403d177a5cSEmil Constantinescu 7413d177a5cSEmil Constantinescu static PetscErrorCode TSGLLEGetAdapt_GLLE(TS ts,TSGLLEAdapt *adapt) 7423d177a5cSEmil Constantinescu { 7433d177a5cSEmil Constantinescu TS_GLLE *gl = (TS_GLLE*)ts->data; 7443d177a5cSEmil Constantinescu 7453d177a5cSEmil Constantinescu PetscFunctionBegin; 7463d177a5cSEmil Constantinescu if (!gl->adapt) { 7479566063dSJacob Faibussowitsch PetscCall(TSGLLEAdaptCreate(PetscObjectComm((PetscObject)ts),&gl->adapt)); 7489566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)gl->adapt,(PetscObject)ts,1)); 7499566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)ts,(PetscObject)gl->adapt)); 7503d177a5cSEmil Constantinescu } 7513d177a5cSEmil Constantinescu *adapt = gl->adapt; 7523d177a5cSEmil Constantinescu PetscFunctionReturn(0); 7533d177a5cSEmil Constantinescu } 7543d177a5cSEmil Constantinescu 7553d177a5cSEmil Constantinescu static PetscErrorCode TSGLLEChooseNextScheme(TS ts,PetscReal h,const PetscReal hmnorm[],PetscInt *next_scheme,PetscReal *next_h,PetscBool *finish) 7563d177a5cSEmil Constantinescu { 7573d177a5cSEmil Constantinescu TS_GLLE *gl = (TS_GLLE*)ts->data; 7583d177a5cSEmil Constantinescu PetscInt i,n,cur_p,cur,next_sc,candidates[64],orders[64]; 7593d177a5cSEmil Constantinescu PetscReal errors[64],costs[64],tleft; 7603d177a5cSEmil Constantinescu 7613d177a5cSEmil Constantinescu PetscFunctionBegin; 7623d177a5cSEmil Constantinescu cur = -1; 7633d177a5cSEmil Constantinescu cur_p = gl->schemes[gl->current_scheme]->p; 7643d177a5cSEmil Constantinescu tleft = ts->max_time - (ts->ptime + ts->time_step); 7653d177a5cSEmil Constantinescu for (i=0,n=0; i<gl->nschemes; i++) { 7663d177a5cSEmil Constantinescu TSGLLEScheme sc = gl->schemes[i]; 7673d177a5cSEmil Constantinescu if (sc->p < gl->min_order || gl->max_order < sc->p) continue; 7683d177a5cSEmil Constantinescu if (sc->p == cur_p - 1) errors[n] = PetscAbsScalar(sc->alpha[0])*hmnorm[0]; 7693d177a5cSEmil Constantinescu else if (sc->p == cur_p) errors[n] = PetscAbsScalar(sc->alpha[0])*hmnorm[1]; 7703d177a5cSEmil Constantinescu else if (sc->p == cur_p+1) errors[n] = PetscAbsScalar(sc->alpha[0])*(hmnorm[2]+hmnorm[3]); 7713d177a5cSEmil Constantinescu else continue; 7723d177a5cSEmil Constantinescu candidates[n] = i; 7733d177a5cSEmil Constantinescu orders[n] = PetscMin(sc->p,sc->q); /* order of global truncation error */ 7743d177a5cSEmil Constantinescu costs[n] = sc->s; /* estimate the cost as the number of stages */ 7753d177a5cSEmil Constantinescu if (i == gl->current_scheme) cur = n; 7763d177a5cSEmil Constantinescu n++; 7773d177a5cSEmil Constantinescu } 778cad9d221SBarry Smith PetscCheck(cur >= 0 && gl->nschemes > cur,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Current scheme not found in scheme list"); 7799566063dSJacob Faibussowitsch PetscCall(TSGLLEAdaptChoose(gl->adapt,n,orders,errors,costs,cur,h,tleft,&next_sc,next_h,finish)); 7803d177a5cSEmil Constantinescu *next_scheme = candidates[next_sc]; 78163a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(ts,"Adapt chose scheme %" PetscInt_FMT " (%" PetscInt_FMT ",%" PetscInt_FMT ",%" PetscInt_FMT ",%" PetscInt_FMT ") with step size %6.2e, finish=%s\n",*next_scheme,gl->schemes[*next_scheme]->p,gl->schemes[*next_scheme]->q,gl->schemes[*next_scheme]->r,gl->schemes[*next_scheme]->s,(double)*next_h,PetscBools[*finish])); 7823d177a5cSEmil Constantinescu PetscFunctionReturn(0); 7833d177a5cSEmil Constantinescu } 7843d177a5cSEmil Constantinescu 7853d177a5cSEmil Constantinescu static PetscErrorCode TSGLLEGetMaxSizes(TS ts,PetscInt *max_r,PetscInt *max_s) 7863d177a5cSEmil Constantinescu { 7873d177a5cSEmil Constantinescu TS_GLLE *gl = (TS_GLLE*)ts->data; 7883d177a5cSEmil Constantinescu 7893d177a5cSEmil Constantinescu PetscFunctionBegin; 7903d177a5cSEmil Constantinescu *max_r = gl->schemes[gl->nschemes-1]->r; 7913d177a5cSEmil Constantinescu *max_s = gl->schemes[gl->nschemes-1]->s; 7923d177a5cSEmil Constantinescu PetscFunctionReturn(0); 7933d177a5cSEmil Constantinescu } 7943d177a5cSEmil Constantinescu 7953d177a5cSEmil Constantinescu static PetscErrorCode TSSolve_GLLE(TS ts) 7963d177a5cSEmil Constantinescu { 7973d177a5cSEmil Constantinescu TS_GLLE *gl = (TS_GLLE*)ts->data; 7983d177a5cSEmil Constantinescu PetscInt i,k,its,lits,max_r,max_s; 7993d177a5cSEmil Constantinescu PetscBool final_step,finish; 8003d177a5cSEmil Constantinescu SNESConvergedReason snesreason; 8013d177a5cSEmil Constantinescu 8023d177a5cSEmil Constantinescu PetscFunctionBegin; 8039566063dSJacob Faibussowitsch PetscCall(TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol)); 8043d177a5cSEmil Constantinescu 8059566063dSJacob Faibussowitsch PetscCall(TSGLLEGetMaxSizes(ts,&max_r,&max_s)); 8069566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol,gl->X[0])); 8073d177a5cSEmil Constantinescu for (i=1; i<max_r; i++) { 8089566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(gl->X[i])); 8093d177a5cSEmil Constantinescu } 8109566063dSJacob Faibussowitsch PetscCall(TSGLLEUpdateWRMS(ts)); 8113d177a5cSEmil Constantinescu 8123d177a5cSEmil Constantinescu if (0) { 8133d177a5cSEmil Constantinescu /* Find consistent initial data for DAE */ 8143d177a5cSEmil Constantinescu gl->stage_time = ts->ptime + ts->time_step; 8153d177a5cSEmil Constantinescu gl->scoeff = 1.; 8163d177a5cSEmil Constantinescu gl->stage = 0; 8173d177a5cSEmil Constantinescu 8189566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol,gl->Z)); 8199566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol,gl->Y)); 8209566063dSJacob Faibussowitsch PetscCall(SNESSolve(ts->snes,NULL,gl->Y)); 8219566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(ts->snes,&its)); 8229566063dSJacob Faibussowitsch PetscCall(SNESGetLinearSolveIterations(ts->snes,&lits)); 8239566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(ts->snes,&snesreason)); 8243d177a5cSEmil Constantinescu 8253d177a5cSEmil Constantinescu ts->snes_its += its; ts->ksp_its += lits; 8263d177a5cSEmil Constantinescu if (snesreason < 0 && ts->max_snes_failures > 0 && ++ts->num_snes_failures >= ts->max_snes_failures) { 8273d177a5cSEmil Constantinescu ts->reason = TS_DIVERGED_NONLINEAR_SOLVE; 82863a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(ts,"Step=%" PetscInt_FMT ", nonlinear solve solve failures %" PetscInt_FMT " greater than current TS allowed, stopping solve\n",ts->steps,ts->num_snes_failures)); 8293d177a5cSEmil Constantinescu PetscFunctionReturn(0); 8303d177a5cSEmil Constantinescu } 8313d177a5cSEmil Constantinescu } 8323d177a5cSEmil Constantinescu 83308401ef6SPierre Jolivet PetscCheck(gl->current_scheme >= 0,PETSC_COMM_SELF,PETSC_ERR_ORDER,"A starting scheme has not been provided"); 8343d177a5cSEmil Constantinescu 8353d177a5cSEmil Constantinescu for (k=0,final_step=PETSC_FALSE,finish=PETSC_FALSE; k<ts->max_steps && !finish; k++) { 8363d177a5cSEmil Constantinescu PetscInt j,r,s,next_scheme = 0,rejections; 8373d177a5cSEmil Constantinescu PetscReal h,hmnorm[4],enorm[3],next_h; 8383d177a5cSEmil Constantinescu PetscBool accept; 8393d177a5cSEmil Constantinescu const PetscScalar *c,*a,*u; 8403d177a5cSEmil Constantinescu Vec *X,*Ydot,Y; 8413d177a5cSEmil Constantinescu TSGLLEScheme scheme = gl->schemes[gl->current_scheme]; 8423d177a5cSEmil Constantinescu 8433d177a5cSEmil Constantinescu r = scheme->r; s = scheme->s; 8443d177a5cSEmil Constantinescu c = scheme->c; 8453d177a5cSEmil Constantinescu a = scheme->a; u = scheme->u; 8463d177a5cSEmil Constantinescu h = ts->time_step; 8473d177a5cSEmil Constantinescu X = gl->X; Ydot = gl->Ydot; Y = gl->Y; 8483d177a5cSEmil Constantinescu 8493d177a5cSEmil Constantinescu if (ts->ptime > ts->max_time) break; 8503d177a5cSEmil Constantinescu 8513d177a5cSEmil Constantinescu /* 8523d177a5cSEmil Constantinescu We only call PreStep at the start of each STEP, not each STAGE. This is because it is 8533d177a5cSEmil Constantinescu possible to fail (have to restart a step) after multiple stages. 8543d177a5cSEmil Constantinescu */ 8559566063dSJacob Faibussowitsch PetscCall(TSPreStep(ts)); 8563d177a5cSEmil Constantinescu 8573d177a5cSEmil Constantinescu rejections = 0; 8583d177a5cSEmil Constantinescu while (1) { 8593d177a5cSEmil Constantinescu for (i=0; i<s; i++) { 8603d177a5cSEmil Constantinescu PetscScalar shift; 8613d177a5cSEmil Constantinescu gl->scoeff = 1./PetscRealPart(a[i*s+i]); 8623d177a5cSEmil Constantinescu shift = gl->scoeff/ts->time_step; 8633d177a5cSEmil Constantinescu gl->stage = i; 8643d177a5cSEmil Constantinescu gl->stage_time = ts->ptime + PetscRealPart(c[i])*h; 8653d177a5cSEmil Constantinescu 8663d177a5cSEmil Constantinescu /* 8673d177a5cSEmil Constantinescu * Stage equation: Y = h A Y' + U X 8683d177a5cSEmil Constantinescu * We assume that A is lower-triangular so that we can solve the stages (Y,Y') sequentially 8693d177a5cSEmil Constantinescu * Build the affine vector z_i = -[1/(h a_ii)](h sum_j a_ij y'_j + sum_j u_ij x_j) 8703d177a5cSEmil Constantinescu * Then y'_i = z + 1/(h a_ii) y_i 8713d177a5cSEmil Constantinescu */ 8729566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(gl->Z)); 8733d177a5cSEmil Constantinescu for (j=0; j<r; j++) { 8749566063dSJacob Faibussowitsch PetscCall(VecAXPY(gl->Z,-shift*u[i*r+j],X[j])); 8753d177a5cSEmil Constantinescu } 8763d177a5cSEmil Constantinescu for (j=0; j<i; j++) { 8779566063dSJacob Faibussowitsch PetscCall(VecAXPY(gl->Z,-shift*h*a[i*s+j],Ydot[j])); 8783d177a5cSEmil Constantinescu } 8793d177a5cSEmil Constantinescu /* Note: Z is used within function evaluation, Ydot = Z + shift*Y */ 8803d177a5cSEmil Constantinescu 8813d177a5cSEmil Constantinescu /* Compute an estimate of Y to start Newton iteration */ 8823d177a5cSEmil Constantinescu if (gl->extrapolate) { 8833d177a5cSEmil Constantinescu if (i==0) { 8843d177a5cSEmil Constantinescu /* Linear extrapolation on the first stage */ 8859566063dSJacob Faibussowitsch PetscCall(VecWAXPY(Y,c[i]*h,X[1],X[0])); 8863d177a5cSEmil Constantinescu } else { 8873d177a5cSEmil Constantinescu /* Linear extrapolation from the last stage */ 8889566063dSJacob Faibussowitsch PetscCall(VecAXPY(Y,(c[i]-c[i-1])*h,Ydot[i-1])); 8893d177a5cSEmil Constantinescu } 8903d177a5cSEmil Constantinescu } else if (i==0) { /* Directly use solution from the last step, otherwise reuse the last stage (do nothing) */ 8919566063dSJacob Faibussowitsch PetscCall(VecCopy(X[0],Y)); 8923d177a5cSEmil Constantinescu } 8933d177a5cSEmil Constantinescu 8943d177a5cSEmil Constantinescu /* Solve this stage (Ydot[i] is computed during function evaluation) */ 8959566063dSJacob Faibussowitsch PetscCall(SNESSolve(ts->snes,NULL,Y)); 8969566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(ts->snes,&its)); 8979566063dSJacob Faibussowitsch PetscCall(SNESGetLinearSolveIterations(ts->snes,&lits)); 8989566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(ts->snes,&snesreason)); 8993d177a5cSEmil Constantinescu ts->snes_its += its; ts->ksp_its += lits; 9003d177a5cSEmil Constantinescu if (snesreason < 0 && ts->max_snes_failures > 0 && ++ts->num_snes_failures >= ts->max_snes_failures) { 9013d177a5cSEmil Constantinescu ts->reason = TS_DIVERGED_NONLINEAR_SOLVE; 90263a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(ts,"Step=%" PetscInt_FMT ", nonlinear solve solve failures %" PetscInt_FMT " greater than current TS allowed, stopping solve\n",ts->steps,ts->num_snes_failures)); 9033d177a5cSEmil Constantinescu PetscFunctionReturn(0); 9043d177a5cSEmil Constantinescu } 9053d177a5cSEmil Constantinescu } 9063d177a5cSEmil Constantinescu 9073d177a5cSEmil Constantinescu gl->stage_time = ts->ptime + ts->time_step; 9083d177a5cSEmil Constantinescu 9099566063dSJacob Faibussowitsch PetscCall((*gl->EstimateHigherMoments)(scheme,h,Ydot,gl->X,gl->himom)); 9103d177a5cSEmil Constantinescu /* hmnorm[i] = h^{p+i}x^{(p+i)} with i=0,1,2; hmnorm[3] = h^{p+2}(dx'/dx) x^{(p+1)} */ 9113d177a5cSEmil Constantinescu for (i=0; i<3; i++) { 9129566063dSJacob Faibussowitsch PetscCall(TSGLLEVecNormWRMS(ts,gl->himom[i],&hmnorm[i+1])); 9133d177a5cSEmil Constantinescu } 9143d177a5cSEmil Constantinescu enorm[0] = PetscRealPart(scheme->alpha[0])*hmnorm[1]; 9153d177a5cSEmil Constantinescu enorm[1] = PetscRealPart(scheme->beta[0]) *hmnorm[2]; 9163d177a5cSEmil Constantinescu enorm[2] = PetscRealPart(scheme->gamma[0])*hmnorm[3]; 9179566063dSJacob Faibussowitsch PetscCall((*gl->Accept)(ts,ts->max_time-gl->stage_time,h,enorm,&accept)); 9183d177a5cSEmil Constantinescu if (accept) goto accepted; 9193d177a5cSEmil Constantinescu rejections++; 92063a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(ts,"Step %" PetscInt_FMT " (t=%g) not accepted, rejections=%" PetscInt_FMT "\n",k,(double)gl->stage_time,rejections)); 9213d177a5cSEmil Constantinescu if (rejections > gl->max_step_rejections) break; 9223d177a5cSEmil Constantinescu /* 9233d177a5cSEmil Constantinescu There are lots of reasons why a step might be rejected, including solvers not converging and other factors that 9243d177a5cSEmil Constantinescu TSGLLEChooseNextScheme does not support. Additionally, the error estimates may be very screwed up, so I'm not 9253d177a5cSEmil Constantinescu convinced that it's safe to just compute a new error estimate using the same interface as the current adaptor 9263d177a5cSEmil Constantinescu (the adaptor interface probably has to change). Here we make an arbitrary and naive choice. This assumes that 9273d177a5cSEmil Constantinescu steps were written in Nordsieck form. The "correct" method would be to re-complete the previous time step with 9283d177a5cSEmil Constantinescu the correct "next" step size. It is unclear to me whether the present ad-hoc method of rescaling X is stable. 9293d177a5cSEmil Constantinescu */ 9303d177a5cSEmil Constantinescu h *= 0.5; 9313d177a5cSEmil Constantinescu for (i=1; i<scheme->r; i++) { 9329566063dSJacob Faibussowitsch PetscCall(VecScale(X[i],PetscPowRealInt(0.5,i))); 9333d177a5cSEmil Constantinescu } 9343d177a5cSEmil Constantinescu } 93563a3b9bcSJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF,PETSC_ERR_CONV_FAILED,"Time step %" PetscInt_FMT " (t=%g) not accepted after %" PetscInt_FMT " failures",k,(double)gl->stage_time,rejections); 9363d177a5cSEmil Constantinescu 9373d177a5cSEmil Constantinescu accepted: 9383d177a5cSEmil Constantinescu /* This term is not error, but it *would* be the leading term for a lower order method */ 9399566063dSJacob Faibussowitsch PetscCall(TSGLLEVecNormWRMS(ts,gl->X[scheme->r-1],&hmnorm[0])); 9403d177a5cSEmil Constantinescu /* Correct scaling so that these are equivalent to norms of the Nordsieck vectors */ 9413d177a5cSEmil Constantinescu 94263a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(ts,"Last moment norm %10.2e, estimated error norms %10.2e %10.2e %10.2e\n",(double)hmnorm[0],(double)enorm[0],(double)enorm[1],(double)enorm[2])); 9433d177a5cSEmil Constantinescu if (!final_step) { 9449566063dSJacob Faibussowitsch PetscCall(TSGLLEChooseNextScheme(ts,h,hmnorm,&next_scheme,&next_h,&final_step)); 9453d177a5cSEmil Constantinescu } else { 9463d177a5cSEmil Constantinescu /* Dummy values to complete the current step in a consistent manner */ 9473d177a5cSEmil Constantinescu next_scheme = gl->current_scheme; 9483d177a5cSEmil Constantinescu next_h = h; 9493d177a5cSEmil Constantinescu finish = PETSC_TRUE; 9503d177a5cSEmil Constantinescu } 9513d177a5cSEmil Constantinescu 9523d177a5cSEmil Constantinescu X = gl->Xold; 9533d177a5cSEmil Constantinescu gl->Xold = gl->X; 9543d177a5cSEmil Constantinescu gl->X = X; 9559566063dSJacob Faibussowitsch PetscCall((*gl->CompleteStep)(scheme,h,gl->schemes[next_scheme],next_h,Ydot,gl->Xold,gl->X)); 9563d177a5cSEmil Constantinescu 9579566063dSJacob Faibussowitsch PetscCall(TSGLLEUpdateWRMS(ts)); 9583d177a5cSEmil Constantinescu 9593d177a5cSEmil Constantinescu /* Post the solution for the user, we could avoid this copy with a small bit of cleverness */ 9609566063dSJacob Faibussowitsch PetscCall(VecCopy(gl->X[0],ts->vec_sol)); 9613d177a5cSEmil Constantinescu ts->ptime += h; 9623d177a5cSEmil Constantinescu ts->steps++; 9633d177a5cSEmil Constantinescu 9649566063dSJacob Faibussowitsch PetscCall(TSPostEvaluate(ts)); 9659566063dSJacob Faibussowitsch PetscCall(TSPostStep(ts)); 9669566063dSJacob Faibussowitsch PetscCall(TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol)); 9673d177a5cSEmil Constantinescu 9683d177a5cSEmil Constantinescu gl->current_scheme = next_scheme; 9693d177a5cSEmil Constantinescu ts->time_step = next_h; 9703d177a5cSEmil Constantinescu } 9713d177a5cSEmil Constantinescu PetscFunctionReturn(0); 9723d177a5cSEmil Constantinescu } 9733d177a5cSEmil Constantinescu 9743d177a5cSEmil Constantinescu /*------------------------------------------------------------*/ 9753d177a5cSEmil Constantinescu 9763d177a5cSEmil Constantinescu static PetscErrorCode TSReset_GLLE(TS ts) 9773d177a5cSEmil Constantinescu { 9783d177a5cSEmil Constantinescu TS_GLLE *gl = (TS_GLLE*)ts->data; 9793d177a5cSEmil Constantinescu PetscInt max_r,max_s; 9803d177a5cSEmil Constantinescu 9813d177a5cSEmil Constantinescu PetscFunctionBegin; 9823d177a5cSEmil Constantinescu if (gl->setupcalled) { 9839566063dSJacob Faibussowitsch PetscCall(TSGLLEGetMaxSizes(ts,&max_r,&max_s)); 9849566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(max_r,&gl->Xold)); 9859566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(max_r,&gl->X)); 9869566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(max_s,&gl->Ydot)); 9879566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(3,&gl->himom)); 9889566063dSJacob Faibussowitsch PetscCall(VecDestroy(&gl->W)); 9899566063dSJacob Faibussowitsch PetscCall(VecDestroy(&gl->Y)); 9909566063dSJacob Faibussowitsch PetscCall(VecDestroy(&gl->Z)); 9913d177a5cSEmil Constantinescu } 9923d177a5cSEmil Constantinescu gl->setupcalled = PETSC_FALSE; 9933d177a5cSEmil Constantinescu PetscFunctionReturn(0); 9943d177a5cSEmil Constantinescu } 9953d177a5cSEmil Constantinescu 9963d177a5cSEmil Constantinescu static PetscErrorCode TSDestroy_GLLE(TS ts) 9973d177a5cSEmil Constantinescu { 9983d177a5cSEmil Constantinescu TS_GLLE *gl = (TS_GLLE*)ts->data; 9993d177a5cSEmil Constantinescu 10003d177a5cSEmil Constantinescu PetscFunctionBegin; 10019566063dSJacob Faibussowitsch PetscCall(TSReset_GLLE(ts)); 1002b3a6b972SJed Brown if (ts->dm) { 10039566063dSJacob Faibussowitsch PetscCall(DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSGLLE,DMRestrictHook_TSGLLE,ts)); 10049566063dSJacob Faibussowitsch PetscCall(DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSGLLE,DMSubDomainRestrictHook_TSGLLE,ts)); 1005b3a6b972SJed Brown } 10069566063dSJacob Faibussowitsch if (gl->adapt) PetscCall(TSGLLEAdaptDestroy(&gl->adapt)); 10079566063dSJacob Faibussowitsch if (gl->Destroy) PetscCall((*gl->Destroy)(gl)); 10089566063dSJacob Faibussowitsch PetscCall(PetscFree(ts->data)); 10099566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSGLLESetType_C",NULL)); 10109566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSGLLESetAcceptType_C",NULL)); 10119566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSGLLEGetAdapt_C",NULL)); 10123d177a5cSEmil Constantinescu PetscFunctionReturn(0); 10133d177a5cSEmil Constantinescu } 10143d177a5cSEmil Constantinescu 10153d177a5cSEmil Constantinescu /* 10163d177a5cSEmil Constantinescu This defines the nonlinear equation that is to be solved with SNES 10173d177a5cSEmil Constantinescu g(x) = f(t,x,z+shift*x) = 0 10183d177a5cSEmil Constantinescu */ 10193d177a5cSEmil Constantinescu static PetscErrorCode SNESTSFormFunction_GLLE(SNES snes,Vec x,Vec f,TS ts) 10203d177a5cSEmil Constantinescu { 10213d177a5cSEmil Constantinescu TS_GLLE *gl = (TS_GLLE*)ts->data; 10223d177a5cSEmil Constantinescu Vec Z,Ydot; 10233d177a5cSEmil Constantinescu DM dm,dmsave; 10243d177a5cSEmil Constantinescu 10253d177a5cSEmil Constantinescu PetscFunctionBegin; 10269566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes,&dm)); 10279566063dSJacob Faibussowitsch PetscCall(TSGLLEGetVecs(ts,dm,&Z,&Ydot)); 10289566063dSJacob Faibussowitsch PetscCall(VecWAXPY(Ydot,gl->scoeff/ts->time_step,x,Z)); 10293d177a5cSEmil Constantinescu dmsave = ts->dm; 10303d177a5cSEmil Constantinescu ts->dm = dm; 10319566063dSJacob Faibussowitsch PetscCall(TSComputeIFunction(ts,gl->stage_time,x,Ydot,f,PETSC_FALSE)); 10323d177a5cSEmil Constantinescu ts->dm = dmsave; 10339566063dSJacob Faibussowitsch PetscCall(TSGLLERestoreVecs(ts,dm,&Z,&Ydot)); 10343d177a5cSEmil Constantinescu PetscFunctionReturn(0); 10353d177a5cSEmil Constantinescu } 10363d177a5cSEmil Constantinescu 10373d177a5cSEmil Constantinescu static PetscErrorCode SNESTSFormJacobian_GLLE(SNES snes,Vec x,Mat A,Mat B,TS ts) 10383d177a5cSEmil Constantinescu { 10393d177a5cSEmil Constantinescu TS_GLLE *gl = (TS_GLLE*)ts->data; 10403d177a5cSEmil Constantinescu Vec Z,Ydot; 10413d177a5cSEmil Constantinescu DM dm,dmsave; 10423d177a5cSEmil Constantinescu 10433d177a5cSEmil Constantinescu PetscFunctionBegin; 10449566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes,&dm)); 10459566063dSJacob Faibussowitsch PetscCall(TSGLLEGetVecs(ts,dm,&Z,&Ydot)); 10463d177a5cSEmil Constantinescu dmsave = ts->dm; 10473d177a5cSEmil Constantinescu ts->dm = dm; 10483d177a5cSEmil Constantinescu /* gl->Xdot will have already been computed in SNESTSFormFunction_GLLE */ 10499566063dSJacob Faibussowitsch PetscCall(TSComputeIJacobian(ts,gl->stage_time,x,gl->Ydot[gl->stage],gl->scoeff/ts->time_step,A,B,PETSC_FALSE)); 10503d177a5cSEmil Constantinescu ts->dm = dmsave; 10519566063dSJacob Faibussowitsch PetscCall(TSGLLERestoreVecs(ts,dm,&Z,&Ydot)); 10523d177a5cSEmil Constantinescu PetscFunctionReturn(0); 10533d177a5cSEmil Constantinescu } 10543d177a5cSEmil Constantinescu 10553d177a5cSEmil Constantinescu static PetscErrorCode TSSetUp_GLLE(TS ts) 10563d177a5cSEmil Constantinescu { 10573d177a5cSEmil Constantinescu TS_GLLE *gl = (TS_GLLE*)ts->data; 10583d177a5cSEmil Constantinescu PetscInt max_r,max_s; 10593d177a5cSEmil Constantinescu DM dm; 10603d177a5cSEmil Constantinescu 10613d177a5cSEmil Constantinescu PetscFunctionBegin; 10623d177a5cSEmil Constantinescu gl->setupcalled = PETSC_TRUE; 10639566063dSJacob Faibussowitsch PetscCall(TSGLLEGetMaxSizes(ts,&max_r,&max_s)); 10649566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vec_sol,max_r,&gl->X)); 10659566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vec_sol,max_r,&gl->Xold)); 10669566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vec_sol,max_s,&gl->Ydot)); 10679566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vec_sol,3,&gl->himom)); 10689566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol,&gl->W)); 10699566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol,&gl->Y)); 10709566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol,&gl->Z)); 10713d177a5cSEmil Constantinescu 10723d177a5cSEmil Constantinescu /* Default acceptance tests and adaptivity */ 10739566063dSJacob Faibussowitsch if (!gl->Accept) PetscCall(TSGLLESetAcceptType(ts,TSGLLEACCEPT_ALWAYS)); 10749566063dSJacob Faibussowitsch if (!gl->adapt) PetscCall(TSGLLEGetAdapt(ts,&gl->adapt)); 10753d177a5cSEmil Constantinescu 10763d177a5cSEmil Constantinescu if (gl->current_scheme < 0) { 10773d177a5cSEmil Constantinescu PetscInt i; 10783d177a5cSEmil Constantinescu for (i=0;; i++) { 10793d177a5cSEmil Constantinescu if (gl->schemes[i]->p == gl->start_order) break; 108063a3b9bcSJacob Faibussowitsch PetscCheck(i+1 != gl->nschemes,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"No schemes available with requested start order %" PetscInt_FMT,i); 10813d177a5cSEmil Constantinescu } 10823d177a5cSEmil Constantinescu gl->current_scheme = i; 10833d177a5cSEmil Constantinescu } 10849566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts,&dm)); 10859566063dSJacob Faibussowitsch PetscCall(DMCoarsenHookAdd(dm,DMCoarsenHook_TSGLLE,DMRestrictHook_TSGLLE,ts)); 10869566063dSJacob Faibussowitsch PetscCall(DMSubDomainHookAdd(dm,DMSubDomainHook_TSGLLE,DMSubDomainRestrictHook_TSGLLE,ts)); 10873d177a5cSEmil Constantinescu PetscFunctionReturn(0); 10883d177a5cSEmil Constantinescu } 10893d177a5cSEmil Constantinescu /*------------------------------------------------------------*/ 10903d177a5cSEmil Constantinescu 1091*dbbe0bcdSBarry Smith static PetscErrorCode TSSetFromOptions_GLLE(TS ts,PetscOptionItems *PetscOptionsObject) 10923d177a5cSEmil Constantinescu { 10933d177a5cSEmil Constantinescu TS_GLLE *gl = (TS_GLLE*)ts->data; 10943d177a5cSEmil Constantinescu char tname[256] = TSGLLE_IRKS,completef[256] = "rescale-and-modify"; 10953d177a5cSEmil Constantinescu 10963d177a5cSEmil Constantinescu PetscFunctionBegin; 1097d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject,"General Linear ODE solver options"); 10983d177a5cSEmil Constantinescu { 10993d177a5cSEmil Constantinescu PetscBool flg; 11009566063dSJacob Faibussowitsch PetscCall(PetscOptionsFList("-ts_gl_type","Type of GL method","TSGLLESetType",TSGLLEList,gl->type_name[0] ? gl->type_name : tname,tname,sizeof(tname),&flg)); 11013d177a5cSEmil Constantinescu if (flg || !gl->type_name[0]) { 11029566063dSJacob Faibussowitsch PetscCall(TSGLLESetType(ts,tname)); 11033d177a5cSEmil Constantinescu } 11049566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-ts_gl_max_step_rejections","Maximum number of times to attempt a step","None",gl->max_step_rejections,&gl->max_step_rejections,NULL)); 11059566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-ts_gl_max_order","Maximum order to try","TSGLLESetMaxOrder",gl->max_order,&gl->max_order,NULL)); 11069566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-ts_gl_min_order","Minimum order to try","TSGLLESetMinOrder",gl->min_order,&gl->min_order,NULL)); 11079566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-ts_gl_start_order","Initial order to try","TSGLLESetMinOrder",gl->start_order,&gl->start_order,NULL)); 11089566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-ts_gl_error_direction","Which direction to look when estimating error","TSGLLESetErrorDirection",TSGLLEErrorDirections,(PetscEnum)gl->error_direction,(PetscEnum*)&gl->error_direction,NULL)); 11099566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-ts_gl_extrapolate","Extrapolate stage solution from previous solution (sometimes unstable)","TSGLLESetExtrapolate",gl->extrapolate,&gl->extrapolate,NULL)); 11109566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_gl_atol","Absolute tolerance","TSGLLESetTolerances",gl->wrms_atol,&gl->wrms_atol,NULL)); 11119566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_gl_rtol","Relative tolerance","TSGLLESetTolerances",gl->wrms_rtol,&gl->wrms_rtol,NULL)); 11129566063dSJacob Faibussowitsch PetscCall(PetscOptionsString("-ts_gl_complete","Method to use for completing the step","none",completef,completef,sizeof(completef),&flg)); 11133d177a5cSEmil Constantinescu if (flg) { 11143d177a5cSEmil Constantinescu PetscBool match1,match2; 11159566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(completef,"rescale",&match1)); 11169566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(completef,"rescale-and-modify",&match2)); 11173d177a5cSEmil Constantinescu if (match1) gl->CompleteStep = TSGLLECompleteStep_Rescale; 11183d177a5cSEmil Constantinescu else if (match2) gl->CompleteStep = TSGLLECompleteStep_RescaleAndModify; 111998921bdaSJacob Faibussowitsch else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"%s",completef); 11203d177a5cSEmil Constantinescu } 11213d177a5cSEmil Constantinescu { 11223d177a5cSEmil Constantinescu char type[256] = TSGLLEACCEPT_ALWAYS; 11239566063dSJacob Faibussowitsch PetscCall(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)); 11243d177a5cSEmil Constantinescu if (flg || !gl->accept_name[0]) { 11259566063dSJacob Faibussowitsch PetscCall(TSGLLESetAcceptType(ts,type)); 11263d177a5cSEmil Constantinescu } 11273d177a5cSEmil Constantinescu } 11283d177a5cSEmil Constantinescu { 11293d177a5cSEmil Constantinescu TSGLLEAdapt adapt; 11309566063dSJacob Faibussowitsch PetscCall(TSGLLEGetAdapt(ts,&adapt)); 1131*dbbe0bcdSBarry Smith PetscCall(TSGLLEAdaptSetFromOptions(adapt,PetscOptionsObject)); 11323d177a5cSEmil Constantinescu } 11333d177a5cSEmil Constantinescu } 1134d0609cedSBarry Smith PetscOptionsHeadEnd(); 11353d177a5cSEmil Constantinescu PetscFunctionReturn(0); 11363d177a5cSEmil Constantinescu } 11373d177a5cSEmil Constantinescu 11383d177a5cSEmil Constantinescu static PetscErrorCode TSView_GLLE(TS ts,PetscViewer viewer) 11393d177a5cSEmil Constantinescu { 11403d177a5cSEmil Constantinescu TS_GLLE *gl = (TS_GLLE*)ts->data; 11413d177a5cSEmil Constantinescu PetscInt i; 11423d177a5cSEmil Constantinescu PetscBool iascii,details; 11433d177a5cSEmil Constantinescu 11443d177a5cSEmil Constantinescu PetscFunctionBegin; 11459566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii)); 11463d177a5cSEmil Constantinescu if (iascii) { 114763a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," min order %" PetscInt_FMT ", max order %" PetscInt_FMT ", current order %" PetscInt_FMT "\n",gl->min_order,gl->max_order,gl->schemes[gl->current_scheme]->p)); 11489566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Error estimation: %s\n",TSGLLEErrorDirections[gl->error_direction])); 11499566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Extrapolation: %s\n",gl->extrapolate ? "yes" : "no")); 11509566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Acceptance test: %s\n",gl->accept_name[0] ? gl->accept_name : "(not yet set)")); 11519566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 11529566063dSJacob Faibussowitsch PetscCall(TSGLLEAdaptView(gl->adapt,viewer)); 11539566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 11549566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," type: %s\n",gl->type_name[0] ? gl->type_name : "(not yet set)")); 115563a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Schemes within family (%" PetscInt_FMT "):\n",gl->nschemes)); 11563d177a5cSEmil Constantinescu details = PETSC_FALSE; 11579566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(((PetscObject)ts)->options,((PetscObject)ts)->prefix,"-ts_gl_view_detailed",&details,NULL)); 11589566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 11593d177a5cSEmil Constantinescu for (i=0; i<gl->nschemes; i++) { 11609566063dSJacob Faibussowitsch PetscCall(TSGLLESchemeView(gl->schemes[i],details,viewer)); 11613d177a5cSEmil Constantinescu } 11621baa6e33SBarry Smith if (gl->View) PetscCall((*gl->View)(gl,viewer)); 11639566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 11643d177a5cSEmil Constantinescu } 11653d177a5cSEmil Constantinescu PetscFunctionReturn(0); 11663d177a5cSEmil Constantinescu } 11673d177a5cSEmil Constantinescu 11683d177a5cSEmil Constantinescu /*@C 11693d177a5cSEmil Constantinescu TSGLLERegister - adds a TSGLLE implementation 11703d177a5cSEmil Constantinescu 11713d177a5cSEmil Constantinescu Not Collective 11723d177a5cSEmil Constantinescu 11733d177a5cSEmil Constantinescu Input Parameters: 11743d177a5cSEmil Constantinescu + name_scheme - name of user-defined general linear scheme 11753d177a5cSEmil Constantinescu - routine_create - routine to create method context 11763d177a5cSEmil Constantinescu 11773d177a5cSEmil Constantinescu Notes: 11783d177a5cSEmil Constantinescu TSGLLERegister() may be called multiple times to add several user-defined families. 11793d177a5cSEmil Constantinescu 11803d177a5cSEmil Constantinescu Sample usage: 11813d177a5cSEmil Constantinescu .vb 11823d177a5cSEmil Constantinescu TSGLLERegister("my_scheme",MySchemeCreate); 11833d177a5cSEmil Constantinescu .ve 11843d177a5cSEmil Constantinescu 11853d177a5cSEmil Constantinescu Then, your scheme can be chosen with the procedural interface via 11863d177a5cSEmil Constantinescu $ TSGLLESetType(ts,"my_scheme") 11873d177a5cSEmil Constantinescu or at runtime via the option 11883d177a5cSEmil Constantinescu $ -ts_gl_type my_scheme 11893d177a5cSEmil Constantinescu 11903d177a5cSEmil Constantinescu Level: advanced 11913d177a5cSEmil Constantinescu 1192db781477SPatrick Sanan .seealso: `TSGLLERegisterAll()` 11933d177a5cSEmil Constantinescu @*/ 11943d177a5cSEmil Constantinescu PetscErrorCode TSGLLERegister(const char sname[],PetscErrorCode (*function)(TS)) 11953d177a5cSEmil Constantinescu { 11963d177a5cSEmil Constantinescu PetscFunctionBegin; 11979566063dSJacob Faibussowitsch PetscCall(TSGLLEInitializePackage()); 11989566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&TSGLLEList,sname,function)); 11993d177a5cSEmil Constantinescu PetscFunctionReturn(0); 12003d177a5cSEmil Constantinescu } 12013d177a5cSEmil Constantinescu 12023d177a5cSEmil Constantinescu /*@C 12033d177a5cSEmil Constantinescu TSGLLEAcceptRegister - adds a TSGLLE acceptance scheme 12043d177a5cSEmil Constantinescu 12053d177a5cSEmil Constantinescu Not Collective 12063d177a5cSEmil Constantinescu 12073d177a5cSEmil Constantinescu Input Parameters: 12083d177a5cSEmil Constantinescu + name_scheme - name of user-defined acceptance scheme 12093d177a5cSEmil Constantinescu - routine_create - routine to create method context 12103d177a5cSEmil Constantinescu 12113d177a5cSEmil Constantinescu Notes: 12123d177a5cSEmil Constantinescu TSGLLEAcceptRegister() may be called multiple times to add several user-defined families. 12133d177a5cSEmil Constantinescu 12143d177a5cSEmil Constantinescu Sample usage: 12153d177a5cSEmil Constantinescu .vb 12163d177a5cSEmil Constantinescu TSGLLEAcceptRegister("my_scheme",MySchemeCreate); 12173d177a5cSEmil Constantinescu .ve 12183d177a5cSEmil Constantinescu 12193d177a5cSEmil Constantinescu Then, your scheme can be chosen with the procedural interface via 12203d177a5cSEmil Constantinescu $ TSGLLESetAcceptType(ts,"my_scheme") 12213d177a5cSEmil Constantinescu or at runtime via the option 12223d177a5cSEmil Constantinescu $ -ts_gl_accept_type my_scheme 12233d177a5cSEmil Constantinescu 12243d177a5cSEmil Constantinescu Level: advanced 12253d177a5cSEmil Constantinescu 1226db781477SPatrick Sanan .seealso: `TSGLLERegisterAll()` 12273d177a5cSEmil Constantinescu @*/ 12283d177a5cSEmil Constantinescu PetscErrorCode TSGLLEAcceptRegister(const char sname[],TSGLLEAcceptFunction function) 12293d177a5cSEmil Constantinescu { 12303d177a5cSEmil Constantinescu PetscFunctionBegin; 12319566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&TSGLLEAcceptList,sname,function)); 12323d177a5cSEmil Constantinescu PetscFunctionReturn(0); 12333d177a5cSEmil Constantinescu } 12343d177a5cSEmil Constantinescu 12353d177a5cSEmil Constantinescu /*@C 12363d177a5cSEmil Constantinescu TSGLLERegisterAll - Registers all of the general linear methods in TSGLLE 12373d177a5cSEmil Constantinescu 12383d177a5cSEmil Constantinescu Not Collective 12393d177a5cSEmil Constantinescu 12403d177a5cSEmil Constantinescu Level: advanced 12413d177a5cSEmil Constantinescu 1242db781477SPatrick Sanan .seealso: `TSGLLERegisterDestroy()` 12433d177a5cSEmil Constantinescu @*/ 12443d177a5cSEmil Constantinescu PetscErrorCode TSGLLERegisterAll(void) 12453d177a5cSEmil Constantinescu { 12463d177a5cSEmil Constantinescu PetscFunctionBegin; 12473d177a5cSEmil Constantinescu if (TSGLLERegisterAllCalled) PetscFunctionReturn(0); 12483d177a5cSEmil Constantinescu TSGLLERegisterAllCalled = PETSC_TRUE; 12493d177a5cSEmil Constantinescu 12509566063dSJacob Faibussowitsch PetscCall(TSGLLERegister(TSGLLE_IRKS, TSGLLECreate_IRKS)); 12519566063dSJacob Faibussowitsch PetscCall(TSGLLEAcceptRegister(TSGLLEACCEPT_ALWAYS,TSGLLEAccept_Always)); 12523d177a5cSEmil Constantinescu PetscFunctionReturn(0); 12533d177a5cSEmil Constantinescu } 12543d177a5cSEmil Constantinescu 12553d177a5cSEmil Constantinescu /*@C 12563d177a5cSEmil Constantinescu TSGLLEInitializePackage - This function initializes everything in the TSGLLE package. It is called 12578a690491SBarry Smith from TSInitializePackage(). 12583d177a5cSEmil Constantinescu 12593d177a5cSEmil Constantinescu Level: developer 12603d177a5cSEmil Constantinescu 1261db781477SPatrick Sanan .seealso: `PetscInitialize()` 12623d177a5cSEmil Constantinescu @*/ 12633d177a5cSEmil Constantinescu PetscErrorCode TSGLLEInitializePackage(void) 12643d177a5cSEmil Constantinescu { 12653d177a5cSEmil Constantinescu PetscFunctionBegin; 12663d177a5cSEmil Constantinescu if (TSGLLEPackageInitialized) PetscFunctionReturn(0); 12673d177a5cSEmil Constantinescu TSGLLEPackageInitialized = PETSC_TRUE; 12689566063dSJacob Faibussowitsch PetscCall(TSGLLERegisterAll()); 12699566063dSJacob Faibussowitsch PetscCall(PetscRegisterFinalize(TSGLLEFinalizePackage)); 12703d177a5cSEmil Constantinescu PetscFunctionReturn(0); 12713d177a5cSEmil Constantinescu } 12723d177a5cSEmil Constantinescu 12733d177a5cSEmil Constantinescu /*@C 12743d177a5cSEmil Constantinescu TSGLLEFinalizePackage - This function destroys everything in the TSGLLE package. It is 12753d177a5cSEmil Constantinescu called from PetscFinalize(). 12763d177a5cSEmil Constantinescu 12773d177a5cSEmil Constantinescu Level: developer 12783d177a5cSEmil Constantinescu 1279db781477SPatrick Sanan .seealso: `PetscFinalize()` 12803d177a5cSEmil Constantinescu @*/ 12813d177a5cSEmil Constantinescu PetscErrorCode TSGLLEFinalizePackage(void) 12823d177a5cSEmil Constantinescu { 12833d177a5cSEmil Constantinescu PetscFunctionBegin; 12849566063dSJacob Faibussowitsch PetscCall(PetscFunctionListDestroy(&TSGLLEList)); 12859566063dSJacob Faibussowitsch PetscCall(PetscFunctionListDestroy(&TSGLLEAcceptList)); 12863d177a5cSEmil Constantinescu TSGLLEPackageInitialized = PETSC_FALSE; 12873d177a5cSEmil Constantinescu TSGLLERegisterAllCalled = PETSC_FALSE; 12883d177a5cSEmil Constantinescu PetscFunctionReturn(0); 12893d177a5cSEmil Constantinescu } 12903d177a5cSEmil Constantinescu 12913d177a5cSEmil Constantinescu /* ------------------------------------------------------------ */ 12923d177a5cSEmil Constantinescu /*MC 12933d177a5cSEmil Constantinescu TSGLLE - DAE solver using implicit General Linear methods 12943d177a5cSEmil Constantinescu 12953d177a5cSEmil Constantinescu These methods contain Runge-Kutta and multistep schemes as special cases. These special cases have some fundamental 12963d177a5cSEmil Constantinescu limitations. For example, diagonally implicit Runge-Kutta cannot have stage order greater than 1 which limits their 12973d177a5cSEmil Constantinescu applicability to very stiff systems. Meanwhile, multistep methods cannot be A-stable for order greater than 2 and BDF 12983d177a5cSEmil Constantinescu are not 0-stable for order greater than 6. GL methods can be A- and L-stable with arbitrarily high stage order and 12993d177a5cSEmil Constantinescu reliable error estimates for both 1 and 2 orders higher to facilitate adaptive step sizes and adaptive order schemes. 13003d177a5cSEmil Constantinescu All this is possible while preserving a singly diagonally implicit structure. 13013d177a5cSEmil Constantinescu 13023d177a5cSEmil Constantinescu Options database keys: 13033d177a5cSEmil Constantinescu + -ts_gl_type <type> - the class of general linear method (irks) 13043d177a5cSEmil Constantinescu . -ts_gl_rtol <tol> - relative error 13053d177a5cSEmil Constantinescu . -ts_gl_atol <tol> - absolute error 13063d177a5cSEmil Constantinescu . -ts_gl_min_order <p> - minimum order method to consider (default=1) 13073d177a5cSEmil Constantinescu . -ts_gl_max_order <p> - maximum order method to consider (default=3) 13083d177a5cSEmil Constantinescu . -ts_gl_start_order <p> - order of starting method (default=1) 13093d177a5cSEmil Constantinescu . -ts_gl_complete <method> - method to use for completing the step (rescale-and-modify or rescale) 13103d177a5cSEmil Constantinescu - -ts_adapt_type <method> - adaptive controller to use (none step both) 13113d177a5cSEmil Constantinescu 13123d177a5cSEmil Constantinescu Notes: 13133d177a5cSEmil Constantinescu This integrator can be applied to DAE. 13143d177a5cSEmil Constantinescu 13153d177a5cSEmil Constantinescu Diagonally implicit general linear (DIGL) methods are a generalization of diagonally implicit Runge-Kutta (DIRK). 13163d177a5cSEmil Constantinescu They are represented by the tableau 13173d177a5cSEmil Constantinescu 13183d177a5cSEmil Constantinescu .vb 13193d177a5cSEmil Constantinescu A | U 13203d177a5cSEmil Constantinescu ------- 13213d177a5cSEmil Constantinescu B | V 13223d177a5cSEmil Constantinescu .ve 13233d177a5cSEmil Constantinescu 13243d177a5cSEmil Constantinescu combined with a vector c of abscissa. "Diagonally implicit" means that A is lower triangular. 13253d177a5cSEmil Constantinescu A step of the general method reads 13263d177a5cSEmil Constantinescu 13273d177a5cSEmil Constantinescu .vb 13283d177a5cSEmil Constantinescu [ Y ] = [A U] [ Y' ] 13293d177a5cSEmil Constantinescu [X^k] = [B V] [X^{k-1}] 13303d177a5cSEmil Constantinescu .ve 13313d177a5cSEmil Constantinescu 13323d177a5cSEmil Constantinescu where Y is the multivector of stage values, Y' is the multivector of stage derivatives, X^k is the Nordsieck vector of 13333d177a5cSEmil Constantinescu the solution at step k. The Nordsieck vector consists of the first r moments of the solution, given by 13343d177a5cSEmil Constantinescu 13353d177a5cSEmil Constantinescu .vb 13363d177a5cSEmil Constantinescu X = [x_0,x_1,...,x_{r-1}] = [x, h x', h^2 x'', ..., h^{r-1} x^{(r-1)} ] 13373d177a5cSEmil Constantinescu .ve 13383d177a5cSEmil Constantinescu 13393d177a5cSEmil Constantinescu If A is lower triangular, we can solve the stages (Y,Y') sequentially 13403d177a5cSEmil Constantinescu 13413d177a5cSEmil Constantinescu .vb 13423d177a5cSEmil 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} 13433d177a5cSEmil Constantinescu .ve 13443d177a5cSEmil Constantinescu 13453d177a5cSEmil Constantinescu and then construct the pieces to carry to the next step 13463d177a5cSEmil Constantinescu 13473d177a5cSEmil Constantinescu .vb 13483d177a5cSEmil 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} 13493d177a5cSEmil Constantinescu .ve 13503d177a5cSEmil Constantinescu 13513d177a5cSEmil Constantinescu Note that when the equations are cast in implicit form, we are using the stage equation to define y'_i 13523d177a5cSEmil Constantinescu in terms of y_i and known stuff (y_j for j<i and x_j for all j). 13533d177a5cSEmil Constantinescu 13543d177a5cSEmil Constantinescu Error estimation 13553d177a5cSEmil Constantinescu 13563d177a5cSEmil Constantinescu At present, the most attractive GL methods for stiff problems are singly diagonally implicit schemes which posses 13573d177a5cSEmil Constantinescu Inherent Runge-Kutta Stability (IRKS). These methods have r=s, the number of items passed between steps is equal to 13583d177a5cSEmil Constantinescu the number of stages. The order and stage-order are one less than the number of stages. We use the error estimates 13593d177a5cSEmil Constantinescu in the 2007 paper which provide the following estimates 13603d177a5cSEmil Constantinescu 13613d177a5cSEmil Constantinescu .vb 13623d177a5cSEmil Constantinescu h^{p+1} X^{(p+1)} = phi_0^T Y' + [0 psi_0^T] Xold 13633d177a5cSEmil Constantinescu h^{p+2} X^{(p+2)} = phi_1^T Y' + [0 psi_1^T] Xold 13643d177a5cSEmil Constantinescu h^{p+2} (dx'/dx) X^{(p+1)} = phi_2^T Y' + [0 psi_2^T] Xold 13653d177a5cSEmil Constantinescu .ve 13663d177a5cSEmil Constantinescu 13673d177a5cSEmil Constantinescu These estimates are accurate to O(h^{p+3}). 13683d177a5cSEmil Constantinescu 13693d177a5cSEmil Constantinescu Changing the step size 13703d177a5cSEmil Constantinescu 13713d177a5cSEmil Constantinescu We use the generalized "rescale and modify" scheme, see equation (4.5) of the 2007 paper. 13723d177a5cSEmil Constantinescu 13733d177a5cSEmil Constantinescu Level: beginner 13743d177a5cSEmil Constantinescu 13753d177a5cSEmil Constantinescu References: 1376606c0280SSatish Balay + * - John Butcher and Z. Jackieweicz and W. Wright, On error propagation in general linear methods for 13773d177a5cSEmil Constantinescu ordinary differential equations, Journal of Complexity, Vol 23, 2007. 1378606c0280SSatish Balay - * - John Butcher, Numerical methods for ordinary differential equations, second edition, Wiley, 2009. 13793d177a5cSEmil Constantinescu 1380db781477SPatrick Sanan .seealso: `TSCreate()`, `TS`, `TSSetType()` 13813d177a5cSEmil Constantinescu 13823d177a5cSEmil Constantinescu M*/ 13833d177a5cSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSCreate_GLLE(TS ts) 13843d177a5cSEmil Constantinescu { 13853d177a5cSEmil Constantinescu TS_GLLE *gl; 13863d177a5cSEmil Constantinescu 13873d177a5cSEmil Constantinescu PetscFunctionBegin; 13889566063dSJacob Faibussowitsch PetscCall(TSGLLEInitializePackage()); 13893d177a5cSEmil Constantinescu 13909566063dSJacob Faibussowitsch PetscCall(PetscNewLog(ts,&gl)); 13913d177a5cSEmil Constantinescu ts->data = (void*)gl; 13923d177a5cSEmil Constantinescu 13933d177a5cSEmil Constantinescu ts->ops->reset = TSReset_GLLE; 13943d177a5cSEmil Constantinescu ts->ops->destroy = TSDestroy_GLLE; 13953d177a5cSEmil Constantinescu ts->ops->view = TSView_GLLE; 13963d177a5cSEmil Constantinescu ts->ops->setup = TSSetUp_GLLE; 13973d177a5cSEmil Constantinescu ts->ops->solve = TSSolve_GLLE; 13983d177a5cSEmil Constantinescu ts->ops->setfromoptions = TSSetFromOptions_GLLE; 13993d177a5cSEmil Constantinescu ts->ops->snesfunction = SNESTSFormFunction_GLLE; 14003d177a5cSEmil Constantinescu ts->ops->snesjacobian = SNESTSFormJacobian_GLLE; 14013d177a5cSEmil Constantinescu 1402efd4aadfSBarry Smith ts->usessnes = PETSC_TRUE; 1403efd4aadfSBarry Smith 14043d177a5cSEmil Constantinescu gl->max_step_rejections = 1; 14053d177a5cSEmil Constantinescu gl->min_order = 1; 14063d177a5cSEmil Constantinescu gl->max_order = 3; 14073d177a5cSEmil Constantinescu gl->start_order = 1; 14083d177a5cSEmil Constantinescu gl->current_scheme = -1; 14093d177a5cSEmil Constantinescu gl->extrapolate = PETSC_FALSE; 14103d177a5cSEmil Constantinescu 14113d177a5cSEmil Constantinescu gl->wrms_atol = 1e-8; 14123d177a5cSEmil Constantinescu gl->wrms_rtol = 1e-5; 14133d177a5cSEmil Constantinescu 14149566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSGLLESetType_C", &TSGLLESetType_GLLE)); 14159566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSGLLESetAcceptType_C",&TSGLLESetAcceptType_GLLE)); 14169566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSGLLEGetAdapt_C", &TSGLLEGetAdapt_GLLE)); 14173d177a5cSEmil Constantinescu PetscFunctionReturn(0); 14183d177a5cSEmil Constantinescu } 1419