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] = {1. / 3., 2. / 3., 1}; 456 const PetscScalar a[3][3] = { 457 {4. / 9., 0, 0 }, 458 {1.03750643704090e+00, 4. / 9., 0 }, 459 {7.67024779410304e-01, -3.81140216918943e-01, 4. / 9.} 460 }; 461 const PetscScalar b[3][3] = { 462 {0.767024779410304, -0.381140216918943, 4. / 9. }, 463 {0.000000000000000, 0.000000000000000, 1.000000000000000}, 464 {-2.075048385225385, 0.621728385225383, 1.277197204924873} 465 }; 466 const PetscScalar u[3][3] = { 467 {1.0000000000000000, -0.1111111111111109, -0.0925925925925922}, 468 {1.0000000000000000, -0.8152842148186744, -0.4199095530877056}, 469 {1.0000000000000000, 0.1696709930641948, 0.0539741070314165 } 470 }; 471 const PetscScalar v[3][3] = { 472 {1.0000000000000000, 0.1696709930641948, 0.0539741070314165}, 473 {0.000000000000000, 0.000000000000000, 0.000000000000000 }, 474 {0.000000000000000, 0.176122795075129, 0.000000000000000 } 475 }; 476 PetscCall(TSGLLESchemeCreate(2, 2, 3, 3, c, *a, *b, *u, *v, &gl->schemes[gl->nschemes++])); 477 } 478 { 479 /* p=q=3, r=s=4: irks(9/40,0,[1:4]/4,[0.3312 1.0050],[0.49541 1;1 0]) */ 480 const PetscScalar c[4] = {0.25, 0.5, 0.75, 1.0}; 481 const PetscScalar a[4][4] = { 482 {9. / 40., 0, 0, 0 }, 483 {2.11286958887701e-01, 9. / 40., 0, 0 }, 484 {9.46338294287584e-01, -3.42942861246094e-01, 9. / 40., 0 }, 485 {0.521490453970721, -0.662474225622980, 0.490476425459734, 9. / 40.} 486 }; 487 const PetscScalar b[4][4] = { 488 {0.521490453970721, -0.662474225622980, 0.490476425459734, 9. / 40. }, 489 {0.000000000000000, 0.000000000000000, 0.000000000000000, 1.000000000000000}, 490 {-0.084677029310348, 1.390757514776085, -1.568157386206001, 2.023192696767826}, 491 {0.465383797936408, 1.478273530625148, -1.930836081010182, 1.644872111193354} 492 }; 493 const PetscScalar u[4][4] = { 494 {1.00000000000000000, 0.02500000000001035, -0.02499999999999053, -0.00442708333332865}, 495 {1.00000000000000000, 0.06371304111232945, -0.04032173972189845, -0.01389438413189452}, 496 {1.00000000000000000, -0.07839543304147778, 0.04738685705116663, 0.02032603595928376 }, 497 {1.00000000000000000, 0.42550734619251651, 0.10800718022400080, -0.01726712647760034} 498 }; 499 const PetscScalar v[4][4] = { 500 {1.00000000000000000, 0.42550734619251651, 0.10800718022400080, -0.01726712647760034}, 501 {0.000000000000000, 0.000000000000000, 0.000000000000000, 0.000000000000000 }, 502 {0.000000000000000, -1.761115796027561, -0.521284157173780, 0.258249384305463 }, 503 {0.000000000000000, -1.657693358744728, -1.052227765232394, 0.521284157173780 } 504 }; 505 PetscCall(TSGLLESchemeCreate(3, 3, 4, 4, c, *a, *b, *u, *v, &gl->schemes[gl->nschemes++])); 506 } 507 { 508 /* p=q=4, r=s=5: 509 irks(3/11,0,[1:5]/5, [0.1715 -0.1238 0.6617],... 510 [ -0.0812 0.4079 1.0000 511 1.0000 0 0 512 0.8270 1.0000 0]) 513 */ 514 const PetscScalar c[5] = {0.2, 0.4, 0.6, 0.8, 1.0}; 515 const PetscScalar a[5][5] = { 516 {2.72727272727352e-01, 0.00000000000000e+00, 0.00000000000000e+00, 0.00000000000000e+00, 0.00000000000000e+00}, 517 {-1.03980153733431e-01, 2.72727272727405e-01, 0.00000000000000e+00, 0.00000000000000e+00, 0.00000000000000e+00}, 518 {-1.58615400341492e+00, 7.44168951881122e-01, 2.72727272727309e-01, 0.00000000000000e+00, 0.00000000000000e+00}, 519 {-8.73658042865628e-01, 5.37884671894595e-01, -1.63298538799523e-01, 2.72727272726996e-01, 0.00000000000000e+00}, 520 {2.95489397443992e-01, -1.18481693910097e+00, -6.68029812659953e-01, 1.00716687860943e+00, 2.72727272727288e-01} 521 }; 522 const PetscScalar b[5][5] = { 523 {2.95489397443992e-01, -1.18481693910097e+00, -6.68029812659953e-01, 1.00716687860943e+00, 2.72727272727288e-01}, 524 {0.00000000000000e+00, 1.11022302462516e-16, -2.22044604925031e-16, 0.00000000000000e+00, 1.00000000000000e+00}, 525 {-4.05882503986005e+00, -4.00924006567769e+00, -1.38930610972481e+00, 4.45223930308488e+00, 6.32331093108427e-01}, 526 {8.35690179937017e+00, -2.26640927349732e+00, 6.86647884973826e+00, -5.22595158025740e+00, 4.50893068837431e+00}, 527 {1.27656267027479e+01, 2.80882153840821e+00, 8.91173096522890e+00, -1.07936444078906e+01, 4.82534148988854e+00} 528 }; 529 const PetscScalar u[5][5] = { 530 {1.00000000000000e+00, -7.27272727273551e-02, -3.45454545454419e-02, -4.12121212119565e-03, -2.96969696964014e-04}, 531 {1.00000000000000e+00, 2.31252881006154e-01, -8.29487834416481e-03, -9.07191207681020e-03, -1.70378403743473e-03}, 532 {1.00000000000000e+00, 1.16925777880663e+00, 3.59268562942635e-02, -4.09013451730615e-02, -1.02411119670164e-02}, 533 {1.00000000000000e+00, 1.02634463704356e+00, 1.59375044913405e-01, 1.89673015035370e-03, -4.89987231897569e-03}, 534 {1.00000000000000e+00, 1.27746320298021e+00, 2.37186008132728e-01, -8.28694373940065e-02, -5.34396510196430e-02} 535 }; 536 const PetscScalar v[5][5] = { 537 {1.00000000000000e+00, 1.27746320298021e+00, 2.37186008132728e-01, -8.28694373940065e-02, -5.34396510196430e-02}, 538 {0.00000000000000e+00, -1.77635683940025e-15, -1.99840144432528e-15, -9.99200722162641e-16, -3.33066907387547e-16}, 539 {0.00000000000000e+00, 4.37280081906924e+00, 5.49221645016377e-02, -8.88913177394943e-02, 1.12879077989154e-01 }, 540 {0.00000000000000e+00, -1.22399504837280e+01, -5.21287338448645e+00, -8.03952325565291e-01, 4.60298678047147e-01 }, 541 {0.00000000000000e+00, -1.85178762883829e+01, -5.21411849862624e+00, -1.04283436528809e+00, 7.49030161063651e-01 } 542 }; 543 PetscCall(TSGLLESchemeCreate(4, 4, 5, 5, c, *a, *b, *u, *v, &gl->schemes[gl->nschemes++])); 544 } 545 { 546 /* p=q=5, r=s=6; 547 irks(1/3,0,[1:6]/6,... 548 [-0.0489 0.4228 -0.8814 0.9021],... 549 [-0.3474 -0.6617 0.6294 0.2129 550 0.0044 -0.4256 -0.1427 -0.8936 551 -0.8267 0.4821 0.1371 -0.2557 552 -0.4426 -0.3855 -0.7514 0.3014]) 553 */ 554 const PetscScalar c[6] = {1. / 6, 2. / 6, 3. / 6, 4. / 6, 5. / 6, 1.}; 555 const PetscScalar a[6][6] = { 556 {3.33333333333940e-01, 0, 0, 0, 0, 0 }, 557 {-8.64423857333350e-02, 3.33333333332888e-01, 0, 0, 0, 0 }, 558 {-2.16850174258252e+00, -2.23619072028839e+00, 3.33333333335204e-01, 0, 0, 0 }, 559 {-4.73160970138997e+00, -3.89265344629268e+00, -2.76318716520933e-01, 3.33333333335759e-01, 0, 0 }, 560 {-6.75187540297338e+00, -7.90756533769377e+00, 7.90245051802259e-01, -4.48352364517632e-01, 3.33333333328483e-01, 0 }, 561 {-4.26488287921548e+00, -1.19320395589302e+01, 3.38924509887755e+00, -2.23969848002481e+00, 6.62807710124007e-01, 3.33333333335440e-01} 562 }; 563 const PetscScalar b[6][6] = { 564 {-4.26488287921548e+00, -1.19320395589302e+01, 3.38924509887755e+00, -2.23969848002481e+00, 6.62807710124007e-01, 3.33333333335440e-01 }, 565 {-8.88178419700125e-16, 4.44089209850063e-16, -1.54737334057131e-15, -8.88178419700125e-16, 0.00000000000000e+00, 1.00000000000001e+00 }, 566 {-2.87780425770651e+01, -1.13520448264971e+01, 2.62002318943161e+01, 2.56943874812797e+01, -3.06702268304488e+01, 6.68067773510103e+00 }, 567 {5.47971245256474e+01, 6.80366875868284e+01, -6.50952588861999e+01, -8.28643975339097e+01, 8.17416943896414e+01, -1.17819043489036e+01}, 568 {-2.33332114788869e+02, 6.12942539462634e+01, -4.91850135865944e+01, 1.82716844135480e+02, -1.29788173979395e+02, 3.09968095651099e+01 }, 569 {-1.72049132343751e+02, 8.60194713593999e+00, 7.98154219170200e-01, 1.50371386053218e+02, -1.18515423962066e+02, 2.50898277784663e+01 } 570 }; 571 const PetscScalar u[6][6] = { 572 {1.00000000000000e+00, -1.66666666666870e-01, -4.16666666664335e-02, -3.85802469124815e-03, -2.25051440302250e-04, -9.64506172339142e-06}, 573 {1.00000000000000e+00, 8.64423857327162e-02, -4.11484912671353e-02, -1.11450903217645e-02, -1.47651050487126e-03, -1.34395070766826e-04}, 574 {1.00000000000000e+00, 4.57135912953434e+00, 1.06514719719137e+00, 1.33517564218007e-01, 1.11365952968659e-02, 6.12382756769504e-04 }, 575 {1.00000000000000e+00, 9.23391519753404e+00, 2.22431212392095e+00, 2.91823807741891e-01, 2.52058456411084e-02, 1.22800542949647e-03 }, 576 {1.00000000000000e+00, 1.48175480533865e+01, 3.73439117461835e+00, 5.14648336541804e-01, 4.76430038853402e-02, 2.56798515502156e-03 }, 577 {1.00000000000000e+00, 1.50512347758335e+01, 4.10099701165164e+00, 5.66039141003603e-01, 3.91213893800891e-02, -2.99136269067853e-03} 578 }; 579 const PetscScalar v[6][6] = { 580 {1.00000000000000e+00, 1.50512347758335e+01, 4.10099701165164e+00, 5.66039141003603e-01, 3.91213893800891e-02, -2.99136269067853e-03}, 581 {0.00000000000000e+00, -4.88498130835069e-15, -6.43929354282591e-15, -3.55271367880050e-15, -1.22124532708767e-15, -3.12250225675825e-16}, 582 {0.00000000000000e+00, 1.22250171233141e+01, -1.77150760606169e+00, 3.54516769879390e-01, 6.22298845883398e-01, 2.31647447450276e-01 }, 583 {0.00000000000000e+00, -4.48339457331040e+01, -3.57363126641880e-01, 5.18750173123425e-01, 6.55727990241799e-02, 1.63175368287079e-01 }, 584 {0.00000000000000e+00, 1.37297394708005e+02, -1.60145272991317e+00, -5.05319555199441e+00, 1.55328940390990e-01, 9.16629423682464e-01 }, 585 {0.00000000000000e+00, 1.05703241119022e+02, -1.16610260983038e+00, -2.99767252773859e+00, -1.13472315553890e-01, 1.09742849254729e+00 } 586 }; 587 PetscCall(TSGLLESchemeCreate(5, 5, 6, 6, c, *a, *b, *u, *v, &gl->schemes[gl->nschemes++])); 588 } 589 PetscFunctionReturn(0); 590 } 591 592 /*@C 593 TSGLLESetType - sets the class of general linear method to use for time-stepping 594 595 Collective on TS 596 597 Input Parameters: 598 + ts - the TS context 599 - type - a method 600 601 Options Database Key: 602 . -ts_gl_type <type> - sets the method, use -help for a list of available method (e.g. irks) 603 604 Notes: 605 See "petsc/include/petscts.h" for available methods (for instance) 606 . TSGLLE_IRKS - Diagonally implicit methods with inherent Runge-Kutta stability (for stiff problems) 607 608 Normally, it is best to use the TSSetFromOptions() command and 609 then set the TSGLLE type from the options database rather than by using 610 this routine. Using the options database provides the user with 611 maximum flexibility in evaluating the many different solvers. 612 The TSGLLESetType() routine is provided for those situations where it 613 is necessary to set the timestepping solver independently of the 614 command line or options database. This might be the case, for example, 615 when the choice of solver changes during the execution of the 616 program, and the user's application is taking responsibility for 617 choosing the appropriate method. 618 619 Level: intermediate 620 621 @*/ 622 PetscErrorCode TSGLLESetType(TS ts, TSGLLEType type) { 623 PetscFunctionBegin; 624 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 625 PetscValidCharPointer(type, 2); 626 PetscTryMethod(ts, "TSGLLESetType_C", (TS, TSGLLEType), (ts, type)); 627 PetscFunctionReturn(0); 628 } 629 630 /*@C 631 TSGLLESetAcceptType - sets the acceptance test 632 633 Time integrators that need to control error must have the option to reject a time step based on local error 634 estimates. This function allows different schemes to be set. 635 636 Logically Collective on TS 637 638 Input Parameters: 639 + ts - the TS context 640 - type - the type 641 642 Options Database Key: 643 . -ts_gl_accept_type <type> - sets the method used to determine whether to accept or reject a step 644 645 Level: intermediate 646 647 .seealso: `TS`, `TSGLLE`, `TSGLLEAcceptRegister()`, `TSGLLEAdapt`, `set` `type` 648 @*/ 649 PetscErrorCode TSGLLESetAcceptType(TS ts, TSGLLEAcceptType type) { 650 PetscFunctionBegin; 651 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 652 PetscValidCharPointer(type, 2); 653 PetscTryMethod(ts, "TSGLLESetAcceptType_C", (TS, TSGLLEAcceptType), (ts, type)); 654 PetscFunctionReturn(0); 655 } 656 657 /*@C 658 TSGLLEGetAdapt - gets the TSGLLEAdapt object from the TS 659 660 Not Collective 661 662 Input Parameter: 663 . ts - the TS context 664 665 Output Parameter: 666 . adapt - the TSGLLEAdapt context 667 668 Notes: 669 This allows the user set options on the TSGLLEAdapt object. Usually it is better to do this using the options 670 database, so this function is rarely needed. 671 672 Level: advanced 673 674 .seealso: `TSGLLEAdapt`, `TSGLLEAdaptRegister()` 675 @*/ 676 PetscErrorCode TSGLLEGetAdapt(TS ts, TSGLLEAdapt *adapt) { 677 PetscFunctionBegin; 678 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 679 PetscValidPointer(adapt, 2); 680 PetscUseMethod(ts, "TSGLLEGetAdapt_C", (TS, TSGLLEAdapt *), (ts, adapt)); 681 PetscFunctionReturn(0); 682 } 683 684 static PetscErrorCode TSGLLEAccept_Always(TS ts, PetscReal tleft, PetscReal h, const PetscReal enorms[], PetscBool *accept) { 685 PetscFunctionBegin; 686 *accept = PETSC_TRUE; 687 PetscFunctionReturn(0); 688 } 689 690 static PetscErrorCode TSGLLEUpdateWRMS(TS ts) { 691 TS_GLLE *gl = (TS_GLLE *)ts->data; 692 PetscScalar *x, *w; 693 PetscInt n, i; 694 695 PetscFunctionBegin; 696 PetscCall(VecGetArray(gl->X[0], &x)); 697 PetscCall(VecGetArray(gl->W, &w)); 698 PetscCall(VecGetLocalSize(gl->W, &n)); 699 for (i = 0; i < n; i++) w[i] = 1. / (gl->wrms_atol + gl->wrms_rtol * PetscAbsScalar(x[i])); 700 PetscCall(VecRestoreArray(gl->X[0], &x)); 701 PetscCall(VecRestoreArray(gl->W, &w)); 702 PetscFunctionReturn(0); 703 } 704 705 static PetscErrorCode TSGLLEVecNormWRMS(TS ts, Vec X, PetscReal *nrm) { 706 TS_GLLE *gl = (TS_GLLE *)ts->data; 707 PetscScalar *x, *w; 708 PetscReal sum = 0.0, gsum; 709 PetscInt n, N, i; 710 711 PetscFunctionBegin; 712 PetscCall(VecGetArray(X, &x)); 713 PetscCall(VecGetArray(gl->W, &w)); 714 PetscCall(VecGetLocalSize(gl->W, &n)); 715 for (i = 0; i < n; i++) sum += PetscAbsScalar(PetscSqr(x[i] * w[i])); 716 PetscCall(VecRestoreArray(X, &x)); 717 PetscCall(VecRestoreArray(gl->W, &w)); 718 PetscCall(MPIU_Allreduce(&sum, &gsum, 1, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)ts))); 719 PetscCall(VecGetSize(gl->W, &N)); 720 *nrm = PetscSqrtReal(gsum / (1. * N)); 721 PetscFunctionReturn(0); 722 } 723 724 static PetscErrorCode TSGLLESetType_GLLE(TS ts, TSGLLEType type) { 725 PetscBool same; 726 TS_GLLE *gl = (TS_GLLE *)ts->data; 727 PetscErrorCode (*r)(TS); 728 729 PetscFunctionBegin; 730 if (gl->type_name[0]) { 731 PetscCall(PetscStrcmp(gl->type_name, type, &same)); 732 if (same) PetscFunctionReturn(0); 733 PetscCall((*gl->Destroy)(gl)); 734 } 735 736 PetscCall(PetscFunctionListFind(TSGLLEList, type, &r)); 737 PetscCheck(r, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown TSGLLE type \"%s\" given", type); 738 PetscCall((*r)(ts)); 739 PetscCall(PetscStrcpy(gl->type_name, type)); 740 PetscFunctionReturn(0); 741 } 742 743 static PetscErrorCode TSGLLESetAcceptType_GLLE(TS ts, TSGLLEAcceptType type) { 744 TSGLLEAcceptFunction r; 745 TS_GLLE *gl = (TS_GLLE *)ts->data; 746 747 PetscFunctionBegin; 748 PetscCall(PetscFunctionListFind(TSGLLEAcceptList, type, &r)); 749 PetscCheck(r, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown TSGLLEAccept type \"%s\" given", type); 750 gl->Accept = r; 751 PetscCall(PetscStrncpy(gl->accept_name, type, sizeof(gl->accept_name))); 752 PetscFunctionReturn(0); 753 } 754 755 static PetscErrorCode TSGLLEGetAdapt_GLLE(TS ts, TSGLLEAdapt *adapt) { 756 TS_GLLE *gl = (TS_GLLE *)ts->data; 757 758 PetscFunctionBegin; 759 if (!gl->adapt) { 760 PetscCall(TSGLLEAdaptCreate(PetscObjectComm((PetscObject)ts), &gl->adapt)); 761 PetscCall(PetscObjectIncrementTabLevel((PetscObject)gl->adapt, (PetscObject)ts, 1)); 762 } 763 *adapt = gl->adapt; 764 PetscFunctionReturn(0); 765 } 766 767 static PetscErrorCode TSGLLEChooseNextScheme(TS ts, PetscReal h, const PetscReal hmnorm[], PetscInt *next_scheme, PetscReal *next_h, PetscBool *finish) { 768 TS_GLLE *gl = (TS_GLLE *)ts->data; 769 PetscInt i, n, cur_p, cur, next_sc, candidates[64], orders[64]; 770 PetscReal errors[64], costs[64], tleft; 771 772 PetscFunctionBegin; 773 cur = -1; 774 cur_p = gl->schemes[gl->current_scheme]->p; 775 tleft = ts->max_time - (ts->ptime + ts->time_step); 776 for (i = 0, n = 0; i < gl->nschemes; i++) { 777 TSGLLEScheme sc = gl->schemes[i]; 778 if (sc->p < gl->min_order || gl->max_order < sc->p) continue; 779 if (sc->p == cur_p - 1) errors[n] = PetscAbsScalar(sc->alpha[0]) * hmnorm[0]; 780 else if (sc->p == cur_p) errors[n] = PetscAbsScalar(sc->alpha[0]) * hmnorm[1]; 781 else if (sc->p == cur_p + 1) errors[n] = PetscAbsScalar(sc->alpha[0]) * (hmnorm[2] + hmnorm[3]); 782 else continue; 783 candidates[n] = i; 784 orders[n] = PetscMin(sc->p, sc->q); /* order of global truncation error */ 785 costs[n] = sc->s; /* estimate the cost as the number of stages */ 786 if (i == gl->current_scheme) cur = n; 787 n++; 788 } 789 PetscCheck(cur >= 0 && gl->nschemes > cur, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Current scheme not found in scheme list"); 790 PetscCall(TSGLLEAdaptChoose(gl->adapt, n, orders, errors, costs, cur, h, tleft, &next_sc, next_h, finish)); 791 *next_scheme = candidates[next_sc]; 792 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, 793 gl->schemes[*next_scheme]->r, gl->schemes[*next_scheme]->s, (double)*next_h, PetscBools[*finish])); 794 PetscFunctionReturn(0); 795 } 796 797 static PetscErrorCode TSGLLEGetMaxSizes(TS ts, PetscInt *max_r, PetscInt *max_s) { 798 TS_GLLE *gl = (TS_GLLE *)ts->data; 799 800 PetscFunctionBegin; 801 *max_r = gl->schemes[gl->nschemes - 1]->r; 802 *max_s = gl->schemes[gl->nschemes - 1]->s; 803 PetscFunctionReturn(0); 804 } 805 806 static PetscErrorCode TSSolve_GLLE(TS ts) { 807 TS_GLLE *gl = (TS_GLLE *)ts->data; 808 PetscInt i, k, its, lits, max_r, max_s; 809 PetscBool final_step, finish; 810 SNESConvergedReason snesreason; 811 812 PetscFunctionBegin; 813 PetscCall(TSMonitor(ts, ts->steps, ts->ptime, ts->vec_sol)); 814 815 PetscCall(TSGLLEGetMaxSizes(ts, &max_r, &max_s)); 816 PetscCall(VecCopy(ts->vec_sol, gl->X[0])); 817 for (i = 1; i < max_r; i++) PetscCall(VecZeroEntries(gl->X[i])); 818 PetscCall(TSGLLEUpdateWRMS(ts)); 819 820 if (0) { 821 /* Find consistent initial data for DAE */ 822 gl->stage_time = ts->ptime + ts->time_step; 823 gl->scoeff = 1.; 824 gl->stage = 0; 825 826 PetscCall(VecCopy(ts->vec_sol, gl->Z)); 827 PetscCall(VecCopy(ts->vec_sol, gl->Y)); 828 PetscCall(SNESSolve(ts->snes, NULL, gl->Y)); 829 PetscCall(SNESGetIterationNumber(ts->snes, &its)); 830 PetscCall(SNESGetLinearSolveIterations(ts->snes, &lits)); 831 PetscCall(SNESGetConvergedReason(ts->snes, &snesreason)); 832 833 ts->snes_its += its; 834 ts->ksp_its += lits; 835 if (snesreason < 0 && ts->max_snes_failures > 0 && ++ts->num_snes_failures >= ts->max_snes_failures) { 836 ts->reason = TS_DIVERGED_NONLINEAR_SOLVE; 837 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)); 838 PetscFunctionReturn(0); 839 } 840 } 841 842 PetscCheck(gl->current_scheme >= 0, PETSC_COMM_SELF, PETSC_ERR_ORDER, "A starting scheme has not been provided"); 843 844 for (k = 0, final_step = PETSC_FALSE, finish = PETSC_FALSE; k < ts->max_steps && !finish; k++) { 845 PetscInt j, r, s, next_scheme = 0, rejections; 846 PetscReal h, hmnorm[4], enorm[3], next_h; 847 PetscBool accept; 848 const PetscScalar *c, *a, *u; 849 Vec *X, *Ydot, Y; 850 TSGLLEScheme scheme = gl->schemes[gl->current_scheme]; 851 852 r = scheme->r; 853 s = scheme->s; 854 c = scheme->c; 855 a = scheme->a; 856 u = scheme->u; 857 h = ts->time_step; 858 X = gl->X; 859 Ydot = gl->Ydot; 860 Y = gl->Y; 861 862 if (ts->ptime > ts->max_time) break; 863 864 /* 865 We only call PreStep at the start of each STEP, not each STAGE. This is because it is 866 possible to fail (have to restart a step) after multiple stages. 867 */ 868 PetscCall(TSPreStep(ts)); 869 870 rejections = 0; 871 while (1) { 872 for (i = 0; i < s; i++) { 873 PetscScalar shift; 874 gl->scoeff = 1. / PetscRealPart(a[i * s + i]); 875 shift = gl->scoeff / ts->time_step; 876 gl->stage = i; 877 gl->stage_time = ts->ptime + PetscRealPart(c[i]) * h; 878 879 /* 880 * Stage equation: Y = h A Y' + U X 881 * We assume that A is lower-triangular so that we can solve the stages (Y,Y') sequentially 882 * Build the affine vector z_i = -[1/(h a_ii)](h sum_j a_ij y'_j + sum_j u_ij x_j) 883 * Then y'_i = z + 1/(h a_ii) y_i 884 */ 885 PetscCall(VecZeroEntries(gl->Z)); 886 for (j = 0; j < r; j++) PetscCall(VecAXPY(gl->Z, -shift * u[i * r + j], X[j])); 887 for (j = 0; j < i; j++) PetscCall(VecAXPY(gl->Z, -shift * h * a[i * s + j], Ydot[j])); 888 /* Note: Z is used within function evaluation, Ydot = Z + shift*Y */ 889 890 /* Compute an estimate of Y to start Newton iteration */ 891 if (gl->extrapolate) { 892 if (i == 0) { 893 /* Linear extrapolation on the first stage */ 894 PetscCall(VecWAXPY(Y, c[i] * h, X[1], X[0])); 895 } else { 896 /* Linear extrapolation from the last stage */ 897 PetscCall(VecAXPY(Y, (c[i] - c[i - 1]) * h, Ydot[i - 1])); 898 } 899 } else if (i == 0) { /* Directly use solution from the last step, otherwise reuse the last stage (do nothing) */ 900 PetscCall(VecCopy(X[0], Y)); 901 } 902 903 /* Solve this stage (Ydot[i] is computed during function evaluation) */ 904 PetscCall(SNESSolve(ts->snes, NULL, Y)); 905 PetscCall(SNESGetIterationNumber(ts->snes, &its)); 906 PetscCall(SNESGetLinearSolveIterations(ts->snes, &lits)); 907 PetscCall(SNESGetConvergedReason(ts->snes, &snesreason)); 908 ts->snes_its += its; 909 ts->ksp_its += lits; 910 if (snesreason < 0 && ts->max_snes_failures > 0 && ++ts->num_snes_failures >= ts->max_snes_failures) { 911 ts->reason = TS_DIVERGED_NONLINEAR_SOLVE; 912 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)); 913 PetscFunctionReturn(0); 914 } 915 } 916 917 gl->stage_time = ts->ptime + ts->time_step; 918 919 PetscCall((*gl->EstimateHigherMoments)(scheme, h, Ydot, gl->X, gl->himom)); 920 /* hmnorm[i] = h^{p+i}x^{(p+i)} with i=0,1,2; hmnorm[3] = h^{p+2}(dx'/dx) x^{(p+1)} */ 921 for (i = 0; i < 3; i++) PetscCall(TSGLLEVecNormWRMS(ts, gl->himom[i], &hmnorm[i + 1])); 922 enorm[0] = PetscRealPart(scheme->alpha[0]) * hmnorm[1]; 923 enorm[1] = PetscRealPart(scheme->beta[0]) * hmnorm[2]; 924 enorm[2] = PetscRealPart(scheme->gamma[0]) * hmnorm[3]; 925 PetscCall((*gl->Accept)(ts, ts->max_time - gl->stage_time, h, enorm, &accept)); 926 if (accept) goto accepted; 927 rejections++; 928 PetscCall(PetscInfo(ts, "Step %" PetscInt_FMT " (t=%g) not accepted, rejections=%" PetscInt_FMT "\n", k, (double)gl->stage_time, rejections)); 929 if (rejections > gl->max_step_rejections) break; 930 /* 931 There are lots of reasons why a step might be rejected, including solvers not converging and other factors that 932 TSGLLEChooseNextScheme does not support. Additionally, the error estimates may be very screwed up, so I'm not 933 convinced that it's safe to just compute a new error estimate using the same interface as the current adaptor 934 (the adaptor interface probably has to change). Here we make an arbitrary and naive choice. This assumes that 935 steps were written in Nordsieck form. The "correct" method would be to re-complete the previous time step with 936 the correct "next" step size. It is unclear to me whether the present ad-hoc method of rescaling X is stable. 937 */ 938 h *= 0.5; 939 for (i = 1; i < scheme->r; i++) PetscCall(VecScale(X[i], PetscPowRealInt(0.5, i))); 940 } 941 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); 942 943 accepted: 944 /* This term is not error, but it *would* be the leading term for a lower order method */ 945 PetscCall(TSGLLEVecNormWRMS(ts, gl->X[scheme->r - 1], &hmnorm[0])); 946 /* Correct scaling so that these are equivalent to norms of the Nordsieck vectors */ 947 948 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])); 949 if (!final_step) { 950 PetscCall(TSGLLEChooseNextScheme(ts, h, hmnorm, &next_scheme, &next_h, &final_step)); 951 } else { 952 /* Dummy values to complete the current step in a consistent manner */ 953 next_scheme = gl->current_scheme; 954 next_h = h; 955 finish = PETSC_TRUE; 956 } 957 958 X = gl->Xold; 959 gl->Xold = gl->X; 960 gl->X = X; 961 PetscCall((*gl->CompleteStep)(scheme, h, gl->schemes[next_scheme], next_h, Ydot, gl->Xold, gl->X)); 962 963 PetscCall(TSGLLEUpdateWRMS(ts)); 964 965 /* Post the solution for the user, we could avoid this copy with a small bit of cleverness */ 966 PetscCall(VecCopy(gl->X[0], ts->vec_sol)); 967 ts->ptime += h; 968 ts->steps++; 969 970 PetscCall(TSPostEvaluate(ts)); 971 PetscCall(TSPostStep(ts)); 972 PetscCall(TSMonitor(ts, ts->steps, ts->ptime, ts->vec_sol)); 973 974 gl->current_scheme = next_scheme; 975 ts->time_step = next_h; 976 } 977 PetscFunctionReturn(0); 978 } 979 980 /*------------------------------------------------------------*/ 981 982 static PetscErrorCode TSReset_GLLE(TS ts) { 983 TS_GLLE *gl = (TS_GLLE *)ts->data; 984 PetscInt max_r, max_s; 985 986 PetscFunctionBegin; 987 if (gl->setupcalled) { 988 PetscCall(TSGLLEGetMaxSizes(ts, &max_r, &max_s)); 989 PetscCall(VecDestroyVecs(max_r, &gl->Xold)); 990 PetscCall(VecDestroyVecs(max_r, &gl->X)); 991 PetscCall(VecDestroyVecs(max_s, &gl->Ydot)); 992 PetscCall(VecDestroyVecs(3, &gl->himom)); 993 PetscCall(VecDestroy(&gl->W)); 994 PetscCall(VecDestroy(&gl->Y)); 995 PetscCall(VecDestroy(&gl->Z)); 996 } 997 gl->setupcalled = PETSC_FALSE; 998 PetscFunctionReturn(0); 999 } 1000 1001 static PetscErrorCode TSDestroy_GLLE(TS ts) { 1002 TS_GLLE *gl = (TS_GLLE *)ts->data; 1003 1004 PetscFunctionBegin; 1005 PetscCall(TSReset_GLLE(ts)); 1006 if (ts->dm) { 1007 PetscCall(DMCoarsenHookRemove(ts->dm, DMCoarsenHook_TSGLLE, DMRestrictHook_TSGLLE, ts)); 1008 PetscCall(DMSubDomainHookRemove(ts->dm, DMSubDomainHook_TSGLLE, DMSubDomainRestrictHook_TSGLLE, ts)); 1009 } 1010 if (gl->adapt) PetscCall(TSGLLEAdaptDestroy(&gl->adapt)); 1011 if (gl->Destroy) PetscCall((*gl->Destroy)(gl)); 1012 PetscCall(PetscFree(ts->data)); 1013 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSGLLESetType_C", NULL)); 1014 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSGLLESetAcceptType_C", NULL)); 1015 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSGLLEGetAdapt_C", NULL)); 1016 PetscFunctionReturn(0); 1017 } 1018 1019 /* 1020 This defines the nonlinear equation that is to be solved with SNES 1021 g(x) = f(t,x,z+shift*x) = 0 1022 */ 1023 static PetscErrorCode SNESTSFormFunction_GLLE(SNES snes, Vec x, Vec f, TS ts) { 1024 TS_GLLE *gl = (TS_GLLE *)ts->data; 1025 Vec Z, Ydot; 1026 DM dm, dmsave; 1027 1028 PetscFunctionBegin; 1029 PetscCall(SNESGetDM(snes, &dm)); 1030 PetscCall(TSGLLEGetVecs(ts, dm, &Z, &Ydot)); 1031 PetscCall(VecWAXPY(Ydot, gl->scoeff / ts->time_step, x, Z)); 1032 dmsave = ts->dm; 1033 ts->dm = dm; 1034 PetscCall(TSComputeIFunction(ts, gl->stage_time, x, Ydot, f, PETSC_FALSE)); 1035 ts->dm = dmsave; 1036 PetscCall(TSGLLERestoreVecs(ts, dm, &Z, &Ydot)); 1037 PetscFunctionReturn(0); 1038 } 1039 1040 static PetscErrorCode SNESTSFormJacobian_GLLE(SNES snes, Vec x, Mat A, Mat B, TS ts) { 1041 TS_GLLE *gl = (TS_GLLE *)ts->data; 1042 Vec Z, Ydot; 1043 DM dm, dmsave; 1044 1045 PetscFunctionBegin; 1046 PetscCall(SNESGetDM(snes, &dm)); 1047 PetscCall(TSGLLEGetVecs(ts, dm, &Z, &Ydot)); 1048 dmsave = ts->dm; 1049 ts->dm = dm; 1050 /* gl->Xdot will have already been computed in SNESTSFormFunction_GLLE */ 1051 PetscCall(TSComputeIJacobian(ts, gl->stage_time, x, gl->Ydot[gl->stage], gl->scoeff / ts->time_step, A, B, PETSC_FALSE)); 1052 ts->dm = dmsave; 1053 PetscCall(TSGLLERestoreVecs(ts, dm, &Z, &Ydot)); 1054 PetscFunctionReturn(0); 1055 } 1056 1057 static PetscErrorCode TSSetUp_GLLE(TS ts) { 1058 TS_GLLE *gl = (TS_GLLE *)ts->data; 1059 PetscInt max_r, max_s; 1060 DM dm; 1061 1062 PetscFunctionBegin; 1063 gl->setupcalled = PETSC_TRUE; 1064 PetscCall(TSGLLEGetMaxSizes(ts, &max_r, &max_s)); 1065 PetscCall(VecDuplicateVecs(ts->vec_sol, max_r, &gl->X)); 1066 PetscCall(VecDuplicateVecs(ts->vec_sol, max_r, &gl->Xold)); 1067 PetscCall(VecDuplicateVecs(ts->vec_sol, max_s, &gl->Ydot)); 1068 PetscCall(VecDuplicateVecs(ts->vec_sol, 3, &gl->himom)); 1069 PetscCall(VecDuplicate(ts->vec_sol, &gl->W)); 1070 PetscCall(VecDuplicate(ts->vec_sol, &gl->Y)); 1071 PetscCall(VecDuplicate(ts->vec_sol, &gl->Z)); 1072 1073 /* Default acceptance tests and adaptivity */ 1074 if (!gl->Accept) PetscCall(TSGLLESetAcceptType(ts, TSGLLEACCEPT_ALWAYS)); 1075 if (!gl->adapt) PetscCall(TSGLLEGetAdapt(ts, &gl->adapt)); 1076 1077 if (gl->current_scheme < 0) { 1078 PetscInt i; 1079 for (i = 0;; i++) { 1080 if (gl->schemes[i]->p == gl->start_order) break; 1081 PetscCheck(i + 1 != gl->nschemes, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No schemes available with requested start order %" PetscInt_FMT, i); 1082 } 1083 gl->current_scheme = i; 1084 } 1085 PetscCall(TSGetDM(ts, &dm)); 1086 PetscCall(DMCoarsenHookAdd(dm, DMCoarsenHook_TSGLLE, DMRestrictHook_TSGLLE, ts)); 1087 PetscCall(DMSubDomainHookAdd(dm, DMSubDomainHook_TSGLLE, DMSubDomainRestrictHook_TSGLLE, ts)); 1088 PetscFunctionReturn(0); 1089 } 1090 /*------------------------------------------------------------*/ 1091 1092 static PetscErrorCode TSSetFromOptions_GLLE(TS ts, PetscOptionItems *PetscOptionsObject) { 1093 TS_GLLE *gl = (TS_GLLE *)ts->data; 1094 char tname[256] = TSGLLE_IRKS, completef[256] = "rescale-and-modify"; 1095 1096 PetscFunctionBegin; 1097 PetscOptionsHeadBegin(PetscOptionsObject, "General Linear ODE solver options"); 1098 { 1099 PetscBool flg; 1100 PetscCall(PetscOptionsFList("-ts_gl_type", "Type of GL method", "TSGLLESetType", TSGLLEList, gl->type_name[0] ? gl->type_name : tname, tname, sizeof(tname), &flg)); 1101 if (flg || !gl->type_name[0]) PetscCall(TSGLLESetType(ts, tname)); 1102 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)); 1103 PetscCall(PetscOptionsInt("-ts_gl_max_order", "Maximum order to try", "TSGLLESetMaxOrder", gl->max_order, &gl->max_order, NULL)); 1104 PetscCall(PetscOptionsInt("-ts_gl_min_order", "Minimum order to try", "TSGLLESetMinOrder", gl->min_order, &gl->min_order, NULL)); 1105 PetscCall(PetscOptionsInt("-ts_gl_start_order", "Initial order to try", "TSGLLESetMinOrder", gl->start_order, &gl->start_order, NULL)); 1106 PetscCall(PetscOptionsEnum("-ts_gl_error_direction", "Which direction to look when estimating error", "TSGLLESetErrorDirection", TSGLLEErrorDirections, (PetscEnum)gl->error_direction, (PetscEnum *)&gl->error_direction, NULL)); 1107 PetscCall(PetscOptionsBool("-ts_gl_extrapolate", "Extrapolate stage solution from previous solution (sometimes unstable)", "TSGLLESetExtrapolate", gl->extrapolate, &gl->extrapolate, NULL)); 1108 PetscCall(PetscOptionsReal("-ts_gl_atol", "Absolute tolerance", "TSGLLESetTolerances", gl->wrms_atol, &gl->wrms_atol, NULL)); 1109 PetscCall(PetscOptionsReal("-ts_gl_rtol", "Relative tolerance", "TSGLLESetTolerances", gl->wrms_rtol, &gl->wrms_rtol, NULL)); 1110 PetscCall(PetscOptionsString("-ts_gl_complete", "Method to use for completing the step", "none", completef, completef, sizeof(completef), &flg)); 1111 if (flg) { 1112 PetscBool match1, match2; 1113 PetscCall(PetscStrcmp(completef, "rescale", &match1)); 1114 PetscCall(PetscStrcmp(completef, "rescale-and-modify", &match2)); 1115 if (match1) gl->CompleteStep = TSGLLECompleteStep_Rescale; 1116 else if (match2) gl->CompleteStep = TSGLLECompleteStep_RescaleAndModify; 1117 else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "%s", completef); 1118 } 1119 { 1120 char type[256] = TSGLLEACCEPT_ALWAYS; 1121 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)); 1122 if (flg || !gl->accept_name[0]) PetscCall(TSGLLESetAcceptType(ts, type)); 1123 } 1124 { 1125 TSGLLEAdapt adapt; 1126 PetscCall(TSGLLEGetAdapt(ts, &adapt)); 1127 PetscCall(TSGLLEAdaptSetFromOptions(adapt, PetscOptionsObject)); 1128 } 1129 } 1130 PetscOptionsHeadEnd(); 1131 PetscFunctionReturn(0); 1132 } 1133 1134 static PetscErrorCode TSView_GLLE(TS ts, PetscViewer viewer) { 1135 TS_GLLE *gl = (TS_GLLE *)ts->data; 1136 PetscInt i; 1137 PetscBool iascii, details; 1138 1139 PetscFunctionBegin; 1140 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 1141 if (iascii) { 1142 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)); 1143 PetscCall(PetscViewerASCIIPrintf(viewer, " Error estimation: %s\n", TSGLLEErrorDirections[gl->error_direction])); 1144 PetscCall(PetscViewerASCIIPrintf(viewer, " Extrapolation: %s\n", gl->extrapolate ? "yes" : "no")); 1145 PetscCall(PetscViewerASCIIPrintf(viewer, " Acceptance test: %s\n", gl->accept_name[0] ? gl->accept_name : "(not yet set)")); 1146 PetscCall(PetscViewerASCIIPushTab(viewer)); 1147 PetscCall(TSGLLEAdaptView(gl->adapt, viewer)); 1148 PetscCall(PetscViewerASCIIPopTab(viewer)); 1149 PetscCall(PetscViewerASCIIPrintf(viewer, " type: %s\n", gl->type_name[0] ? gl->type_name : "(not yet set)")); 1150 PetscCall(PetscViewerASCIIPrintf(viewer, "Schemes within family (%" PetscInt_FMT "):\n", gl->nschemes)); 1151 details = PETSC_FALSE; 1152 PetscCall(PetscOptionsGetBool(((PetscObject)ts)->options, ((PetscObject)ts)->prefix, "-ts_gl_view_detailed", &details, NULL)); 1153 PetscCall(PetscViewerASCIIPushTab(viewer)); 1154 for (i = 0; i < gl->nschemes; i++) PetscCall(TSGLLESchemeView(gl->schemes[i], details, viewer)); 1155 if (gl->View) PetscCall((*gl->View)(gl, viewer)); 1156 PetscCall(PetscViewerASCIIPopTab(viewer)); 1157 } 1158 PetscFunctionReturn(0); 1159 } 1160 1161 /*@C 1162 TSGLLERegister - adds a TSGLLE implementation 1163 1164 Not Collective 1165 1166 Input Parameters: 1167 + name_scheme - name of user-defined general linear scheme 1168 - routine_create - routine to create method context 1169 1170 Notes: 1171 TSGLLERegister() may be called multiple times to add several user-defined families. 1172 1173 Sample usage: 1174 .vb 1175 TSGLLERegister("my_scheme",MySchemeCreate); 1176 .ve 1177 1178 Then, your scheme can be chosen with the procedural interface via 1179 $ TSGLLESetType(ts,"my_scheme") 1180 or at runtime via the option 1181 $ -ts_gl_type my_scheme 1182 1183 Level: advanced 1184 1185 .seealso: `TSGLLERegisterAll()` 1186 @*/ 1187 PetscErrorCode TSGLLERegister(const char sname[], PetscErrorCode (*function)(TS)) { 1188 PetscFunctionBegin; 1189 PetscCall(TSGLLEInitializePackage()); 1190 PetscCall(PetscFunctionListAdd(&TSGLLEList, sname, function)); 1191 PetscFunctionReturn(0); 1192 } 1193 1194 /*@C 1195 TSGLLEAcceptRegister - adds a TSGLLE acceptance scheme 1196 1197 Not Collective 1198 1199 Input Parameters: 1200 + name_scheme - name of user-defined acceptance scheme 1201 - routine_create - routine to create method context 1202 1203 Notes: 1204 TSGLLEAcceptRegister() may be called multiple times to add several user-defined families. 1205 1206 Sample usage: 1207 .vb 1208 TSGLLEAcceptRegister("my_scheme",MySchemeCreate); 1209 .ve 1210 1211 Then, your scheme can be chosen with the procedural interface via 1212 $ TSGLLESetAcceptType(ts,"my_scheme") 1213 or at runtime via the option 1214 $ -ts_gl_accept_type my_scheme 1215 1216 Level: advanced 1217 1218 .seealso: `TSGLLERegisterAll()` 1219 @*/ 1220 PetscErrorCode TSGLLEAcceptRegister(const char sname[], TSGLLEAcceptFunction function) { 1221 PetscFunctionBegin; 1222 PetscCall(PetscFunctionListAdd(&TSGLLEAcceptList, sname, function)); 1223 PetscFunctionReturn(0); 1224 } 1225 1226 /*@C 1227 TSGLLERegisterAll - Registers all of the general linear methods in TSGLLE 1228 1229 Not Collective 1230 1231 Level: advanced 1232 1233 .seealso: `TSGLLERegisterDestroy()` 1234 @*/ 1235 PetscErrorCode TSGLLERegisterAll(void) { 1236 PetscFunctionBegin; 1237 if (TSGLLERegisterAllCalled) PetscFunctionReturn(0); 1238 TSGLLERegisterAllCalled = PETSC_TRUE; 1239 1240 PetscCall(TSGLLERegister(TSGLLE_IRKS, TSGLLECreate_IRKS)); 1241 PetscCall(TSGLLEAcceptRegister(TSGLLEACCEPT_ALWAYS, TSGLLEAccept_Always)); 1242 PetscFunctionReturn(0); 1243 } 1244 1245 /*@C 1246 TSGLLEInitializePackage - This function initializes everything in the TSGLLE package. It is called 1247 from TSInitializePackage(). 1248 1249 Level: developer 1250 1251 .seealso: `PetscInitialize()` 1252 @*/ 1253 PetscErrorCode TSGLLEInitializePackage(void) { 1254 PetscFunctionBegin; 1255 if (TSGLLEPackageInitialized) PetscFunctionReturn(0); 1256 TSGLLEPackageInitialized = PETSC_TRUE; 1257 PetscCall(TSGLLERegisterAll()); 1258 PetscCall(PetscRegisterFinalize(TSGLLEFinalizePackage)); 1259 PetscFunctionReturn(0); 1260 } 1261 1262 /*@C 1263 TSGLLEFinalizePackage - This function destroys everything in the TSGLLE package. It is 1264 called from PetscFinalize(). 1265 1266 Level: developer 1267 1268 .seealso: `PetscFinalize()` 1269 @*/ 1270 PetscErrorCode TSGLLEFinalizePackage(void) { 1271 PetscFunctionBegin; 1272 PetscCall(PetscFunctionListDestroy(&TSGLLEList)); 1273 PetscCall(PetscFunctionListDestroy(&TSGLLEAcceptList)); 1274 TSGLLEPackageInitialized = PETSC_FALSE; 1275 TSGLLERegisterAllCalled = PETSC_FALSE; 1276 PetscFunctionReturn(0); 1277 } 1278 1279 /* ------------------------------------------------------------ */ 1280 /*MC 1281 TSGLLE - DAE solver using implicit General Linear methods 1282 1283 These methods contain Runge-Kutta and multistep schemes as special cases. These special cases have some fundamental 1284 limitations. For example, diagonally implicit Runge-Kutta cannot have stage order greater than 1 which limits their 1285 applicability to very stiff systems. Meanwhile, multistep methods cannot be A-stable for order greater than 2 and BDF 1286 are not 0-stable for order greater than 6. GL methods can be A- and L-stable with arbitrarily high stage order and 1287 reliable error estimates for both 1 and 2 orders higher to facilitate adaptive step sizes and adaptive order schemes. 1288 All this is possible while preserving a singly diagonally implicit structure. 1289 1290 Options database keys: 1291 + -ts_gl_type <type> - the class of general linear method (irks) 1292 . -ts_gl_rtol <tol> - relative error 1293 . -ts_gl_atol <tol> - absolute error 1294 . -ts_gl_min_order <p> - minimum order method to consider (default=1) 1295 . -ts_gl_max_order <p> - maximum order method to consider (default=3) 1296 . -ts_gl_start_order <p> - order of starting method (default=1) 1297 . -ts_gl_complete <method> - method to use for completing the step (rescale-and-modify or rescale) 1298 - -ts_adapt_type <method> - adaptive controller to use (none step both) 1299 1300 Notes: 1301 This integrator can be applied to DAE. 1302 1303 Diagonally implicit general linear (DIGL) methods are a generalization of diagonally implicit Runge-Kutta (DIRK). 1304 They are represented by the tableau 1305 1306 .vb 1307 A | U 1308 ------- 1309 B | V 1310 .ve 1311 1312 combined with a vector c of abscissa. "Diagonally implicit" means that A is lower triangular. 1313 A step of the general method reads 1314 1315 .vb 1316 [ Y ] = [A U] [ Y' ] 1317 [X^k] = [B V] [X^{k-1}] 1318 .ve 1319 1320 where Y is the multivector of stage values, Y' is the multivector of stage derivatives, X^k is the Nordsieck vector of 1321 the solution at step k. The Nordsieck vector consists of the first r moments of the solution, given by 1322 1323 .vb 1324 X = [x_0,x_1,...,x_{r-1}] = [x, h x', h^2 x'', ..., h^{r-1} x^{(r-1)} ] 1325 .ve 1326 1327 If A is lower triangular, we can solve the stages (Y,Y') sequentially 1328 1329 .vb 1330 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} 1331 .ve 1332 1333 and then construct the pieces to carry to the next step 1334 1335 .vb 1336 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} 1337 .ve 1338 1339 Note that when the equations are cast in implicit form, we are using the stage equation to define y'_i 1340 in terms of y_i and known stuff (y_j for j<i and x_j for all j). 1341 1342 Error estimation 1343 1344 At present, the most attractive GL methods for stiff problems are singly diagonally implicit schemes which posses 1345 Inherent Runge-Kutta Stability (IRKS). These methods have r=s, the number of items passed between steps is equal to 1346 the number of stages. The order and stage-order are one less than the number of stages. We use the error estimates 1347 in the 2007 paper which provide the following estimates 1348 1349 .vb 1350 h^{p+1} X^{(p+1)} = phi_0^T Y' + [0 psi_0^T] Xold 1351 h^{p+2} X^{(p+2)} = phi_1^T Y' + [0 psi_1^T] Xold 1352 h^{p+2} (dx'/dx) X^{(p+1)} = phi_2^T Y' + [0 psi_2^T] Xold 1353 .ve 1354 1355 These estimates are accurate to O(h^{p+3}). 1356 1357 Changing the step size 1358 1359 We use the generalized "rescale and modify" scheme, see equation (4.5) of the 2007 paper. 1360 1361 Level: beginner 1362 1363 References: 1364 + * - John Butcher and Z. Jackieweicz and W. Wright, On error propagation in general linear methods for 1365 ordinary differential equations, Journal of Complexity, Vol 23, 2007. 1366 - * - John Butcher, Numerical methods for ordinary differential equations, second edition, Wiley, 2009. 1367 1368 .seealso: `TSCreate()`, `TS`, `TSSetType()` 1369 1370 M*/ 1371 PETSC_EXTERN PetscErrorCode TSCreate_GLLE(TS ts) { 1372 TS_GLLE *gl; 1373 1374 PetscFunctionBegin; 1375 PetscCall(TSGLLEInitializePackage()); 1376 1377 PetscCall(PetscNew(&gl)); 1378 ts->data = (void *)gl; 1379 1380 ts->ops->reset = TSReset_GLLE; 1381 ts->ops->destroy = TSDestroy_GLLE; 1382 ts->ops->view = TSView_GLLE; 1383 ts->ops->setup = TSSetUp_GLLE; 1384 ts->ops->solve = TSSolve_GLLE; 1385 ts->ops->setfromoptions = TSSetFromOptions_GLLE; 1386 ts->ops->snesfunction = SNESTSFormFunction_GLLE; 1387 ts->ops->snesjacobian = SNESTSFormJacobian_GLLE; 1388 1389 ts->usessnes = PETSC_TRUE; 1390 1391 gl->max_step_rejections = 1; 1392 gl->min_order = 1; 1393 gl->max_order = 3; 1394 gl->start_order = 1; 1395 gl->current_scheme = -1; 1396 gl->extrapolate = PETSC_FALSE; 1397 1398 gl->wrms_atol = 1e-8; 1399 gl->wrms_rtol = 1e-5; 1400 1401 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSGLLESetType_C", &TSGLLESetType_GLLE)); 1402 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSGLLESetAcceptType_C", &TSGLLESetAcceptType_GLLE)); 1403 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSGLLEGetAdapt_C", &TSGLLEGetAdapt_GLLE)); 1404 PetscFunctionReturn(0); 1405 } 1406