1 2 #include <../src/ts/impls/implicit/glle/glle.h> /*I "petscts.h" I*/ 3 #include <petscdm.h> 4 #include <petscblaslapack.h> 5 6 static const char *TSGLLEErrorDirections[] = {"FORWARD", "BACKWARD", "TSGLLEErrorDirection", "TSGLLEERROR_", NULL}; 7 static PetscFunctionList TSGLLEList; 8 static PetscFunctionList TSGLLEAcceptList; 9 static PetscBool TSGLLEPackageInitialized; 10 static PetscBool TSGLLERegisterAllCalled; 11 12 /* This function is pure */ 13 static PetscScalar Factorial(PetscInt n) { 14 PetscInt i; 15 if (n < 12) { /* Can compute with 32-bit integers */ 16 PetscInt f = 1; 17 for (i = 2; i <= n; i++) f *= i; 18 return (PetscScalar)f; 19 } else { 20 PetscScalar f = 1.; 21 for (i = 2; i <= n; i++) f *= (PetscScalar)i; 22 return f; 23 } 24 } 25 26 /* This function is pure */ 27 static PetscScalar CPowF(PetscScalar c, PetscInt p) { 28 return PetscPowRealInt(PetscRealPart(c), p) / Factorial(p); 29 } 30 31 static PetscErrorCode TSGLLEGetVecs(TS ts, DM dm, Vec *Z, Vec *Ydotstage) { 32 TS_GLLE *gl = (TS_GLLE *)ts->data; 33 34 PetscFunctionBegin; 35 if (Z) { 36 if (dm && dm != ts->dm) { 37 PetscCall(DMGetNamedGlobalVector(dm, "TSGLLE_Z", Z)); 38 } else *Z = gl->Z; 39 } 40 if (Ydotstage) { 41 if (dm && dm != ts->dm) { 42 PetscCall(DMGetNamedGlobalVector(dm, "TSGLLE_Ydot", Ydotstage)); 43 } else *Ydotstage = gl->Ydot[gl->stage]; 44 } 45 PetscFunctionReturn(0); 46 } 47 48 static PetscErrorCode TSGLLERestoreVecs(TS ts, DM dm, Vec *Z, Vec *Ydotstage) { 49 PetscFunctionBegin; 50 if (Z) { 51 if (dm && dm != ts->dm) { PetscCall(DMRestoreNamedGlobalVector(dm, "TSGLLE_Z", Z)); } 52 } 53 if (Ydotstage) { 54 if (dm && dm != ts->dm) { PetscCall(DMRestoreNamedGlobalVector(dm, "TSGLLE_Ydot", Ydotstage)); } 55 } 56 PetscFunctionReturn(0); 57 } 58 59 static PetscErrorCode DMCoarsenHook_TSGLLE(DM fine, DM coarse, void *ctx) { 60 PetscFunctionBegin; 61 PetscFunctionReturn(0); 62 } 63 64 static PetscErrorCode DMRestrictHook_TSGLLE(DM fine, Mat restrct, Vec rscale, Mat inject, DM coarse, void *ctx) { 65 TS ts = (TS)ctx; 66 Vec Ydot, Ydot_c; 67 68 PetscFunctionBegin; 69 PetscCall(TSGLLEGetVecs(ts, fine, NULL, &Ydot)); 70 PetscCall(TSGLLEGetVecs(ts, coarse, NULL, &Ydot_c)); 71 PetscCall(MatRestrict(restrct, Ydot, Ydot_c)); 72 PetscCall(VecPointwiseMult(Ydot_c, rscale, Ydot_c)); 73 PetscCall(TSGLLERestoreVecs(ts, fine, NULL, &Ydot)); 74 PetscCall(TSGLLERestoreVecs(ts, coarse, NULL, &Ydot_c)); 75 PetscFunctionReturn(0); 76 } 77 78 static PetscErrorCode DMSubDomainHook_TSGLLE(DM dm, DM subdm, void *ctx) { 79 PetscFunctionBegin; 80 PetscFunctionReturn(0); 81 } 82 83 static PetscErrorCode DMSubDomainRestrictHook_TSGLLE(DM dm, VecScatter gscat, VecScatter lscat, DM subdm, void *ctx) { 84 TS ts = (TS)ctx; 85 Vec Ydot, Ydot_s; 86 87 PetscFunctionBegin; 88 PetscCall(TSGLLEGetVecs(ts, dm, NULL, &Ydot)); 89 PetscCall(TSGLLEGetVecs(ts, subdm, NULL, &Ydot_s)); 90 91 PetscCall(VecScatterBegin(gscat, Ydot, Ydot_s, INSERT_VALUES, SCATTER_FORWARD)); 92 PetscCall(VecScatterEnd(gscat, Ydot, Ydot_s, INSERT_VALUES, SCATTER_FORWARD)); 93 94 PetscCall(TSGLLERestoreVecs(ts, dm, NULL, &Ydot)); 95 PetscCall(TSGLLERestoreVecs(ts, subdm, NULL, &Ydot_s)); 96 PetscFunctionReturn(0); 97 } 98 99 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) { 100 TSGLLEScheme scheme; 101 PetscInt j; 102 103 PetscFunctionBegin; 104 PetscCheck(p >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Scheme order must be positive"); 105 PetscCheck(r >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "At least one item must be carried between steps"); 106 PetscCheck(s >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "At least one stage is required"); 107 PetscValidPointer(inscheme, 10); 108 *inscheme = NULL; 109 PetscCall(PetscNew(&scheme)); 110 scheme->p = p; 111 scheme->q = q; 112 scheme->r = r; 113 scheme->s = s; 114 115 PetscCall(PetscMalloc5(s, &scheme->c, s * s, &scheme->a, r * s, &scheme->b, r * s, &scheme->u, r * r, &scheme->v)); 116 PetscCall(PetscArraycpy(scheme->c, c, s)); 117 for (j = 0; j < s * s; j++) scheme->a[j] = (PetscAbsScalar(a[j]) < 1e-12) ? 0 : a[j]; 118 for (j = 0; j < r * s; j++) scheme->b[j] = (PetscAbsScalar(b[j]) < 1e-12) ? 0 : b[j]; 119 for (j = 0; j < s * r; j++) scheme->u[j] = (PetscAbsScalar(u[j]) < 1e-12) ? 0 : u[j]; 120 for (j = 0; j < r * r; j++) scheme->v[j] = (PetscAbsScalar(v[j]) < 1e-12) ? 0 : v[j]; 121 122 PetscCall(PetscMalloc6(r, &scheme->alpha, r, &scheme->beta, r, &scheme->gamma, 3 * s, &scheme->phi, 3 * r, &scheme->psi, r, &scheme->stage_error)); 123 { 124 PetscInt i, j, k, ss = s + 2; 125 PetscBLASInt m, n, one = 1, *ipiv, lwork = 4 * ((s + 3) * 3 + 3), info, ldb; 126 PetscReal rcond, *sing, *workreal; 127 PetscScalar *ImV, *H, *bmat, *workscalar, *c = scheme->c, *a = scheme->a, *b = scheme->b, *u = scheme->u, *v = scheme->v; 128 PetscBLASInt rank; 129 PetscCall(PetscMalloc7(PetscSqr(r), &ImV, 3 * s, &H, 3 * ss, &bmat, lwork, &workscalar, 5 * (3 + r), &workreal, r + s, &sing, r + s, &ipiv)); 130 131 /* column-major input */ 132 for (i = 0; i < r - 1; i++) { 133 for (j = 0; j < r - 1; j++) ImV[i + j * r] = 1.0 * (i == j) - v[(i + 1) * r + j + 1]; 134 } 135 /* Build right hand side for alpha (tp - glm.B(2:end,:)*(glm.c.^(p)./factorial(p))) */ 136 for (i = 1; i < r; i++) { 137 scheme->alpha[i] = 1. / Factorial(p + 1 - i); 138 for (j = 0; j < s; j++) scheme->alpha[i] -= b[i * s + j] * CPowF(c[j], p); 139 } 140 PetscCall(PetscBLASIntCast(r - 1, &m)); 141 PetscCall(PetscBLASIntCast(r, &n)); 142 PetscCallBLAS("LAPACKgesv", LAPACKgesv_(&m, &one, ImV, &n, ipiv, scheme->alpha + 1, &n, &info)); 143 PetscCheck(info >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "Bad argument to GESV"); 144 PetscCheck(info <= 0, PETSC_COMM_SELF, PETSC_ERR_MAT_LU_ZRPVT, "Bad LU factorization"); 145 146 /* Build right hand side for beta (tp1 - glm.B(2:end,:)*(glm.c.^(p+1)./factorial(p+1)) - e.alpha) */ 147 for (i = 1; i < r; i++) { 148 scheme->beta[i] = 1. / Factorial(p + 2 - i) - scheme->alpha[i]; 149 for (j = 0; j < s; j++) scheme->beta[i] -= b[i * s + j] * CPowF(c[j], p + 1); 150 } 151 PetscCallBLAS("LAPACKgetrs", LAPACKgetrs_("No transpose", &m, &one, ImV, &n, ipiv, scheme->beta + 1, &n, &info)); 152 PetscCheck(info >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "Bad argument to GETRS"); 153 PetscCheck(info <= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "Should not happen"); 154 155 /* Build stage_error vector 156 xi = glm.c.^(p+1)/factorial(p+1) - glm.A*glm.c.^p/factorial(p) + glm.U(:,2:end)*e.alpha; 157 */ 158 for (i = 0; i < s; i++) { 159 scheme->stage_error[i] = CPowF(c[i], p + 1); 160 for (j = 0; j < s; j++) scheme->stage_error[i] -= a[i * s + j] * CPowF(c[j], p); 161 for (j = 1; j < r; j++) scheme->stage_error[i] += u[i * r + j] * scheme->alpha[j]; 162 } 163 164 /* alpha[0] (epsilon in B,J,W 2007) 165 epsilon = 1/factorial(p+1) - B(1,:)*c.^p/factorial(p) + V(1,2:end)*e.alpha; 166 */ 167 scheme->alpha[0] = 1. / Factorial(p + 1); 168 for (j = 0; j < s; j++) scheme->alpha[0] -= b[0 * s + j] * CPowF(c[j], p); 169 for (j = 1; j < r; j++) scheme->alpha[0] += v[0 * r + j] * scheme->alpha[j]; 170 171 /* right hand side for gamma (glm.B(2:end,:)*e.xi - e.epsilon*eye(s-1,1)) */ 172 for (i = 1; i < r; i++) { 173 scheme->gamma[i] = (i == 1 ? -1. : 0) * scheme->alpha[0]; 174 for (j = 0; j < s; j++) scheme->gamma[i] += b[i * s + j] * scheme->stage_error[j]; 175 } 176 PetscCallBLAS("LAPACKgetrs", LAPACKgetrs_("No transpose", &m, &one, ImV, &n, ipiv, scheme->gamma + 1, &n, &info)); 177 PetscCheck(info >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "Bad argument to GETRS"); 178 PetscCheck(info <= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "Should not happen"); 179 180 /* beta[0] (rho in B,J,W 2007) 181 e.rho = 1/factorial(p+2) - glm.B(1,:)*glm.c.^(p+1)/factorial(p+1) ... 182 + glm.V(1,2:end)*e.beta;% - e.epsilon; 183 % Note: The paper (B,J,W 2007) includes the last term in their definition 184 * */ 185 scheme->beta[0] = 1. / Factorial(p + 2); 186 for (j = 0; j < s; j++) scheme->beta[0] -= b[0 * s + j] * CPowF(c[j], p + 1); 187 for (j = 1; j < r; j++) scheme->beta[0] += v[0 * r + j] * scheme->beta[j]; 188 189 /* gamma[0] (sigma in B,J,W 2007) 190 * e.sigma = glm.B(1,:)*e.xi + glm.V(1,2:end)*e.gamma; 191 * */ 192 scheme->gamma[0] = 0.0; 193 for (j = 0; j < s; j++) scheme->gamma[0] += b[0 * s + j] * scheme->stage_error[j]; 194 for (j = 1; j < r; j++) scheme->gamma[0] += v[0 * s + j] * scheme->gamma[j]; 195 196 /* Assemble H 197 * % " PetscInt_FMT "etermine the error estimators phi 198 H = [[cpow(glm.c,p) + C*e.alpha] [cpow(glm.c,p+1) + C*e.beta] ... 199 [e.xi - C*(e.gamma + 0*e.epsilon*eye(s-1,1))]]'; 200 % Paper has formula above without the 0, but that term must be left 201 % out to satisfy the conditions they propose and to make the 202 % example schemes work 203 e.H = H; 204 e.phi = (H \ [1 0 0;1 1 0;0 0 -1])'; 205 e.psi = -e.phi*C; 206 * */ 207 for (j = 0; j < s; j++) { 208 H[0 + j * 3] = CPowF(c[j], p); 209 H[1 + j * 3] = CPowF(c[j], p + 1); 210 H[2 + j * 3] = scheme->stage_error[j]; 211 for (k = 1; k < r; k++) { 212 H[0 + j * 3] += CPowF(c[j], k - 1) * scheme->alpha[k]; 213 H[1 + j * 3] += CPowF(c[j], k - 1) * scheme->beta[k]; 214 H[2 + j * 3] -= CPowF(c[j], k - 1) * scheme->gamma[k]; 215 } 216 } 217 bmat[0 + 0 * ss] = 1.; 218 bmat[0 + 1 * ss] = 0.; 219 bmat[0 + 2 * ss] = 0.; 220 bmat[1 + 0 * ss] = 1.; 221 bmat[1 + 1 * ss] = 1.; 222 bmat[1 + 2 * ss] = 0.; 223 bmat[2 + 0 * ss] = 0.; 224 bmat[2 + 1 * ss] = 0.; 225 bmat[2 + 2 * ss] = -1.; 226 m = 3; 227 PetscCall(PetscBLASIntCast(s, &n)); 228 PetscCall(PetscBLASIntCast(ss, &ldb)); 229 rcond = 1e-12; 230 #if defined(PETSC_USE_COMPLEX) 231 /* ZGELSS( M, N, NRHS, A, LDA, B, LDB, S, RCOND, RANK, WORK, LWORK, RWORK, INFO) */ 232 PetscCallBLAS("LAPACKgelss", LAPACKgelss_(&m, &n, &m, H, &m, bmat, &ldb, sing, &rcond, &rank, workscalar, &lwork, workreal, &info)); 233 #else 234 /* DGELSS( M, N, NRHS, A, LDA, B, LDB, S, RCOND, RANK, WORK, LWORK, INFO) */ 235 PetscCallBLAS("LAPACKgelss", LAPACKgelss_(&m, &n, &m, H, &m, bmat, &ldb, sing, &rcond, &rank, workscalar, &lwork, &info)); 236 #endif 237 PetscCheck(info >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "Bad argument to GELSS"); 238 PetscCheck(info <= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "SVD failed to converge"); 239 240 for (j = 0; j < 3; j++) { 241 for (k = 0; k < s; k++) scheme->phi[k + j * s] = bmat[k + j * ss]; 242 } 243 244 /* the other part of the error estimator, psi in B,J,W 2007 */ 245 scheme->psi[0 * r + 0] = 0.; 246 scheme->psi[1 * r + 0] = 0.; 247 scheme->psi[2 * r + 0] = 0.; 248 for (j = 1; j < r; j++) { 249 scheme->psi[0 * r + j] = 0.; 250 scheme->psi[1 * r + j] = 0.; 251 scheme->psi[2 * r + j] = 0.; 252 for (k = 0; k < s; k++) { 253 scheme->psi[0 * r + j] -= CPowF(c[k], j - 1) * scheme->phi[0 * s + k]; 254 scheme->psi[1 * r + j] -= CPowF(c[k], j - 1) * scheme->phi[1 * s + k]; 255 scheme->psi[2 * r + j] -= CPowF(c[k], j - 1) * scheme->phi[2 * s + k]; 256 } 257 } 258 PetscCall(PetscFree7(ImV, H, bmat, workscalar, workreal, sing, ipiv)); 259 } 260 /* Check which properties are satisfied */ 261 scheme->stiffly_accurate = PETSC_TRUE; 262 if (scheme->c[s - 1] != 1.) scheme->stiffly_accurate = PETSC_FALSE; 263 for (j = 0; j < s; j++) 264 if (a[(s - 1) * s + j] != b[j]) scheme->stiffly_accurate = PETSC_FALSE; 265 for (j = 0; j < r; j++) 266 if (u[(s - 1) * r + j] != v[j]) scheme->stiffly_accurate = PETSC_FALSE; 267 scheme->fsal = scheme->stiffly_accurate; /* FSAL is stronger */ 268 for (j = 0; j < s - 1; j++) 269 if (r > 1 && b[1 * s + j] != 0.) scheme->fsal = PETSC_FALSE; 270 if (b[1 * s + r - 1] != 1.) scheme->fsal = PETSC_FALSE; 271 for (j = 0; j < r; j++) 272 if (r > 1 && v[1 * r + j] != 0.) scheme->fsal = PETSC_FALSE; 273 274 *inscheme = scheme; 275 PetscFunctionReturn(0); 276 } 277 278 static PetscErrorCode TSGLLESchemeDestroy(TSGLLEScheme sc) { 279 PetscFunctionBegin; 280 PetscCall(PetscFree5(sc->c, sc->a, sc->b, sc->u, sc->v)); 281 PetscCall(PetscFree6(sc->alpha, sc->beta, sc->gamma, sc->phi, sc->psi, sc->stage_error)); 282 PetscCall(PetscFree(sc)); 283 PetscFunctionReturn(0); 284 } 285 286 static PetscErrorCode TSGLLEDestroy_Default(TS_GLLE *gl) { 287 PetscInt i; 288 289 PetscFunctionBegin; 290 for (i = 0; i < gl->nschemes; i++) { 291 if (gl->schemes[i]) PetscCall(TSGLLESchemeDestroy(gl->schemes[i])); 292 } 293 PetscCall(PetscFree(gl->schemes)); 294 gl->nschemes = 0; 295 PetscCall(PetscMemzero(gl->type_name, sizeof(gl->type_name))); 296 PetscFunctionReturn(0); 297 } 298 299 static PetscErrorCode TSGLLEViewTable_Private(PetscViewer viewer, PetscInt m, PetscInt n, const PetscScalar a[], const char name[]) { 300 PetscBool iascii; 301 PetscInt i, j; 302 303 PetscFunctionBegin; 304 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 305 if (iascii) { 306 PetscCall(PetscViewerASCIIPrintf(viewer, "%30s = [", name)); 307 for (i = 0; i < m; i++) { 308 if (i) PetscCall(PetscViewerASCIIPrintf(viewer, "%30s [", "")); 309 PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE)); 310 for (j = 0; j < n; j++) { PetscCall(PetscViewerASCIIPrintf(viewer, " %12.8g", (double)PetscRealPart(a[i * n + j]))); } 311 PetscCall(PetscViewerASCIIPrintf(viewer, "]\n")); 312 PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE)); 313 } 314 } 315 PetscFunctionReturn(0); 316 } 317 318 static PetscErrorCode TSGLLESchemeView(TSGLLEScheme sc, PetscBool view_details, PetscViewer viewer) { 319 PetscBool iascii; 320 321 PetscFunctionBegin; 322 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 323 if (iascii) { 324 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)); 325 PetscCall(PetscViewerASCIIPushTab(viewer)); 326 PetscCall(PetscViewerASCIIPrintf(viewer, "Stiffly accurate: %s, FSAL: %s\n", sc->stiffly_accurate ? "yes" : "no", sc->fsal ? "yes" : "no")); 327 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]))); 328 PetscCall(TSGLLEViewTable_Private(viewer, 1, sc->s, sc->c, "Abscissas c")); 329 if (view_details) { 330 PetscCall(TSGLLEViewTable_Private(viewer, sc->s, sc->s, sc->a, "A")); 331 PetscCall(TSGLLEViewTable_Private(viewer, sc->r, sc->s, sc->b, "B")); 332 PetscCall(TSGLLEViewTable_Private(viewer, sc->s, sc->r, sc->u, "U")); 333 PetscCall(TSGLLEViewTable_Private(viewer, sc->r, sc->r, sc->v, "V")); 334 335 PetscCall(TSGLLEViewTable_Private(viewer, 3, sc->s, sc->phi, "Error estimate phi")); 336 PetscCall(TSGLLEViewTable_Private(viewer, 3, sc->r, sc->psi, "Error estimate psi")); 337 PetscCall(TSGLLEViewTable_Private(viewer, 1, sc->r, sc->alpha, "Modify alpha")); 338 PetscCall(TSGLLEViewTable_Private(viewer, 1, sc->r, sc->beta, "Modify beta")); 339 PetscCall(TSGLLEViewTable_Private(viewer, 1, sc->r, sc->gamma, "Modify gamma")); 340 PetscCall(TSGLLEViewTable_Private(viewer, 1, sc->s, sc->stage_error, "Stage error xi")); 341 } 342 PetscCall(PetscViewerASCIIPopTab(viewer)); 343 } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Viewer type %s not supported", ((PetscObject)viewer)->type_name); 344 PetscFunctionReturn(0); 345 } 346 347 static PetscErrorCode TSGLLEEstimateHigherMoments_Default(TSGLLEScheme sc, PetscReal h, Vec Ydot[], Vec Xold[], Vec hm[]) { 348 PetscInt i; 349 350 PetscFunctionBegin; 351 PetscCheck(sc->r <= 64 && sc->s <= 64, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Ridiculous number of stages or items passed between stages"); 352 /* build error vectors*/ 353 for (i = 0; i < 3; i++) { 354 PetscScalar phih[64]; 355 PetscInt j; 356 for (j = 0; j < sc->s; j++) phih[j] = sc->phi[i * sc->s + j] * h; 357 PetscCall(VecZeroEntries(hm[i])); 358 PetscCall(VecMAXPY(hm[i], sc->s, phih, Ydot)); 359 PetscCall(VecMAXPY(hm[i], sc->r, &sc->psi[i * sc->r], Xold)); 360 } 361 PetscFunctionReturn(0); 362 } 363 364 static PetscErrorCode TSGLLECompleteStep_Rescale(TSGLLEScheme sc, PetscReal h, TSGLLEScheme next_sc, PetscReal next_h, Vec Ydot[], Vec Xold[], Vec X[]) { 365 PetscScalar brow[32], vrow[32]; 366 PetscInt i, j, r, s; 367 368 PetscFunctionBegin; 369 /* Build the new solution from (X,Ydot) */ 370 r = sc->r; 371 s = sc->s; 372 for (i = 0; i < r; i++) { 373 PetscCall(VecZeroEntries(X[i])); 374 for (j = 0; j < s; j++) brow[j] = h * sc->b[i * s + j]; 375 PetscCall(VecMAXPY(X[i], s, brow, Ydot)); 376 for (j = 0; j < r; j++) vrow[j] = sc->v[i * r + j]; 377 PetscCall(VecMAXPY(X[i], r, vrow, Xold)); 378 } 379 PetscFunctionReturn(0); 380 } 381 382 static PetscErrorCode TSGLLECompleteStep_RescaleAndModify(TSGLLEScheme sc, PetscReal h, TSGLLEScheme next_sc, PetscReal next_h, Vec Ydot[], Vec Xold[], Vec X[]) { 383 PetscScalar brow[32], vrow[32]; 384 PetscReal ratio; 385 PetscInt i, j, p, r, s; 386 387 PetscFunctionBegin; 388 /* Build the new solution from (X,Ydot) */ 389 p = sc->p; 390 r = sc->r; 391 s = sc->s; 392 ratio = next_h / h; 393 for (i = 0; i < r; i++) { 394 PetscCall(VecZeroEntries(X[i])); 395 for (j = 0; j < s; j++) { 396 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])); 397 } 398 PetscCall(VecMAXPY(X[i], s, brow, Ydot)); 399 for (j = 0; j < r; j++) { 400 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])); 401 } 402 PetscCall(VecMAXPY(X[i], r, vrow, Xold)); 403 } 404 if (r < next_sc->r) { 405 PetscCheck(r + 1 == next_sc->r, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Cannot accommodate jump in r greater than 1"); 406 PetscCall(VecZeroEntries(X[r])); 407 for (j = 0; j < s; j++) brow[j] = h * PetscPowRealInt(ratio, p + 1) * sc->phi[0 * s + j]; 408 PetscCall(VecMAXPY(X[r], s, brow, Ydot)); 409 for (j = 0; j < r; j++) vrow[j] = PetscPowRealInt(ratio, p + 1) * sc->psi[0 * r + j]; 410 PetscCall(VecMAXPY(X[r], r, vrow, Xold)); 411 } 412 PetscFunctionReturn(0); 413 } 414 415 static PetscErrorCode TSGLLECreate_IRKS(TS ts) { 416 TS_GLLE *gl = (TS_GLLE *)ts->data; 417 418 PetscFunctionBegin; 419 gl->Destroy = TSGLLEDestroy_Default; 420 gl->EstimateHigherMoments = TSGLLEEstimateHigherMoments_Default; 421 gl->CompleteStep = TSGLLECompleteStep_RescaleAndModify; 422 PetscCall(PetscMalloc1(10, &gl->schemes)); 423 gl->nschemes = 0; 424 425 { 426 /* p=1,q=1, r=s=2, A- and L-stable with error estimates of order 2 and 3 427 * Listed in Butcher & Podhaisky 2006. On error estimation in general linear methods for stiff ODE. 428 * irks(0.3,0,[.3,1],[1],1) 429 * Note: can be made to have classical order (not stage order) 2 by replacing 0.3 with 1-sqrt(1/2) 430 * but doing so would sacrifice the error estimator. 431 */ 432 const PetscScalar c[2] = {3. / 10., 1.}; 433 const PetscScalar a[2][2] = { 434 {3. / 10., 0 }, 435 {7. / 10., 3. / 10.} 436 }; 437 const PetscScalar b[2][2] = { 438 {7. / 10., 3. / 10.}, 439 {0, 1 } 440 }; 441 const PetscScalar u[2][2] = { 442 {1, 0}, 443 {1, 0} 444 }; 445 const PetscScalar v[2][2] = { 446 {1, 0}, 447 {0, 0} 448 }; 449 PetscCall(TSGLLESchemeCreate(1, 1, 2, 2, c, *a, *b, *u, *v, &gl->schemes[gl->nschemes++])); 450 } 451 452 { 453 /* p=q=2, r=s=3: irks(4/9,0,[1:3]/3,[0.33852],1) */ 454 /* http://www.math.auckland.ac.nz/~hpod/atlas/i2a.html */ 455 const PetscScalar c[3] = 456 { 457 1. / 3., 2. / 3., 1 458 }, 459 a[3][3] = {{4. / 9., 0, 0}, {1.03750643704090e+00, 4. / 9., 0}, {7.67024779410304e-01, -3.81140216918943e-01, 4. / 9.}}, b[3][3] = {{0.767024779410304, -0.381140216918943, 4. / 9.}, {0.000000000000000, 0.000000000000000, 1.000000000000000}, {-2.075048385225385, 0.621728385225383, 1.277197204924873}}, u[3][3] = {{1.0000000000000000, -0.1111111111111109, -0.0925925925925922}, {1.0000000000000000, -0.8152842148186744, -0.4199095530877056}, {1.0000000000000000, 0.1696709930641948, 0.0539741070314165}}, v[3][3] = {{1.0000000000000000, 0.1696709930641948, 0.0539741070314165}, {0.000000000000000, 0.000000000000000, 0.000000000000000}, {0.000000000000000, 0.176122795075129, 0.000000000000000}}; 460 PetscCall(TSGLLESchemeCreate(2, 2, 3, 3, c, *a, *b, *u, *v, &gl->schemes[gl->nschemes++])); 461 } 462 { 463 /* p=q=3, r=s=4: irks(9/40,0,[1:4]/4,[0.3312 1.0050],[0.49541 1;1 0]) */ 464 const PetscScalar c[4] = 465 { 466 0.25, 0.5, 0.75, 1.0 467 }, 468 a[4][4] = {{9. / 40., 0, 0, 0}, {2.11286958887701e-01, 9. / 40., 0, 0}, {9.46338294287584e-01, -3.42942861246094e-01, 9. / 40., 0}, {0.521490453970721, -0.662474225622980, 0.490476425459734, 9. / 40.}}, b[4][4] = {{0.521490453970721, -0.662474225622980, 0.490476425459734, 9. / 40.}, {0.000000000000000, 0.000000000000000, 0.000000000000000, 1.000000000000000}, {-0.084677029310348, 1.390757514776085, -1.568157386206001, 2.023192696767826}, {0.465383797936408, 1.478273530625148, -1.930836081010182, 1.644872111193354}}, u[4][4] = {{1.00000000000000000, 0.02500000000001035, -0.02499999999999053, -0.00442708333332865}, {1.00000000000000000, 0.06371304111232945, -0.04032173972189845, -0.01389438413189452}, {1.00000000000000000, -0.07839543304147778, 0.04738685705116663, 0.02032603595928376}, {1.00000000000000000, 0.42550734619251651, 0.10800718022400080, -0.01726712647760034}}, v[4][4] = {{1.00000000000000000, 0.42550734619251651, 0.10800718022400080, -0.01726712647760034}, {0.000000000000000, 0.000000000000000, 0.000000000000000, 0.000000000000000}, {0.000000000000000, -1.761115796027561, -0.521284157173780, 0.258249384305463}, {0.000000000000000, -1.657693358744728, -1.052227765232394, 0.521284157173780}}; 469 PetscCall(TSGLLESchemeCreate(3, 3, 4, 4, c, *a, *b, *u, *v, &gl->schemes[gl->nschemes++])); 470 } 471 { 472 /* p=q=4, r=s=5: 473 irks(3/11,0,[1:5]/5, [0.1715 -0.1238 0.6617],... 474 [ -0.0812 0.4079 1.0000 475 1.0000 0 0 476 0.8270 1.0000 0]) 477 */ 478 const PetscScalar c[5] = 479 { 480 0.2, 0.4, 0.6, 0.8, 1.0 481 }, 482 a[5][5] = {{2.72727272727352e-01, 0.00000000000000e+00, 0.00000000000000e+00, 0.00000000000000e+00, 0.00000000000000e+00}, {-1.03980153733431e-01, 2.72727272727405e-01, 0.00000000000000e+00, 0.00000000000000e+00, 0.00000000000000e+00}, {-1.58615400341492e+00, 7.44168951881122e-01, 2.72727272727309e-01, 0.00000000000000e+00, 0.00000000000000e+00}, {-8.73658042865628e-01, 5.37884671894595e-01, -1.63298538799523e-01, 2.72727272726996e-01, 0.00000000000000e+00}, {2.95489397443992e-01, -1.18481693910097e+00, -6.68029812659953e-01, 1.00716687860943e+00, 2.72727272727288e-01}}, b[5][5] = {{2.95489397443992e-01, -1.18481693910097e+00, -6.68029812659953e-01, 1.00716687860943e+00, 2.72727272727288e-01}, {0.00000000000000e+00, 1.11022302462516e-16, -2.22044604925031e-16, 0.00000000000000e+00, 1.00000000000000e+00}, {-4.05882503986005e+00, -4.00924006567769e+00, -1.38930610972481e+00, 4.45223930308488e+00, 6.32331093108427e-01}, {8.35690179937017e+00, -2.26640927349732e+00, 6.86647884973826e+00, -5.22595158025740e+00, 4.50893068837431e+00}, {1.27656267027479e+01, 2.80882153840821e+00, 8.91173096522890e+00, -1.07936444078906e+01, 4.82534148988854e+00}}, u[5][5] = {{1.00000000000000e+00, -7.27272727273551e-02, -3.45454545454419e-02, -4.12121212119565e-03, -2.96969696964014e-04}, {1.00000000000000e+00, 2.31252881006154e-01, -8.29487834416481e-03, -9.07191207681020e-03, -1.70378403743473e-03}, {1.00000000000000e+00, 1.16925777880663e+00, 3.59268562942635e-02, -4.09013451730615e-02, -1.02411119670164e-02}, {1.00000000000000e+00, 1.02634463704356e+00, 1.59375044913405e-01, 1.89673015035370e-03, -4.89987231897569e-03}, {1.00000000000000e+00, 1.27746320298021e+00, 2.37186008132728e-01, -8.28694373940065e-02, -5.34396510196430e-02}}, v[5][5] = {{1.00000000000000e+00, 1.27746320298021e+00, 2.37186008132728e-01, -8.28694373940065e-02, -5.34396510196430e-02}, {0.00000000000000e+00, -1.77635683940025e-15, -1.99840144432528e-15, -9.99200722162641e-16, -3.33066907387547e-16}, {0.00000000000000e+00, 4.37280081906924e+00, 5.49221645016377e-02, -8.88913177394943e-02, 1.12879077989154e-01}, {0.00000000000000e+00, -1.22399504837280e+01, -5.21287338448645e+00, -8.03952325565291e-01, 4.60298678047147e-01}, {0.00000000000000e+00, -1.85178762883829e+01, -5.21411849862624e+00, -1.04283436528809e+00, 7.49030161063651e-01}}; 483 PetscCall(TSGLLESchemeCreate(4, 4, 5, 5, c, *a, *b, *u, *v, &gl->schemes[gl->nschemes++])); 484 } 485 { 486 /* p=q=5, r=s=6; 487 irks(1/3,0,[1:6]/6,... 488 [-0.0489 0.4228 -0.8814 0.9021],... 489 [-0.3474 -0.6617 0.6294 0.2129 490 0.0044 -0.4256 -0.1427 -0.8936 491 -0.8267 0.4821 0.1371 -0.2557 492 -0.4426 -0.3855 -0.7514 0.3014]) 493 */ 494 const PetscScalar c[6] = 495 { 496 1. / 6, 2. / 6, 3. / 6, 4. / 6, 5. / 6, 1. 497 }, 498 a[6][6] = {{3.33333333333940e-01, 0, 0, 0, 0, 0}, {-8.64423857333350e-02, 3.33333333332888e-01, 0, 0, 0, 0}, {-2.16850174258252e+00, -2.23619072028839e+00, 3.33333333335204e-01, 0, 0, 0}, {-4.73160970138997e+00, -3.89265344629268e+00, -2.76318716520933e-01, 3.33333333335759e-01, 0, 0}, {-6.75187540297338e+00, -7.90756533769377e+00, 7.90245051802259e-01, -4.48352364517632e-01, 3.33333333328483e-01, 0}, {-4.26488287921548e+00, -1.19320395589302e+01, 3.38924509887755e+00, -2.23969848002481e+00, 6.62807710124007e-01, 3.33333333335440e-01}}, b[6][6] = {{-4.26488287921548e+00, -1.19320395589302e+01, 3.38924509887755e+00, -2.23969848002481e+00, 6.62807710124007e-01, 3.33333333335440e-01}, {-8.88178419700125e-16, 4.44089209850063e-16, -1.54737334057131e-15, -8.88178419700125e-16, 0.00000000000000e+00, 1.00000000000001e+00}, {-2.87780425770651e+01, -1.13520448264971e+01, 2.62002318943161e+01, 2.56943874812797e+01, -3.06702268304488e+01, 6.68067773510103e+00}, {5.47971245256474e+01, 6.80366875868284e+01, -6.50952588861999e+01, -8.28643975339097e+01, 8.17416943896414e+01, -1.17819043489036e+01}, {-2.33332114788869e+02, 6.12942539462634e+01, -4.91850135865944e+01, 1.82716844135480e+02, -1.29788173979395e+02, 3.09968095651099e+01}, {-1.72049132343751e+02, 8.60194713593999e+00, 7.98154219170200e-01, 1.50371386053218e+02, -1.18515423962066e+02, 2.50898277784663e+01}}, u[6][6] = {{1.00000000000000e+00, -1.66666666666870e-01, -4.16666666664335e-02, -3.85802469124815e-03, -2.25051440302250e-04, -9.64506172339142e-06}, {1.00000000000000e+00, 8.64423857327162e-02, -4.11484912671353e-02, -1.11450903217645e-02, -1.47651050487126e-03, -1.34395070766826e-04}, {1.00000000000000e+00, 4.57135912953434e+00, 1.06514719719137e+00, 1.33517564218007e-01, 1.11365952968659e-02, 6.12382756769504e-04}, {1.00000000000000e+00, 9.23391519753404e+00, 2.22431212392095e+00, 2.91823807741891e-01, 2.52058456411084e-02, 1.22800542949647e-03}, {1.00000000000000e+00, 1.48175480533865e+01, 3.73439117461835e+00, 5.14648336541804e-01, 4.76430038853402e-02, 2.56798515502156e-03}, {1.00000000000000e+00, 1.50512347758335e+01, 4.10099701165164e+00, 5.66039141003603e-01, 3.91213893800891e-02, -2.99136269067853e-03}}, v[6][6] = {{1.00000000000000e+00, 1.50512347758335e+01, 4.10099701165164e+00, 5.66039141003603e-01, 3.91213893800891e-02, -2.99136269067853e-03}, {0.00000000000000e+00, -4.88498130835069e-15, -6.43929354282591e-15, -3.55271367880050e-15, -1.22124532708767e-15, -3.12250225675825e-16}, {0.00000000000000e+00, 1.22250171233141e+01, -1.77150760606169e+00, 3.54516769879390e-01, 6.22298845883398e-01, 2.31647447450276e-01}, {0.00000000000000e+00, -4.48339457331040e+01, -3.57363126641880e-01, 5.18750173123425e-01, 6.55727990241799e-02, 1.63175368287079e-01}, {0.00000000000000e+00, 1.37297394708005e+02, -1.60145272991317e+00, -5.05319555199441e+00, 1.55328940390990e-01, 9.16629423682464e-01}, {0.00000000000000e+00, 1.05703241119022e+02, -1.16610260983038e+00, -2.99767252773859e+00, -1.13472315553890e-01, 1.09742849254729e+00}}; 499 PetscCall(TSGLLESchemeCreate(5, 5, 6, 6, c, *a, *b, *u, *v, &gl->schemes[gl->nschemes++])); 500 } 501 PetscFunctionReturn(0); 502 } 503 504 /*@C 505 TSGLLESetType - sets the class of general linear method to use for time-stepping 506 507 Collective on TS 508 509 Input Parameters: 510 + ts - the TS context 511 - type - a method 512 513 Options Database Key: 514 . -ts_gl_type <type> - sets the method, use -help for a list of available method (e.g. irks) 515 516 Notes: 517 See "petsc/include/petscts.h" for available methods (for instance) 518 . TSGLLE_IRKS - Diagonally implicit methods with inherent Runge-Kutta stability (for stiff problems) 519 520 Normally, it is best to use the TSSetFromOptions() command and 521 then set the TSGLLE type from the options database rather than by using 522 this routine. Using the options database provides the user with 523 maximum flexibility in evaluating the many different solvers. 524 The TSGLLESetType() routine is provided for those situations where it 525 is necessary to set the timestepping solver independently of the 526 command line or options database. This might be the case, for example, 527 when the choice of solver changes during the execution of the 528 program, and the user's application is taking responsibility for 529 choosing the appropriate method. 530 531 Level: intermediate 532 533 @*/ 534 PetscErrorCode TSGLLESetType(TS ts, TSGLLEType type) { 535 PetscFunctionBegin; 536 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 537 PetscValidCharPointer(type, 2); 538 PetscTryMethod(ts, "TSGLLESetType_C", (TS, TSGLLEType), (ts, type)); 539 PetscFunctionReturn(0); 540 } 541 542 /*@C 543 TSGLLESetAcceptType - sets the acceptance test 544 545 Time integrators that need to control error must have the option to reject a time step based on local error 546 estimates. This function allows different schemes to be set. 547 548 Logically Collective on TS 549 550 Input Parameters: 551 + ts - the TS context 552 - type - the type 553 554 Options Database Key: 555 . -ts_gl_accept_type <type> - sets the method used to determine whether to accept or reject a step 556 557 Level: intermediate 558 559 .seealso: `TS`, `TSGLLE`, `TSGLLEAcceptRegister()`, `TSGLLEAdapt`, `set` `type` 560 @*/ 561 PetscErrorCode TSGLLESetAcceptType(TS ts, TSGLLEAcceptType type) { 562 PetscFunctionBegin; 563 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 564 PetscValidCharPointer(type, 2); 565 PetscTryMethod(ts, "TSGLLESetAcceptType_C", (TS, TSGLLEAcceptType), (ts, type)); 566 PetscFunctionReturn(0); 567 } 568 569 /*@C 570 TSGLLEGetAdapt - gets the TSGLLEAdapt object from the TS 571 572 Not Collective 573 574 Input Parameter: 575 . ts - the TS context 576 577 Output Parameter: 578 . adapt - the TSGLLEAdapt context 579 580 Notes: 581 This allows the user set options on the TSGLLEAdapt object. Usually it is better to do this using the options 582 database, so this function is rarely needed. 583 584 Level: advanced 585 586 .seealso: `TSGLLEAdapt`, `TSGLLEAdaptRegister()` 587 @*/ 588 PetscErrorCode TSGLLEGetAdapt(TS ts, TSGLLEAdapt *adapt) { 589 PetscFunctionBegin; 590 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 591 PetscValidPointer(adapt, 2); 592 PetscUseMethod(ts, "TSGLLEGetAdapt_C", (TS, TSGLLEAdapt *), (ts, adapt)); 593 PetscFunctionReturn(0); 594 } 595 596 static PetscErrorCode TSGLLEAccept_Always(TS ts, PetscReal tleft, PetscReal h, const PetscReal enorms[], PetscBool *accept) { 597 PetscFunctionBegin; 598 *accept = PETSC_TRUE; 599 PetscFunctionReturn(0); 600 } 601 602 static PetscErrorCode TSGLLEUpdateWRMS(TS ts) { 603 TS_GLLE *gl = (TS_GLLE *)ts->data; 604 PetscScalar *x, *w; 605 PetscInt n, i; 606 607 PetscFunctionBegin; 608 PetscCall(VecGetArray(gl->X[0], &x)); 609 PetscCall(VecGetArray(gl->W, &w)); 610 PetscCall(VecGetLocalSize(gl->W, &n)); 611 for (i = 0; i < n; i++) w[i] = 1. / (gl->wrms_atol + gl->wrms_rtol * PetscAbsScalar(x[i])); 612 PetscCall(VecRestoreArray(gl->X[0], &x)); 613 PetscCall(VecRestoreArray(gl->W, &w)); 614 PetscFunctionReturn(0); 615 } 616 617 static PetscErrorCode TSGLLEVecNormWRMS(TS ts, Vec X, PetscReal *nrm) { 618 TS_GLLE *gl = (TS_GLLE *)ts->data; 619 PetscScalar *x, *w; 620 PetscReal sum = 0.0, gsum; 621 PetscInt n, N, i; 622 623 PetscFunctionBegin; 624 PetscCall(VecGetArray(X, &x)); 625 PetscCall(VecGetArray(gl->W, &w)); 626 PetscCall(VecGetLocalSize(gl->W, &n)); 627 for (i = 0; i < n; i++) sum += PetscAbsScalar(PetscSqr(x[i] * w[i])); 628 PetscCall(VecRestoreArray(X, &x)); 629 PetscCall(VecRestoreArray(gl->W, &w)); 630 PetscCall(MPIU_Allreduce(&sum, &gsum, 1, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)ts))); 631 PetscCall(VecGetSize(gl->W, &N)); 632 *nrm = PetscSqrtReal(gsum / (1. * N)); 633 PetscFunctionReturn(0); 634 } 635 636 static PetscErrorCode TSGLLESetType_GLLE(TS ts, TSGLLEType type) { 637 PetscBool same; 638 TS_GLLE *gl = (TS_GLLE *)ts->data; 639 PetscErrorCode (*r)(TS); 640 641 PetscFunctionBegin; 642 if (gl->type_name[0]) { 643 PetscCall(PetscStrcmp(gl->type_name, type, &same)); 644 if (same) PetscFunctionReturn(0); 645 PetscCall((*gl->Destroy)(gl)); 646 } 647 648 PetscCall(PetscFunctionListFind(TSGLLEList, type, &r)); 649 PetscCheck(r, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown TSGLLE type \"%s\" given", type); 650 PetscCall((*r)(ts)); 651 PetscCall(PetscStrcpy(gl->type_name, type)); 652 PetscFunctionReturn(0); 653 } 654 655 static PetscErrorCode TSGLLESetAcceptType_GLLE(TS ts, TSGLLEAcceptType type) { 656 TSGLLEAcceptFunction r; 657 TS_GLLE *gl = (TS_GLLE *)ts->data; 658 659 PetscFunctionBegin; 660 PetscCall(PetscFunctionListFind(TSGLLEAcceptList, type, &r)); 661 PetscCheck(r, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown TSGLLEAccept type \"%s\" given", type); 662 gl->Accept = r; 663 PetscCall(PetscStrncpy(gl->accept_name, type, sizeof(gl->accept_name))); 664 PetscFunctionReturn(0); 665 } 666 667 static PetscErrorCode TSGLLEGetAdapt_GLLE(TS ts, TSGLLEAdapt *adapt) { 668 TS_GLLE *gl = (TS_GLLE *)ts->data; 669 670 PetscFunctionBegin; 671 if (!gl->adapt) { 672 PetscCall(TSGLLEAdaptCreate(PetscObjectComm((PetscObject)ts), &gl->adapt)); 673 PetscCall(PetscObjectIncrementTabLevel((PetscObject)gl->adapt, (PetscObject)ts, 1)); 674 PetscCall(PetscLogObjectParent((PetscObject)ts, (PetscObject)gl->adapt)); 675 } 676 *adapt = gl->adapt; 677 PetscFunctionReturn(0); 678 } 679 680 static PetscErrorCode TSGLLEChooseNextScheme(TS ts, PetscReal h, const PetscReal hmnorm[], PetscInt *next_scheme, PetscReal *next_h, PetscBool *finish) { 681 TS_GLLE *gl = (TS_GLLE *)ts->data; 682 PetscInt i, n, cur_p, cur, next_sc, candidates[64], orders[64]; 683 PetscReal errors[64], costs[64], tleft; 684 685 PetscFunctionBegin; 686 cur = -1; 687 cur_p = gl->schemes[gl->current_scheme]->p; 688 tleft = ts->max_time - (ts->ptime + ts->time_step); 689 for (i = 0, n = 0; i < gl->nschemes; i++) { 690 TSGLLEScheme sc = gl->schemes[i]; 691 if (sc->p < gl->min_order || gl->max_order < sc->p) continue; 692 if (sc->p == cur_p - 1) errors[n] = PetscAbsScalar(sc->alpha[0]) * hmnorm[0]; 693 else if (sc->p == cur_p) errors[n] = PetscAbsScalar(sc->alpha[0]) * hmnorm[1]; 694 else if (sc->p == cur_p + 1) errors[n] = PetscAbsScalar(sc->alpha[0]) * (hmnorm[2] + hmnorm[3]); 695 else continue; 696 candidates[n] = i; 697 orders[n] = PetscMin(sc->p, sc->q); /* order of global truncation error */ 698 costs[n] = sc->s; /* estimate the cost as the number of stages */ 699 if (i == gl->current_scheme) cur = n; 700 n++; 701 } 702 PetscCheck(cur >= 0 && gl->nschemes > cur, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Current scheme not found in scheme list"); 703 PetscCall(TSGLLEAdaptChoose(gl->adapt, n, orders, errors, costs, cur, h, tleft, &next_sc, next_h, finish)); 704 *next_scheme = candidates[next_sc]; 705 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, 706 gl->schemes[*next_scheme]->r, gl->schemes[*next_scheme]->s, (double)*next_h, PetscBools[*finish])); 707 PetscFunctionReturn(0); 708 } 709 710 static PetscErrorCode TSGLLEGetMaxSizes(TS ts, PetscInt *max_r, PetscInt *max_s) { 711 TS_GLLE *gl = (TS_GLLE *)ts->data; 712 713 PetscFunctionBegin; 714 *max_r = gl->schemes[gl->nschemes - 1]->r; 715 *max_s = gl->schemes[gl->nschemes - 1]->s; 716 PetscFunctionReturn(0); 717 } 718 719 static PetscErrorCode TSSolve_GLLE(TS ts) { 720 TS_GLLE *gl = (TS_GLLE *)ts->data; 721 PetscInt i, k, its, lits, max_r, max_s; 722 PetscBool final_step, finish; 723 SNESConvergedReason snesreason; 724 725 PetscFunctionBegin; 726 PetscCall(TSMonitor(ts, ts->steps, ts->ptime, ts->vec_sol)); 727 728 PetscCall(TSGLLEGetMaxSizes(ts, &max_r, &max_s)); 729 PetscCall(VecCopy(ts->vec_sol, gl->X[0])); 730 for (i = 1; i < max_r; i++) { PetscCall(VecZeroEntries(gl->X[i])); } 731 PetscCall(TSGLLEUpdateWRMS(ts)); 732 733 if (0) { 734 /* Find consistent initial data for DAE */ 735 gl->stage_time = ts->ptime + ts->time_step; 736 gl->scoeff = 1.; 737 gl->stage = 0; 738 739 PetscCall(VecCopy(ts->vec_sol, gl->Z)); 740 PetscCall(VecCopy(ts->vec_sol, gl->Y)); 741 PetscCall(SNESSolve(ts->snes, NULL, gl->Y)); 742 PetscCall(SNESGetIterationNumber(ts->snes, &its)); 743 PetscCall(SNESGetLinearSolveIterations(ts->snes, &lits)); 744 PetscCall(SNESGetConvergedReason(ts->snes, &snesreason)); 745 746 ts->snes_its += its; 747 ts->ksp_its += lits; 748 if (snesreason < 0 && ts->max_snes_failures > 0 && ++ts->num_snes_failures >= ts->max_snes_failures) { 749 ts->reason = TS_DIVERGED_NONLINEAR_SOLVE; 750 PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", nonlinear solve solve failures %" PetscInt_FMT " greater than current TS allowed, stopping solve\n", ts->steps, ts->num_snes_failures)); 751 PetscFunctionReturn(0); 752 } 753 } 754 755 PetscCheck(gl->current_scheme >= 0, PETSC_COMM_SELF, PETSC_ERR_ORDER, "A starting scheme has not been provided"); 756 757 for (k = 0, final_step = PETSC_FALSE, finish = PETSC_FALSE; k < ts->max_steps && !finish; k++) { 758 PetscInt j, r, s, next_scheme = 0, rejections; 759 PetscReal h, hmnorm[4], enorm[3], next_h; 760 PetscBool accept; 761 const PetscScalar *c, *a, *u; 762 Vec *X, *Ydot, Y; 763 TSGLLEScheme scheme = gl->schemes[gl->current_scheme]; 764 765 r = scheme->r; 766 s = scheme->s; 767 c = scheme->c; 768 a = scheme->a; 769 u = scheme->u; 770 h = ts->time_step; 771 X = gl->X; 772 Ydot = gl->Ydot; 773 Y = gl->Y; 774 775 if (ts->ptime > ts->max_time) break; 776 777 /* 778 We only call PreStep at the start of each STEP, not each STAGE. This is because it is 779 possible to fail (have to restart a step) after multiple stages. 780 */ 781 PetscCall(TSPreStep(ts)); 782 783 rejections = 0; 784 while (1) { 785 for (i = 0; i < s; i++) { 786 PetscScalar shift; 787 gl->scoeff = 1. / PetscRealPart(a[i * s + i]); 788 shift = gl->scoeff / ts->time_step; 789 gl->stage = i; 790 gl->stage_time = ts->ptime + PetscRealPart(c[i]) * h; 791 792 /* 793 * Stage equation: Y = h A Y' + U X 794 * We assume that A is lower-triangular so that we can solve the stages (Y,Y') sequentially 795 * Build the affine vector z_i = -[1/(h a_ii)](h sum_j a_ij y'_j + sum_j u_ij x_j) 796 * Then y'_i = z + 1/(h a_ii) y_i 797 */ 798 PetscCall(VecZeroEntries(gl->Z)); 799 for (j = 0; j < r; j++) { PetscCall(VecAXPY(gl->Z, -shift * u[i * r + j], X[j])); } 800 for (j = 0; j < i; j++) { PetscCall(VecAXPY(gl->Z, -shift * h * a[i * s + j], Ydot[j])); } 801 /* Note: Z is used within function evaluation, Ydot = Z + shift*Y */ 802 803 /* Compute an estimate of Y to start Newton iteration */ 804 if (gl->extrapolate) { 805 if (i == 0) { 806 /* Linear extrapolation on the first stage */ 807 PetscCall(VecWAXPY(Y, c[i] * h, X[1], X[0])); 808 } else { 809 /* Linear extrapolation from the last stage */ 810 PetscCall(VecAXPY(Y, (c[i] - c[i - 1]) * h, Ydot[i - 1])); 811 } 812 } else if (i == 0) { /* Directly use solution from the last step, otherwise reuse the last stage (do nothing) */ 813 PetscCall(VecCopy(X[0], Y)); 814 } 815 816 /* Solve this stage (Ydot[i] is computed during function evaluation) */ 817 PetscCall(SNESSolve(ts->snes, NULL, Y)); 818 PetscCall(SNESGetIterationNumber(ts->snes, &its)); 819 PetscCall(SNESGetLinearSolveIterations(ts->snes, &lits)); 820 PetscCall(SNESGetConvergedReason(ts->snes, &snesreason)); 821 ts->snes_its += its; 822 ts->ksp_its += lits; 823 if (snesreason < 0 && ts->max_snes_failures > 0 && ++ts->num_snes_failures >= ts->max_snes_failures) { 824 ts->reason = TS_DIVERGED_NONLINEAR_SOLVE; 825 PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", nonlinear solve solve failures %" PetscInt_FMT " greater than current TS allowed, stopping solve\n", ts->steps, ts->num_snes_failures)); 826 PetscFunctionReturn(0); 827 } 828 } 829 830 gl->stage_time = ts->ptime + ts->time_step; 831 832 PetscCall((*gl->EstimateHigherMoments)(scheme, h, Ydot, gl->X, gl->himom)); 833 /* hmnorm[i] = h^{p+i}x^{(p+i)} with i=0,1,2; hmnorm[3] = h^{p+2}(dx'/dx) x^{(p+1)} */ 834 for (i = 0; i < 3; i++) { PetscCall(TSGLLEVecNormWRMS(ts, gl->himom[i], &hmnorm[i + 1])); } 835 enorm[0] = PetscRealPart(scheme->alpha[0]) * hmnorm[1]; 836 enorm[1] = PetscRealPart(scheme->beta[0]) * hmnorm[2]; 837 enorm[2] = PetscRealPart(scheme->gamma[0]) * hmnorm[3]; 838 PetscCall((*gl->Accept)(ts, ts->max_time - gl->stage_time, h, enorm, &accept)); 839 if (accept) goto accepted; 840 rejections++; 841 PetscCall(PetscInfo(ts, "Step %" PetscInt_FMT " (t=%g) not accepted, rejections=%" PetscInt_FMT "\n", k, (double)gl->stage_time, rejections)); 842 if (rejections > gl->max_step_rejections) break; 843 /* 844 There are lots of reasons why a step might be rejected, including solvers not converging and other factors that 845 TSGLLEChooseNextScheme does not support. Additionally, the error estimates may be very screwed up, so I'm not 846 convinced that it's safe to just compute a new error estimate using the same interface as the current adaptor 847 (the adaptor interface probably has to change). Here we make an arbitrary and naive choice. This assumes that 848 steps were written in Nordsieck form. The "correct" method would be to re-complete the previous time step with 849 the correct "next" step size. It is unclear to me whether the present ad-hoc method of rescaling X is stable. 850 */ 851 h *= 0.5; 852 for (i = 1; i < scheme->r; i++) { PetscCall(VecScale(X[i], PetscPowRealInt(0.5, i))); } 853 } 854 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); 855 856 accepted: 857 /* This term is not error, but it *would* be the leading term for a lower order method */ 858 PetscCall(TSGLLEVecNormWRMS(ts, gl->X[scheme->r - 1], &hmnorm[0])); 859 /* Correct scaling so that these are equivalent to norms of the Nordsieck vectors */ 860 861 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])); 862 if (!final_step) { 863 PetscCall(TSGLLEChooseNextScheme(ts, h, hmnorm, &next_scheme, &next_h, &final_step)); 864 } else { 865 /* Dummy values to complete the current step in a consistent manner */ 866 next_scheme = gl->current_scheme; 867 next_h = h; 868 finish = PETSC_TRUE; 869 } 870 871 X = gl->Xold; 872 gl->Xold = gl->X; 873 gl->X = X; 874 PetscCall((*gl->CompleteStep)(scheme, h, gl->schemes[next_scheme], next_h, Ydot, gl->Xold, gl->X)); 875 876 PetscCall(TSGLLEUpdateWRMS(ts)); 877 878 /* Post the solution for the user, we could avoid this copy with a small bit of cleverness */ 879 PetscCall(VecCopy(gl->X[0], ts->vec_sol)); 880 ts->ptime += h; 881 ts->steps++; 882 883 PetscCall(TSPostEvaluate(ts)); 884 PetscCall(TSPostStep(ts)); 885 PetscCall(TSMonitor(ts, ts->steps, ts->ptime, ts->vec_sol)); 886 887 gl->current_scheme = next_scheme; 888 ts->time_step = next_h; 889 } 890 PetscFunctionReturn(0); 891 } 892 893 /*------------------------------------------------------------*/ 894 895 static PetscErrorCode TSReset_GLLE(TS ts) { 896 TS_GLLE *gl = (TS_GLLE *)ts->data; 897 PetscInt max_r, max_s; 898 899 PetscFunctionBegin; 900 if (gl->setupcalled) { 901 PetscCall(TSGLLEGetMaxSizes(ts, &max_r, &max_s)); 902 PetscCall(VecDestroyVecs(max_r, &gl->Xold)); 903 PetscCall(VecDestroyVecs(max_r, &gl->X)); 904 PetscCall(VecDestroyVecs(max_s, &gl->Ydot)); 905 PetscCall(VecDestroyVecs(3, &gl->himom)); 906 PetscCall(VecDestroy(&gl->W)); 907 PetscCall(VecDestroy(&gl->Y)); 908 PetscCall(VecDestroy(&gl->Z)); 909 } 910 gl->setupcalled = PETSC_FALSE; 911 PetscFunctionReturn(0); 912 } 913 914 static PetscErrorCode TSDestroy_GLLE(TS ts) { 915 TS_GLLE *gl = (TS_GLLE *)ts->data; 916 917 PetscFunctionBegin; 918 PetscCall(TSReset_GLLE(ts)); 919 if (ts->dm) { 920 PetscCall(DMCoarsenHookRemove(ts->dm, DMCoarsenHook_TSGLLE, DMRestrictHook_TSGLLE, ts)); 921 PetscCall(DMSubDomainHookRemove(ts->dm, DMSubDomainHook_TSGLLE, DMSubDomainRestrictHook_TSGLLE, ts)); 922 } 923 if (gl->adapt) PetscCall(TSGLLEAdaptDestroy(&gl->adapt)); 924 if (gl->Destroy) PetscCall((*gl->Destroy)(gl)); 925 PetscCall(PetscFree(ts->data)); 926 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSGLLESetType_C", NULL)); 927 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSGLLESetAcceptType_C", NULL)); 928 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSGLLEGetAdapt_C", NULL)); 929 PetscFunctionReturn(0); 930 } 931 932 /* 933 This defines the nonlinear equation that is to be solved with SNES 934 g(x) = f(t,x,z+shift*x) = 0 935 */ 936 static PetscErrorCode SNESTSFormFunction_GLLE(SNES snes, Vec x, Vec f, TS ts) { 937 TS_GLLE *gl = (TS_GLLE *)ts->data; 938 Vec Z, Ydot; 939 DM dm, dmsave; 940 941 PetscFunctionBegin; 942 PetscCall(SNESGetDM(snes, &dm)); 943 PetscCall(TSGLLEGetVecs(ts, dm, &Z, &Ydot)); 944 PetscCall(VecWAXPY(Ydot, gl->scoeff / ts->time_step, x, Z)); 945 dmsave = ts->dm; 946 ts->dm = dm; 947 PetscCall(TSComputeIFunction(ts, gl->stage_time, x, Ydot, f, PETSC_FALSE)); 948 ts->dm = dmsave; 949 PetscCall(TSGLLERestoreVecs(ts, dm, &Z, &Ydot)); 950 PetscFunctionReturn(0); 951 } 952 953 static PetscErrorCode SNESTSFormJacobian_GLLE(SNES snes, Vec x, Mat A, Mat B, TS ts) { 954 TS_GLLE *gl = (TS_GLLE *)ts->data; 955 Vec Z, Ydot; 956 DM dm, dmsave; 957 958 PetscFunctionBegin; 959 PetscCall(SNESGetDM(snes, &dm)); 960 PetscCall(TSGLLEGetVecs(ts, dm, &Z, &Ydot)); 961 dmsave = ts->dm; 962 ts->dm = dm; 963 /* gl->Xdot will have already been computed in SNESTSFormFunction_GLLE */ 964 PetscCall(TSComputeIJacobian(ts, gl->stage_time, x, gl->Ydot[gl->stage], gl->scoeff / ts->time_step, A, B, PETSC_FALSE)); 965 ts->dm = dmsave; 966 PetscCall(TSGLLERestoreVecs(ts, dm, &Z, &Ydot)); 967 PetscFunctionReturn(0); 968 } 969 970 static PetscErrorCode TSSetUp_GLLE(TS ts) { 971 TS_GLLE *gl = (TS_GLLE *)ts->data; 972 PetscInt max_r, max_s; 973 DM dm; 974 975 PetscFunctionBegin; 976 gl->setupcalled = PETSC_TRUE; 977 PetscCall(TSGLLEGetMaxSizes(ts, &max_r, &max_s)); 978 PetscCall(VecDuplicateVecs(ts->vec_sol, max_r, &gl->X)); 979 PetscCall(VecDuplicateVecs(ts->vec_sol, max_r, &gl->Xold)); 980 PetscCall(VecDuplicateVecs(ts->vec_sol, max_s, &gl->Ydot)); 981 PetscCall(VecDuplicateVecs(ts->vec_sol, 3, &gl->himom)); 982 PetscCall(VecDuplicate(ts->vec_sol, &gl->W)); 983 PetscCall(VecDuplicate(ts->vec_sol, &gl->Y)); 984 PetscCall(VecDuplicate(ts->vec_sol, &gl->Z)); 985 986 /* Default acceptance tests and adaptivity */ 987 if (!gl->Accept) PetscCall(TSGLLESetAcceptType(ts, TSGLLEACCEPT_ALWAYS)); 988 if (!gl->adapt) PetscCall(TSGLLEGetAdapt(ts, &gl->adapt)); 989 990 if (gl->current_scheme < 0) { 991 PetscInt i; 992 for (i = 0;; i++) { 993 if (gl->schemes[i]->p == gl->start_order) break; 994 PetscCheck(i + 1 != gl->nschemes, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No schemes available with requested start order %" PetscInt_FMT, i); 995 } 996 gl->current_scheme = i; 997 } 998 PetscCall(TSGetDM(ts, &dm)); 999 PetscCall(DMCoarsenHookAdd(dm, DMCoarsenHook_TSGLLE, DMRestrictHook_TSGLLE, ts)); 1000 PetscCall(DMSubDomainHookAdd(dm, DMSubDomainHook_TSGLLE, DMSubDomainRestrictHook_TSGLLE, ts)); 1001 PetscFunctionReturn(0); 1002 } 1003 /*------------------------------------------------------------*/ 1004 1005 static PetscErrorCode TSSetFromOptions_GLLE(TS ts, PetscOptionItems *PetscOptionsObject) { 1006 TS_GLLE *gl = (TS_GLLE *)ts->data; 1007 char tname[256] = TSGLLE_IRKS, completef[256] = "rescale-and-modify"; 1008 1009 PetscFunctionBegin; 1010 PetscOptionsHeadBegin(PetscOptionsObject, "General Linear ODE solver options"); 1011 { 1012 PetscBool flg; 1013 PetscCall(PetscOptionsFList("-ts_gl_type", "Type of GL method", "TSGLLESetType", TSGLLEList, gl->type_name[0] ? gl->type_name : tname, tname, sizeof(tname), &flg)); 1014 if (flg || !gl->type_name[0]) { PetscCall(TSGLLESetType(ts, tname)); } 1015 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)); 1016 PetscCall(PetscOptionsInt("-ts_gl_max_order", "Maximum order to try", "TSGLLESetMaxOrder", gl->max_order, &gl->max_order, NULL)); 1017 PetscCall(PetscOptionsInt("-ts_gl_min_order", "Minimum order to try", "TSGLLESetMinOrder", gl->min_order, &gl->min_order, NULL)); 1018 PetscCall(PetscOptionsInt("-ts_gl_start_order", "Initial order to try", "TSGLLESetMinOrder", gl->start_order, &gl->start_order, NULL)); 1019 PetscCall(PetscOptionsEnum("-ts_gl_error_direction", "Which direction to look when estimating error", "TSGLLESetErrorDirection", TSGLLEErrorDirections, (PetscEnum)gl->error_direction, (PetscEnum *)&gl->error_direction, NULL)); 1020 PetscCall(PetscOptionsBool("-ts_gl_extrapolate", "Extrapolate stage solution from previous solution (sometimes unstable)", "TSGLLESetExtrapolate", gl->extrapolate, &gl->extrapolate, NULL)); 1021 PetscCall(PetscOptionsReal("-ts_gl_atol", "Absolute tolerance", "TSGLLESetTolerances", gl->wrms_atol, &gl->wrms_atol, NULL)); 1022 PetscCall(PetscOptionsReal("-ts_gl_rtol", "Relative tolerance", "TSGLLESetTolerances", gl->wrms_rtol, &gl->wrms_rtol, NULL)); 1023 PetscCall(PetscOptionsString("-ts_gl_complete", "Method to use for completing the step", "none", completef, completef, sizeof(completef), &flg)); 1024 if (flg) { 1025 PetscBool match1, match2; 1026 PetscCall(PetscStrcmp(completef, "rescale", &match1)); 1027 PetscCall(PetscStrcmp(completef, "rescale-and-modify", &match2)); 1028 if (match1) gl->CompleteStep = TSGLLECompleteStep_Rescale; 1029 else if (match2) gl->CompleteStep = TSGLLECompleteStep_RescaleAndModify; 1030 else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "%s", completef); 1031 } 1032 { 1033 char type[256] = TSGLLEACCEPT_ALWAYS; 1034 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)); 1035 if (flg || !gl->accept_name[0]) { PetscCall(TSGLLESetAcceptType(ts, type)); } 1036 } 1037 { 1038 TSGLLEAdapt adapt; 1039 PetscCall(TSGLLEGetAdapt(ts, &adapt)); 1040 PetscCall(TSGLLEAdaptSetFromOptions(adapt, PetscOptionsObject)); 1041 } 1042 } 1043 PetscOptionsHeadEnd(); 1044 PetscFunctionReturn(0); 1045 } 1046 1047 static PetscErrorCode TSView_GLLE(TS ts, PetscViewer viewer) { 1048 TS_GLLE *gl = (TS_GLLE *)ts->data; 1049 PetscInt i; 1050 PetscBool iascii, details; 1051 1052 PetscFunctionBegin; 1053 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 1054 if (iascii) { 1055 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)); 1056 PetscCall(PetscViewerASCIIPrintf(viewer, " Error estimation: %s\n", TSGLLEErrorDirections[gl->error_direction])); 1057 PetscCall(PetscViewerASCIIPrintf(viewer, " Extrapolation: %s\n", gl->extrapolate ? "yes" : "no")); 1058 PetscCall(PetscViewerASCIIPrintf(viewer, " Acceptance test: %s\n", gl->accept_name[0] ? gl->accept_name : "(not yet set)")); 1059 PetscCall(PetscViewerASCIIPushTab(viewer)); 1060 PetscCall(TSGLLEAdaptView(gl->adapt, viewer)); 1061 PetscCall(PetscViewerASCIIPopTab(viewer)); 1062 PetscCall(PetscViewerASCIIPrintf(viewer, " type: %s\n", gl->type_name[0] ? gl->type_name : "(not yet set)")); 1063 PetscCall(PetscViewerASCIIPrintf(viewer, "Schemes within family (%" PetscInt_FMT "):\n", gl->nschemes)); 1064 details = PETSC_FALSE; 1065 PetscCall(PetscOptionsGetBool(((PetscObject)ts)->options, ((PetscObject)ts)->prefix, "-ts_gl_view_detailed", &details, NULL)); 1066 PetscCall(PetscViewerASCIIPushTab(viewer)); 1067 for (i = 0; i < gl->nschemes; i++) { PetscCall(TSGLLESchemeView(gl->schemes[i], details, viewer)); } 1068 if (gl->View) PetscCall((*gl->View)(gl, viewer)); 1069 PetscCall(PetscViewerASCIIPopTab(viewer)); 1070 } 1071 PetscFunctionReturn(0); 1072 } 1073 1074 /*@C 1075 TSGLLERegister - adds a TSGLLE implementation 1076 1077 Not Collective 1078 1079 Input Parameters: 1080 + name_scheme - name of user-defined general linear scheme 1081 - routine_create - routine to create method context 1082 1083 Notes: 1084 TSGLLERegister() may be called multiple times to add several user-defined families. 1085 1086 Sample usage: 1087 .vb 1088 TSGLLERegister("my_scheme",MySchemeCreate); 1089 .ve 1090 1091 Then, your scheme can be chosen with the procedural interface via 1092 $ TSGLLESetType(ts,"my_scheme") 1093 or at runtime via the option 1094 $ -ts_gl_type my_scheme 1095 1096 Level: advanced 1097 1098 .seealso: `TSGLLERegisterAll()` 1099 @*/ 1100 PetscErrorCode TSGLLERegister(const char sname[], PetscErrorCode (*function)(TS)) { 1101 PetscFunctionBegin; 1102 PetscCall(TSGLLEInitializePackage()); 1103 PetscCall(PetscFunctionListAdd(&TSGLLEList, sname, function)); 1104 PetscFunctionReturn(0); 1105 } 1106 1107 /*@C 1108 TSGLLEAcceptRegister - adds a TSGLLE acceptance scheme 1109 1110 Not Collective 1111 1112 Input Parameters: 1113 + name_scheme - name of user-defined acceptance scheme 1114 - routine_create - routine to create method context 1115 1116 Notes: 1117 TSGLLEAcceptRegister() may be called multiple times to add several user-defined families. 1118 1119 Sample usage: 1120 .vb 1121 TSGLLEAcceptRegister("my_scheme",MySchemeCreate); 1122 .ve 1123 1124 Then, your scheme can be chosen with the procedural interface via 1125 $ TSGLLESetAcceptType(ts,"my_scheme") 1126 or at runtime via the option 1127 $ -ts_gl_accept_type my_scheme 1128 1129 Level: advanced 1130 1131 .seealso: `TSGLLERegisterAll()` 1132 @*/ 1133 PetscErrorCode TSGLLEAcceptRegister(const char sname[], TSGLLEAcceptFunction function) { 1134 PetscFunctionBegin; 1135 PetscCall(PetscFunctionListAdd(&TSGLLEAcceptList, sname, function)); 1136 PetscFunctionReturn(0); 1137 } 1138 1139 /*@C 1140 TSGLLERegisterAll - Registers all of the general linear methods in TSGLLE 1141 1142 Not Collective 1143 1144 Level: advanced 1145 1146 .seealso: `TSGLLERegisterDestroy()` 1147 @*/ 1148 PetscErrorCode TSGLLERegisterAll(void) { 1149 PetscFunctionBegin; 1150 if (TSGLLERegisterAllCalled) PetscFunctionReturn(0); 1151 TSGLLERegisterAllCalled = PETSC_TRUE; 1152 1153 PetscCall(TSGLLERegister(TSGLLE_IRKS, TSGLLECreate_IRKS)); 1154 PetscCall(TSGLLEAcceptRegister(TSGLLEACCEPT_ALWAYS, TSGLLEAccept_Always)); 1155 PetscFunctionReturn(0); 1156 } 1157 1158 /*@C 1159 TSGLLEInitializePackage - This function initializes everything in the TSGLLE package. It is called 1160 from TSInitializePackage(). 1161 1162 Level: developer 1163 1164 .seealso: `PetscInitialize()` 1165 @*/ 1166 PetscErrorCode TSGLLEInitializePackage(void) { 1167 PetscFunctionBegin; 1168 if (TSGLLEPackageInitialized) PetscFunctionReturn(0); 1169 TSGLLEPackageInitialized = PETSC_TRUE; 1170 PetscCall(TSGLLERegisterAll()); 1171 PetscCall(PetscRegisterFinalize(TSGLLEFinalizePackage)); 1172 PetscFunctionReturn(0); 1173 } 1174 1175 /*@C 1176 TSGLLEFinalizePackage - This function destroys everything in the TSGLLE package. It is 1177 called from PetscFinalize(). 1178 1179 Level: developer 1180 1181 .seealso: `PetscFinalize()` 1182 @*/ 1183 PetscErrorCode TSGLLEFinalizePackage(void) { 1184 PetscFunctionBegin; 1185 PetscCall(PetscFunctionListDestroy(&TSGLLEList)); 1186 PetscCall(PetscFunctionListDestroy(&TSGLLEAcceptList)); 1187 TSGLLEPackageInitialized = PETSC_FALSE; 1188 TSGLLERegisterAllCalled = PETSC_FALSE; 1189 PetscFunctionReturn(0); 1190 } 1191 1192 /* ------------------------------------------------------------ */ 1193 /*MC 1194 TSGLLE - DAE solver using implicit General Linear methods 1195 1196 These methods contain Runge-Kutta and multistep schemes as special cases. These special cases have some fundamental 1197 limitations. For example, diagonally implicit Runge-Kutta cannot have stage order greater than 1 which limits their 1198 applicability to very stiff systems. Meanwhile, multistep methods cannot be A-stable for order greater than 2 and BDF 1199 are not 0-stable for order greater than 6. GL methods can be A- and L-stable with arbitrarily high stage order and 1200 reliable error estimates for both 1 and 2 orders higher to facilitate adaptive step sizes and adaptive order schemes. 1201 All this is possible while preserving a singly diagonally implicit structure. 1202 1203 Options database keys: 1204 + -ts_gl_type <type> - the class of general linear method (irks) 1205 . -ts_gl_rtol <tol> - relative error 1206 . -ts_gl_atol <tol> - absolute error 1207 . -ts_gl_min_order <p> - minimum order method to consider (default=1) 1208 . -ts_gl_max_order <p> - maximum order method to consider (default=3) 1209 . -ts_gl_start_order <p> - order of starting method (default=1) 1210 . -ts_gl_complete <method> - method to use for completing the step (rescale-and-modify or rescale) 1211 - -ts_adapt_type <method> - adaptive controller to use (none step both) 1212 1213 Notes: 1214 This integrator can be applied to DAE. 1215 1216 Diagonally implicit general linear (DIGL) methods are a generalization of diagonally implicit Runge-Kutta (DIRK). 1217 They are represented by the tableau 1218 1219 .vb 1220 A | U 1221 ------- 1222 B | V 1223 .ve 1224 1225 combined with a vector c of abscissa. "Diagonally implicit" means that A is lower triangular. 1226 A step of the general method reads 1227 1228 .vb 1229 [ Y ] = [A U] [ Y' ] 1230 [X^k] = [B V] [X^{k-1}] 1231 .ve 1232 1233 where Y is the multivector of stage values, Y' is the multivector of stage derivatives, X^k is the Nordsieck vector of 1234 the solution at step k. The Nordsieck vector consists of the first r moments of the solution, given by 1235 1236 .vb 1237 X = [x_0,x_1,...,x_{r-1}] = [x, h x', h^2 x'', ..., h^{r-1} x^{(r-1)} ] 1238 .ve 1239 1240 If A is lower triangular, we can solve the stages (Y,Y') sequentially 1241 1242 .vb 1243 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} 1244 .ve 1245 1246 and then construct the pieces to carry to the next step 1247 1248 .vb 1249 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} 1250 .ve 1251 1252 Note that when the equations are cast in implicit form, we are using the stage equation to define y'_i 1253 in terms of y_i and known stuff (y_j for j<i and x_j for all j). 1254 1255 Error estimation 1256 1257 At present, the most attractive GL methods for stiff problems are singly diagonally implicit schemes which posses 1258 Inherent Runge-Kutta Stability (IRKS). These methods have r=s, the number of items passed between steps is equal to 1259 the number of stages. The order and stage-order are one less than the number of stages. We use the error estimates 1260 in the 2007 paper which provide the following estimates 1261 1262 .vb 1263 h^{p+1} X^{(p+1)} = phi_0^T Y' + [0 psi_0^T] Xold 1264 h^{p+2} X^{(p+2)} = phi_1^T Y' + [0 psi_1^T] Xold 1265 h^{p+2} (dx'/dx) X^{(p+1)} = phi_2^T Y' + [0 psi_2^T] Xold 1266 .ve 1267 1268 These estimates are accurate to O(h^{p+3}). 1269 1270 Changing the step size 1271 1272 We use the generalized "rescale and modify" scheme, see equation (4.5) of the 2007 paper. 1273 1274 Level: beginner 1275 1276 References: 1277 + * - John Butcher and Z. Jackieweicz and W. Wright, On error propagation in general linear methods for 1278 ordinary differential equations, Journal of Complexity, Vol 23, 2007. 1279 - * - John Butcher, Numerical methods for ordinary differential equations, second edition, Wiley, 2009. 1280 1281 .seealso: `TSCreate()`, `TS`, `TSSetType()` 1282 1283 M*/ 1284 PETSC_EXTERN PetscErrorCode TSCreate_GLLE(TS ts) { 1285 TS_GLLE *gl; 1286 1287 PetscFunctionBegin; 1288 PetscCall(TSGLLEInitializePackage()); 1289 1290 PetscCall(PetscNewLog(ts, &gl)); 1291 ts->data = (void *)gl; 1292 1293 ts->ops->reset = TSReset_GLLE; 1294 ts->ops->destroy = TSDestroy_GLLE; 1295 ts->ops->view = TSView_GLLE; 1296 ts->ops->setup = TSSetUp_GLLE; 1297 ts->ops->solve = TSSolve_GLLE; 1298 ts->ops->setfromoptions = TSSetFromOptions_GLLE; 1299 ts->ops->snesfunction = SNESTSFormFunction_GLLE; 1300 ts->ops->snesjacobian = SNESTSFormJacobian_GLLE; 1301 1302 ts->usessnes = PETSC_TRUE; 1303 1304 gl->max_step_rejections = 1; 1305 gl->min_order = 1; 1306 gl->max_order = 3; 1307 gl->start_order = 1; 1308 gl->current_scheme = -1; 1309 gl->extrapolate = PETSC_FALSE; 1310 1311 gl->wrms_atol = 1e-8; 1312 gl->wrms_rtol = 1e-5; 1313 1314 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSGLLESetType_C", &TSGLLESetType_GLLE)); 1315 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSGLLESetAcceptType_C", &TSGLLESetAcceptType_GLLE)); 1316 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSGLLEGetAdapt_C", &TSGLLEGetAdapt_GLLE)); 1317 PetscFunctionReturn(0); 1318 } 1319