1211a84d6SLisandro Dalcin /* 2211a84d6SLisandro Dalcin Code for timestepping with BDF methods 3211a84d6SLisandro Dalcin */ 4211a84d6SLisandro Dalcin #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 51117961dSLisandro Dalcin #include <petscdm.h> 6211a84d6SLisandro Dalcin 7211a84d6SLisandro Dalcin static PetscBool cited = PETSC_FALSE; 89371c9d4SSatish Balay static const char citation[] = "@book{Brenan1995,\n" 9211a84d6SLisandro Dalcin " title = {Numerical Solution of Initial-Value Problems in Differential-Algebraic Equations},\n" 10211a84d6SLisandro Dalcin " author = {Brenan, K. and Campbell, S. and Petzold, L.},\n" 11211a84d6SLisandro Dalcin " publisher = {Society for Industrial and Applied Mathematics},\n" 12211a84d6SLisandro Dalcin " year = {1995},\n" 13211a84d6SLisandro Dalcin " doi = {10.1137/1.9781611971224},\n}\n"; 14211a84d6SLisandro Dalcin 15211a84d6SLisandro Dalcin typedef struct { 16211a84d6SLisandro Dalcin PetscInt k, n; 17211a84d6SLisandro Dalcin PetscReal time[6 + 2]; 18211a84d6SLisandro Dalcin Vec work[6 + 2]; 19e3c11fc1SJed Brown Vec tvwork[6 + 2]; 20211a84d6SLisandro Dalcin PetscReal shift; 21e3c11fc1SJed Brown Vec vec_dot; /* Xdot when !transientvar, else Cdot where C(X) is the transient variable. */ 221117961dSLisandro Dalcin Vec vec_wrk; 23211a84d6SLisandro Dalcin Vec vec_lte; 24211a84d6SLisandro Dalcin 25e3c11fc1SJed Brown PetscBool transientvar; 2699e85243SStefano Zampini PetscBool extrapolate; 27211a84d6SLisandro Dalcin PetscInt order; 28211a84d6SLisandro Dalcin TSStepStatus status; 29211a84d6SLisandro Dalcin } TS_BDF; 30211a84d6SLisandro Dalcin 31e3c11fc1SJed Brown /* Compute Lagrange polynomials on T[:n] evaluated at t. 32e3c11fc1SJed Brown * If one has data (T[i], Y[i]), then the interpolation/extrapolation is f(t) = \sum_i L[i]*Y[i]. 33e3c11fc1SJed Brown */ 34d71ae5a4SJacob Faibussowitsch static inline void LagrangeBasisVals(PetscInt n, PetscReal t, const PetscReal T[], PetscScalar L[]) 35d71ae5a4SJacob Faibussowitsch { 36211a84d6SLisandro Dalcin PetscInt k, j; 37211a84d6SLisandro Dalcin for (k = 0; k < n; k++) 38211a84d6SLisandro Dalcin for (L[k] = 1, j = 0; j < n; j++) 399371c9d4SSatish Balay if (j != k) L[k] *= (t - T[j]) / (T[k] - T[j]); 40211a84d6SLisandro Dalcin } 41211a84d6SLisandro Dalcin 42d71ae5a4SJacob Faibussowitsch static inline void LagrangeBasisDers(PetscInt n, PetscReal t, const PetscReal T[], PetscScalar dL[]) 43d71ae5a4SJacob Faibussowitsch { 44211a84d6SLisandro Dalcin PetscInt k, j, i; 45211a84d6SLisandro Dalcin for (k = 0; k < n; k++) 46211a84d6SLisandro Dalcin for (dL[k] = 0, j = 0; j < n; j++) 47211a84d6SLisandro Dalcin if (j != k) { 48211a84d6SLisandro Dalcin PetscReal L = 1 / (T[k] - T[j]); 49211a84d6SLisandro Dalcin for (i = 0; i < n; i++) 509371c9d4SSatish Balay if (i != j && i != k) L *= (t - T[i]) / (T[k] - T[i]); 51211a84d6SLisandro Dalcin dL[k] += L; 52211a84d6SLisandro Dalcin } 53211a84d6SLisandro Dalcin } 54211a84d6SLisandro Dalcin 55d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSBDF_GetVecs(TS ts, DM dm, Vec *Xdot, Vec *Ydot) 56d71ae5a4SJacob Faibussowitsch { 571117961dSLisandro Dalcin TS_BDF *bdf = (TS_BDF *)ts->data; 581117961dSLisandro Dalcin 591117961dSLisandro Dalcin PetscFunctionBegin; 601117961dSLisandro Dalcin if (dm && dm != ts->dm) { 619566063dSJacob Faibussowitsch PetscCall(DMGetNamedGlobalVector(dm, "TSBDF_Vec_Xdot", Xdot)); 629566063dSJacob Faibussowitsch PetscCall(DMGetNamedGlobalVector(dm, "TSBDF_Vec_Ydot", Ydot)); 631117961dSLisandro Dalcin } else { 641117961dSLisandro Dalcin *Xdot = bdf->vec_dot; 651117961dSLisandro Dalcin *Ydot = bdf->vec_wrk; 661117961dSLisandro Dalcin } 673ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 681117961dSLisandro Dalcin } 691117961dSLisandro Dalcin 70d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSBDF_RestoreVecs(TS ts, DM dm, Vec *Xdot, Vec *Ydot) 71d71ae5a4SJacob Faibussowitsch { 721117961dSLisandro Dalcin TS_BDF *bdf = (TS_BDF *)ts->data; 731117961dSLisandro Dalcin 741117961dSLisandro Dalcin PetscFunctionBegin; 751117961dSLisandro Dalcin if (dm && dm != ts->dm) { 769566063dSJacob Faibussowitsch PetscCall(DMRestoreNamedGlobalVector(dm, "TSBDF_Vec_Xdot", Xdot)); 779566063dSJacob Faibussowitsch PetscCall(DMRestoreNamedGlobalVector(dm, "TSBDF_Vec_Ydot", Ydot)); 781117961dSLisandro Dalcin } else { 793c633725SBarry Smith PetscCheck(*Xdot == bdf->vec_dot, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_INCOMP, "Vec does not match the cache"); 803c633725SBarry Smith PetscCheck(*Ydot == bdf->vec_wrk, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_INCOMP, "Vec does not match the cache"); 811117961dSLisandro Dalcin *Xdot = NULL; 821117961dSLisandro Dalcin *Ydot = NULL; 831117961dSLisandro Dalcin } 843ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 851117961dSLisandro Dalcin } 861117961dSLisandro Dalcin 87*2a8381b2SBarry Smith static PetscErrorCode DMCoarsenHook_TSBDF(DM fine, DM coarse, PetscCtx ctx) 88d71ae5a4SJacob Faibussowitsch { 891117961dSLisandro Dalcin PetscFunctionBegin; 903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 911117961dSLisandro Dalcin } 921117961dSLisandro Dalcin 93*2a8381b2SBarry Smith static PetscErrorCode DMRestrictHook_TSBDF(DM fine, Mat restrct, Vec rscale, Mat inject, DM coarse, PetscCtx ctx) 94d71ae5a4SJacob Faibussowitsch { 951117961dSLisandro Dalcin TS ts = (TS)ctx; 961117961dSLisandro Dalcin Vec Ydot, Ydot_c; 971117961dSLisandro Dalcin Vec Xdot, Xdot_c; 981117961dSLisandro Dalcin 991117961dSLisandro Dalcin PetscFunctionBegin; 1009566063dSJacob Faibussowitsch PetscCall(TSBDF_GetVecs(ts, fine, &Xdot, &Ydot)); 1019566063dSJacob Faibussowitsch PetscCall(TSBDF_GetVecs(ts, coarse, &Xdot_c, &Ydot_c)); 1021117961dSLisandro Dalcin 1039566063dSJacob Faibussowitsch PetscCall(MatRestrict(restrct, Ydot, Ydot_c)); 1049566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(Ydot_c, rscale, Ydot_c)); 1051117961dSLisandro Dalcin 1069566063dSJacob Faibussowitsch PetscCall(TSBDF_RestoreVecs(ts, fine, &Xdot, &Ydot)); 1079566063dSJacob Faibussowitsch PetscCall(TSBDF_RestoreVecs(ts, coarse, &Xdot_c, &Ydot_c)); 1083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1091117961dSLisandro Dalcin } 1101117961dSLisandro Dalcin 111d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSBDF_Advance(TS ts, PetscReal t, Vec X) 112d71ae5a4SJacob Faibussowitsch { 113211a84d6SLisandro Dalcin TS_BDF *bdf = (TS_BDF *)ts->data; 114400f6f02SBarry Smith PetscInt i, n = PETSC_STATIC_ARRAY_LENGTH(bdf->work); 115e3c11fc1SJed Brown Vec tail = bdf->work[n - 1], tvtail = bdf->tvwork[n - 1]; 116211a84d6SLisandro Dalcin 117211a84d6SLisandro Dalcin PetscFunctionBegin; 118211a84d6SLisandro Dalcin for (i = n - 1; i >= 2; i--) { 119211a84d6SLisandro Dalcin bdf->time[i] = bdf->time[i - 1]; 120211a84d6SLisandro Dalcin bdf->work[i] = bdf->work[i - 1]; 121e3c11fc1SJed Brown bdf->tvwork[i] = bdf->tvwork[i - 1]; 122211a84d6SLisandro Dalcin } 123211a84d6SLisandro Dalcin bdf->n = PetscMin(bdf->n + 1, n - 1); 124211a84d6SLisandro Dalcin bdf->time[1] = t; 125211a84d6SLisandro Dalcin bdf->work[1] = tail; 126e3c11fc1SJed Brown bdf->tvwork[1] = tvtail; 1279566063dSJacob Faibussowitsch PetscCall(VecCopy(X, tail)); 1289566063dSJacob Faibussowitsch PetscCall(TSComputeTransientVariable(ts, tail, tvtail)); 1293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 130211a84d6SLisandro Dalcin } 131211a84d6SLisandro Dalcin 132d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSBDF_VecLTE(TS ts, PetscInt order, Vec lte) 133d71ae5a4SJacob Faibussowitsch { 134211a84d6SLisandro Dalcin TS_BDF *bdf = (TS_BDF *)ts->data; 135211a84d6SLisandro Dalcin PetscInt i, n = order + 1; 136211a84d6SLisandro Dalcin PetscReal *time = bdf->time; 137211a84d6SLisandro Dalcin Vec *vecs = bdf->work; 138211a84d6SLisandro Dalcin PetscScalar a[8], b[8], alpha[8]; 139211a84d6SLisandro Dalcin 140211a84d6SLisandro Dalcin PetscFunctionBegin; 1419371c9d4SSatish Balay LagrangeBasisDers(n + 0, time[0], time, a); 1429371c9d4SSatish Balay a[n] = 0; 143211a84d6SLisandro Dalcin LagrangeBasisDers(n + 1, time[0], time, b); 144211a84d6SLisandro Dalcin for (i = 0; i < n + 1; i++) alpha[i] = (a[i] - b[i]) / a[0]; 1459566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(lte)); 1469566063dSJacob Faibussowitsch PetscCall(VecMAXPY(lte, n + 1, alpha, vecs)); 1473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 148211a84d6SLisandro Dalcin } 149211a84d6SLisandro Dalcin 150d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSBDF_Extrapolate(TS ts, PetscInt order, PetscReal t, Vec X) 151d71ae5a4SJacob Faibussowitsch { 152211a84d6SLisandro Dalcin TS_BDF *bdf = (TS_BDF *)ts->data; 153211a84d6SLisandro Dalcin PetscInt n = order + 1; 154211a84d6SLisandro Dalcin PetscReal *time = bdf->time + 1; 155211a84d6SLisandro Dalcin Vec *vecs = bdf->work + 1; 156211a84d6SLisandro Dalcin PetscScalar alpha[7]; 157211a84d6SLisandro Dalcin 158211a84d6SLisandro Dalcin PetscFunctionBegin; 159211a84d6SLisandro Dalcin n = PetscMin(n, bdf->n); 160211a84d6SLisandro Dalcin LagrangeBasisVals(n, t, time, alpha); 1619566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(X)); 1629566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X, n, alpha, vecs)); 1633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 164211a84d6SLisandro Dalcin } 165211a84d6SLisandro Dalcin 166d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSBDF_Interpolate(TS ts, PetscInt order, PetscReal t, Vec X) 167d71ae5a4SJacob Faibussowitsch { 168211a84d6SLisandro Dalcin TS_BDF *bdf = (TS_BDF *)ts->data; 169211a84d6SLisandro Dalcin PetscInt n = order + 1; 170211a84d6SLisandro Dalcin PetscReal *time = bdf->time; 171211a84d6SLisandro Dalcin Vec *vecs = bdf->work; 172211a84d6SLisandro Dalcin PetscScalar alpha[7]; 173211a84d6SLisandro Dalcin 174211a84d6SLisandro Dalcin PetscFunctionBegin; 175211a84d6SLisandro Dalcin LagrangeBasisVals(n, t, time, alpha); 1769566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(X)); 1779566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X, n, alpha, vecs)); 1783ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 179211a84d6SLisandro Dalcin } 180211a84d6SLisandro Dalcin 181e3c11fc1SJed Brown /* Compute the affine term V0 such that Xdot = shift*X + V0. 182e3c11fc1SJed Brown * 183e3c11fc1SJed Brown * When using transient variables, we're computing Cdot = shift*C(X) + V0, and thus choose a linear combination of tvwork. 184e3c11fc1SJed Brown */ 185d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSBDF_PreSolve(TS ts) 186d71ae5a4SJacob Faibussowitsch { 1871117961dSLisandro Dalcin TS_BDF *bdf = (TS_BDF *)ts->data; 1881117961dSLisandro Dalcin PetscInt i, n = PetscMax(bdf->k, 1) + 1; 1891117961dSLisandro Dalcin Vec V, V0; 1901117961dSLisandro Dalcin Vec vecs[7]; 1911117961dSLisandro Dalcin PetscScalar alpha[7]; 1921117961dSLisandro Dalcin 1931117961dSLisandro Dalcin PetscFunctionBegin; 1949566063dSJacob Faibussowitsch PetscCall(TSBDF_GetVecs(ts, NULL, &V, &V0)); 1951117961dSLisandro Dalcin LagrangeBasisDers(n, bdf->time[0], bdf->time, alpha); 196ad540459SPierre Jolivet for (i = 1; i < n; i++) vecs[i] = bdf->transientvar ? bdf->tvwork[i] : bdf->work[i]; 1979566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(V0)); 1989566063dSJacob Faibussowitsch PetscCall(VecMAXPY(V0, n - 1, alpha + 1, vecs + 1)); 1991117961dSLisandro Dalcin bdf->shift = PetscRealPart(alpha[0]); 2009566063dSJacob Faibussowitsch PetscCall(TSBDF_RestoreVecs(ts, NULL, &V, &V0)); 2013ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2021117961dSLisandro Dalcin } 2031117961dSLisandro Dalcin 204d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSBDF_SNESSolve(TS ts, Vec b, Vec x) 205d71ae5a4SJacob Faibussowitsch { 206211a84d6SLisandro Dalcin PetscInt nits, lits; 207211a84d6SLisandro Dalcin 208211a84d6SLisandro Dalcin PetscFunctionBegin; 2099566063dSJacob Faibussowitsch PetscCall(TSBDF_PreSolve(ts)); 2109566063dSJacob Faibussowitsch PetscCall(SNESSolve(ts->snes, b, x)); 2119566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(ts->snes, &nits)); 2129566063dSJacob Faibussowitsch PetscCall(SNESGetLinearSolveIterations(ts->snes, &lits)); 2139371c9d4SSatish Balay ts->snes_its += nits; 2149371c9d4SSatish Balay ts->ksp_its += lits; 2153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 216211a84d6SLisandro Dalcin } 217211a84d6SLisandro Dalcin 218d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSBDF_Restart(TS ts, PetscBool *accept) 219d71ae5a4SJacob Faibussowitsch { 220211a84d6SLisandro Dalcin TS_BDF *bdf = (TS_BDF *)ts->data; 221211a84d6SLisandro Dalcin 222211a84d6SLisandro Dalcin PetscFunctionBegin; 2239371c9d4SSatish Balay bdf->k = 1; 2249371c9d4SSatish Balay bdf->n = 0; 2259566063dSJacob Faibussowitsch PetscCall(TSBDF_Advance(ts, ts->ptime, ts->vec_sol)); 2263db8ccc7SStefano Zampini if (bdf->order == 1) { 2273db8ccc7SStefano Zampini *accept = PETSC_TRUE; 2283db8ccc7SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 2293db8ccc7SStefano Zampini } 230211a84d6SLisandro Dalcin bdf->time[0] = ts->ptime + ts->time_step / 2; 2319566063dSJacob Faibussowitsch PetscCall(VecCopy(bdf->work[1], bdf->work[0])); 2329566063dSJacob Faibussowitsch PetscCall(TSPreStage(ts, bdf->time[0])); 2339566063dSJacob Faibussowitsch PetscCall(TSBDF_SNESSolve(ts, NULL, bdf->work[0])); 2349566063dSJacob Faibussowitsch PetscCall(TSPostStage(ts, bdf->time[0], 0, &bdf->work[0])); 2359566063dSJacob Faibussowitsch PetscCall(TSAdaptCheckStage(ts->adapt, ts, bdf->time[0], bdf->work[0], accept)); 2363ba16761SJacob Faibussowitsch if (!*accept) PetscFunctionReturn(PETSC_SUCCESS); 237211a84d6SLisandro Dalcin 2389371c9d4SSatish Balay bdf->k = PetscMin(2, bdf->order); 2399371c9d4SSatish Balay bdf->n++; 2409566063dSJacob Faibussowitsch PetscCall(VecCopy(bdf->work[0], bdf->work[2])); 241211a84d6SLisandro Dalcin bdf->time[2] = bdf->time[0]; 2429566063dSJacob Faibussowitsch PetscCall(TSComputeTransientVariable(ts, bdf->work[2], bdf->tvwork[2])); 2433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 244211a84d6SLisandro Dalcin } 245211a84d6SLisandro Dalcin 246211a84d6SLisandro Dalcin static const char *const BDF_SchemeName[] = {"", "1", "2", "3", "4", "5", "6"}; 247211a84d6SLisandro Dalcin 248d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSStep_BDF(TS ts) 249d71ae5a4SJacob Faibussowitsch { 250211a84d6SLisandro Dalcin TS_BDF *bdf = (TS_BDF *)ts->data; 251211a84d6SLisandro Dalcin PetscInt rejections = 0; 252211a84d6SLisandro Dalcin PetscBool stageok, accept = PETSC_TRUE; 253211a84d6SLisandro Dalcin PetscReal next_time_step = ts->time_step; 254211a84d6SLisandro Dalcin 255211a84d6SLisandro Dalcin PetscFunctionBegin; 2569566063dSJacob Faibussowitsch PetscCall(PetscCitationsRegister(citation, &cited)); 257211a84d6SLisandro Dalcin 258211a84d6SLisandro Dalcin if (!ts->steprollback && !ts->steprestart) { 259211a84d6SLisandro Dalcin bdf->k = PetscMin(bdf->k + 1, bdf->order); 2609566063dSJacob Faibussowitsch PetscCall(TSBDF_Advance(ts, ts->ptime, ts->vec_sol)); 261211a84d6SLisandro Dalcin } 262211a84d6SLisandro Dalcin 263211a84d6SLisandro Dalcin bdf->status = TS_STEP_INCOMPLETE; 264211a84d6SLisandro Dalcin while (!ts->reason && bdf->status != TS_STEP_COMPLETE) { 265211a84d6SLisandro Dalcin if (ts->steprestart) { 2669566063dSJacob Faibussowitsch PetscCall(TSBDF_Restart(ts, &stageok)); 267211a84d6SLisandro Dalcin if (!stageok) goto reject_step; 268211a84d6SLisandro Dalcin } 269211a84d6SLisandro Dalcin 270211a84d6SLisandro Dalcin bdf->time[0] = ts->ptime + ts->time_step; 27199e85243SStefano Zampini if (bdf->extrapolate) PetscCall(TSBDF_Extrapolate(ts, bdf->k - (accept ? 0 : 1), bdf->time[0], bdf->work[0])); 272e4b02790SStefano Zampini else if (!accept) PetscCall(VecCopy(ts->vec_sol, bdf->work[0])); 2739566063dSJacob Faibussowitsch PetscCall(TSPreStage(ts, bdf->time[0])); 2749566063dSJacob Faibussowitsch PetscCall(TSBDF_SNESSolve(ts, NULL, bdf->work[0])); 2759566063dSJacob Faibussowitsch PetscCall(TSPostStage(ts, bdf->time[0], 0, &bdf->work[0])); 2769566063dSJacob Faibussowitsch PetscCall(TSAdaptCheckStage(ts->adapt, ts, bdf->time[0], bdf->work[0], &stageok)); 277211a84d6SLisandro Dalcin if (!stageok) goto reject_step; 278211a84d6SLisandro Dalcin 279211a84d6SLisandro Dalcin bdf->status = TS_STEP_PENDING; 2809566063dSJacob Faibussowitsch PetscCall(TSAdaptCandidatesClear(ts->adapt)); 2819566063dSJacob Faibussowitsch PetscCall(TSAdaptCandidateAdd(ts->adapt, BDF_SchemeName[bdf->k], bdf->k, 1, 1.0, 1.0, PETSC_TRUE)); 2829566063dSJacob Faibussowitsch PetscCall(TSAdaptChoose(ts->adapt, ts, ts->time_step, NULL, &next_time_step, &accept)); 283211a84d6SLisandro Dalcin bdf->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE; 2849371c9d4SSatish Balay if (!accept) { 2859371c9d4SSatish Balay ts->time_step = next_time_step; 2869371c9d4SSatish Balay goto reject_step; 2879371c9d4SSatish Balay } 288211a84d6SLisandro Dalcin 2899566063dSJacob Faibussowitsch PetscCall(VecCopy(bdf->work[0], ts->vec_sol)); 290211a84d6SLisandro Dalcin ts->ptime += ts->time_step; 291211a84d6SLisandro Dalcin ts->time_step = next_time_step; 292211a84d6SLisandro Dalcin break; 293211a84d6SLisandro Dalcin 294211a84d6SLisandro Dalcin reject_step: 2959371c9d4SSatish Balay ts->reject++; 2969371c9d4SSatish Balay accept = PETSC_FALSE; 297211a84d6SLisandro Dalcin if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) { 29863a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", step rejections %" PetscInt_FMT " greater than current TS allowed, stopping solve\n", ts->steps, rejections)); 299211a84d6SLisandro Dalcin ts->reason = TS_DIVERGED_STEP_REJECTED; 300211a84d6SLisandro Dalcin } 301211a84d6SLisandro Dalcin } 3023ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 303211a84d6SLisandro Dalcin } 304211a84d6SLisandro Dalcin 305d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSInterpolate_BDF(TS ts, PetscReal t, Vec X) 306d71ae5a4SJacob Faibussowitsch { 307211a84d6SLisandro Dalcin TS_BDF *bdf = (TS_BDF *)ts->data; 308211a84d6SLisandro Dalcin 309211a84d6SLisandro Dalcin PetscFunctionBegin; 3109566063dSJacob Faibussowitsch PetscCall(TSBDF_Interpolate(ts, bdf->k, t, X)); 3113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 312211a84d6SLisandro Dalcin } 313211a84d6SLisandro Dalcin 314d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSEvaluateWLTE_BDF(TS ts, NormType wnormtype, PetscInt *order, PetscReal *wlte) 315d71ae5a4SJacob Faibussowitsch { 316211a84d6SLisandro Dalcin TS_BDF *bdf = (TS_BDF *)ts->data; 317211a84d6SLisandro Dalcin PetscInt k = bdf->k; 3187453f775SEmil Constantinescu PetscReal wltea, wlter; 319211a84d6SLisandro Dalcin Vec X = bdf->work[0], Y = bdf->vec_lte; 320211a84d6SLisandro Dalcin 321211a84d6SLisandro Dalcin PetscFunctionBegin; 322211a84d6SLisandro Dalcin k = PetscMin(k, bdf->n - 1); 3233db8ccc7SStefano Zampini if (k == 0) { 3243db8ccc7SStefano Zampini *wlte = -1; 3253db8ccc7SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 3263db8ccc7SStefano Zampini } 3279566063dSJacob Faibussowitsch PetscCall(TSBDF_VecLTE(ts, k, Y)); 3289566063dSJacob Faibussowitsch PetscCall(VecAXPY(Y, 1, X)); 3299566063dSJacob Faibussowitsch PetscCall(TSErrorWeightedNorm(ts, X, Y, wnormtype, wlte, &wltea, &wlter)); 330211a84d6SLisandro Dalcin if (order) *order = k + 1; 3313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 332211a84d6SLisandro Dalcin } 333211a84d6SLisandro Dalcin 33436de76e7SStefano Zampini static PetscErrorCode TSResizeRegister_BDF(TS ts, PetscBool reg) 33536de76e7SStefano Zampini { 33636de76e7SStefano Zampini TS_BDF *bdf = (TS_BDF *)ts->data; 337ecc87898SStefano Zampini const char *names[] = {"", "ts:bdf:1", "ts:bdf:2", "ts:bdf:3", "ts:bdf:4", "ts:bdf:5", "ts:bdf:6", "ts:bdf:7"}; 338ecc87898SStefano Zampini PetscInt i, maxn = PETSC_STATIC_ARRAY_LENGTH(bdf->work); 33936de76e7SStefano Zampini 34036de76e7SStefano Zampini PetscFunctionBegin; 341ecc87898SStefano Zampini PetscAssert(maxn == 8, PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "names need to be redefined"); 34236de76e7SStefano Zampini if (reg) { 3433a7d0413SPierre Jolivet for (i = 1; i < PetscMin(bdf->n + 1, maxn); i++) PetscCall(TSResizeRegisterVec(ts, names[i], bdf->work[i])); 34436de76e7SStefano Zampini } else { 34536de76e7SStefano Zampini for (i = 1; i < maxn; i++) { 34636de76e7SStefano Zampini PetscCall(TSResizeRetrieveVec(ts, names[i], &bdf->work[i])); 34736de76e7SStefano Zampini if (!bdf->work[i]) break; 34836de76e7SStefano Zampini PetscCall(PetscObjectReference((PetscObject)bdf->work[i])); 34936de76e7SStefano Zampini if (bdf->transientvar) { 35036de76e7SStefano Zampini PetscCall(VecDuplicate(bdf->work[i], &bdf->tvwork[i])); 35136de76e7SStefano Zampini PetscCall(TSComputeTransientVariable(ts, bdf->work[i], bdf->tvwork[i])); 35236de76e7SStefano Zampini } 35336de76e7SStefano Zampini } 35436de76e7SStefano Zampini } 35536de76e7SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 35636de76e7SStefano Zampini } 35736de76e7SStefano Zampini 358d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormFunction_BDF(SNES snes, Vec X, Vec F, TS ts) 359d71ae5a4SJacob Faibussowitsch { 360211a84d6SLisandro Dalcin TS_BDF *bdf = (TS_BDF *)ts->data; 3611117961dSLisandro Dalcin DM dm, dmsave = ts->dm; 362211a84d6SLisandro Dalcin PetscReal t = bdf->time[0]; 3631117961dSLisandro Dalcin PetscReal shift = bdf->shift; 3641117961dSLisandro Dalcin Vec V, V0; 365211a84d6SLisandro Dalcin 366211a84d6SLisandro Dalcin PetscFunctionBegin; 3679566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 3689566063dSJacob Faibussowitsch PetscCall(TSBDF_GetVecs(ts, dm, &V, &V0)); 369e3c11fc1SJed Brown if (bdf->transientvar) { /* shift*C(X) + V0 */ 3709566063dSJacob Faibussowitsch PetscCall(TSComputeTransientVariable(ts, X, V)); 3719566063dSJacob Faibussowitsch PetscCall(VecAYPX(V, shift, V0)); 372e3c11fc1SJed Brown } else { /* shift*X + V0 */ 3739566063dSJacob Faibussowitsch PetscCall(VecWAXPY(V, shift, X, V0)); 374e3c11fc1SJed Brown } 3751117961dSLisandro Dalcin 376211a84d6SLisandro Dalcin /* F = Function(t,X,V) */ 3771117961dSLisandro Dalcin ts->dm = dm; 3789566063dSJacob Faibussowitsch PetscCall(TSComputeIFunction(ts, t, X, V, F, PETSC_FALSE)); 3791117961dSLisandro Dalcin ts->dm = dmsave; 3801117961dSLisandro Dalcin 3819566063dSJacob Faibussowitsch PetscCall(TSBDF_RestoreVecs(ts, dm, &V, &V0)); 3823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 383211a84d6SLisandro Dalcin } 384211a84d6SLisandro Dalcin 385d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormJacobian_BDF(SNES snes, Vec X, Mat J, Mat P, TS ts) 386d71ae5a4SJacob Faibussowitsch { 387211a84d6SLisandro Dalcin TS_BDF *bdf = (TS_BDF *)ts->data; 3881117961dSLisandro Dalcin DM dm, dmsave = ts->dm; 389211a84d6SLisandro Dalcin PetscReal t = bdf->time[0]; 3901117961dSLisandro Dalcin PetscReal shift = bdf->shift; 3911117961dSLisandro Dalcin Vec V, V0; 392211a84d6SLisandro Dalcin 393211a84d6SLisandro Dalcin PetscFunctionBegin; 3949566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 3959566063dSJacob Faibussowitsch PetscCall(TSBDF_GetVecs(ts, dm, &V, &V0)); 3961117961dSLisandro Dalcin 397211a84d6SLisandro Dalcin /* J,P = Jacobian(t,X,V) */ 3981117961dSLisandro Dalcin ts->dm = dm; 3999566063dSJacob Faibussowitsch PetscCall(TSComputeIJacobian(ts, t, X, V, shift, J, P, PETSC_FALSE)); 4001117961dSLisandro Dalcin ts->dm = dmsave; 4011117961dSLisandro Dalcin 4029566063dSJacob Faibussowitsch PetscCall(TSBDF_RestoreVecs(ts, dm, &V, &V0)); 4033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 404211a84d6SLisandro Dalcin } 405211a84d6SLisandro Dalcin 406d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSReset_BDF(TS ts) 407d71ae5a4SJacob Faibussowitsch { 408211a84d6SLisandro Dalcin TS_BDF *bdf = (TS_BDF *)ts->data; 409400f6f02SBarry Smith size_t i, n = PETSC_STATIC_ARRAY_LENGTH(bdf->work); 410211a84d6SLisandro Dalcin 411211a84d6SLisandro Dalcin PetscFunctionBegin; 412e3c11fc1SJed Brown for (i = 0; i < n; i++) { 4139566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bdf->work[i])); 4149566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bdf->tvwork[i])); 415e3c11fc1SJed Brown } 4169566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bdf->vec_dot)); 4179566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bdf->vec_wrk)); 4189566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bdf->vec_lte)); 4199566063dSJacob Faibussowitsch if (ts->dm) PetscCall(DMCoarsenHookRemove(ts->dm, DMCoarsenHook_TSBDF, DMRestrictHook_TSBDF, ts)); 4203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 421211a84d6SLisandro Dalcin } 422211a84d6SLisandro Dalcin 423d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSDestroy_BDF(TS ts) 424d71ae5a4SJacob Faibussowitsch { 425211a84d6SLisandro Dalcin PetscFunctionBegin; 4269566063dSJacob Faibussowitsch PetscCall(TSReset_BDF(ts)); 4279566063dSJacob Faibussowitsch PetscCall(PetscFree(ts->data)); 4289566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSBDFSetOrder_C", NULL)); 4299566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSBDFGetOrder_C", NULL)); 4303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 431211a84d6SLisandro Dalcin } 432211a84d6SLisandro Dalcin 433d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetUp_BDF(TS ts) 434d71ae5a4SJacob Faibussowitsch { 435211a84d6SLisandro Dalcin TS_BDF *bdf = (TS_BDF *)ts->data; 436400f6f02SBarry Smith size_t n = PETSC_STATIC_ARRAY_LENGTH(bdf->work); 4372ffb9264SLisandro Dalcin PetscReal low, high, two = 2; 43836de76e7SStefano Zampini PetscInt cnt = 0; 439211a84d6SLisandro Dalcin 440211a84d6SLisandro Dalcin PetscFunctionBegin; 4419566063dSJacob Faibussowitsch PetscCall(TSHasTransientVariable(ts, &bdf->transientvar)); 44236de76e7SStefano Zampini for (size_t i = 0; i < n; i++) { 44336de76e7SStefano Zampini if (!bdf->work[i]) PetscCall(VecDuplicate(ts->vec_sol, &bdf->work[i])); 44436de76e7SStefano Zampini else cnt++; 44536de76e7SStefano Zampini if (i && bdf->transientvar && !bdf->tvwork[i]) PetscCall(VecDuplicate(ts->vec_sol, &bdf->tvwork[i])); 446e3c11fc1SJed Brown } 44736de76e7SStefano Zampini if (!cnt) bdf->k = bdf->n = 0; 4489566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &bdf->vec_dot)); 4499566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &bdf->vec_wrk)); 4509566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &bdf->vec_lte)); 4519566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &ts->dm)); 4529566063dSJacob Faibussowitsch PetscCall(DMCoarsenHookAdd(ts->dm, DMCoarsenHook_TSBDF, DMRestrictHook_TSBDF, ts)); 453211a84d6SLisandro Dalcin 4549566063dSJacob Faibussowitsch PetscCall(TSGetAdapt(ts, &ts->adapt)); 4559566063dSJacob Faibussowitsch PetscCall(TSAdaptCandidatesClear(ts->adapt)); 4569566063dSJacob Faibussowitsch PetscCall(TSAdaptGetClip(ts->adapt, &low, &high)); 4579566063dSJacob Faibussowitsch PetscCall(TSAdaptSetClip(ts->adapt, low, PetscMin(high, two))); 458211a84d6SLisandro Dalcin 4599566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts, &ts->snes)); 4603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 461211a84d6SLisandro Dalcin } 462211a84d6SLisandro Dalcin 463ce78bad3SBarry Smith static PetscErrorCode TSSetFromOptions_BDF(TS ts, PetscOptionItems PetscOptionsObject) 464d71ae5a4SJacob Faibussowitsch { 46599e85243SStefano Zampini TS_BDF *bdf = (TS_BDF *)ts->data; 46699e85243SStefano Zampini 467211a84d6SLisandro Dalcin PetscFunctionBegin; 468d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "BDF ODE solver options"); 469211a84d6SLisandro Dalcin { 470211a84d6SLisandro Dalcin PetscBool flg; 471e5b8ffdfSLisandro Dalcin PetscInt order; 4729566063dSJacob Faibussowitsch PetscCall(TSBDFGetOrder(ts, &order)); 4739566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-ts_bdf_order", "Order of the BDF method", "TSBDFSetOrder", order, &order, &flg)); 4749566063dSJacob Faibussowitsch if (flg) PetscCall(TSBDFSetOrder(ts, order)); 47599e85243SStefano Zampini PetscCall(PetscOptionsBool("-ts_bdf_initial_guess_extrapolate", "Extrapolate the initial guess of the nonlinear solve from previous time steps", "", bdf->extrapolate, &bdf->extrapolate, NULL)); 476211a84d6SLisandro Dalcin } 477d0609cedSBarry Smith PetscOptionsHeadEnd(); 4783ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 479211a84d6SLisandro Dalcin } 480211a84d6SLisandro Dalcin 481d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSView_BDF(TS ts, PetscViewer viewer) 482d71ae5a4SJacob Faibussowitsch { 483211a84d6SLisandro Dalcin TS_BDF *bdf = (TS_BDF *)ts->data; 4849f196a02SMartin Diehl PetscBool isascii; 485211a84d6SLisandro Dalcin 486211a84d6SLisandro Dalcin PetscFunctionBegin; 4879f196a02SMartin Diehl PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii)); 4889f196a02SMartin Diehl if (isascii) PetscCall(PetscViewerASCIIPrintf(viewer, " Order=%" PetscInt_FMT "\n", bdf->order)); 4893ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 490211a84d6SLisandro Dalcin } 491211a84d6SLisandro Dalcin 492d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSBDFSetOrder_BDF(TS ts, PetscInt order) 493d71ae5a4SJacob Faibussowitsch { 494211a84d6SLisandro Dalcin TS_BDF *bdf = (TS_BDF *)ts->data; 495211a84d6SLisandro Dalcin 496211a84d6SLisandro Dalcin PetscFunctionBegin; 4973ba16761SJacob Faibussowitsch if (order == bdf->order) PetscFunctionReturn(PETSC_SUCCESS); 49863a3b9bcSJacob Faibussowitsch PetscCheck(order >= 1 && order <= 6, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_OUTOFRANGE, "BDF Order %" PetscInt_FMT " not implemented", order); 499211a84d6SLisandro Dalcin bdf->order = order; 5003ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 501211a84d6SLisandro Dalcin } 502211a84d6SLisandro Dalcin 503d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSBDFGetOrder_BDF(TS ts, PetscInt *order) 504d71ae5a4SJacob Faibussowitsch { 505211a84d6SLisandro Dalcin TS_BDF *bdf = (TS_BDF *)ts->data; 506211a84d6SLisandro Dalcin 507211a84d6SLisandro Dalcin PetscFunctionBegin; 508211a84d6SLisandro Dalcin *order = bdf->order; 5093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 510211a84d6SLisandro Dalcin } 511211a84d6SLisandro Dalcin 512211a84d6SLisandro Dalcin /*MC 513efa39862SBarry Smith TSBDF - DAE solver using implicit backward differentiation formula (BDF) methods suitable for stiff ODEs. 514211a84d6SLisandro Dalcin 515211a84d6SLisandro Dalcin Level: beginner 516211a84d6SLisandro Dalcin 517efa39862SBarry Smith Options Database Keys: 518efa39862SBarry Smith + -ts_bdf_order <n> - Order of the BDF method 519efa39862SBarry Smith - -ts_bdf_initial_guess_extrapolate <true, false> - Extrapolate the initial guess of the nonlinear solve from previous time steps, defaults to true 520efa39862SBarry Smith 521efa39862SBarry Smith .seealso: [](ch_ts), `TS`, `TSCreate()`, `TSSetType()`, `TSType`, `TSBDFSetOrder()` 522211a84d6SLisandro Dalcin M*/ 523d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSCreate_BDF(TS ts) 524d71ae5a4SJacob Faibussowitsch { 525211a84d6SLisandro Dalcin TS_BDF *bdf; 526211a84d6SLisandro Dalcin 527211a84d6SLisandro Dalcin PetscFunctionBegin; 528211a84d6SLisandro Dalcin ts->ops->reset = TSReset_BDF; 529211a84d6SLisandro Dalcin ts->ops->destroy = TSDestroy_BDF; 530211a84d6SLisandro Dalcin ts->ops->view = TSView_BDF; 531211a84d6SLisandro Dalcin ts->ops->setup = TSSetUp_BDF; 532211a84d6SLisandro Dalcin ts->ops->setfromoptions = TSSetFromOptions_BDF; 533211a84d6SLisandro Dalcin ts->ops->step = TSStep_BDF; 534211a84d6SLisandro Dalcin ts->ops->evaluatewlte = TSEvaluateWLTE_BDF; 535211a84d6SLisandro Dalcin ts->ops->interpolate = TSInterpolate_BDF; 53636de76e7SStefano Zampini ts->ops->resizeregister = TSResizeRegister_BDF; 537211a84d6SLisandro Dalcin ts->ops->snesfunction = SNESTSFormFunction_BDF; 538211a84d6SLisandro Dalcin ts->ops->snesjacobian = SNESTSFormJacobian_BDF; 5392ffb9264SLisandro Dalcin ts->default_adapt_type = TSADAPTBASIC; 540211a84d6SLisandro Dalcin 541efd4aadfSBarry Smith ts->usessnes = PETSC_TRUE; 542efd4aadfSBarry Smith 5434dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&bdf)); 544211a84d6SLisandro Dalcin ts->data = (void *)bdf; 545211a84d6SLisandro Dalcin 54699e85243SStefano Zampini bdf->extrapolate = PETSC_TRUE; 547211a84d6SLisandro Dalcin bdf->status = TS_STEP_COMPLETE; 548ac530a7eSPierre Jolivet for (size_t i = 0; i < PETSC_STATIC_ARRAY_LENGTH(bdf->work); i++) bdf->work[i] = bdf->tvwork[i] = NULL; 549211a84d6SLisandro Dalcin 5509566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSBDFSetOrder_C", TSBDFSetOrder_BDF)); 5519566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSBDFGetOrder_C", TSBDFGetOrder_BDF)); 5529566063dSJacob Faibussowitsch PetscCall(TSBDFSetOrder(ts, 2)); 5533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 554211a84d6SLisandro Dalcin } 555211a84d6SLisandro Dalcin 556211a84d6SLisandro Dalcin /*@ 557bcf0153eSBarry Smith TSBDFSetOrder - Set the order of the `TSBDF` method 558211a84d6SLisandro Dalcin 559c3339decSBarry Smith Logically Collective 560211a84d6SLisandro Dalcin 561d8d19677SJose E. Roman Input Parameters: 562211a84d6SLisandro Dalcin + ts - timestepping context 563211a84d6SLisandro Dalcin - order - order of the method 564211a84d6SLisandro Dalcin 565bcf0153eSBarry Smith Options Database Key: 566147403d9SBarry Smith . -ts_bdf_order <order> - select the order 567211a84d6SLisandro Dalcin 568211a84d6SLisandro Dalcin Level: intermediate 569211a84d6SLisandro Dalcin 570bcf0153eSBarry Smith .seealso: `TSBDFGetOrder()`, `TS`, `TSBDF` 571211a84d6SLisandro Dalcin @*/ 572d71ae5a4SJacob Faibussowitsch PetscErrorCode TSBDFSetOrder(TS ts, PetscInt order) 573d71ae5a4SJacob Faibussowitsch { 574211a84d6SLisandro Dalcin PetscFunctionBegin; 575211a84d6SLisandro Dalcin PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 576211a84d6SLisandro Dalcin PetscValidLogicalCollectiveInt(ts, order, 2); 577cac4c232SBarry Smith PetscTryMethod(ts, "TSBDFSetOrder_C", (TS, PetscInt), (ts, order)); 5783ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 579211a84d6SLisandro Dalcin } 580211a84d6SLisandro Dalcin 581211a84d6SLisandro Dalcin /*@ 582bcf0153eSBarry Smith TSBDFGetOrder - Get the order of the `TSBDF` method 583211a84d6SLisandro Dalcin 584211a84d6SLisandro Dalcin Not Collective 585211a84d6SLisandro Dalcin 586211a84d6SLisandro Dalcin Input Parameter: 587211a84d6SLisandro Dalcin . ts - timestepping context 588211a84d6SLisandro Dalcin 589211a84d6SLisandro Dalcin Output Parameter: 590211a84d6SLisandro Dalcin . order - order of the method 591211a84d6SLisandro Dalcin 592211a84d6SLisandro Dalcin Level: intermediate 593211a84d6SLisandro Dalcin 594bcf0153eSBarry Smith .seealso: `TSBDFSetOrder()`, `TS`, `TSBDF` 595211a84d6SLisandro Dalcin @*/ 596d71ae5a4SJacob Faibussowitsch PetscErrorCode TSBDFGetOrder(TS ts, PetscInt *order) 597d71ae5a4SJacob Faibussowitsch { 598211a84d6SLisandro Dalcin PetscFunctionBegin; 599211a84d6SLisandro Dalcin PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 6004f572ea9SToby Isaac PetscAssertPointer(order, 2); 601cac4c232SBarry Smith PetscUseMethod(ts, "TSBDFGetOrder_C", (TS, PetscInt *), (ts, order)); 6023ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 603211a84d6SLisandro Dalcin } 604