13d177a5cSEmil Constantinescu #include <../src/ts/impls/implicit/glle/glle.h> /*I "petscts.h" I*/ 23d177a5cSEmil Constantinescu #include <petscdm.h> 33d177a5cSEmil Constantinescu #include <petscblaslapack.h> 43d177a5cSEmil Constantinescu 5c793f718SLisandro Dalcin static const char *TSGLLEErrorDirections[] = {"FORWARD", "BACKWARD", "TSGLLEErrorDirection", "TSGLLEERROR_", NULL}; 63d177a5cSEmil Constantinescu static PetscFunctionList TSGLLEList; 73d177a5cSEmil Constantinescu static PetscFunctionList TSGLLEAcceptList; 83d177a5cSEmil Constantinescu static PetscBool TSGLLEPackageInitialized; 93d177a5cSEmil Constantinescu static PetscBool TSGLLERegisterAllCalled; 103d177a5cSEmil Constantinescu 113d177a5cSEmil Constantinescu /* This function is pure */ 12d71ae5a4SJacob Faibussowitsch static PetscScalar Factorial(PetscInt n) 13d71ae5a4SJacob Faibussowitsch { 143d177a5cSEmil Constantinescu PetscInt i; 153d177a5cSEmil Constantinescu if (n < 12) { /* Can compute with 32-bit integers */ 163d177a5cSEmil Constantinescu PetscInt f = 1; 173d177a5cSEmil Constantinescu for (i = 2; i <= n; i++) f *= i; 183d177a5cSEmil Constantinescu return (PetscScalar)f; 193d177a5cSEmil Constantinescu } else { 203d177a5cSEmil Constantinescu PetscScalar f = 1.; 213d177a5cSEmil Constantinescu for (i = 2; i <= n; i++) f *= (PetscScalar)i; 223d177a5cSEmil Constantinescu return f; 233d177a5cSEmil Constantinescu } 243d177a5cSEmil Constantinescu } 253d177a5cSEmil Constantinescu 263d177a5cSEmil Constantinescu /* This function is pure */ 27d71ae5a4SJacob Faibussowitsch static PetscScalar CPowF(PetscScalar c, PetscInt p) 28d71ae5a4SJacob Faibussowitsch { 293d177a5cSEmil Constantinescu return PetscPowRealInt(PetscRealPart(c), p) / Factorial(p); 303d177a5cSEmil Constantinescu } 313d177a5cSEmil Constantinescu 32d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSGLLEGetVecs(TS ts, DM dm, Vec *Z, Vec *Ydotstage) 33d71ae5a4SJacob Faibussowitsch { 343d177a5cSEmil Constantinescu TS_GLLE *gl = (TS_GLLE *)ts->data; 353d177a5cSEmil Constantinescu 363d177a5cSEmil Constantinescu PetscFunctionBegin; 373d177a5cSEmil Constantinescu if (Z) { 383d177a5cSEmil Constantinescu if (dm && dm != ts->dm) { 399566063dSJacob Faibussowitsch PetscCall(DMGetNamedGlobalVector(dm, "TSGLLE_Z", Z)); 403d177a5cSEmil Constantinescu } else *Z = gl->Z; 413d177a5cSEmil Constantinescu } 423d177a5cSEmil Constantinescu if (Ydotstage) { 433d177a5cSEmil Constantinescu if (dm && dm != ts->dm) { 449566063dSJacob Faibussowitsch PetscCall(DMGetNamedGlobalVector(dm, "TSGLLE_Ydot", Ydotstage)); 453d177a5cSEmil Constantinescu } else *Ydotstage = gl->Ydot[gl->stage]; 463d177a5cSEmil Constantinescu } 473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 483d177a5cSEmil Constantinescu } 493d177a5cSEmil Constantinescu 50d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSGLLERestoreVecs(TS ts, DM dm, Vec *Z, Vec *Ydotstage) 51d71ae5a4SJacob Faibussowitsch { 523d177a5cSEmil Constantinescu PetscFunctionBegin; 533d177a5cSEmil Constantinescu if (Z) { 5448a46eb9SPierre Jolivet if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSGLLE_Z", Z)); 553d177a5cSEmil Constantinescu } 563d177a5cSEmil Constantinescu if (Ydotstage) { 5748a46eb9SPierre Jolivet if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSGLLE_Ydot", Ydotstage)); 583d177a5cSEmil Constantinescu } 593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 603d177a5cSEmil Constantinescu } 613d177a5cSEmil Constantinescu 62d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMCoarsenHook_TSGLLE(DM fine, DM coarse, void *ctx) 63d71ae5a4SJacob Faibussowitsch { 643d177a5cSEmil Constantinescu PetscFunctionBegin; 653ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 663d177a5cSEmil Constantinescu } 673d177a5cSEmil Constantinescu 68d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMRestrictHook_TSGLLE(DM fine, Mat restrct, Vec rscale, Mat inject, DM coarse, void *ctx) 69d71ae5a4SJacob Faibussowitsch { 703d177a5cSEmil Constantinescu TS ts = (TS)ctx; 713d177a5cSEmil Constantinescu Vec Ydot, Ydot_c; 723d177a5cSEmil Constantinescu 733d177a5cSEmil Constantinescu PetscFunctionBegin; 749566063dSJacob Faibussowitsch PetscCall(TSGLLEGetVecs(ts, fine, NULL, &Ydot)); 759566063dSJacob Faibussowitsch PetscCall(TSGLLEGetVecs(ts, coarse, NULL, &Ydot_c)); 769566063dSJacob Faibussowitsch PetscCall(MatRestrict(restrct, Ydot, Ydot_c)); 779566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(Ydot_c, rscale, Ydot_c)); 789566063dSJacob Faibussowitsch PetscCall(TSGLLERestoreVecs(ts, fine, NULL, &Ydot)); 799566063dSJacob Faibussowitsch PetscCall(TSGLLERestoreVecs(ts, coarse, NULL, &Ydot_c)); 803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 813d177a5cSEmil Constantinescu } 823d177a5cSEmil Constantinescu 83d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSubDomainHook_TSGLLE(DM dm, DM subdm, void *ctx) 84d71ae5a4SJacob Faibussowitsch { 853d177a5cSEmil Constantinescu PetscFunctionBegin; 863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 873d177a5cSEmil Constantinescu } 883d177a5cSEmil Constantinescu 89d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSubDomainRestrictHook_TSGLLE(DM dm, VecScatter gscat, VecScatter lscat, DM subdm, void *ctx) 90d71ae5a4SJacob Faibussowitsch { 913d177a5cSEmil Constantinescu TS ts = (TS)ctx; 923d177a5cSEmil Constantinescu Vec Ydot, Ydot_s; 933d177a5cSEmil Constantinescu 943d177a5cSEmil Constantinescu PetscFunctionBegin; 959566063dSJacob Faibussowitsch PetscCall(TSGLLEGetVecs(ts, dm, NULL, &Ydot)); 969566063dSJacob Faibussowitsch PetscCall(TSGLLEGetVecs(ts, subdm, NULL, &Ydot_s)); 973d177a5cSEmil Constantinescu 989566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(gscat, Ydot, Ydot_s, INSERT_VALUES, SCATTER_FORWARD)); 999566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(gscat, Ydot, Ydot_s, INSERT_VALUES, SCATTER_FORWARD)); 1003d177a5cSEmil Constantinescu 1019566063dSJacob Faibussowitsch PetscCall(TSGLLERestoreVecs(ts, dm, NULL, &Ydot)); 1029566063dSJacob Faibussowitsch PetscCall(TSGLLERestoreVecs(ts, subdm, NULL, &Ydot_s)); 1033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1043d177a5cSEmil Constantinescu } 1053d177a5cSEmil Constantinescu 106d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSGLLESchemeCreate(PetscInt p, PetscInt q, PetscInt r, PetscInt s, const PetscScalar *c, const PetscScalar *a, const PetscScalar *b, const PetscScalar *u, const PetscScalar *v, TSGLLEScheme *inscheme) 107d71ae5a4SJacob Faibussowitsch { 1083d177a5cSEmil Constantinescu TSGLLEScheme scheme; 1093d177a5cSEmil Constantinescu PetscInt j; 1103d177a5cSEmil Constantinescu 1113d177a5cSEmil Constantinescu PetscFunctionBegin; 11208401ef6SPierre Jolivet PetscCheck(p >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Scheme order must be positive"); 11308401ef6SPierre Jolivet PetscCheck(r >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "At least one item must be carried between steps"); 11408401ef6SPierre Jolivet PetscCheck(s >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "At least one stage is required"); 1154f572ea9SToby Isaac PetscAssertPointer(inscheme, 10); 116c793f718SLisandro Dalcin *inscheme = NULL; 1179566063dSJacob Faibussowitsch PetscCall(PetscNew(&scheme)); 1183d177a5cSEmil Constantinescu scheme->p = p; 1193d177a5cSEmil Constantinescu scheme->q = q; 1203d177a5cSEmil Constantinescu scheme->r = r; 1213d177a5cSEmil Constantinescu scheme->s = s; 1223d177a5cSEmil Constantinescu 1239566063dSJacob Faibussowitsch PetscCall(PetscMalloc5(s, &scheme->c, s * s, &scheme->a, r * s, &scheme->b, r * s, &scheme->u, r * r, &scheme->v)); 1249566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(scheme->c, c, s)); 1253d177a5cSEmil Constantinescu for (j = 0; j < s * s; j++) scheme->a[j] = (PetscAbsScalar(a[j]) < 1e-12) ? 0 : a[j]; 1263d177a5cSEmil Constantinescu for (j = 0; j < r * s; j++) scheme->b[j] = (PetscAbsScalar(b[j]) < 1e-12) ? 0 : b[j]; 1273d177a5cSEmil Constantinescu for (j = 0; j < s * r; j++) scheme->u[j] = (PetscAbsScalar(u[j]) < 1e-12) ? 0 : u[j]; 1283d177a5cSEmil Constantinescu for (j = 0; j < r * r; j++) scheme->v[j] = (PetscAbsScalar(v[j]) < 1e-12) ? 0 : v[j]; 1293d177a5cSEmil Constantinescu 1309566063dSJacob Faibussowitsch PetscCall(PetscMalloc6(r, &scheme->alpha, r, &scheme->beta, r, &scheme->gamma, 3 * s, &scheme->phi, 3 * r, &scheme->psi, r, &scheme->stage_error)); 1313d177a5cSEmil Constantinescu { 1323d177a5cSEmil Constantinescu PetscInt i, j, k, ss = s + 2; 1336497c311SBarry Smith PetscBLASInt m, n, one = 1, *ipiv, lwork, info, ldb; 1343d177a5cSEmil Constantinescu PetscReal rcond, *sing, *workreal; 1353d177a5cSEmil Constantinescu PetscScalar *ImV, *H, *bmat, *workscalar, *c = scheme->c, *a = scheme->a, *b = scheme->b, *u = scheme->u, *v = scheme->v; 1363d177a5cSEmil Constantinescu PetscBLASInt rank; 1376497c311SBarry Smith 1386497c311SBarry Smith PetscCall(PetscBLASIntCast(4 * ((s + 3) * 3 + 3), &lwork)); 1399566063dSJacob Faibussowitsch PetscCall(PetscMalloc7(PetscSqr(r), &ImV, 3 * s, &H, 3 * ss, &bmat, lwork, &workscalar, 5 * (3 + r), &workreal, r + s, &sing, r + s, &ipiv)); 1403d177a5cSEmil Constantinescu 1413d177a5cSEmil Constantinescu /* column-major input */ 1423d177a5cSEmil Constantinescu for (i = 0; i < r - 1; i++) { 1433d177a5cSEmil Constantinescu for (j = 0; j < r - 1; j++) ImV[i + j * r] = 1.0 * (i == j) - v[(i + 1) * r + j + 1]; 1443d177a5cSEmil Constantinescu } 145dd8e379bSPierre Jolivet /* Build right-hand side for alpha (tp - glm.B(2:end,:)*(glm.c.^(p)./factorial(p))) */ 1463d177a5cSEmil Constantinescu for (i = 1; i < r; i++) { 1473d177a5cSEmil Constantinescu scheme->alpha[i] = 1. / Factorial(p + 1 - i); 1483d177a5cSEmil Constantinescu for (j = 0; j < s; j++) scheme->alpha[i] -= b[i * s + j] * CPowF(c[j], p); 1493d177a5cSEmil Constantinescu } 1509566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(r - 1, &m)); 1519566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(r, &n)); 152792fecdfSBarry Smith PetscCallBLAS("LAPACKgesv", LAPACKgesv_(&m, &one, ImV, &n, ipiv, scheme->alpha + 1, &n, &info)); 15308401ef6SPierre Jolivet PetscCheck(info >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "Bad argument to GESV"); 15408401ef6SPierre Jolivet PetscCheck(info <= 0, PETSC_COMM_SELF, PETSC_ERR_MAT_LU_ZRPVT, "Bad LU factorization"); 1553d177a5cSEmil Constantinescu 156dd8e379bSPierre Jolivet /* Build right-hand side for beta (tp1 - glm.B(2:end,:)*(glm.c.^(p+1)./factorial(p+1)) - e.alpha) */ 1573d177a5cSEmil Constantinescu for (i = 1; i < r; i++) { 1583d177a5cSEmil Constantinescu scheme->beta[i] = 1. / Factorial(p + 2 - i) - scheme->alpha[i]; 1593d177a5cSEmil Constantinescu for (j = 0; j < s; j++) scheme->beta[i] -= b[i * s + j] * CPowF(c[j], p + 1); 1603d177a5cSEmil Constantinescu } 161792fecdfSBarry Smith PetscCallBLAS("LAPACKgetrs", LAPACKgetrs_("No transpose", &m, &one, ImV, &n, ipiv, scheme->beta + 1, &n, &info)); 16208401ef6SPierre Jolivet PetscCheck(info >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "Bad argument to GETRS"); 16308401ef6SPierre Jolivet PetscCheck(info <= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "Should not happen"); 1643d177a5cSEmil Constantinescu 1653d177a5cSEmil Constantinescu /* Build stage_error vector 1663d177a5cSEmil Constantinescu xi = glm.c.^(p+1)/factorial(p+1) - glm.A*glm.c.^p/factorial(p) + glm.U(:,2:end)*e.alpha; 1673d177a5cSEmil Constantinescu */ 1683d177a5cSEmil Constantinescu for (i = 0; i < s; i++) { 1693d177a5cSEmil Constantinescu scheme->stage_error[i] = CPowF(c[i], p + 1); 1703d177a5cSEmil Constantinescu for (j = 0; j < s; j++) scheme->stage_error[i] -= a[i * s + j] * CPowF(c[j], p); 1713d177a5cSEmil Constantinescu for (j = 1; j < r; j++) scheme->stage_error[i] += u[i * r + j] * scheme->alpha[j]; 1723d177a5cSEmil Constantinescu } 1733d177a5cSEmil Constantinescu 1743d177a5cSEmil Constantinescu /* alpha[0] (epsilon in B,J,W 2007) 1753d177a5cSEmil Constantinescu epsilon = 1/factorial(p+1) - B(1,:)*c.^p/factorial(p) + V(1,2:end)*e.alpha; 1763d177a5cSEmil Constantinescu */ 1773d177a5cSEmil Constantinescu scheme->alpha[0] = 1. / Factorial(p + 1); 1783d177a5cSEmil Constantinescu for (j = 0; j < s; j++) scheme->alpha[0] -= b[0 * s + j] * CPowF(c[j], p); 1793d177a5cSEmil Constantinescu for (j = 1; j < r; j++) scheme->alpha[0] += v[0 * r + j] * scheme->alpha[j]; 1803d177a5cSEmil Constantinescu 181dd8e379bSPierre Jolivet /* right-hand side for gamma (glm.B(2:end,:)*e.xi - e.epsilon*eye(s-1,1)) */ 1823d177a5cSEmil Constantinescu for (i = 1; i < r; i++) { 1833d177a5cSEmil Constantinescu scheme->gamma[i] = (i == 1 ? -1. : 0) * scheme->alpha[0]; 1843d177a5cSEmil Constantinescu for (j = 0; j < s; j++) scheme->gamma[i] += b[i * s + j] * scheme->stage_error[j]; 1853d177a5cSEmil Constantinescu } 186792fecdfSBarry Smith PetscCallBLAS("LAPACKgetrs", LAPACKgetrs_("No transpose", &m, &one, ImV, &n, ipiv, scheme->gamma + 1, &n, &info)); 18708401ef6SPierre Jolivet PetscCheck(info >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "Bad argument to GETRS"); 18808401ef6SPierre Jolivet PetscCheck(info <= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "Should not happen"); 1893d177a5cSEmil Constantinescu 1903d177a5cSEmil Constantinescu /* beta[0] (rho in B,J,W 2007) 1913d177a5cSEmil Constantinescu e.rho = 1/factorial(p+2) - glm.B(1,:)*glm.c.^(p+1)/factorial(p+1) ... 1923d177a5cSEmil Constantinescu + glm.V(1,2:end)*e.beta;% - e.epsilon; 1933d177a5cSEmil Constantinescu % Note: The paper (B,J,W 2007) includes the last term in their definition 1943d177a5cSEmil Constantinescu * */ 1953d177a5cSEmil Constantinescu scheme->beta[0] = 1. / Factorial(p + 2); 1963d177a5cSEmil Constantinescu for (j = 0; j < s; j++) scheme->beta[0] -= b[0 * s + j] * CPowF(c[j], p + 1); 1973d177a5cSEmil Constantinescu for (j = 1; j < r; j++) scheme->beta[0] += v[0 * r + j] * scheme->beta[j]; 1983d177a5cSEmil Constantinescu 1993d177a5cSEmil Constantinescu /* gamma[0] (sigma in B,J,W 2007) 2003d177a5cSEmil Constantinescu * e.sigma = glm.B(1,:)*e.xi + glm.V(1,2:end)*e.gamma; 2013d177a5cSEmil Constantinescu * */ 2023d177a5cSEmil Constantinescu scheme->gamma[0] = 0.0; 2033d177a5cSEmil Constantinescu for (j = 0; j < s; j++) scheme->gamma[0] += b[0 * s + j] * scheme->stage_error[j]; 2043d177a5cSEmil Constantinescu for (j = 1; j < r; j++) scheme->gamma[0] += v[0 * s + j] * scheme->gamma[j]; 2053d177a5cSEmil Constantinescu 2063d177a5cSEmil Constantinescu /* Assemble H 20763a3b9bcSJacob Faibussowitsch * % " PetscInt_FMT "etermine the error estimators phi 2083d177a5cSEmil Constantinescu H = [[cpow(glm.c,p) + C*e.alpha] [cpow(glm.c,p+1) + C*e.beta] ... 2093d177a5cSEmil Constantinescu [e.xi - C*(e.gamma + 0*e.epsilon*eye(s-1,1))]]'; 2103d177a5cSEmil Constantinescu % Paper has formula above without the 0, but that term must be left 2113d177a5cSEmil Constantinescu % out to satisfy the conditions they propose and to make the 2123d177a5cSEmil Constantinescu % example schemes work 2133d177a5cSEmil Constantinescu e.H = H; 2143d177a5cSEmil Constantinescu e.phi = (H \ [1 0 0;1 1 0;0 0 -1])'; 2153d177a5cSEmil Constantinescu e.psi = -e.phi*C; 2163d177a5cSEmil Constantinescu * */ 2173d177a5cSEmil Constantinescu for (j = 0; j < s; j++) { 2183d177a5cSEmil Constantinescu H[0 + j * 3] = CPowF(c[j], p); 2193d177a5cSEmil Constantinescu H[1 + j * 3] = CPowF(c[j], p + 1); 2203d177a5cSEmil Constantinescu H[2 + j * 3] = scheme->stage_error[j]; 2213d177a5cSEmil Constantinescu for (k = 1; k < r; k++) { 2223d177a5cSEmil Constantinescu H[0 + j * 3] += CPowF(c[j], k - 1) * scheme->alpha[k]; 2233d177a5cSEmil Constantinescu H[1 + j * 3] += CPowF(c[j], k - 1) * scheme->beta[k]; 2243d177a5cSEmil Constantinescu H[2 + j * 3] -= CPowF(c[j], k - 1) * scheme->gamma[k]; 2253d177a5cSEmil Constantinescu } 2263d177a5cSEmil Constantinescu } 2279371c9d4SSatish Balay bmat[0 + 0 * ss] = 1.; 2289371c9d4SSatish Balay bmat[0 + 1 * ss] = 0.; 2299371c9d4SSatish Balay bmat[0 + 2 * ss] = 0.; 2309371c9d4SSatish Balay bmat[1 + 0 * ss] = 1.; 2319371c9d4SSatish Balay bmat[1 + 1 * ss] = 1.; 2329371c9d4SSatish Balay bmat[1 + 2 * ss] = 0.; 2339371c9d4SSatish Balay bmat[2 + 0 * ss] = 0.; 2349371c9d4SSatish Balay bmat[2 + 1 * ss] = 0.; 2359371c9d4SSatish Balay bmat[2 + 2 * ss] = -1.; 2363d177a5cSEmil Constantinescu m = 3; 2379566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(s, &n)); 2389566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(ss, &ldb)); 2393d177a5cSEmil Constantinescu rcond = 1e-12; 2403d177a5cSEmil Constantinescu #if defined(PETSC_USE_COMPLEX) 2413d177a5cSEmil Constantinescu /* ZGELSS( M, N, NRHS, A, LDA, B, LDB, S, RCOND, RANK, WORK, LWORK, RWORK, INFO) */ 242792fecdfSBarry Smith PetscCallBLAS("LAPACKgelss", LAPACKgelss_(&m, &n, &m, H, &m, bmat, &ldb, sing, &rcond, &rank, workscalar, &lwork, workreal, &info)); 2433d177a5cSEmil Constantinescu #else 2443d177a5cSEmil Constantinescu /* DGELSS( M, N, NRHS, A, LDA, B, LDB, S, RCOND, RANK, WORK, LWORK, INFO) */ 245792fecdfSBarry Smith PetscCallBLAS("LAPACKgelss", LAPACKgelss_(&m, &n, &m, H, &m, bmat, &ldb, sing, &rcond, &rank, workscalar, &lwork, &info)); 2463d177a5cSEmil Constantinescu #endif 24708401ef6SPierre Jolivet PetscCheck(info >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "Bad argument to GELSS"); 24808401ef6SPierre Jolivet PetscCheck(info <= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "SVD failed to converge"); 2493d177a5cSEmil Constantinescu 2503d177a5cSEmil Constantinescu for (j = 0; j < 3; j++) { 2513d177a5cSEmil Constantinescu for (k = 0; k < s; k++) scheme->phi[k + j * s] = bmat[k + j * ss]; 2523d177a5cSEmil Constantinescu } 2533d177a5cSEmil Constantinescu 2543d177a5cSEmil Constantinescu /* the other part of the error estimator, psi in B,J,W 2007 */ 2553d177a5cSEmil Constantinescu scheme->psi[0 * r + 0] = 0.; 2563d177a5cSEmil Constantinescu scheme->psi[1 * r + 0] = 0.; 2573d177a5cSEmil Constantinescu scheme->psi[2 * r + 0] = 0.; 2583d177a5cSEmil Constantinescu for (j = 1; j < r; j++) { 2593d177a5cSEmil Constantinescu scheme->psi[0 * r + j] = 0.; 2603d177a5cSEmil Constantinescu scheme->psi[1 * r + j] = 0.; 2613d177a5cSEmil Constantinescu scheme->psi[2 * r + j] = 0.; 2623d177a5cSEmil Constantinescu for (k = 0; k < s; k++) { 2633d177a5cSEmil Constantinescu scheme->psi[0 * r + j] -= CPowF(c[k], j - 1) * scheme->phi[0 * s + k]; 2643d177a5cSEmil Constantinescu scheme->psi[1 * r + j] -= CPowF(c[k], j - 1) * scheme->phi[1 * s + k]; 2653d177a5cSEmil Constantinescu scheme->psi[2 * r + j] -= CPowF(c[k], j - 1) * scheme->phi[2 * s + k]; 2663d177a5cSEmil Constantinescu } 2673d177a5cSEmil Constantinescu } 2689566063dSJacob Faibussowitsch PetscCall(PetscFree7(ImV, H, bmat, workscalar, workreal, sing, ipiv)); 2693d177a5cSEmil Constantinescu } 2703d177a5cSEmil Constantinescu /* Check which properties are satisfied */ 2713d177a5cSEmil Constantinescu scheme->stiffly_accurate = PETSC_TRUE; 2723d177a5cSEmil Constantinescu if (scheme->c[s - 1] != 1.) scheme->stiffly_accurate = PETSC_FALSE; 2739371c9d4SSatish Balay for (j = 0; j < s; j++) 2749371c9d4SSatish Balay if (a[(s - 1) * s + j] != b[j]) scheme->stiffly_accurate = PETSC_FALSE; 2759371c9d4SSatish Balay for (j = 0; j < r; j++) 2769371c9d4SSatish Balay if (u[(s - 1) * r + j] != v[j]) scheme->stiffly_accurate = PETSC_FALSE; 2773d177a5cSEmil Constantinescu scheme->fsal = scheme->stiffly_accurate; /* FSAL is stronger */ 2789371c9d4SSatish Balay for (j = 0; j < s - 1; j++) 2799371c9d4SSatish Balay if (r > 1 && b[1 * s + j] != 0.) scheme->fsal = PETSC_FALSE; 2803d177a5cSEmil Constantinescu if (b[1 * s + r - 1] != 1.) scheme->fsal = PETSC_FALSE; 2819371c9d4SSatish Balay for (j = 0; j < r; j++) 2829371c9d4SSatish Balay if (r > 1 && v[1 * r + j] != 0.) scheme->fsal = PETSC_FALSE; 2833d177a5cSEmil Constantinescu 2843d177a5cSEmil Constantinescu *inscheme = scheme; 2853ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2863d177a5cSEmil Constantinescu } 2873d177a5cSEmil Constantinescu 288d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSGLLESchemeDestroy(TSGLLEScheme sc) 289d71ae5a4SJacob Faibussowitsch { 2903d177a5cSEmil Constantinescu PetscFunctionBegin; 2919566063dSJacob Faibussowitsch PetscCall(PetscFree5(sc->c, sc->a, sc->b, sc->u, sc->v)); 2929566063dSJacob Faibussowitsch PetscCall(PetscFree6(sc->alpha, sc->beta, sc->gamma, sc->phi, sc->psi, sc->stage_error)); 2939566063dSJacob Faibussowitsch PetscCall(PetscFree(sc)); 2943ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2953d177a5cSEmil Constantinescu } 2963d177a5cSEmil Constantinescu 297d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSGLLEDestroy_Default(TS_GLLE *gl) 298d71ae5a4SJacob Faibussowitsch { 2993d177a5cSEmil Constantinescu PetscInt i; 3003d177a5cSEmil Constantinescu 3013d177a5cSEmil Constantinescu PetscFunctionBegin; 3023d177a5cSEmil Constantinescu for (i = 0; i < gl->nschemes; i++) { 3039566063dSJacob Faibussowitsch if (gl->schemes[i]) PetscCall(TSGLLESchemeDestroy(gl->schemes[i])); 3043d177a5cSEmil Constantinescu } 3059566063dSJacob Faibussowitsch PetscCall(PetscFree(gl->schemes)); 3063d177a5cSEmil Constantinescu gl->nschemes = 0; 3079566063dSJacob Faibussowitsch PetscCall(PetscMemzero(gl->type_name, sizeof(gl->type_name))); 3083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3093d177a5cSEmil Constantinescu } 3103d177a5cSEmil Constantinescu 311d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSGLLEViewTable_Private(PetscViewer viewer, PetscInt m, PetscInt n, const PetscScalar a[], const char name[]) 312d71ae5a4SJacob Faibussowitsch { 313*9f196a02SMartin Diehl PetscBool isascii; 3143d177a5cSEmil Constantinescu PetscInt i, j; 3153d177a5cSEmil Constantinescu 3163d177a5cSEmil Constantinescu PetscFunctionBegin; 317*9f196a02SMartin Diehl PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii)); 318*9f196a02SMartin Diehl if (isascii) { 3199566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "%30s = [", name)); 3203d177a5cSEmil Constantinescu for (i = 0; i < m; i++) { 3219566063dSJacob Faibussowitsch if (i) PetscCall(PetscViewerASCIIPrintf(viewer, "%30s [", "")); 3229566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE)); 32348a46eb9SPierre Jolivet for (j = 0; j < n; j++) PetscCall(PetscViewerASCIIPrintf(viewer, " %12.8g", (double)PetscRealPart(a[i * n + j]))); 3249566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "]\n")); 3259566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE)); 3263d177a5cSEmil Constantinescu } 3273d177a5cSEmil Constantinescu } 3283ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3293d177a5cSEmil Constantinescu } 3303d177a5cSEmil Constantinescu 331d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSGLLESchemeView(TSGLLEScheme sc, PetscBool view_details, PetscViewer viewer) 332d71ae5a4SJacob Faibussowitsch { 333*9f196a02SMartin Diehl PetscBool isascii; 3343d177a5cSEmil Constantinescu 3353d177a5cSEmil Constantinescu PetscFunctionBegin; 336*9f196a02SMartin Diehl PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii)); 337*9f196a02SMartin Diehl if (isascii) { 33863a3b9bcSJacob 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)); 3399566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 3409566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Stiffly accurate: %s, FSAL: %s\n", sc->stiffly_accurate ? "yes" : "no", sc->fsal ? "yes" : "no")); 3419371c9d4SSatish Balay PetscCall(PetscViewerASCIIPrintf(viewer, "Leading error constants: %10.3e %10.3e %10.3e\n", (double)PetscRealPart(sc->alpha[0]), (double)PetscRealPart(sc->beta[0]), (double)PetscRealPart(sc->gamma[0]))); 3429566063dSJacob Faibussowitsch PetscCall(TSGLLEViewTable_Private(viewer, 1, sc->s, sc->c, "Abscissas c")); 3433d177a5cSEmil Constantinescu if (view_details) { 3449566063dSJacob Faibussowitsch PetscCall(TSGLLEViewTable_Private(viewer, sc->s, sc->s, sc->a, "A")); 3459566063dSJacob Faibussowitsch PetscCall(TSGLLEViewTable_Private(viewer, sc->r, sc->s, sc->b, "B")); 3469566063dSJacob Faibussowitsch PetscCall(TSGLLEViewTable_Private(viewer, sc->s, sc->r, sc->u, "U")); 3479566063dSJacob Faibussowitsch PetscCall(TSGLLEViewTable_Private(viewer, sc->r, sc->r, sc->v, "V")); 3483d177a5cSEmil Constantinescu 3499566063dSJacob Faibussowitsch PetscCall(TSGLLEViewTable_Private(viewer, 3, sc->s, sc->phi, "Error estimate phi")); 3509566063dSJacob Faibussowitsch PetscCall(TSGLLEViewTable_Private(viewer, 3, sc->r, sc->psi, "Error estimate psi")); 3519566063dSJacob Faibussowitsch PetscCall(TSGLLEViewTable_Private(viewer, 1, sc->r, sc->alpha, "Modify alpha")); 3529566063dSJacob Faibussowitsch PetscCall(TSGLLEViewTable_Private(viewer, 1, sc->r, sc->beta, "Modify beta")); 3539566063dSJacob Faibussowitsch PetscCall(TSGLLEViewTable_Private(viewer, 1, sc->r, sc->gamma, "Modify gamma")); 3549566063dSJacob Faibussowitsch PetscCall(TSGLLEViewTable_Private(viewer, 1, sc->s, sc->stage_error, "Stage error xi")); 3553d177a5cSEmil Constantinescu } 3569566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 35798921bdaSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Viewer type %s not supported", ((PetscObject)viewer)->type_name); 3583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3593d177a5cSEmil Constantinescu } 3603d177a5cSEmil Constantinescu 361d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSGLLEEstimateHigherMoments_Default(TSGLLEScheme sc, PetscReal h, Vec Ydot[], Vec Xold[], Vec hm[]) 362d71ae5a4SJacob Faibussowitsch { 3633d177a5cSEmil Constantinescu PetscInt i; 3643d177a5cSEmil Constantinescu 3653d177a5cSEmil Constantinescu PetscFunctionBegin; 366cad9d221SBarry Smith PetscCheck(sc->r <= 64 && sc->s <= 64, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Ridiculous number of stages or items passed between stages"); 3673d177a5cSEmil Constantinescu /* build error vectors*/ 3683d177a5cSEmil Constantinescu for (i = 0; i < 3; i++) { 3693d177a5cSEmil Constantinescu PetscScalar phih[64]; 3703d177a5cSEmil Constantinescu PetscInt j; 3713d177a5cSEmil Constantinescu for (j = 0; j < sc->s; j++) phih[j] = sc->phi[i * sc->s + j] * h; 3729566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(hm[i])); 3739566063dSJacob Faibussowitsch PetscCall(VecMAXPY(hm[i], sc->s, phih, Ydot)); 3749566063dSJacob Faibussowitsch PetscCall(VecMAXPY(hm[i], sc->r, &sc->psi[i * sc->r], Xold)); 3753d177a5cSEmil Constantinescu } 3763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3773d177a5cSEmil Constantinescu } 3783d177a5cSEmil Constantinescu 379d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSGLLECompleteStep_Rescale(TSGLLEScheme sc, PetscReal h, TSGLLEScheme next_sc, PetscReal next_h, Vec Ydot[], Vec Xold[], Vec X[]) 380d71ae5a4SJacob Faibussowitsch { 3813d177a5cSEmil Constantinescu PetscScalar brow[32], vrow[32]; 3823d177a5cSEmil Constantinescu PetscInt i, j, r, s; 3833d177a5cSEmil Constantinescu 3843d177a5cSEmil Constantinescu PetscFunctionBegin; 3853d177a5cSEmil Constantinescu /* Build the new solution from (X,Ydot) */ 3863d177a5cSEmil Constantinescu r = sc->r; 3873d177a5cSEmil Constantinescu s = sc->s; 3883d177a5cSEmil Constantinescu for (i = 0; i < r; i++) { 3899566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(X[i])); 3903d177a5cSEmil Constantinescu for (j = 0; j < s; j++) brow[j] = h * sc->b[i * s + j]; 3919566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X[i], s, brow, Ydot)); 3923d177a5cSEmil Constantinescu for (j = 0; j < r; j++) vrow[j] = sc->v[i * r + j]; 3939566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X[i], r, vrow, Xold)); 3943d177a5cSEmil Constantinescu } 3953ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3963d177a5cSEmil Constantinescu } 3973d177a5cSEmil Constantinescu 398d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSGLLECompleteStep_RescaleAndModify(TSGLLEScheme sc, PetscReal h, TSGLLEScheme next_sc, PetscReal next_h, Vec Ydot[], Vec Xold[], Vec X[]) 399d71ae5a4SJacob Faibussowitsch { 4003d177a5cSEmil Constantinescu PetscScalar brow[32], vrow[32]; 4013d177a5cSEmil Constantinescu PetscReal ratio; 4023d177a5cSEmil Constantinescu PetscInt i, j, p, r, s; 4033d177a5cSEmil Constantinescu 4043d177a5cSEmil Constantinescu PetscFunctionBegin; 4053d177a5cSEmil Constantinescu /* Build the new solution from (X,Ydot) */ 4063d177a5cSEmil Constantinescu p = sc->p; 4073d177a5cSEmil Constantinescu r = sc->r; 4083d177a5cSEmil Constantinescu s = sc->s; 4093d177a5cSEmil Constantinescu ratio = next_h / h; 4103d177a5cSEmil Constantinescu for (i = 0; i < r; i++) { 4119566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(X[i])); 4123d177a5cSEmil Constantinescu for (j = 0; j < s; j++) { 4139371c9d4SSatish Balay brow[j] = h * (PetscPowRealInt(ratio, i) * sc->b[i * s + j] + (PetscPowRealInt(ratio, i) - PetscPowRealInt(ratio, p + 1)) * (+sc->alpha[i] * sc->phi[0 * s + j]) + (PetscPowRealInt(ratio, i) - PetscPowRealInt(ratio, p + 2)) * (+sc->beta[i] * sc->phi[1 * s + j] + sc->gamma[i] * sc->phi[2 * s + j])); 4143d177a5cSEmil Constantinescu } 4159566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X[i], s, brow, Ydot)); 4163d177a5cSEmil Constantinescu for (j = 0; j < r; j++) { 4179371c9d4SSatish Balay vrow[j] = (PetscPowRealInt(ratio, i) * sc->v[i * r + j] + (PetscPowRealInt(ratio, i) - PetscPowRealInt(ratio, p + 1)) * (+sc->alpha[i] * sc->psi[0 * r + j]) + (PetscPowRealInt(ratio, i) - PetscPowRealInt(ratio, p + 2)) * (+sc->beta[i] * sc->psi[1 * r + j] + sc->gamma[i] * sc->psi[2 * r + j])); 4183d177a5cSEmil Constantinescu } 4199566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X[i], r, vrow, Xold)); 4203d177a5cSEmil Constantinescu } 4213d177a5cSEmil Constantinescu if (r < next_sc->r) { 42208401ef6SPierre Jolivet PetscCheck(r + 1 == next_sc->r, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Cannot accommodate jump in r greater than 1"); 4239566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(X[r])); 4243d177a5cSEmil Constantinescu for (j = 0; j < s; j++) brow[j] = h * PetscPowRealInt(ratio, p + 1) * sc->phi[0 * s + j]; 4259566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X[r], s, brow, Ydot)); 4263d177a5cSEmil Constantinescu for (j = 0; j < r; j++) vrow[j] = PetscPowRealInt(ratio, p + 1) * sc->psi[0 * r + j]; 4279566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X[r], r, vrow, Xold)); 4283d177a5cSEmil Constantinescu } 4293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4303d177a5cSEmil Constantinescu } 4313d177a5cSEmil Constantinescu 432d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSGLLECreate_IRKS(TS ts) 433d71ae5a4SJacob Faibussowitsch { 4343d177a5cSEmil Constantinescu TS_GLLE *gl = (TS_GLLE *)ts->data; 4353d177a5cSEmil Constantinescu 4363d177a5cSEmil Constantinescu PetscFunctionBegin; 4373d177a5cSEmil Constantinescu gl->Destroy = TSGLLEDestroy_Default; 4383d177a5cSEmil Constantinescu gl->EstimateHigherMoments = TSGLLEEstimateHigherMoments_Default; 4393d177a5cSEmil Constantinescu gl->CompleteStep = TSGLLECompleteStep_RescaleAndModify; 4409566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(10, &gl->schemes)); 4413d177a5cSEmil Constantinescu gl->nschemes = 0; 4423d177a5cSEmil Constantinescu 4433d177a5cSEmil Constantinescu { 4443d177a5cSEmil Constantinescu /* p=1,q=1, r=s=2, A- and L-stable with error estimates of order 2 and 3 4453d177a5cSEmil Constantinescu * Listed in Butcher & Podhaisky 2006. On error estimation in general linear methods for stiff ODE. 4463d177a5cSEmil Constantinescu * irks(0.3,0,[.3,1],[1],1) 4473d177a5cSEmil Constantinescu * Note: can be made to have classical order (not stage order) 2 by replacing 0.3 with 1-sqrt(1/2) 4483d177a5cSEmil Constantinescu * but doing so would sacrifice the error estimator. 4493d177a5cSEmil Constantinescu */ 4503d177a5cSEmil Constantinescu const PetscScalar c[2] = {3. / 10., 1.}; 4519371c9d4SSatish Balay const PetscScalar a[2][2] = { 4529371c9d4SSatish Balay {3. / 10., 0 }, 4539371c9d4SSatish Balay {7. / 10., 3. / 10.} 4549371c9d4SSatish Balay }; 4559371c9d4SSatish Balay const PetscScalar b[2][2] = { 4569371c9d4SSatish Balay {7. / 10., 3. / 10.}, 4579371c9d4SSatish Balay {0, 1 } 4589371c9d4SSatish Balay }; 4599371c9d4SSatish Balay const PetscScalar u[2][2] = { 4609371c9d4SSatish Balay {1, 0}, 4619371c9d4SSatish Balay {1, 0} 4629371c9d4SSatish Balay }; 4639371c9d4SSatish Balay const PetscScalar v[2][2] = { 4649371c9d4SSatish Balay {1, 0}, 4659371c9d4SSatish Balay {0, 0} 4669371c9d4SSatish Balay }; 4679566063dSJacob Faibussowitsch PetscCall(TSGLLESchemeCreate(1, 1, 2, 2, c, *a, *b, *u, *v, &gl->schemes[gl->nschemes++])); 4683d177a5cSEmil Constantinescu } 4693d177a5cSEmil Constantinescu 4703d177a5cSEmil Constantinescu { 4713d177a5cSEmil Constantinescu /* p=q=2, r=s=3: irks(4/9,0,[1:3]/3,[0.33852],1) */ 4723d177a5cSEmil Constantinescu /* http://www.math.auckland.ac.nz/~hpod/atlas/i2a.html */ 47397a1619fSSatish Balay const PetscScalar c[3] = {1. / 3., 2. / 3., 1}; 47497a1619fSSatish Balay const PetscScalar a[3][3] = { 47597a1619fSSatish Balay {4. / 9., 0, 0 }, 47697a1619fSSatish Balay {1.03750643704090e+00, 4. / 9., 0 }, 47797a1619fSSatish Balay {7.67024779410304e-01, -3.81140216918943e-01, 4. / 9.} 47897a1619fSSatish Balay }; 47997a1619fSSatish Balay const PetscScalar b[3][3] = { 48097a1619fSSatish Balay {0.767024779410304, -0.381140216918943, 4. / 9. }, 48197a1619fSSatish Balay {0.000000000000000, 0.000000000000000, 1.000000000000000}, 48297a1619fSSatish Balay {-2.075048385225385, 0.621728385225383, 1.277197204924873} 48397a1619fSSatish Balay }; 48497a1619fSSatish Balay const PetscScalar u[3][3] = { 48597a1619fSSatish Balay {1.0000000000000000, -0.1111111111111109, -0.0925925925925922}, 48697a1619fSSatish Balay {1.0000000000000000, -0.8152842148186744, -0.4199095530877056}, 48797a1619fSSatish Balay {1.0000000000000000, 0.1696709930641948, 0.0539741070314165 } 48897a1619fSSatish Balay }; 48997a1619fSSatish Balay const PetscScalar v[3][3] = { 49097a1619fSSatish Balay {1.0000000000000000, 0.1696709930641948, 0.0539741070314165}, 49197a1619fSSatish Balay {0.000000000000000, 0.000000000000000, 0.000000000000000 }, 49297a1619fSSatish Balay {0.000000000000000, 0.176122795075129, 0.000000000000000 } 49397a1619fSSatish Balay }; 4949566063dSJacob Faibussowitsch PetscCall(TSGLLESchemeCreate(2, 2, 3, 3, c, *a, *b, *u, *v, &gl->schemes[gl->nschemes++])); 4953d177a5cSEmil Constantinescu } 4963d177a5cSEmil Constantinescu { 4973d177a5cSEmil Constantinescu /* p=q=3, r=s=4: irks(9/40,0,[1:4]/4,[0.3312 1.0050],[0.49541 1;1 0]) */ 49897a1619fSSatish Balay const PetscScalar c[4] = {0.25, 0.5, 0.75, 1.0}; 49997a1619fSSatish Balay const PetscScalar a[4][4] = { 50097a1619fSSatish Balay {9. / 40., 0, 0, 0 }, 50197a1619fSSatish Balay {2.11286958887701e-01, 9. / 40., 0, 0 }, 50297a1619fSSatish Balay {9.46338294287584e-01, -3.42942861246094e-01, 9. / 40., 0 }, 50397a1619fSSatish Balay {0.521490453970721, -0.662474225622980, 0.490476425459734, 9. / 40.} 50497a1619fSSatish Balay }; 50597a1619fSSatish Balay const PetscScalar b[4][4] = { 50697a1619fSSatish Balay {0.521490453970721, -0.662474225622980, 0.490476425459734, 9. / 40. }, 50797a1619fSSatish Balay {0.000000000000000, 0.000000000000000, 0.000000000000000, 1.000000000000000}, 50897a1619fSSatish Balay {-0.084677029310348, 1.390757514776085, -1.568157386206001, 2.023192696767826}, 50997a1619fSSatish Balay {0.465383797936408, 1.478273530625148, -1.930836081010182, 1.644872111193354} 51097a1619fSSatish Balay }; 51197a1619fSSatish Balay const PetscScalar u[4][4] = { 51297a1619fSSatish Balay {1.00000000000000000, 0.02500000000001035, -0.02499999999999053, -0.00442708333332865}, 51397a1619fSSatish Balay {1.00000000000000000, 0.06371304111232945, -0.04032173972189845, -0.01389438413189452}, 51497a1619fSSatish Balay {1.00000000000000000, -0.07839543304147778, 0.04738685705116663, 0.02032603595928376 }, 51597a1619fSSatish Balay {1.00000000000000000, 0.42550734619251651, 0.10800718022400080, -0.01726712647760034} 51697a1619fSSatish Balay }; 51797a1619fSSatish Balay const PetscScalar v[4][4] = { 51897a1619fSSatish Balay {1.00000000000000000, 0.42550734619251651, 0.10800718022400080, -0.01726712647760034}, 51997a1619fSSatish Balay {0.000000000000000, 0.000000000000000, 0.000000000000000, 0.000000000000000 }, 52097a1619fSSatish Balay {0.000000000000000, -1.761115796027561, -0.521284157173780, 0.258249384305463 }, 52197a1619fSSatish Balay {0.000000000000000, -1.657693358744728, -1.052227765232394, 0.521284157173780 } 52297a1619fSSatish Balay }; 5239566063dSJacob Faibussowitsch PetscCall(TSGLLESchemeCreate(3, 3, 4, 4, c, *a, *b, *u, *v, &gl->schemes[gl->nschemes++])); 5243d177a5cSEmil Constantinescu } 5253d177a5cSEmil Constantinescu { 5263d177a5cSEmil Constantinescu /* p=q=4, r=s=5: 5273d177a5cSEmil Constantinescu irks(3/11,0,[1:5]/5, [0.1715 -0.1238 0.6617],... 5283d177a5cSEmil Constantinescu [ -0.0812 0.4079 1.0000 5293d177a5cSEmil Constantinescu 1.0000 0 0 5303d177a5cSEmil Constantinescu 0.8270 1.0000 0]) 5313d177a5cSEmil Constantinescu */ 53297a1619fSSatish Balay const PetscScalar c[5] = {0.2, 0.4, 0.6, 0.8, 1.0}; 53397a1619fSSatish Balay const PetscScalar a[5][5] = { 53497a1619fSSatish Balay {2.72727272727352e-01, 0.00000000000000e+00, 0.00000000000000e+00, 0.00000000000000e+00, 0.00000000000000e+00}, 53597a1619fSSatish Balay {-1.03980153733431e-01, 2.72727272727405e-01, 0.00000000000000e+00, 0.00000000000000e+00, 0.00000000000000e+00}, 53697a1619fSSatish Balay {-1.58615400341492e+00, 7.44168951881122e-01, 2.72727272727309e-01, 0.00000000000000e+00, 0.00000000000000e+00}, 53797a1619fSSatish Balay {-8.73658042865628e-01, 5.37884671894595e-01, -1.63298538799523e-01, 2.72727272726996e-01, 0.00000000000000e+00}, 53897a1619fSSatish Balay {2.95489397443992e-01, -1.18481693910097e+00, -6.68029812659953e-01, 1.00716687860943e+00, 2.72727272727288e-01} 53997a1619fSSatish Balay }; 54097a1619fSSatish Balay const PetscScalar b[5][5] = { 54197a1619fSSatish Balay {2.95489397443992e-01, -1.18481693910097e+00, -6.68029812659953e-01, 1.00716687860943e+00, 2.72727272727288e-01}, 54297a1619fSSatish Balay {0.00000000000000e+00, 1.11022302462516e-16, -2.22044604925031e-16, 0.00000000000000e+00, 1.00000000000000e+00}, 54397a1619fSSatish Balay {-4.05882503986005e+00, -4.00924006567769e+00, -1.38930610972481e+00, 4.45223930308488e+00, 6.32331093108427e-01}, 54497a1619fSSatish Balay {8.35690179937017e+00, -2.26640927349732e+00, 6.86647884973826e+00, -5.22595158025740e+00, 4.50893068837431e+00}, 54597a1619fSSatish Balay {1.27656267027479e+01, 2.80882153840821e+00, 8.91173096522890e+00, -1.07936444078906e+01, 4.82534148988854e+00} 54697a1619fSSatish Balay }; 54797a1619fSSatish Balay const PetscScalar u[5][5] = { 54897a1619fSSatish Balay {1.00000000000000e+00, -7.27272727273551e-02, -3.45454545454419e-02, -4.12121212119565e-03, -2.96969696964014e-04}, 54997a1619fSSatish Balay {1.00000000000000e+00, 2.31252881006154e-01, -8.29487834416481e-03, -9.07191207681020e-03, -1.70378403743473e-03}, 55097a1619fSSatish Balay {1.00000000000000e+00, 1.16925777880663e+00, 3.59268562942635e-02, -4.09013451730615e-02, -1.02411119670164e-02}, 55197a1619fSSatish Balay {1.00000000000000e+00, 1.02634463704356e+00, 1.59375044913405e-01, 1.89673015035370e-03, -4.89987231897569e-03}, 55297a1619fSSatish Balay {1.00000000000000e+00, 1.27746320298021e+00, 2.37186008132728e-01, -8.28694373940065e-02, -5.34396510196430e-02} 55397a1619fSSatish Balay }; 55497a1619fSSatish Balay const PetscScalar v[5][5] = { 55597a1619fSSatish Balay {1.00000000000000e+00, 1.27746320298021e+00, 2.37186008132728e-01, -8.28694373940065e-02, -5.34396510196430e-02}, 55697a1619fSSatish Balay {0.00000000000000e+00, -1.77635683940025e-15, -1.99840144432528e-15, -9.99200722162641e-16, -3.33066907387547e-16}, 55797a1619fSSatish Balay {0.00000000000000e+00, 4.37280081906924e+00, 5.49221645016377e-02, -8.88913177394943e-02, 1.12879077989154e-01 }, 55897a1619fSSatish Balay {0.00000000000000e+00, -1.22399504837280e+01, -5.21287338448645e+00, -8.03952325565291e-01, 4.60298678047147e-01 }, 55997a1619fSSatish Balay {0.00000000000000e+00, -1.85178762883829e+01, -5.21411849862624e+00, -1.04283436528809e+00, 7.49030161063651e-01 } 56097a1619fSSatish Balay }; 5619566063dSJacob Faibussowitsch PetscCall(TSGLLESchemeCreate(4, 4, 5, 5, c, *a, *b, *u, *v, &gl->schemes[gl->nschemes++])); 5623d177a5cSEmil Constantinescu } 5633d177a5cSEmil Constantinescu { 5643d177a5cSEmil Constantinescu /* p=q=5, r=s=6; 5653d177a5cSEmil Constantinescu irks(1/3,0,[1:6]/6,... 5663d177a5cSEmil Constantinescu [-0.0489 0.4228 -0.8814 0.9021],... 5673d177a5cSEmil Constantinescu [-0.3474 -0.6617 0.6294 0.2129 5683d177a5cSEmil Constantinescu 0.0044 -0.4256 -0.1427 -0.8936 5693d177a5cSEmil Constantinescu -0.8267 0.4821 0.1371 -0.2557 5703d177a5cSEmil Constantinescu -0.4426 -0.3855 -0.7514 0.3014]) 5713d177a5cSEmil Constantinescu */ 57297a1619fSSatish Balay const PetscScalar c[6] = {1. / 6, 2. / 6, 3. / 6, 4. / 6, 5. / 6, 1.}; 57397a1619fSSatish Balay const PetscScalar a[6][6] = { 57497a1619fSSatish Balay {3.33333333333940e-01, 0, 0, 0, 0, 0 }, 57597a1619fSSatish Balay {-8.64423857333350e-02, 3.33333333332888e-01, 0, 0, 0, 0 }, 57697a1619fSSatish Balay {-2.16850174258252e+00, -2.23619072028839e+00, 3.33333333335204e-01, 0, 0, 0 }, 57797a1619fSSatish Balay {-4.73160970138997e+00, -3.89265344629268e+00, -2.76318716520933e-01, 3.33333333335759e-01, 0, 0 }, 57897a1619fSSatish Balay {-6.75187540297338e+00, -7.90756533769377e+00, 7.90245051802259e-01, -4.48352364517632e-01, 3.33333333328483e-01, 0 }, 57997a1619fSSatish Balay {-4.26488287921548e+00, -1.19320395589302e+01, 3.38924509887755e+00, -2.23969848002481e+00, 6.62807710124007e-01, 3.33333333335440e-01} 58097a1619fSSatish Balay }; 58197a1619fSSatish Balay const PetscScalar b[6][6] = { 58297a1619fSSatish Balay {-4.26488287921548e+00, -1.19320395589302e+01, 3.38924509887755e+00, -2.23969848002481e+00, 6.62807710124007e-01, 3.33333333335440e-01 }, 58397a1619fSSatish Balay {-8.88178419700125e-16, 4.44089209850063e-16, -1.54737334057131e-15, -8.88178419700125e-16, 0.00000000000000e+00, 1.00000000000001e+00 }, 58497a1619fSSatish Balay {-2.87780425770651e+01, -1.13520448264971e+01, 2.62002318943161e+01, 2.56943874812797e+01, -3.06702268304488e+01, 6.68067773510103e+00 }, 58597a1619fSSatish Balay {5.47971245256474e+01, 6.80366875868284e+01, -6.50952588861999e+01, -8.28643975339097e+01, 8.17416943896414e+01, -1.17819043489036e+01}, 58697a1619fSSatish Balay {-2.33332114788869e+02, 6.12942539462634e+01, -4.91850135865944e+01, 1.82716844135480e+02, -1.29788173979395e+02, 3.09968095651099e+01 }, 58797a1619fSSatish Balay {-1.72049132343751e+02, 8.60194713593999e+00, 7.98154219170200e-01, 1.50371386053218e+02, -1.18515423962066e+02, 2.50898277784663e+01 } 58897a1619fSSatish Balay }; 58997a1619fSSatish Balay const PetscScalar u[6][6] = { 59097a1619fSSatish Balay {1.00000000000000e+00, -1.66666666666870e-01, -4.16666666664335e-02, -3.85802469124815e-03, -2.25051440302250e-04, -9.64506172339142e-06}, 59197a1619fSSatish Balay {1.00000000000000e+00, 8.64423857327162e-02, -4.11484912671353e-02, -1.11450903217645e-02, -1.47651050487126e-03, -1.34395070766826e-04}, 59297a1619fSSatish Balay {1.00000000000000e+00, 4.57135912953434e+00, 1.06514719719137e+00, 1.33517564218007e-01, 1.11365952968659e-02, 6.12382756769504e-04 }, 59397a1619fSSatish Balay {1.00000000000000e+00, 9.23391519753404e+00, 2.22431212392095e+00, 2.91823807741891e-01, 2.52058456411084e-02, 1.22800542949647e-03 }, 59497a1619fSSatish Balay {1.00000000000000e+00, 1.48175480533865e+01, 3.73439117461835e+00, 5.14648336541804e-01, 4.76430038853402e-02, 2.56798515502156e-03 }, 59597a1619fSSatish Balay {1.00000000000000e+00, 1.50512347758335e+01, 4.10099701165164e+00, 5.66039141003603e-01, 3.91213893800891e-02, -2.99136269067853e-03} 59697a1619fSSatish Balay }; 59797a1619fSSatish Balay const PetscScalar v[6][6] = { 59897a1619fSSatish Balay {1.00000000000000e+00, 1.50512347758335e+01, 4.10099701165164e+00, 5.66039141003603e-01, 3.91213893800891e-02, -2.99136269067853e-03}, 59997a1619fSSatish Balay {0.00000000000000e+00, -4.88498130835069e-15, -6.43929354282591e-15, -3.55271367880050e-15, -1.22124532708767e-15, -3.12250225675825e-16}, 60097a1619fSSatish Balay {0.00000000000000e+00, 1.22250171233141e+01, -1.77150760606169e+00, 3.54516769879390e-01, 6.22298845883398e-01, 2.31647447450276e-01 }, 60197a1619fSSatish Balay {0.00000000000000e+00, -4.48339457331040e+01, -3.57363126641880e-01, 5.18750173123425e-01, 6.55727990241799e-02, 1.63175368287079e-01 }, 60297a1619fSSatish Balay {0.00000000000000e+00, 1.37297394708005e+02, -1.60145272991317e+00, -5.05319555199441e+00, 1.55328940390990e-01, 9.16629423682464e-01 }, 60397a1619fSSatish Balay {0.00000000000000e+00, 1.05703241119022e+02, -1.16610260983038e+00, -2.99767252773859e+00, -1.13472315553890e-01, 1.09742849254729e+00 } 60497a1619fSSatish Balay }; 6059566063dSJacob Faibussowitsch PetscCall(TSGLLESchemeCreate(5, 5, 6, 6, c, *a, *b, *u, *v, &gl->schemes[gl->nschemes++])); 6063d177a5cSEmil Constantinescu } 6073ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6083d177a5cSEmil Constantinescu } 6093d177a5cSEmil Constantinescu 610cc4c1da9SBarry Smith /*@ 611bcf0153eSBarry Smith TSGLLESetType - sets the class of general linear method, `TSGLLE` to use for time-stepping 6123d177a5cSEmil Constantinescu 613c3339decSBarry Smith Collective 6143d177a5cSEmil Constantinescu 6153d177a5cSEmil Constantinescu Input Parameters: 616bcf0153eSBarry Smith + ts - the `TS` context 6173d177a5cSEmil Constantinescu - type - a method 6183d177a5cSEmil Constantinescu 6193d177a5cSEmil Constantinescu Options Database Key: 6203d177a5cSEmil Constantinescu . -ts_gl_type <type> - sets the method, use -help for a list of available method (e.g. irks) 6213d177a5cSEmil Constantinescu 622bcf0153eSBarry Smith Level: intermediate 623bcf0153eSBarry Smith 6243d177a5cSEmil Constantinescu Notes: 6253d177a5cSEmil Constantinescu See "petsc/include/petscts.h" for available methods (for instance) 6263d177a5cSEmil Constantinescu . TSGLLE_IRKS - Diagonally implicit methods with inherent Runge-Kutta stability (for stiff problems) 6273d177a5cSEmil Constantinescu 62814d0ab18SJacob Faibussowitsch Normally, it is best to use the `TSSetFromOptions()` command and then set the `TSGLLE` type 62914d0ab18SJacob Faibussowitsch from the options database rather than by using this routine. Using the options database 63014d0ab18SJacob Faibussowitsch provides the user with maximum flexibility in evaluating the many different solvers. The 63114d0ab18SJacob Faibussowitsch `TSGLLESetType()` routine is provided for those situations where it is necessary to set the 63214d0ab18SJacob Faibussowitsch timestepping solver independently of the command line or options database. This might be the 63314d0ab18SJacob Faibussowitsch case, for example, when the choice of solver changes during the execution of the program, and 63414d0ab18SJacob Faibussowitsch the user's application is taking responsibility for choosing the appropriate method. 6353d177a5cSEmil Constantinescu 6361cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSGLLEType`, `TSGLLE` 6373d177a5cSEmil Constantinescu @*/ 638d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGLLESetType(TS ts, TSGLLEType type) 639d71ae5a4SJacob Faibussowitsch { 6403d177a5cSEmil Constantinescu PetscFunctionBegin; 6413d177a5cSEmil Constantinescu PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 6424f572ea9SToby Isaac PetscAssertPointer(type, 2); 643cac4c232SBarry Smith PetscTryMethod(ts, "TSGLLESetType_C", (TS, TSGLLEType), (ts, type)); 6443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6453d177a5cSEmil Constantinescu } 6463d177a5cSEmil Constantinescu 6473d177a5cSEmil Constantinescu /*@C 648bcf0153eSBarry Smith TSGLLESetAcceptType - sets the acceptance test for `TSGLLE` 6493d177a5cSEmil Constantinescu 650c3339decSBarry Smith Logically Collective 6513d177a5cSEmil Constantinescu 6523d177a5cSEmil Constantinescu Input Parameters: 653bcf0153eSBarry Smith + ts - the `TS` context 6543d177a5cSEmil Constantinescu - type - the type 6553d177a5cSEmil Constantinescu 6563d177a5cSEmil Constantinescu Options Database Key: 6573d177a5cSEmil Constantinescu . -ts_gl_accept_type <type> - sets the method used to determine whether to accept or reject a step 6583d177a5cSEmil Constantinescu 6593d177a5cSEmil Constantinescu Level: intermediate 6603d177a5cSEmil Constantinescu 66114d0ab18SJacob Faibussowitsch Notes: 66214d0ab18SJacob Faibussowitsch Time integrators that need to control error must have the option to reject a time step based 66314d0ab18SJacob Faibussowitsch on local error estimates. This function allows different schemes to be set. 66414d0ab18SJacob Faibussowitsch 6651cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSGLLE`, `TSGLLEAcceptRegister()`, `TSGLLEAdapt` 6663d177a5cSEmil Constantinescu @*/ 667d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGLLESetAcceptType(TS ts, TSGLLEAcceptType type) 668d71ae5a4SJacob Faibussowitsch { 6693d177a5cSEmil Constantinescu PetscFunctionBegin; 6703d177a5cSEmil Constantinescu PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 6714f572ea9SToby Isaac PetscAssertPointer(type, 2); 672cac4c232SBarry Smith PetscTryMethod(ts, "TSGLLESetAcceptType_C", (TS, TSGLLEAcceptType), (ts, type)); 6733ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6743d177a5cSEmil Constantinescu } 6753d177a5cSEmil Constantinescu 676cc4c1da9SBarry Smith /*@ 677bcf0153eSBarry Smith TSGLLEGetAdapt - gets the `TSGLLEAdapt` object from the `TS` 6783d177a5cSEmil Constantinescu 6793d177a5cSEmil Constantinescu Not Collective 6803d177a5cSEmil Constantinescu 6813d177a5cSEmil Constantinescu Input Parameter: 682bcf0153eSBarry Smith . ts - the `TS` context 6833d177a5cSEmil Constantinescu 6843d177a5cSEmil Constantinescu Output Parameter: 685bcf0153eSBarry Smith . adapt - the `TSGLLEAdapt` context 6863d177a5cSEmil Constantinescu 6873d177a5cSEmil Constantinescu Level: advanced 6883d177a5cSEmil Constantinescu 689bcf0153eSBarry Smith Note: 69014d0ab18SJacob Faibussowitsch This allows the user set options on the `TSGLLEAdapt` object. Usually it is better to do this 69114d0ab18SJacob Faibussowitsch using the options database, so this function is rarely needed. 692bcf0153eSBarry Smith 6931cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSGLLE`, `TSGLLEAdapt`, `TSGLLEAdaptRegister()` 6943d177a5cSEmil Constantinescu @*/ 695d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGLLEGetAdapt(TS ts, TSGLLEAdapt *adapt) 696d71ae5a4SJacob Faibussowitsch { 6973d177a5cSEmil Constantinescu PetscFunctionBegin; 6983d177a5cSEmil Constantinescu PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 6994f572ea9SToby Isaac PetscAssertPointer(adapt, 2); 700cac4c232SBarry Smith PetscUseMethod(ts, "TSGLLEGetAdapt_C", (TS, TSGLLEAdapt *), (ts, adapt)); 7013ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7023d177a5cSEmil Constantinescu } 7033d177a5cSEmil Constantinescu 704d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSGLLEAccept_Always(TS ts, PetscReal tleft, PetscReal h, const PetscReal enorms[], PetscBool *accept) 705d71ae5a4SJacob Faibussowitsch { 7063d177a5cSEmil Constantinescu PetscFunctionBegin; 7073d177a5cSEmil Constantinescu *accept = PETSC_TRUE; 7083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7093d177a5cSEmil Constantinescu } 7103d177a5cSEmil Constantinescu 711d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSGLLEUpdateWRMS(TS ts) 712d71ae5a4SJacob Faibussowitsch { 7133d177a5cSEmil Constantinescu TS_GLLE *gl = (TS_GLLE *)ts->data; 7143d177a5cSEmil Constantinescu PetscScalar *x, *w; 7153d177a5cSEmil Constantinescu PetscInt n, i; 7163d177a5cSEmil Constantinescu 7173d177a5cSEmil Constantinescu PetscFunctionBegin; 7189566063dSJacob Faibussowitsch PetscCall(VecGetArray(gl->X[0], &x)); 7199566063dSJacob Faibussowitsch PetscCall(VecGetArray(gl->W, &w)); 7209566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(gl->W, &n)); 7213d177a5cSEmil Constantinescu for (i = 0; i < n; i++) w[i] = 1. / (gl->wrms_atol + gl->wrms_rtol * PetscAbsScalar(x[i])); 7229566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(gl->X[0], &x)); 7239566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(gl->W, &w)); 7243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7253d177a5cSEmil Constantinescu } 7263d177a5cSEmil Constantinescu 727d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSGLLEVecNormWRMS(TS ts, Vec X, PetscReal *nrm) 728d71ae5a4SJacob Faibussowitsch { 7293d177a5cSEmil Constantinescu TS_GLLE *gl = (TS_GLLE *)ts->data; 7303d177a5cSEmil Constantinescu PetscScalar *x, *w; 7313d177a5cSEmil Constantinescu PetscReal sum = 0.0, gsum; 7323d177a5cSEmil Constantinescu PetscInt n, N, i; 7333d177a5cSEmil Constantinescu 7343d177a5cSEmil Constantinescu PetscFunctionBegin; 7359566063dSJacob Faibussowitsch PetscCall(VecGetArray(X, &x)); 7369566063dSJacob Faibussowitsch PetscCall(VecGetArray(gl->W, &w)); 7379566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(gl->W, &n)); 7383d177a5cSEmil Constantinescu for (i = 0; i < n; i++) sum += PetscAbsScalar(PetscSqr(x[i] * w[i])); 7399566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(X, &x)); 7409566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(gl->W, &w)); 741462c564dSBarry Smith PetscCallMPI(MPIU_Allreduce(&sum, &gsum, 1, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)ts))); 7429566063dSJacob Faibussowitsch PetscCall(VecGetSize(gl->W, &N)); 7433d177a5cSEmil Constantinescu *nrm = PetscSqrtReal(gsum / (1. * N)); 7443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7453d177a5cSEmil Constantinescu } 7463d177a5cSEmil Constantinescu 747d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSGLLESetType_GLLE(TS ts, TSGLLEType type) 748d71ae5a4SJacob Faibussowitsch { 7493d177a5cSEmil Constantinescu PetscBool same; 7503d177a5cSEmil Constantinescu TS_GLLE *gl = (TS_GLLE *)ts->data; 7515f80ce2aSJacob Faibussowitsch PetscErrorCode (*r)(TS); 7523d177a5cSEmil Constantinescu 7533d177a5cSEmil Constantinescu PetscFunctionBegin; 7543d177a5cSEmil Constantinescu if (gl->type_name[0]) { 7559566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(gl->type_name, type, &same)); 7563ba16761SJacob Faibussowitsch if (same) PetscFunctionReturn(PETSC_SUCCESS); 7579566063dSJacob Faibussowitsch PetscCall((*gl->Destroy)(gl)); 7583d177a5cSEmil Constantinescu } 7593d177a5cSEmil Constantinescu 7609566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(TSGLLEList, type, &r)); 7616adde796SStefano Zampini PetscCheck(r, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown TSGLLE type \"%s\" given", type); 7629566063dSJacob Faibussowitsch PetscCall((*r)(ts)); 763c6a7a370SJeremy L Thompson PetscCall(PetscStrncpy(gl->type_name, type, sizeof(gl->type_name))); 7643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7653d177a5cSEmil Constantinescu } 7663d177a5cSEmil Constantinescu 767d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSGLLESetAcceptType_GLLE(TS ts, TSGLLEAcceptType type) 768d71ae5a4SJacob Faibussowitsch { 7698434afd1SBarry Smith TSGLLEAcceptFn *r; 7703d177a5cSEmil Constantinescu TS_GLLE *gl = (TS_GLLE *)ts->data; 7713d177a5cSEmil Constantinescu 7723d177a5cSEmil Constantinescu PetscFunctionBegin; 7739566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(TSGLLEAcceptList, type, &r)); 7746adde796SStefano Zampini PetscCheck(r, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown TSGLLEAccept type \"%s\" given", type); 7753d177a5cSEmil Constantinescu gl->Accept = r; 7769566063dSJacob Faibussowitsch PetscCall(PetscStrncpy(gl->accept_name, type, sizeof(gl->accept_name))); 7773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7783d177a5cSEmil Constantinescu } 7793d177a5cSEmil Constantinescu 780d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSGLLEGetAdapt_GLLE(TS ts, TSGLLEAdapt *adapt) 781d71ae5a4SJacob Faibussowitsch { 7823d177a5cSEmil Constantinescu TS_GLLE *gl = (TS_GLLE *)ts->data; 7833d177a5cSEmil Constantinescu 7843d177a5cSEmil Constantinescu PetscFunctionBegin; 7853d177a5cSEmil Constantinescu if (!gl->adapt) { 7869566063dSJacob Faibussowitsch PetscCall(TSGLLEAdaptCreate(PetscObjectComm((PetscObject)ts), &gl->adapt)); 7879566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)gl->adapt, (PetscObject)ts, 1)); 7883d177a5cSEmil Constantinescu } 7893d177a5cSEmil Constantinescu *adapt = gl->adapt; 7903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7913d177a5cSEmil Constantinescu } 7923d177a5cSEmil Constantinescu 793d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSGLLEChooseNextScheme(TS ts, PetscReal h, const PetscReal hmnorm[], PetscInt *next_scheme, PetscReal *next_h, PetscBool *finish) 794d71ae5a4SJacob Faibussowitsch { 7953d177a5cSEmil Constantinescu TS_GLLE *gl = (TS_GLLE *)ts->data; 7963d177a5cSEmil Constantinescu PetscInt i, n, cur_p, cur, next_sc, candidates[64], orders[64]; 7973d177a5cSEmil Constantinescu PetscReal errors[64], costs[64], tleft; 7983d177a5cSEmil Constantinescu 7993d177a5cSEmil Constantinescu PetscFunctionBegin; 8003d177a5cSEmil Constantinescu cur = -1; 8013d177a5cSEmil Constantinescu cur_p = gl->schemes[gl->current_scheme]->p; 8023d177a5cSEmil Constantinescu tleft = ts->max_time - (ts->ptime + ts->time_step); 8033d177a5cSEmil Constantinescu for (i = 0, n = 0; i < gl->nschemes; i++) { 8043d177a5cSEmil Constantinescu TSGLLEScheme sc = gl->schemes[i]; 8053d177a5cSEmil Constantinescu if (sc->p < gl->min_order || gl->max_order < sc->p) continue; 8063d177a5cSEmil Constantinescu if (sc->p == cur_p - 1) errors[n] = PetscAbsScalar(sc->alpha[0]) * hmnorm[0]; 8073d177a5cSEmil Constantinescu else if (sc->p == cur_p) errors[n] = PetscAbsScalar(sc->alpha[0]) * hmnorm[1]; 8083d177a5cSEmil Constantinescu else if (sc->p == cur_p + 1) errors[n] = PetscAbsScalar(sc->alpha[0]) * (hmnorm[2] + hmnorm[3]); 8093d177a5cSEmil Constantinescu else continue; 8103d177a5cSEmil Constantinescu candidates[n] = i; 8113d177a5cSEmil Constantinescu orders[n] = PetscMin(sc->p, sc->q); /* order of global truncation error */ 8123d177a5cSEmil Constantinescu costs[n] = sc->s; /* estimate the cost as the number of stages */ 8133d177a5cSEmil Constantinescu if (i == gl->current_scheme) cur = n; 8143d177a5cSEmil Constantinescu n++; 8153d177a5cSEmil Constantinescu } 816cad9d221SBarry Smith PetscCheck(cur >= 0 && gl->nschemes > cur, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Current scheme not found in scheme list"); 8179566063dSJacob Faibussowitsch PetscCall(TSGLLEAdaptChoose(gl->adapt, n, orders, errors, costs, cur, h, tleft, &next_sc, next_h, finish)); 8183d177a5cSEmil Constantinescu *next_scheme = candidates[next_sc]; 8199371c9d4SSatish Balay 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, 8209371c9d4SSatish Balay gl->schemes[*next_scheme]->r, gl->schemes[*next_scheme]->s, (double)*next_h, PetscBools[*finish])); 8213ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8223d177a5cSEmil Constantinescu } 8233d177a5cSEmil Constantinescu 824d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSGLLEGetMaxSizes(TS ts, PetscInt *max_r, PetscInt *max_s) 825d71ae5a4SJacob Faibussowitsch { 8263d177a5cSEmil Constantinescu TS_GLLE *gl = (TS_GLLE *)ts->data; 8273d177a5cSEmil Constantinescu 8283d177a5cSEmil Constantinescu PetscFunctionBegin; 8293d177a5cSEmil Constantinescu *max_r = gl->schemes[gl->nschemes - 1]->r; 8303d177a5cSEmil Constantinescu *max_s = gl->schemes[gl->nschemes - 1]->s; 8313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8323d177a5cSEmil Constantinescu } 8333d177a5cSEmil Constantinescu 834d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSolve_GLLE(TS ts) 835d71ae5a4SJacob Faibussowitsch { 8363d177a5cSEmil Constantinescu TS_GLLE *gl = (TS_GLLE *)ts->data; 8373d177a5cSEmil Constantinescu PetscInt i, k, its, lits, max_r, max_s; 8383d177a5cSEmil Constantinescu PetscBool final_step, finish; 8393d177a5cSEmil Constantinescu SNESConvergedReason snesreason; 8403d177a5cSEmil Constantinescu 8413d177a5cSEmil Constantinescu PetscFunctionBegin; 8429566063dSJacob Faibussowitsch PetscCall(TSMonitor(ts, ts->steps, ts->ptime, ts->vec_sol)); 8433d177a5cSEmil Constantinescu 8449566063dSJacob Faibussowitsch PetscCall(TSGLLEGetMaxSizes(ts, &max_r, &max_s)); 8459566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol, gl->X[0])); 84648a46eb9SPierre Jolivet for (i = 1; i < max_r; i++) PetscCall(VecZeroEntries(gl->X[i])); 8479566063dSJacob Faibussowitsch PetscCall(TSGLLEUpdateWRMS(ts)); 8483d177a5cSEmil Constantinescu 8493d177a5cSEmil Constantinescu if (0) { 8503d177a5cSEmil Constantinescu /* Find consistent initial data for DAE */ 8513d177a5cSEmil Constantinescu gl->stage_time = ts->ptime + ts->time_step; 8523d177a5cSEmil Constantinescu gl->scoeff = 1.; 8533d177a5cSEmil Constantinescu gl->stage = 0; 8543d177a5cSEmil Constantinescu 8559566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol, gl->Z)); 8569566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol, gl->Y)); 8579566063dSJacob Faibussowitsch PetscCall(SNESSolve(ts->snes, NULL, gl->Y)); 8589566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(ts->snes, &its)); 8599566063dSJacob Faibussowitsch PetscCall(SNESGetLinearSolveIterations(ts->snes, &lits)); 8609566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(ts->snes, &snesreason)); 8613d177a5cSEmil Constantinescu 8629371c9d4SSatish Balay ts->snes_its += its; 8639371c9d4SSatish Balay ts->ksp_its += lits; 86409cb0f53SBarry Smith if (snesreason < 0 && ts->max_snes_failures != PETSC_UNLIMITED && ++ts->num_snes_failures >= ts->max_snes_failures) { 8653d177a5cSEmil Constantinescu ts->reason = TS_DIVERGED_NONLINEAR_SOLVE; 86615229ffcSPierre Jolivet PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", nonlinear solve failures %" PetscInt_FMT " greater than current TS allowed, stopping solve\n", ts->steps, ts->num_snes_failures)); 8673ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8683d177a5cSEmil Constantinescu } 8693d177a5cSEmil Constantinescu } 8703d177a5cSEmil Constantinescu 87108401ef6SPierre Jolivet PetscCheck(gl->current_scheme >= 0, PETSC_COMM_SELF, PETSC_ERR_ORDER, "A starting scheme has not been provided"); 8723d177a5cSEmil Constantinescu 8733d177a5cSEmil Constantinescu for (k = 0, final_step = PETSC_FALSE, finish = PETSC_FALSE; k < ts->max_steps && !finish; k++) { 8743d177a5cSEmil Constantinescu PetscInt j, r, s, next_scheme = 0, rejections; 8753d177a5cSEmil Constantinescu PetscReal h, hmnorm[4], enorm[3], next_h; 8763d177a5cSEmil Constantinescu PetscBool accept; 8773d177a5cSEmil Constantinescu const PetscScalar *c, *a, *u; 8783d177a5cSEmil Constantinescu Vec *X, *Ydot, Y; 8793d177a5cSEmil Constantinescu TSGLLEScheme scheme = gl->schemes[gl->current_scheme]; 8803d177a5cSEmil Constantinescu 8819371c9d4SSatish Balay r = scheme->r; 8829371c9d4SSatish Balay s = scheme->s; 8833d177a5cSEmil Constantinescu c = scheme->c; 8849371c9d4SSatish Balay a = scheme->a; 8859371c9d4SSatish Balay u = scheme->u; 8863d177a5cSEmil Constantinescu h = ts->time_step; 8879371c9d4SSatish Balay X = gl->X; 8889371c9d4SSatish Balay Ydot = gl->Ydot; 8899371c9d4SSatish Balay Y = gl->Y; 8903d177a5cSEmil Constantinescu 8913d177a5cSEmil Constantinescu if (ts->ptime > ts->max_time) break; 8923d177a5cSEmil Constantinescu 8933d177a5cSEmil Constantinescu /* 8943d177a5cSEmil Constantinescu We only call PreStep at the start of each STEP, not each STAGE. This is because it is 8953d177a5cSEmil Constantinescu possible to fail (have to restart a step) after multiple stages. 8963d177a5cSEmil Constantinescu */ 8979566063dSJacob Faibussowitsch PetscCall(TSPreStep(ts)); 8983d177a5cSEmil Constantinescu 8993d177a5cSEmil Constantinescu rejections = 0; 9003d177a5cSEmil Constantinescu while (1) { 9013d177a5cSEmil Constantinescu for (i = 0; i < s; i++) { 9023d177a5cSEmil Constantinescu PetscScalar shift; 9033d177a5cSEmil Constantinescu gl->scoeff = 1. / PetscRealPart(a[i * s + i]); 9043d177a5cSEmil Constantinescu shift = gl->scoeff / ts->time_step; 9053d177a5cSEmil Constantinescu gl->stage = i; 9063d177a5cSEmil Constantinescu gl->stage_time = ts->ptime + PetscRealPart(c[i]) * h; 9073d177a5cSEmil Constantinescu 9083d177a5cSEmil Constantinescu /* 9093d177a5cSEmil Constantinescu * Stage equation: Y = h A Y' + U X 9103d177a5cSEmil Constantinescu * We assume that A is lower-triangular so that we can solve the stages (Y,Y') sequentially 9113d177a5cSEmil Constantinescu * Build the affine vector z_i = -[1/(h a_ii)](h sum_j a_ij y'_j + sum_j u_ij x_j) 9123d177a5cSEmil Constantinescu * Then y'_i = z + 1/(h a_ii) y_i 9133d177a5cSEmil Constantinescu */ 9149566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(gl->Z)); 91548a46eb9SPierre Jolivet for (j = 0; j < r; j++) PetscCall(VecAXPY(gl->Z, -shift * u[i * r + j], X[j])); 91648a46eb9SPierre Jolivet for (j = 0; j < i; j++) PetscCall(VecAXPY(gl->Z, -shift * h * a[i * s + j], Ydot[j])); 9173d177a5cSEmil Constantinescu /* Note: Z is used within function evaluation, Ydot = Z + shift*Y */ 9183d177a5cSEmil Constantinescu 9193d177a5cSEmil Constantinescu /* Compute an estimate of Y to start Newton iteration */ 9203d177a5cSEmil Constantinescu if (gl->extrapolate) { 9213d177a5cSEmil Constantinescu if (i == 0) { 9223d177a5cSEmil Constantinescu /* Linear extrapolation on the first stage */ 9239566063dSJacob Faibussowitsch PetscCall(VecWAXPY(Y, c[i] * h, X[1], X[0])); 9243d177a5cSEmil Constantinescu } else { 9253d177a5cSEmil Constantinescu /* Linear extrapolation from the last stage */ 9269566063dSJacob Faibussowitsch PetscCall(VecAXPY(Y, (c[i] - c[i - 1]) * h, Ydot[i - 1])); 9273d177a5cSEmil Constantinescu } 9283d177a5cSEmil Constantinescu } else if (i == 0) { /* Directly use solution from the last step, otherwise reuse the last stage (do nothing) */ 9299566063dSJacob Faibussowitsch PetscCall(VecCopy(X[0], Y)); 9303d177a5cSEmil Constantinescu } 9313d177a5cSEmil Constantinescu 9323d177a5cSEmil Constantinescu /* Solve this stage (Ydot[i] is computed during function evaluation) */ 9339566063dSJacob Faibussowitsch PetscCall(SNESSolve(ts->snes, NULL, Y)); 9349566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(ts->snes, &its)); 9359566063dSJacob Faibussowitsch PetscCall(SNESGetLinearSolveIterations(ts->snes, &lits)); 9369566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(ts->snes, &snesreason)); 9379371c9d4SSatish Balay ts->snes_its += its; 9389371c9d4SSatish Balay ts->ksp_its += lits; 93909cb0f53SBarry Smith if (snesreason < 0 && ts->max_snes_failures != PETSC_UNLIMITED && ++ts->num_snes_failures >= ts->max_snes_failures) { 9403d177a5cSEmil Constantinescu ts->reason = TS_DIVERGED_NONLINEAR_SOLVE; 94115229ffcSPierre Jolivet PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", nonlinear solve failures %" PetscInt_FMT " greater than current TS allowed, stopping solve\n", ts->steps, ts->num_snes_failures)); 9423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9433d177a5cSEmil Constantinescu } 9443d177a5cSEmil Constantinescu } 9453d177a5cSEmil Constantinescu 9463d177a5cSEmil Constantinescu gl->stage_time = ts->ptime + ts->time_step; 9473d177a5cSEmil Constantinescu 9489566063dSJacob Faibussowitsch PetscCall((*gl->EstimateHigherMoments)(scheme, h, Ydot, gl->X, gl->himom)); 9493d177a5cSEmil Constantinescu /* hmnorm[i] = h^{p+i}x^{(p+i)} with i=0,1,2; hmnorm[3] = h^{p+2}(dx'/dx) x^{(p+1)} */ 95048a46eb9SPierre Jolivet for (i = 0; i < 3; i++) PetscCall(TSGLLEVecNormWRMS(ts, gl->himom[i], &hmnorm[i + 1])); 9513d177a5cSEmil Constantinescu enorm[0] = PetscRealPart(scheme->alpha[0]) * hmnorm[1]; 9523d177a5cSEmil Constantinescu enorm[1] = PetscRealPart(scheme->beta[0]) * hmnorm[2]; 9533d177a5cSEmil Constantinescu enorm[2] = PetscRealPart(scheme->gamma[0]) * hmnorm[3]; 9549566063dSJacob Faibussowitsch PetscCall((*gl->Accept)(ts, ts->max_time - gl->stage_time, h, enorm, &accept)); 9553d177a5cSEmil Constantinescu if (accept) goto accepted; 9563d177a5cSEmil Constantinescu rejections++; 95763a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(ts, "Step %" PetscInt_FMT " (t=%g) not accepted, rejections=%" PetscInt_FMT "\n", k, (double)gl->stage_time, rejections)); 9583d177a5cSEmil Constantinescu if (rejections > gl->max_step_rejections) break; 9593d177a5cSEmil Constantinescu /* 9603d177a5cSEmil Constantinescu There are lots of reasons why a step might be rejected, including solvers not converging and other factors that 9613d177a5cSEmil Constantinescu TSGLLEChooseNextScheme does not support. Additionally, the error estimates may be very screwed up, so I'm not 9623d177a5cSEmil Constantinescu convinced that it's safe to just compute a new error estimate using the same interface as the current adaptor 9633d177a5cSEmil Constantinescu (the adaptor interface probably has to change). Here we make an arbitrary and naive choice. This assumes that 9643d177a5cSEmil Constantinescu steps were written in Nordsieck form. The "correct" method would be to re-complete the previous time step with 9653d177a5cSEmil Constantinescu the correct "next" step size. It is unclear to me whether the present ad-hoc method of rescaling X is stable. 9663d177a5cSEmil Constantinescu */ 9673d177a5cSEmil Constantinescu h *= 0.5; 96848a46eb9SPierre Jolivet for (i = 1; i < scheme->r; i++) PetscCall(VecScale(X[i], PetscPowRealInt(0.5, i))); 9693d177a5cSEmil Constantinescu } 97063a3b9bcSJacob 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); 9713d177a5cSEmil Constantinescu 9723d177a5cSEmil Constantinescu accepted: 9733d177a5cSEmil Constantinescu /* This term is not error, but it *would* be the leading term for a lower order method */ 9749566063dSJacob Faibussowitsch PetscCall(TSGLLEVecNormWRMS(ts, gl->X[scheme->r - 1], &hmnorm[0])); 9753d177a5cSEmil Constantinescu /* Correct scaling so that these are equivalent to norms of the Nordsieck vectors */ 9763d177a5cSEmil Constantinescu 97763a3b9bcSJacob 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])); 9783d177a5cSEmil Constantinescu if (!final_step) { 9799566063dSJacob Faibussowitsch PetscCall(TSGLLEChooseNextScheme(ts, h, hmnorm, &next_scheme, &next_h, &final_step)); 9803d177a5cSEmil Constantinescu } else { 9813d177a5cSEmil Constantinescu /* Dummy values to complete the current step in a consistent manner */ 9823d177a5cSEmil Constantinescu next_scheme = gl->current_scheme; 9833d177a5cSEmil Constantinescu next_h = h; 9843d177a5cSEmil Constantinescu finish = PETSC_TRUE; 9853d177a5cSEmil Constantinescu } 9863d177a5cSEmil Constantinescu 9873d177a5cSEmil Constantinescu X = gl->Xold; 9883d177a5cSEmil Constantinescu gl->Xold = gl->X; 9893d177a5cSEmil Constantinescu gl->X = X; 9909566063dSJacob Faibussowitsch PetscCall((*gl->CompleteStep)(scheme, h, gl->schemes[next_scheme], next_h, Ydot, gl->Xold, gl->X)); 9913d177a5cSEmil Constantinescu 9929566063dSJacob Faibussowitsch PetscCall(TSGLLEUpdateWRMS(ts)); 9933d177a5cSEmil Constantinescu 9943d177a5cSEmil Constantinescu /* Post the solution for the user, we could avoid this copy with a small bit of cleverness */ 9959566063dSJacob Faibussowitsch PetscCall(VecCopy(gl->X[0], ts->vec_sol)); 9963d177a5cSEmil Constantinescu ts->ptime += h; 9973d177a5cSEmil Constantinescu ts->steps++; 9983d177a5cSEmil Constantinescu 9999566063dSJacob Faibussowitsch PetscCall(TSPostEvaluate(ts)); 10009566063dSJacob Faibussowitsch PetscCall(TSPostStep(ts)); 10019566063dSJacob Faibussowitsch PetscCall(TSMonitor(ts, ts->steps, ts->ptime, ts->vec_sol)); 10023d177a5cSEmil Constantinescu 10033d177a5cSEmil Constantinescu gl->current_scheme = next_scheme; 10043d177a5cSEmil Constantinescu ts->time_step = next_h; 10053d177a5cSEmil Constantinescu } 10063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 10073d177a5cSEmil Constantinescu } 10083d177a5cSEmil Constantinescu 10093d177a5cSEmil Constantinescu /*------------------------------------------------------------*/ 10103d177a5cSEmil Constantinescu 1011d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSReset_GLLE(TS ts) 1012d71ae5a4SJacob Faibussowitsch { 10133d177a5cSEmil Constantinescu TS_GLLE *gl = (TS_GLLE *)ts->data; 10143d177a5cSEmil Constantinescu PetscInt max_r, max_s; 10153d177a5cSEmil Constantinescu 10163d177a5cSEmil Constantinescu PetscFunctionBegin; 10173d177a5cSEmil Constantinescu if (gl->setupcalled) { 10189566063dSJacob Faibussowitsch PetscCall(TSGLLEGetMaxSizes(ts, &max_r, &max_s)); 10199566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(max_r, &gl->Xold)); 10209566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(max_r, &gl->X)); 10219566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(max_s, &gl->Ydot)); 10229566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(3, &gl->himom)); 10239566063dSJacob Faibussowitsch PetscCall(VecDestroy(&gl->W)); 10249566063dSJacob Faibussowitsch PetscCall(VecDestroy(&gl->Y)); 10259566063dSJacob Faibussowitsch PetscCall(VecDestroy(&gl->Z)); 10263d177a5cSEmil Constantinescu } 10273d177a5cSEmil Constantinescu gl->setupcalled = PETSC_FALSE; 10283ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 10293d177a5cSEmil Constantinescu } 10303d177a5cSEmil Constantinescu 1031d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSDestroy_GLLE(TS ts) 1032d71ae5a4SJacob Faibussowitsch { 10333d177a5cSEmil Constantinescu TS_GLLE *gl = (TS_GLLE *)ts->data; 10343d177a5cSEmil Constantinescu 10353d177a5cSEmil Constantinescu PetscFunctionBegin; 10369566063dSJacob Faibussowitsch PetscCall(TSReset_GLLE(ts)); 1037b3a6b972SJed Brown if (ts->dm) { 10389566063dSJacob Faibussowitsch PetscCall(DMCoarsenHookRemove(ts->dm, DMCoarsenHook_TSGLLE, DMRestrictHook_TSGLLE, ts)); 10399566063dSJacob Faibussowitsch PetscCall(DMSubDomainHookRemove(ts->dm, DMSubDomainHook_TSGLLE, DMSubDomainRestrictHook_TSGLLE, ts)); 1040b3a6b972SJed Brown } 10419566063dSJacob Faibussowitsch if (gl->adapt) PetscCall(TSGLLEAdaptDestroy(&gl->adapt)); 10429566063dSJacob Faibussowitsch if (gl->Destroy) PetscCall((*gl->Destroy)(gl)); 10439566063dSJacob Faibussowitsch PetscCall(PetscFree(ts->data)); 10449566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSGLLESetType_C", NULL)); 10459566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSGLLESetAcceptType_C", NULL)); 10469566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSGLLEGetAdapt_C", NULL)); 10473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 10483d177a5cSEmil Constantinescu } 10493d177a5cSEmil Constantinescu 10503d177a5cSEmil Constantinescu /* 10513d177a5cSEmil Constantinescu This defines the nonlinear equation that is to be solved with SNES 10523d177a5cSEmil Constantinescu g(x) = f(t,x,z+shift*x) = 0 10533d177a5cSEmil Constantinescu */ 1054d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormFunction_GLLE(SNES snes, Vec x, Vec f, TS ts) 1055d71ae5a4SJacob Faibussowitsch { 10563d177a5cSEmil Constantinescu TS_GLLE *gl = (TS_GLLE *)ts->data; 10573d177a5cSEmil Constantinescu Vec Z, Ydot; 10583d177a5cSEmil Constantinescu DM dm, dmsave; 10593d177a5cSEmil Constantinescu 10603d177a5cSEmil Constantinescu PetscFunctionBegin; 10619566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 10629566063dSJacob Faibussowitsch PetscCall(TSGLLEGetVecs(ts, dm, &Z, &Ydot)); 10639566063dSJacob Faibussowitsch PetscCall(VecWAXPY(Ydot, gl->scoeff / ts->time_step, x, Z)); 10643d177a5cSEmil Constantinescu dmsave = ts->dm; 10653d177a5cSEmil Constantinescu ts->dm = dm; 10669566063dSJacob Faibussowitsch PetscCall(TSComputeIFunction(ts, gl->stage_time, x, Ydot, f, PETSC_FALSE)); 10673d177a5cSEmil Constantinescu ts->dm = dmsave; 10689566063dSJacob Faibussowitsch PetscCall(TSGLLERestoreVecs(ts, dm, &Z, &Ydot)); 10693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 10703d177a5cSEmil Constantinescu } 10713d177a5cSEmil Constantinescu 1072d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormJacobian_GLLE(SNES snes, Vec x, Mat A, Mat B, TS ts) 1073d71ae5a4SJacob Faibussowitsch { 10743d177a5cSEmil Constantinescu TS_GLLE *gl = (TS_GLLE *)ts->data; 10753d177a5cSEmil Constantinescu Vec Z, Ydot; 10763d177a5cSEmil Constantinescu DM dm, dmsave; 10773d177a5cSEmil Constantinescu 10783d177a5cSEmil Constantinescu PetscFunctionBegin; 10799566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 10809566063dSJacob Faibussowitsch PetscCall(TSGLLEGetVecs(ts, dm, &Z, &Ydot)); 10813d177a5cSEmil Constantinescu dmsave = ts->dm; 10823d177a5cSEmil Constantinescu ts->dm = dm; 10833d177a5cSEmil Constantinescu /* gl->Xdot will have already been computed in SNESTSFormFunction_GLLE */ 10849566063dSJacob Faibussowitsch PetscCall(TSComputeIJacobian(ts, gl->stage_time, x, gl->Ydot[gl->stage], gl->scoeff / ts->time_step, A, B, PETSC_FALSE)); 10853d177a5cSEmil Constantinescu ts->dm = dmsave; 10869566063dSJacob Faibussowitsch PetscCall(TSGLLERestoreVecs(ts, dm, &Z, &Ydot)); 10873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 10883d177a5cSEmil Constantinescu } 10893d177a5cSEmil Constantinescu 1090d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetUp_GLLE(TS ts) 1091d71ae5a4SJacob Faibussowitsch { 10923d177a5cSEmil Constantinescu TS_GLLE *gl = (TS_GLLE *)ts->data; 10933d177a5cSEmil Constantinescu PetscInt max_r, max_s; 10943d177a5cSEmil Constantinescu DM dm; 10953d177a5cSEmil Constantinescu 10963d177a5cSEmil Constantinescu PetscFunctionBegin; 10973d177a5cSEmil Constantinescu gl->setupcalled = PETSC_TRUE; 10989566063dSJacob Faibussowitsch PetscCall(TSGLLEGetMaxSizes(ts, &max_r, &max_s)); 10999566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vec_sol, max_r, &gl->X)); 11009566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vec_sol, max_r, &gl->Xold)); 11019566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vec_sol, max_s, &gl->Ydot)); 11029566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vec_sol, 3, &gl->himom)); 11039566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &gl->W)); 11049566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &gl->Y)); 11059566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &gl->Z)); 11063d177a5cSEmil Constantinescu 11073d177a5cSEmil Constantinescu /* Default acceptance tests and adaptivity */ 11089566063dSJacob Faibussowitsch if (!gl->Accept) PetscCall(TSGLLESetAcceptType(ts, TSGLLEACCEPT_ALWAYS)); 11099566063dSJacob Faibussowitsch if (!gl->adapt) PetscCall(TSGLLEGetAdapt(ts, &gl->adapt)); 11103d177a5cSEmil Constantinescu 11113d177a5cSEmil Constantinescu if (gl->current_scheme < 0) { 11123d177a5cSEmil Constantinescu PetscInt i; 11133d177a5cSEmil Constantinescu for (i = 0;; i++) { 11143d177a5cSEmil Constantinescu if (gl->schemes[i]->p == gl->start_order) break; 111563a3b9bcSJacob Faibussowitsch PetscCheck(i + 1 != gl->nschemes, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No schemes available with requested start order %" PetscInt_FMT, i); 11163d177a5cSEmil Constantinescu } 11173d177a5cSEmil Constantinescu gl->current_scheme = i; 11183d177a5cSEmil Constantinescu } 11199566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 11209566063dSJacob Faibussowitsch PetscCall(DMCoarsenHookAdd(dm, DMCoarsenHook_TSGLLE, DMRestrictHook_TSGLLE, ts)); 11219566063dSJacob Faibussowitsch PetscCall(DMSubDomainHookAdd(dm, DMSubDomainHook_TSGLLE, DMSubDomainRestrictHook_TSGLLE, ts)); 11223ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11233d177a5cSEmil Constantinescu } 11243d177a5cSEmil Constantinescu /*------------------------------------------------------------*/ 11253d177a5cSEmil Constantinescu 1126ce78bad3SBarry Smith static PetscErrorCode TSSetFromOptions_GLLE(TS ts, PetscOptionItems PetscOptionsObject) 1127d71ae5a4SJacob Faibussowitsch { 11283d177a5cSEmil Constantinescu TS_GLLE *gl = (TS_GLLE *)ts->data; 11293d177a5cSEmil Constantinescu char tname[256] = TSGLLE_IRKS, completef[256] = "rescale-and-modify"; 11303d177a5cSEmil Constantinescu 11313d177a5cSEmil Constantinescu PetscFunctionBegin; 1132d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "General Linear ODE solver options"); 11333d177a5cSEmil Constantinescu { 11343d177a5cSEmil Constantinescu PetscBool flg; 11359566063dSJacob Faibussowitsch PetscCall(PetscOptionsFList("-ts_gl_type", "Type of GL method", "TSGLLESetType", TSGLLEList, gl->type_name[0] ? gl->type_name : tname, tname, sizeof(tname), &flg)); 113648a46eb9SPierre Jolivet if (flg || !gl->type_name[0]) PetscCall(TSGLLESetType(ts, tname)); 11379566063dSJacob 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)); 11389566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-ts_gl_max_order", "Maximum order to try", "TSGLLESetMaxOrder", gl->max_order, &gl->max_order, NULL)); 11399566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-ts_gl_min_order", "Minimum order to try", "TSGLLESetMinOrder", gl->min_order, &gl->min_order, NULL)); 11409566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-ts_gl_start_order", "Initial order to try", "TSGLLESetMinOrder", gl->start_order, &gl->start_order, NULL)); 11419566063dSJacob 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)); 11429566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-ts_gl_extrapolate", "Extrapolate stage solution from previous solution (sometimes unstable)", "TSGLLESetExtrapolate", gl->extrapolate, &gl->extrapolate, NULL)); 11439566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_gl_atol", "Absolute tolerance", "TSGLLESetTolerances", gl->wrms_atol, &gl->wrms_atol, NULL)); 11449566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_gl_rtol", "Relative tolerance", "TSGLLESetTolerances", gl->wrms_rtol, &gl->wrms_rtol, NULL)); 11459566063dSJacob Faibussowitsch PetscCall(PetscOptionsString("-ts_gl_complete", "Method to use for completing the step", "none", completef, completef, sizeof(completef), &flg)); 11463d177a5cSEmil Constantinescu if (flg) { 11473d177a5cSEmil Constantinescu PetscBool match1, match2; 11489566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(completef, "rescale", &match1)); 11499566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(completef, "rescale-and-modify", &match2)); 11503d177a5cSEmil Constantinescu if (match1) gl->CompleteStep = TSGLLECompleteStep_Rescale; 11513d177a5cSEmil Constantinescu else if (match2) gl->CompleteStep = TSGLLECompleteStep_RescaleAndModify; 11526adde796SStefano Zampini else SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_UNKNOWN_TYPE, "%s", completef); 11533d177a5cSEmil Constantinescu } 11543d177a5cSEmil Constantinescu { 11553d177a5cSEmil Constantinescu char type[256] = TSGLLEACCEPT_ALWAYS; 11569566063dSJacob 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)); 115748a46eb9SPierre Jolivet if (flg || !gl->accept_name[0]) PetscCall(TSGLLESetAcceptType(ts, type)); 11583d177a5cSEmil Constantinescu } 11593d177a5cSEmil Constantinescu { 11603d177a5cSEmil Constantinescu TSGLLEAdapt adapt; 11619566063dSJacob Faibussowitsch PetscCall(TSGLLEGetAdapt(ts, &adapt)); 1162dbbe0bcdSBarry Smith PetscCall(TSGLLEAdaptSetFromOptions(adapt, PetscOptionsObject)); 11633d177a5cSEmil Constantinescu } 11643d177a5cSEmil Constantinescu } 1165d0609cedSBarry Smith PetscOptionsHeadEnd(); 11663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11673d177a5cSEmil Constantinescu } 11683d177a5cSEmil Constantinescu 1169d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSView_GLLE(TS ts, PetscViewer viewer) 1170d71ae5a4SJacob Faibussowitsch { 11713d177a5cSEmil Constantinescu TS_GLLE *gl = (TS_GLLE *)ts->data; 11723d177a5cSEmil Constantinescu PetscInt i; 1173*9f196a02SMartin Diehl PetscBool isascii, details; 11743d177a5cSEmil Constantinescu 11753d177a5cSEmil Constantinescu PetscFunctionBegin; 1176*9f196a02SMartin Diehl PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii)); 1177*9f196a02SMartin Diehl if (isascii) { 117863a3b9bcSJacob 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)); 11799566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Error estimation: %s\n", TSGLLEErrorDirections[gl->error_direction])); 11809566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Extrapolation: %s\n", gl->extrapolate ? "yes" : "no")); 11819566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Acceptance test: %s\n", gl->accept_name[0] ? gl->accept_name : "(not yet set)")); 11829566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 11839566063dSJacob Faibussowitsch PetscCall(TSGLLEAdaptView(gl->adapt, viewer)); 11849566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 11859566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " type: %s\n", gl->type_name[0] ? gl->type_name : "(not yet set)")); 118663a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Schemes within family (%" PetscInt_FMT "):\n", gl->nschemes)); 11873d177a5cSEmil Constantinescu details = PETSC_FALSE; 11889566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(((PetscObject)ts)->options, ((PetscObject)ts)->prefix, "-ts_gl_view_detailed", &details, NULL)); 11899566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 119048a46eb9SPierre Jolivet for (i = 0; i < gl->nschemes; i++) PetscCall(TSGLLESchemeView(gl->schemes[i], details, viewer)); 11911baa6e33SBarry Smith if (gl->View) PetscCall((*gl->View)(gl, viewer)); 11929566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 11933d177a5cSEmil Constantinescu } 11943ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11953d177a5cSEmil Constantinescu } 11963d177a5cSEmil Constantinescu 11973d177a5cSEmil Constantinescu /*@C 1198bcf0153eSBarry Smith TSGLLERegister - adds a `TSGLLE` implementation 11993d177a5cSEmil Constantinescu 1200cc4c1da9SBarry Smith Not Collective, No Fortran Support 12013d177a5cSEmil Constantinescu 12023d177a5cSEmil Constantinescu Input Parameters: 120320f4b53cSBarry Smith + sname - name of user-defined general linear scheme 120420f4b53cSBarry Smith - function - routine to create method context 12053d177a5cSEmil Constantinescu 1206bcf0153eSBarry Smith Level: advanced 1207bcf0153eSBarry Smith 1208bcf0153eSBarry Smith Note: 1209bcf0153eSBarry Smith `TSGLLERegister()` may be called multiple times to add several user-defined families. 12103d177a5cSEmil Constantinescu 1211b43aa488SJacob Faibussowitsch Example Usage: 12123d177a5cSEmil Constantinescu .vb 12133d177a5cSEmil Constantinescu TSGLLERegister("my_scheme", MySchemeCreate); 12143d177a5cSEmil Constantinescu .ve 12153d177a5cSEmil Constantinescu 12163d177a5cSEmil Constantinescu Then, your scheme can be chosen with the procedural interface via 1217b44f4de4SBarry Smith .vb 1218b44f4de4SBarry Smith TSGLLESetType(ts, "my_scheme") 1219b44f4de4SBarry Smith .ve 12203d177a5cSEmil Constantinescu or at runtime via the option 1221b44f4de4SBarry Smith .vb 1222b44f4de4SBarry Smith -ts_gl_type my_scheme 1223b44f4de4SBarry Smith .ve 12243d177a5cSEmil Constantinescu 12251cc06b55SBarry Smith .seealso: [](ch_ts), `TSGLLE`, `TSGLLEType`, `TSGLLERegisterAll()` 12263d177a5cSEmil Constantinescu @*/ 1227d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGLLERegister(const char sname[], PetscErrorCode (*function)(TS)) 1228d71ae5a4SJacob Faibussowitsch { 12293d177a5cSEmil Constantinescu PetscFunctionBegin; 12309566063dSJacob Faibussowitsch PetscCall(TSGLLEInitializePackage()); 12319566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&TSGLLEList, sname, function)); 12323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 12333d177a5cSEmil Constantinescu } 12343d177a5cSEmil Constantinescu 12353d177a5cSEmil Constantinescu /*@C 1236bcf0153eSBarry Smith TSGLLEAcceptRegister - adds a `TSGLLE` acceptance scheme 12373d177a5cSEmil Constantinescu 12383d177a5cSEmil Constantinescu Not Collective 12393d177a5cSEmil Constantinescu 12403d177a5cSEmil Constantinescu Input Parameters: 124120f4b53cSBarry Smith + sname - name of user-defined acceptance scheme 12428434afd1SBarry Smith - function - routine to create method context, see `TSGLLEAcceptFn` for the calling sequence 12433d177a5cSEmil Constantinescu 1244bcf0153eSBarry Smith Level: advanced 1245bcf0153eSBarry Smith 1246bcf0153eSBarry Smith Note: 1247bcf0153eSBarry Smith `TSGLLEAcceptRegister()` may be called multiple times to add several user-defined families. 12483d177a5cSEmil Constantinescu 1249b43aa488SJacob Faibussowitsch Example Usage: 12503d177a5cSEmil Constantinescu .vb 12513d177a5cSEmil Constantinescu TSGLLEAcceptRegister("my_scheme", MySchemeCreate); 12523d177a5cSEmil Constantinescu .ve 12533d177a5cSEmil Constantinescu 12543d177a5cSEmil Constantinescu Then, your scheme can be chosen with the procedural interface via 1255d1c5d1fcSBarry Smith .vb 1256d1c5d1fcSBarry Smith TSGLLESetAcceptType(ts, "my_scheme") 1257d1c5d1fcSBarry Smith .ve 1258d1c5d1fcSBarry Smith or at runtime via the option `-ts_gl_accept_type my_scheme` 12593d177a5cSEmil Constantinescu 12608434afd1SBarry Smith .seealso: [](ch_ts), `TSGLLE`, `TSGLLEType`, `TSGLLERegisterAll()`, `TSGLLEAcceptFn` 12613d177a5cSEmil Constantinescu @*/ 12628434afd1SBarry Smith PetscErrorCode TSGLLEAcceptRegister(const char sname[], TSGLLEAcceptFn *function) 1263d71ae5a4SJacob Faibussowitsch { 12643d177a5cSEmil Constantinescu PetscFunctionBegin; 12659566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&TSGLLEAcceptList, sname, function)); 12663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 12673d177a5cSEmil Constantinescu } 12683d177a5cSEmil Constantinescu 12693d177a5cSEmil Constantinescu /*@C 1270bcf0153eSBarry Smith TSGLLERegisterAll - Registers all of the general linear methods in `TSGLLE` 12713d177a5cSEmil Constantinescu 12723d177a5cSEmil Constantinescu Not Collective 12733d177a5cSEmil Constantinescu 12743d177a5cSEmil Constantinescu Level: advanced 12753d177a5cSEmil Constantinescu 12761cc06b55SBarry Smith .seealso: [](ch_ts), `TSGLLE`, `TSGLLERegisterDestroy()` 12773d177a5cSEmil Constantinescu @*/ 1278d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGLLERegisterAll(void) 1279d71ae5a4SJacob Faibussowitsch { 12803d177a5cSEmil Constantinescu PetscFunctionBegin; 12813ba16761SJacob Faibussowitsch if (TSGLLERegisterAllCalled) PetscFunctionReturn(PETSC_SUCCESS); 12823d177a5cSEmil Constantinescu TSGLLERegisterAllCalled = PETSC_TRUE; 12833d177a5cSEmil Constantinescu 12849566063dSJacob Faibussowitsch PetscCall(TSGLLERegister(TSGLLE_IRKS, TSGLLECreate_IRKS)); 12859566063dSJacob Faibussowitsch PetscCall(TSGLLEAcceptRegister(TSGLLEACCEPT_ALWAYS, TSGLLEAccept_Always)); 12863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 12873d177a5cSEmil Constantinescu } 12883d177a5cSEmil Constantinescu 12893d177a5cSEmil Constantinescu /*@C 1290bcf0153eSBarry Smith TSGLLEInitializePackage - This function initializes everything in the `TSGLLE` package. It is called 1291bcf0153eSBarry Smith from `TSInitializePackage()`. 12923d177a5cSEmil Constantinescu 12933d177a5cSEmil Constantinescu Level: developer 12943d177a5cSEmil Constantinescu 12951cc06b55SBarry Smith .seealso: [](ch_ts), `PetscInitialize()`, `TSInitializePackage()`, `TSGLLEFinalizePackage()` 12963d177a5cSEmil Constantinescu @*/ 1297d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGLLEInitializePackage(void) 1298d71ae5a4SJacob Faibussowitsch { 12993d177a5cSEmil Constantinescu PetscFunctionBegin; 13003ba16761SJacob Faibussowitsch if (TSGLLEPackageInitialized) PetscFunctionReturn(PETSC_SUCCESS); 13013d177a5cSEmil Constantinescu TSGLLEPackageInitialized = PETSC_TRUE; 13029566063dSJacob Faibussowitsch PetscCall(TSGLLERegisterAll()); 13039566063dSJacob Faibussowitsch PetscCall(PetscRegisterFinalize(TSGLLEFinalizePackage)); 13043ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13053d177a5cSEmil Constantinescu } 13063d177a5cSEmil Constantinescu 13073d177a5cSEmil Constantinescu /*@C 1308bcf0153eSBarry Smith TSGLLEFinalizePackage - This function destroys everything in the `TSGLLE` package. It is 1309bcf0153eSBarry Smith called from `PetscFinalize()`. 13103d177a5cSEmil Constantinescu 13113d177a5cSEmil Constantinescu Level: developer 13123d177a5cSEmil Constantinescu 13131cc06b55SBarry Smith .seealso: [](ch_ts), `PetscFinalize()`, `TSGLLEInitializePackage()`, `TSInitializePackage()` 13143d177a5cSEmil Constantinescu @*/ 1315d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGLLEFinalizePackage(void) 1316d71ae5a4SJacob Faibussowitsch { 13173d177a5cSEmil Constantinescu PetscFunctionBegin; 13189566063dSJacob Faibussowitsch PetscCall(PetscFunctionListDestroy(&TSGLLEList)); 13199566063dSJacob Faibussowitsch PetscCall(PetscFunctionListDestroy(&TSGLLEAcceptList)); 13203d177a5cSEmil Constantinescu TSGLLEPackageInitialized = PETSC_FALSE; 13213d177a5cSEmil Constantinescu TSGLLERegisterAllCalled = PETSC_FALSE; 13223ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13233d177a5cSEmil Constantinescu } 13243d177a5cSEmil Constantinescu 13253d177a5cSEmil Constantinescu /* ------------------------------------------------------------ */ 13263d177a5cSEmil Constantinescu /*MC 13271d27aa22SBarry Smith TSGLLE - DAE solver using implicit General Linear methods {cite}`butcher_2007` {cite}`butcher2016numerical` 13283d177a5cSEmil Constantinescu 1329bcf0153eSBarry Smith Options Database Keys: 13303d177a5cSEmil Constantinescu + -ts_gl_type <type> - the class of general linear method (irks) 13313d177a5cSEmil Constantinescu . -ts_gl_rtol <tol> - relative error 13323d177a5cSEmil Constantinescu . -ts_gl_atol <tol> - absolute error 13333d177a5cSEmil Constantinescu . -ts_gl_min_order <p> - minimum order method to consider (default=1) 13343d177a5cSEmil Constantinescu . -ts_gl_max_order <p> - maximum order method to consider (default=3) 13353d177a5cSEmil Constantinescu . -ts_gl_start_order <p> - order of starting method (default=1) 13363d177a5cSEmil Constantinescu . -ts_gl_complete <method> - method to use for completing the step (rescale-and-modify or rescale) 13373d177a5cSEmil Constantinescu - -ts_adapt_type <method> - adaptive controller to use (none step both) 13383d177a5cSEmil Constantinescu 1339bcf0153eSBarry Smith Level: beginner 1340bcf0153eSBarry Smith 13413d177a5cSEmil Constantinescu Notes: 134214d0ab18SJacob Faibussowitsch These methods contain Runge-Kutta and multistep schemes as special cases. These special cases 134314d0ab18SJacob Faibussowitsch have some fundamental limitations. For example, diagonally implicit Runge-Kutta cannot have 134414d0ab18SJacob Faibussowitsch stage order greater than 1 which limits their applicability to very stiff systems. 134514d0ab18SJacob Faibussowitsch Meanwhile, multistep methods cannot be A-stable for order greater than 2 and BDF are not 134614d0ab18SJacob Faibussowitsch 0-stable for order greater than 6. GL methods can be A- and L-stable with arbitrarily high 134714d0ab18SJacob Faibussowitsch stage order and reliable error estimates for both 1 and 2 orders higher to facilitate 134814d0ab18SJacob Faibussowitsch adaptive step sizes and adaptive order schemes. All this is possible while preserving a 134914d0ab18SJacob Faibussowitsch singly diagonally implicit structure. 135014d0ab18SJacob Faibussowitsch 13513d177a5cSEmil Constantinescu This integrator can be applied to DAE. 13523d177a5cSEmil Constantinescu 135314d0ab18SJacob Faibussowitsch Diagonally implicit general linear (DIGL) methods are a generalization of diagonally implicit 135414d0ab18SJacob Faibussowitsch Runge-Kutta (DIRK). They are represented by the tableau 13553d177a5cSEmil Constantinescu 13563d177a5cSEmil Constantinescu .vb 13573d177a5cSEmil Constantinescu A | U 13583d177a5cSEmil Constantinescu ------- 13593d177a5cSEmil Constantinescu B | V 13603d177a5cSEmil Constantinescu .ve 13613d177a5cSEmil Constantinescu 1362562efe2eSBarry Smith combined with a vector c of abscissa. "Diagonally implicit" means that $A$ is lower 136314d0ab18SJacob Faibussowitsch triangular. A step of the general method reads 13643d177a5cSEmil Constantinescu 136514d0ab18SJacob Faibussowitsch $$ 1366562efe2eSBarry Smith \begin{align*} 1367562efe2eSBarry Smith [ Y ] = [A U] [ Y' ] \\ 13683d177a5cSEmil Constantinescu [X^k] = [B V] [X^{k-1}] 1369562efe2eSBarry Smith \end{align*} 137014d0ab18SJacob Faibussowitsch $$ 13713d177a5cSEmil Constantinescu 1372562efe2eSBarry Smith where Y is the multivector of stage values, $Y'$ is the multivector of stage derivatives, $X^k$ 1373562efe2eSBarry Smith is the Nordsieck vector of the solution at step $k$. The Nordsieck vector consists of the first 1374562efe2eSBarry Smith $r$ moments of the solution, given by 13753d177a5cSEmil Constantinescu 137614d0ab18SJacob Faibussowitsch $$ 13773d177a5cSEmil Constantinescu X = [x_0,x_1,...,x_{r-1}] = [x, h x', h^2 x'', ..., h^{r-1} x^{(r-1)} ] 137814d0ab18SJacob Faibussowitsch $$ 13793d177a5cSEmil Constantinescu 1380562efe2eSBarry Smith If $A$ is lower triangular, we can solve the stages $(Y, Y')$ sequentially 13813d177a5cSEmil Constantinescu 138214d0ab18SJacob Faibussowitsch $$ 1383562efe2eSBarry Smith 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} 138414d0ab18SJacob Faibussowitsch $$ 13853d177a5cSEmil Constantinescu 13863d177a5cSEmil Constantinescu and then construct the pieces to carry to the next step 13873d177a5cSEmil Constantinescu 138814d0ab18SJacob Faibussowitsch $$ 1389562efe2eSBarry Smith 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} 139014d0ab18SJacob Faibussowitsch $$ 13913d177a5cSEmil Constantinescu 139214d0ab18SJacob Faibussowitsch Note that when the equations are cast in implicit form, we are using the stage equation to 1393562efe2eSBarry Smith define $y'_i$ in terms of $y_i$ and known stuff ($y_j$ for $j<i$ and $x_j$ for all $j$). 13943d177a5cSEmil Constantinescu 13953d177a5cSEmil Constantinescu Error estimation 13963d177a5cSEmil Constantinescu 139714d0ab18SJacob Faibussowitsch At present, the most attractive GL methods for stiff problems are singly diagonally implicit 1398562efe2eSBarry Smith schemes which posses Inherent Runge-Kutta Stability (`TSIRKS`). These methods have $r=s$, the 139914d0ab18SJacob Faibussowitsch number of items passed between steps is equal to the number of stages. The order and 140014d0ab18SJacob Faibussowitsch stage-order are one less than the number of stages. We use the error estimates in the 2007 140114d0ab18SJacob Faibussowitsch paper which provide the following estimates 14023d177a5cSEmil Constantinescu 140314d0ab18SJacob Faibussowitsch $$ 1404562efe2eSBarry Smith \begin{align*} 1405562efe2eSBarry Smith h^{p+1} X^{(p+1)} = \phi_0^T Y' + [0 \psi_0^T] Xold \\ 1406562efe2eSBarry Smith h^{p+2} X^{(p+2)} = \phi_1^T Y' + [0 \psi_1^T] Xold \\ 1407562efe2eSBarry Smith h^{p+2} (dx'/dx) X^{(p+1)} = \phi_2^T Y' + [0 \psi_2^T] Xold 1408562efe2eSBarry Smith \end{align*} 140914d0ab18SJacob Faibussowitsch $$ 14103d177a5cSEmil Constantinescu 1411562efe2eSBarry Smith These estimates are accurate to $ O(h^{p+3})$. 14123d177a5cSEmil Constantinescu 14133d177a5cSEmil Constantinescu Changing the step size 14143d177a5cSEmil Constantinescu 14151d27aa22SBarry Smith Uses the generalized "rescale and modify" scheme, see equation (4.5) of {cite}`butcher_2007`. 14163d177a5cSEmil Constantinescu 14171cc06b55SBarry Smith .seealso: [](ch_ts), `TSCreate()`, `TS`, `TSSetType()`, `TSType` 14183d177a5cSEmil Constantinescu M*/ 1419d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSCreate_GLLE(TS ts) 1420d71ae5a4SJacob Faibussowitsch { 14213d177a5cSEmil Constantinescu TS_GLLE *gl; 14223d177a5cSEmil Constantinescu 14233d177a5cSEmil Constantinescu PetscFunctionBegin; 14249566063dSJacob Faibussowitsch PetscCall(TSGLLEInitializePackage()); 14253d177a5cSEmil Constantinescu 14264dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&gl)); 14273d177a5cSEmil Constantinescu ts->data = (void *)gl; 14283d177a5cSEmil Constantinescu 14293d177a5cSEmil Constantinescu ts->ops->reset = TSReset_GLLE; 14303d177a5cSEmil Constantinescu ts->ops->destroy = TSDestroy_GLLE; 14313d177a5cSEmil Constantinescu ts->ops->view = TSView_GLLE; 14323d177a5cSEmil Constantinescu ts->ops->setup = TSSetUp_GLLE; 14333d177a5cSEmil Constantinescu ts->ops->solve = TSSolve_GLLE; 14343d177a5cSEmil Constantinescu ts->ops->setfromoptions = TSSetFromOptions_GLLE; 14353d177a5cSEmil Constantinescu ts->ops->snesfunction = SNESTSFormFunction_GLLE; 14363d177a5cSEmil Constantinescu ts->ops->snesjacobian = SNESTSFormJacobian_GLLE; 14373d177a5cSEmil Constantinescu 1438efd4aadfSBarry Smith ts->usessnes = PETSC_TRUE; 1439efd4aadfSBarry Smith 14403d177a5cSEmil Constantinescu gl->max_step_rejections = 1; 14413d177a5cSEmil Constantinescu gl->min_order = 1; 14423d177a5cSEmil Constantinescu gl->max_order = 3; 14433d177a5cSEmil Constantinescu gl->start_order = 1; 14443d177a5cSEmil Constantinescu gl->current_scheme = -1; 14453d177a5cSEmil Constantinescu gl->extrapolate = PETSC_FALSE; 14463d177a5cSEmil Constantinescu 14473d177a5cSEmil Constantinescu gl->wrms_atol = 1e-8; 14483d177a5cSEmil Constantinescu gl->wrms_rtol = 1e-5; 14493d177a5cSEmil Constantinescu 14509566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSGLLESetType_C", &TSGLLESetType_GLLE)); 14519566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSGLLESetAcceptType_C", &TSGLLESetAcceptType_GLLE)); 14529566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSGLLEGetAdapt_C", &TSGLLEGetAdapt_GLLE)); 14533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14543d177a5cSEmil Constantinescu } 1455