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; 26211a84d6SLisandro Dalcin PetscInt order; 27211a84d6SLisandro Dalcin TSStepStatus status; 28211a84d6SLisandro Dalcin } TS_BDF; 29211a84d6SLisandro Dalcin 30e3c11fc1SJed Brown /* Compute Lagrange polynomials on T[:n] evaluated at t. 31e3c11fc1SJed Brown * If one has data (T[i], Y[i]), then the interpolation/extrapolation is f(t) = \sum_i L[i]*Y[i]. 32e3c11fc1SJed Brown */ 33d71ae5a4SJacob Faibussowitsch static inline void LagrangeBasisVals(PetscInt n, PetscReal t, const PetscReal T[], PetscScalar L[]) 34d71ae5a4SJacob Faibussowitsch { 35211a84d6SLisandro Dalcin PetscInt k, j; 36211a84d6SLisandro Dalcin for (k = 0; k < n; k++) 37211a84d6SLisandro Dalcin for (L[k] = 1, j = 0; j < n; j++) 389371c9d4SSatish Balay if (j != k) L[k] *= (t - T[j]) / (T[k] - T[j]); 39211a84d6SLisandro Dalcin } 40211a84d6SLisandro Dalcin 41d71ae5a4SJacob Faibussowitsch static inline void LagrangeBasisDers(PetscInt n, PetscReal t, const PetscReal T[], PetscScalar dL[]) 42d71ae5a4SJacob Faibussowitsch { 43211a84d6SLisandro Dalcin PetscInt k, j, i; 44211a84d6SLisandro Dalcin for (k = 0; k < n; k++) 45211a84d6SLisandro Dalcin for (dL[k] = 0, j = 0; j < n; j++) 46211a84d6SLisandro Dalcin if (j != k) { 47211a84d6SLisandro Dalcin PetscReal L = 1 / (T[k] - T[j]); 48211a84d6SLisandro Dalcin for (i = 0; i < n; i++) 499371c9d4SSatish Balay if (i != j && i != k) L *= (t - T[i]) / (T[k] - T[i]); 50211a84d6SLisandro Dalcin dL[k] += L; 51211a84d6SLisandro Dalcin } 52211a84d6SLisandro Dalcin } 53211a84d6SLisandro Dalcin 54d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSBDF_GetVecs(TS ts, DM dm, Vec *Xdot, Vec *Ydot) 55d71ae5a4SJacob Faibussowitsch { 561117961dSLisandro Dalcin TS_BDF *bdf = (TS_BDF *)ts->data; 571117961dSLisandro Dalcin 581117961dSLisandro Dalcin PetscFunctionBegin; 591117961dSLisandro Dalcin if (dm && dm != ts->dm) { 609566063dSJacob Faibussowitsch PetscCall(DMGetNamedGlobalVector(dm, "TSBDF_Vec_Xdot", Xdot)); 619566063dSJacob Faibussowitsch PetscCall(DMGetNamedGlobalVector(dm, "TSBDF_Vec_Ydot", Ydot)); 621117961dSLisandro Dalcin } else { 631117961dSLisandro Dalcin *Xdot = bdf->vec_dot; 641117961dSLisandro Dalcin *Ydot = bdf->vec_wrk; 651117961dSLisandro Dalcin } 663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 671117961dSLisandro Dalcin } 681117961dSLisandro Dalcin 69d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSBDF_RestoreVecs(TS ts, DM dm, Vec *Xdot, Vec *Ydot) 70d71ae5a4SJacob Faibussowitsch { 711117961dSLisandro Dalcin TS_BDF *bdf = (TS_BDF *)ts->data; 721117961dSLisandro Dalcin 731117961dSLisandro Dalcin PetscFunctionBegin; 741117961dSLisandro Dalcin if (dm && dm != ts->dm) { 759566063dSJacob Faibussowitsch PetscCall(DMRestoreNamedGlobalVector(dm, "TSBDF_Vec_Xdot", Xdot)); 769566063dSJacob Faibussowitsch PetscCall(DMRestoreNamedGlobalVector(dm, "TSBDF_Vec_Ydot", Ydot)); 771117961dSLisandro Dalcin } else { 783c633725SBarry Smith PetscCheck(*Xdot == bdf->vec_dot, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_INCOMP, "Vec does not match the cache"); 793c633725SBarry Smith PetscCheck(*Ydot == bdf->vec_wrk, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_INCOMP, "Vec does not match the cache"); 801117961dSLisandro Dalcin *Xdot = NULL; 811117961dSLisandro Dalcin *Ydot = NULL; 821117961dSLisandro Dalcin } 833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 841117961dSLisandro Dalcin } 851117961dSLisandro Dalcin 86d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMCoarsenHook_TSBDF(DM fine, DM coarse, void *ctx) 87d71ae5a4SJacob Faibussowitsch { 881117961dSLisandro Dalcin PetscFunctionBegin; 893ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 901117961dSLisandro Dalcin } 911117961dSLisandro Dalcin 92d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMRestrictHook_TSBDF(DM fine, Mat restrct, Vec rscale, Mat inject, DM coarse, void *ctx) 93d71ae5a4SJacob Faibussowitsch { 941117961dSLisandro Dalcin TS ts = (TS)ctx; 951117961dSLisandro Dalcin Vec Ydot, Ydot_c; 961117961dSLisandro Dalcin Vec Xdot, Xdot_c; 971117961dSLisandro Dalcin 981117961dSLisandro Dalcin PetscFunctionBegin; 999566063dSJacob Faibussowitsch PetscCall(TSBDF_GetVecs(ts, fine, &Xdot, &Ydot)); 1009566063dSJacob Faibussowitsch PetscCall(TSBDF_GetVecs(ts, coarse, &Xdot_c, &Ydot_c)); 1011117961dSLisandro Dalcin 1029566063dSJacob Faibussowitsch PetscCall(MatRestrict(restrct, Ydot, Ydot_c)); 1039566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(Ydot_c, rscale, Ydot_c)); 1041117961dSLisandro Dalcin 1059566063dSJacob Faibussowitsch PetscCall(TSBDF_RestoreVecs(ts, fine, &Xdot, &Ydot)); 1069566063dSJacob Faibussowitsch PetscCall(TSBDF_RestoreVecs(ts, coarse, &Xdot_c, &Ydot_c)); 1073ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1081117961dSLisandro Dalcin } 1091117961dSLisandro Dalcin 110d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSBDF_Advance(TS ts, PetscReal t, Vec X) 111d71ae5a4SJacob Faibussowitsch { 112211a84d6SLisandro Dalcin TS_BDF *bdf = (TS_BDF *)ts->data; 113400f6f02SBarry Smith PetscInt i, n = PETSC_STATIC_ARRAY_LENGTH(bdf->work); 114e3c11fc1SJed Brown Vec tail = bdf->work[n - 1], tvtail = bdf->tvwork[n - 1]; 115211a84d6SLisandro Dalcin 116211a84d6SLisandro Dalcin PetscFunctionBegin; 117211a84d6SLisandro Dalcin for (i = n - 1; i >= 2; i--) { 118211a84d6SLisandro Dalcin bdf->time[i] = bdf->time[i - 1]; 119211a84d6SLisandro Dalcin bdf->work[i] = bdf->work[i - 1]; 120e3c11fc1SJed Brown bdf->tvwork[i] = bdf->tvwork[i - 1]; 121211a84d6SLisandro Dalcin } 122211a84d6SLisandro Dalcin bdf->n = PetscMin(bdf->n + 1, n - 1); 123211a84d6SLisandro Dalcin bdf->time[1] = t; 124211a84d6SLisandro Dalcin bdf->work[1] = tail; 125e3c11fc1SJed Brown bdf->tvwork[1] = tvtail; 1269566063dSJacob Faibussowitsch PetscCall(VecCopy(X, tail)); 1279566063dSJacob Faibussowitsch PetscCall(TSComputeTransientVariable(ts, tail, tvtail)); 1283ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 129211a84d6SLisandro Dalcin } 130211a84d6SLisandro Dalcin 131d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSBDF_VecLTE(TS ts, PetscInt order, Vec lte) 132d71ae5a4SJacob Faibussowitsch { 133211a84d6SLisandro Dalcin TS_BDF *bdf = (TS_BDF *)ts->data; 134211a84d6SLisandro Dalcin PetscInt i, n = order + 1; 135211a84d6SLisandro Dalcin PetscReal *time = bdf->time; 136211a84d6SLisandro Dalcin Vec *vecs = bdf->work; 137211a84d6SLisandro Dalcin PetscScalar a[8], b[8], alpha[8]; 138211a84d6SLisandro Dalcin 139211a84d6SLisandro Dalcin PetscFunctionBegin; 1409371c9d4SSatish Balay LagrangeBasisDers(n + 0, time[0], time, a); 1419371c9d4SSatish Balay a[n] = 0; 142211a84d6SLisandro Dalcin LagrangeBasisDers(n + 1, time[0], time, b); 143211a84d6SLisandro Dalcin for (i = 0; i < n + 1; i++) alpha[i] = (a[i] - b[i]) / a[0]; 1449566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(lte)); 1459566063dSJacob Faibussowitsch PetscCall(VecMAXPY(lte, n + 1, alpha, vecs)); 1463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 147211a84d6SLisandro Dalcin } 148211a84d6SLisandro Dalcin 149d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSBDF_Extrapolate(TS ts, PetscInt order, PetscReal t, Vec X) 150d71ae5a4SJacob Faibussowitsch { 151211a84d6SLisandro Dalcin TS_BDF *bdf = (TS_BDF *)ts->data; 152211a84d6SLisandro Dalcin PetscInt n = order + 1; 153211a84d6SLisandro Dalcin PetscReal *time = bdf->time + 1; 154211a84d6SLisandro Dalcin Vec *vecs = bdf->work + 1; 155211a84d6SLisandro Dalcin PetscScalar alpha[7]; 156211a84d6SLisandro Dalcin 157211a84d6SLisandro Dalcin PetscFunctionBegin; 158211a84d6SLisandro Dalcin n = PetscMin(n, bdf->n); 159211a84d6SLisandro Dalcin LagrangeBasisVals(n, t, time, alpha); 1609566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(X)); 1619566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X, n, alpha, vecs)); 1623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 163211a84d6SLisandro Dalcin } 164211a84d6SLisandro Dalcin 165d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSBDF_Interpolate(TS ts, PetscInt order, PetscReal t, Vec X) 166d71ae5a4SJacob Faibussowitsch { 167211a84d6SLisandro Dalcin TS_BDF *bdf = (TS_BDF *)ts->data; 168211a84d6SLisandro Dalcin PetscInt n = order + 1; 169211a84d6SLisandro Dalcin PetscReal *time = bdf->time; 170211a84d6SLisandro Dalcin Vec *vecs = bdf->work; 171211a84d6SLisandro Dalcin PetscScalar alpha[7]; 172211a84d6SLisandro Dalcin 173211a84d6SLisandro Dalcin PetscFunctionBegin; 174211a84d6SLisandro Dalcin LagrangeBasisVals(n, t, time, alpha); 1759566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(X)); 1769566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X, n, alpha, vecs)); 1773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 178211a84d6SLisandro Dalcin } 179211a84d6SLisandro Dalcin 180e3c11fc1SJed Brown /* Compute the affine term V0 such that Xdot = shift*X + V0. 181e3c11fc1SJed Brown * 182e3c11fc1SJed Brown * When using transient variables, we're computing Cdot = shift*C(X) + V0, and thus choose a linear combination of tvwork. 183e3c11fc1SJed Brown */ 184d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSBDF_PreSolve(TS ts) 185d71ae5a4SJacob Faibussowitsch { 1861117961dSLisandro Dalcin TS_BDF *bdf = (TS_BDF *)ts->data; 1871117961dSLisandro Dalcin PetscInt i, n = PetscMax(bdf->k, 1) + 1; 1881117961dSLisandro Dalcin Vec V, V0; 1891117961dSLisandro Dalcin Vec vecs[7]; 1901117961dSLisandro Dalcin PetscScalar alpha[7]; 1911117961dSLisandro Dalcin 1921117961dSLisandro Dalcin PetscFunctionBegin; 1939566063dSJacob Faibussowitsch PetscCall(TSBDF_GetVecs(ts, NULL, &V, &V0)); 1941117961dSLisandro Dalcin LagrangeBasisDers(n, bdf->time[0], bdf->time, alpha); 195ad540459SPierre Jolivet for (i = 1; i < n; i++) vecs[i] = bdf->transientvar ? bdf->tvwork[i] : bdf->work[i]; 1969566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(V0)); 1979566063dSJacob Faibussowitsch PetscCall(VecMAXPY(V0, n - 1, alpha + 1, vecs + 1)); 1981117961dSLisandro Dalcin bdf->shift = PetscRealPart(alpha[0]); 1999566063dSJacob Faibussowitsch PetscCall(TSBDF_RestoreVecs(ts, NULL, &V, &V0)); 2003ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2011117961dSLisandro Dalcin } 2021117961dSLisandro Dalcin 203d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSBDF_SNESSolve(TS ts, Vec b, Vec x) 204d71ae5a4SJacob Faibussowitsch { 205211a84d6SLisandro Dalcin PetscInt nits, lits; 206211a84d6SLisandro Dalcin 207211a84d6SLisandro Dalcin PetscFunctionBegin; 2089566063dSJacob Faibussowitsch PetscCall(TSBDF_PreSolve(ts)); 2099566063dSJacob Faibussowitsch PetscCall(SNESSolve(ts->snes, b, x)); 2109566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(ts->snes, &nits)); 2119566063dSJacob Faibussowitsch PetscCall(SNESGetLinearSolveIterations(ts->snes, &lits)); 2129371c9d4SSatish Balay ts->snes_its += nits; 2139371c9d4SSatish Balay ts->ksp_its += lits; 2143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 215211a84d6SLisandro Dalcin } 216211a84d6SLisandro Dalcin 217d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSBDF_Restart(TS ts, PetscBool *accept) 218d71ae5a4SJacob Faibussowitsch { 219211a84d6SLisandro Dalcin TS_BDF *bdf = (TS_BDF *)ts->data; 220211a84d6SLisandro Dalcin 221211a84d6SLisandro Dalcin PetscFunctionBegin; 2229371c9d4SSatish Balay bdf->k = 1; 2239371c9d4SSatish Balay bdf->n = 0; 2249566063dSJacob Faibussowitsch PetscCall(TSBDF_Advance(ts, ts->ptime, ts->vec_sol)); 225*3db8ccc7SStefano Zampini if (bdf->order == 1) { 226*3db8ccc7SStefano Zampini *accept = PETSC_TRUE; 227*3db8ccc7SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 228*3db8ccc7SStefano Zampini } 229211a84d6SLisandro Dalcin bdf->time[0] = ts->ptime + ts->time_step / 2; 2309566063dSJacob Faibussowitsch PetscCall(VecCopy(bdf->work[1], bdf->work[0])); 2319566063dSJacob Faibussowitsch PetscCall(TSPreStage(ts, bdf->time[0])); 2329566063dSJacob Faibussowitsch PetscCall(TSBDF_SNESSolve(ts, NULL, bdf->work[0])); 2339566063dSJacob Faibussowitsch PetscCall(TSPostStage(ts, bdf->time[0], 0, &bdf->work[0])); 2349566063dSJacob Faibussowitsch PetscCall(TSAdaptCheckStage(ts->adapt, ts, bdf->time[0], bdf->work[0], accept)); 2353ba16761SJacob Faibussowitsch if (!*accept) PetscFunctionReturn(PETSC_SUCCESS); 236211a84d6SLisandro Dalcin 2379371c9d4SSatish Balay bdf->k = PetscMin(2, bdf->order); 2389371c9d4SSatish Balay bdf->n++; 2399566063dSJacob Faibussowitsch PetscCall(VecCopy(bdf->work[0], bdf->work[2])); 240211a84d6SLisandro Dalcin bdf->time[2] = bdf->time[0]; 2419566063dSJacob Faibussowitsch PetscCall(TSComputeTransientVariable(ts, bdf->work[2], bdf->tvwork[2])); 2423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 243211a84d6SLisandro Dalcin } 244211a84d6SLisandro Dalcin 245211a84d6SLisandro Dalcin static const char *const BDF_SchemeName[] = {"", "1", "2", "3", "4", "5", "6"}; 246211a84d6SLisandro Dalcin 247d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSStep_BDF(TS ts) 248d71ae5a4SJacob Faibussowitsch { 249211a84d6SLisandro Dalcin TS_BDF *bdf = (TS_BDF *)ts->data; 250211a84d6SLisandro Dalcin PetscInt rejections = 0; 251211a84d6SLisandro Dalcin PetscBool stageok, accept = PETSC_TRUE; 252211a84d6SLisandro Dalcin PetscReal next_time_step = ts->time_step; 253211a84d6SLisandro Dalcin 254211a84d6SLisandro Dalcin PetscFunctionBegin; 2559566063dSJacob Faibussowitsch PetscCall(PetscCitationsRegister(citation, &cited)); 256211a84d6SLisandro Dalcin 257211a84d6SLisandro Dalcin if (!ts->steprollback && !ts->steprestart) { 258211a84d6SLisandro Dalcin bdf->k = PetscMin(bdf->k + 1, bdf->order); 2599566063dSJacob Faibussowitsch PetscCall(TSBDF_Advance(ts, ts->ptime, ts->vec_sol)); 260211a84d6SLisandro Dalcin } 261211a84d6SLisandro Dalcin 262211a84d6SLisandro Dalcin bdf->status = TS_STEP_INCOMPLETE; 263211a84d6SLisandro Dalcin while (!ts->reason && bdf->status != TS_STEP_COMPLETE) { 264211a84d6SLisandro Dalcin if (ts->steprestart) { 2659566063dSJacob Faibussowitsch PetscCall(TSBDF_Restart(ts, &stageok)); 266211a84d6SLisandro Dalcin if (!stageok) goto reject_step; 267211a84d6SLisandro Dalcin } 268211a84d6SLisandro Dalcin 269211a84d6SLisandro Dalcin bdf->time[0] = ts->ptime + ts->time_step; 2709566063dSJacob Faibussowitsch PetscCall(TSBDF_Extrapolate(ts, bdf->k - (accept ? 0 : 1), bdf->time[0], bdf->work[0])); 2719566063dSJacob Faibussowitsch PetscCall(TSPreStage(ts, bdf->time[0])); 2729566063dSJacob Faibussowitsch PetscCall(TSBDF_SNESSolve(ts, NULL, bdf->work[0])); 2739566063dSJacob Faibussowitsch PetscCall(TSPostStage(ts, bdf->time[0], 0, &bdf->work[0])); 2749566063dSJacob Faibussowitsch PetscCall(TSAdaptCheckStage(ts->adapt, ts, bdf->time[0], bdf->work[0], &stageok)); 275211a84d6SLisandro Dalcin if (!stageok) goto reject_step; 276211a84d6SLisandro Dalcin 277211a84d6SLisandro Dalcin bdf->status = TS_STEP_PENDING; 2789566063dSJacob Faibussowitsch PetscCall(TSAdaptCandidatesClear(ts->adapt)); 2799566063dSJacob Faibussowitsch PetscCall(TSAdaptCandidateAdd(ts->adapt, BDF_SchemeName[bdf->k], bdf->k, 1, 1.0, 1.0, PETSC_TRUE)); 2809566063dSJacob Faibussowitsch PetscCall(TSAdaptChoose(ts->adapt, ts, ts->time_step, NULL, &next_time_step, &accept)); 281211a84d6SLisandro Dalcin bdf->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE; 2829371c9d4SSatish Balay if (!accept) { 2839371c9d4SSatish Balay ts->time_step = next_time_step; 2849371c9d4SSatish Balay goto reject_step; 2859371c9d4SSatish Balay } 286211a84d6SLisandro Dalcin 2879566063dSJacob Faibussowitsch PetscCall(VecCopy(bdf->work[0], ts->vec_sol)); 288211a84d6SLisandro Dalcin ts->ptime += ts->time_step; 289211a84d6SLisandro Dalcin ts->time_step = next_time_step; 290211a84d6SLisandro Dalcin break; 291211a84d6SLisandro Dalcin 292211a84d6SLisandro Dalcin reject_step: 2939371c9d4SSatish Balay ts->reject++; 2949371c9d4SSatish Balay accept = PETSC_FALSE; 295211a84d6SLisandro Dalcin if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) { 29663a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", step rejections %" PetscInt_FMT " greater than current TS allowed, stopping solve\n", ts->steps, rejections)); 297211a84d6SLisandro Dalcin ts->reason = TS_DIVERGED_STEP_REJECTED; 298211a84d6SLisandro Dalcin } 299211a84d6SLisandro Dalcin } 3003ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 301211a84d6SLisandro Dalcin } 302211a84d6SLisandro Dalcin 303d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSInterpolate_BDF(TS ts, PetscReal t, Vec X) 304d71ae5a4SJacob Faibussowitsch { 305211a84d6SLisandro Dalcin TS_BDF *bdf = (TS_BDF *)ts->data; 306211a84d6SLisandro Dalcin 307211a84d6SLisandro Dalcin PetscFunctionBegin; 3089566063dSJacob Faibussowitsch PetscCall(TSBDF_Interpolate(ts, bdf->k, t, X)); 3093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 310211a84d6SLisandro Dalcin } 311211a84d6SLisandro Dalcin 312d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSEvaluateWLTE_BDF(TS ts, NormType wnormtype, PetscInt *order, PetscReal *wlte) 313d71ae5a4SJacob Faibussowitsch { 314211a84d6SLisandro Dalcin TS_BDF *bdf = (TS_BDF *)ts->data; 315211a84d6SLisandro Dalcin PetscInt k = bdf->k; 3167453f775SEmil Constantinescu PetscReal wltea, wlter; 317211a84d6SLisandro Dalcin Vec X = bdf->work[0], Y = bdf->vec_lte; 318211a84d6SLisandro Dalcin 319211a84d6SLisandro Dalcin PetscFunctionBegin; 320211a84d6SLisandro Dalcin k = PetscMin(k, bdf->n - 1); 321*3db8ccc7SStefano Zampini if (k == 0) { 322*3db8ccc7SStefano Zampini *wlte = -1; 323*3db8ccc7SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 324*3db8ccc7SStefano Zampini } 3259566063dSJacob Faibussowitsch PetscCall(TSBDF_VecLTE(ts, k, Y)); 3269566063dSJacob Faibussowitsch PetscCall(VecAXPY(Y, 1, X)); 3279566063dSJacob Faibussowitsch PetscCall(TSErrorWeightedNorm(ts, X, Y, wnormtype, wlte, &wltea, &wlter)); 328211a84d6SLisandro Dalcin if (order) *order = k + 1; 3293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 330211a84d6SLisandro Dalcin } 331211a84d6SLisandro Dalcin 33236de76e7SStefano Zampini static PetscErrorCode TSResizeRegister_BDF(TS ts, PetscBool reg) 33336de76e7SStefano Zampini { 33436de76e7SStefano Zampini TS_BDF *bdf = (TS_BDF *)ts->data; 335ecc87898SStefano 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"}; 336ecc87898SStefano Zampini PetscInt i, maxn = PETSC_STATIC_ARRAY_LENGTH(bdf->work); 33736de76e7SStefano Zampini 33836de76e7SStefano Zampini PetscFunctionBegin; 339ecc87898SStefano Zampini PetscAssert(maxn == 8, PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "names need to be redefined"); 34036de76e7SStefano Zampini if (reg) { 34136de76e7SStefano Zampini for (i = 1; i < PetscMin(bdf->n + 1, maxn); i++) { PetscCall(TSResizeRegisterVec(ts, names[i], bdf->work[i])); } 34236de76e7SStefano Zampini } else { 34336de76e7SStefano Zampini for (i = 1; i < maxn; i++) { 34436de76e7SStefano Zampini PetscCall(TSResizeRetrieveVec(ts, names[i], &bdf->work[i])); 34536de76e7SStefano Zampini if (!bdf->work[i]) break; 34636de76e7SStefano Zampini PetscCall(PetscObjectReference((PetscObject)bdf->work[i])); 34736de76e7SStefano Zampini if (bdf->transientvar) { 34836de76e7SStefano Zampini PetscCall(VecDuplicate(bdf->work[i], &bdf->tvwork[i])); 34936de76e7SStefano Zampini PetscCall(TSComputeTransientVariable(ts, bdf->work[i], bdf->tvwork[i])); 35036de76e7SStefano Zampini } 35136de76e7SStefano Zampini } 35236de76e7SStefano Zampini } 35336de76e7SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 35436de76e7SStefano Zampini } 35536de76e7SStefano Zampini 356d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormFunction_BDF(SNES snes, Vec X, Vec F, TS ts) 357d71ae5a4SJacob Faibussowitsch { 358211a84d6SLisandro Dalcin TS_BDF *bdf = (TS_BDF *)ts->data; 3591117961dSLisandro Dalcin DM dm, dmsave = ts->dm; 360211a84d6SLisandro Dalcin PetscReal t = bdf->time[0]; 3611117961dSLisandro Dalcin PetscReal shift = bdf->shift; 3621117961dSLisandro Dalcin Vec V, V0; 363211a84d6SLisandro Dalcin 364211a84d6SLisandro Dalcin PetscFunctionBegin; 3659566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 3669566063dSJacob Faibussowitsch PetscCall(TSBDF_GetVecs(ts, dm, &V, &V0)); 367e3c11fc1SJed Brown if (bdf->transientvar) { /* shift*C(X) + V0 */ 3689566063dSJacob Faibussowitsch PetscCall(TSComputeTransientVariable(ts, X, V)); 3699566063dSJacob Faibussowitsch PetscCall(VecAYPX(V, shift, V0)); 370e3c11fc1SJed Brown } else { /* shift*X + V0 */ 3719566063dSJacob Faibussowitsch PetscCall(VecWAXPY(V, shift, X, V0)); 372e3c11fc1SJed Brown } 3731117961dSLisandro Dalcin 374211a84d6SLisandro Dalcin /* F = Function(t,X,V) */ 3751117961dSLisandro Dalcin ts->dm = dm; 3769566063dSJacob Faibussowitsch PetscCall(TSComputeIFunction(ts, t, X, V, F, PETSC_FALSE)); 3771117961dSLisandro Dalcin ts->dm = dmsave; 3781117961dSLisandro Dalcin 3799566063dSJacob Faibussowitsch PetscCall(TSBDF_RestoreVecs(ts, dm, &V, &V0)); 3803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 381211a84d6SLisandro Dalcin } 382211a84d6SLisandro Dalcin 383d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormJacobian_BDF(SNES snes, Vec X, Mat J, Mat P, TS ts) 384d71ae5a4SJacob Faibussowitsch { 385211a84d6SLisandro Dalcin TS_BDF *bdf = (TS_BDF *)ts->data; 3861117961dSLisandro Dalcin DM dm, dmsave = ts->dm; 387211a84d6SLisandro Dalcin PetscReal t = bdf->time[0]; 3881117961dSLisandro Dalcin PetscReal shift = bdf->shift; 3891117961dSLisandro Dalcin Vec V, V0; 390211a84d6SLisandro Dalcin 391211a84d6SLisandro Dalcin PetscFunctionBegin; 3929566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 3939566063dSJacob Faibussowitsch PetscCall(TSBDF_GetVecs(ts, dm, &V, &V0)); 3941117961dSLisandro Dalcin 395211a84d6SLisandro Dalcin /* J,P = Jacobian(t,X,V) */ 3961117961dSLisandro Dalcin ts->dm = dm; 3979566063dSJacob Faibussowitsch PetscCall(TSComputeIJacobian(ts, t, X, V, shift, J, P, PETSC_FALSE)); 3981117961dSLisandro Dalcin ts->dm = dmsave; 3991117961dSLisandro Dalcin 4009566063dSJacob Faibussowitsch PetscCall(TSBDF_RestoreVecs(ts, dm, &V, &V0)); 4013ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 402211a84d6SLisandro Dalcin } 403211a84d6SLisandro Dalcin 404d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSReset_BDF(TS ts) 405d71ae5a4SJacob Faibussowitsch { 406211a84d6SLisandro Dalcin TS_BDF *bdf = (TS_BDF *)ts->data; 407400f6f02SBarry Smith size_t i, n = PETSC_STATIC_ARRAY_LENGTH(bdf->work); 408211a84d6SLisandro Dalcin 409211a84d6SLisandro Dalcin PetscFunctionBegin; 410e3c11fc1SJed Brown for (i = 0; i < n; i++) { 4119566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bdf->work[i])); 4129566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bdf->tvwork[i])); 413e3c11fc1SJed Brown } 4149566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bdf->vec_dot)); 4159566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bdf->vec_wrk)); 4169566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bdf->vec_lte)); 4179566063dSJacob Faibussowitsch if (ts->dm) PetscCall(DMCoarsenHookRemove(ts->dm, DMCoarsenHook_TSBDF, DMRestrictHook_TSBDF, ts)); 4183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 419211a84d6SLisandro Dalcin } 420211a84d6SLisandro Dalcin 421d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSDestroy_BDF(TS ts) 422d71ae5a4SJacob Faibussowitsch { 423211a84d6SLisandro Dalcin PetscFunctionBegin; 4249566063dSJacob Faibussowitsch PetscCall(TSReset_BDF(ts)); 4259566063dSJacob Faibussowitsch PetscCall(PetscFree(ts->data)); 4269566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSBDFSetOrder_C", NULL)); 4279566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSBDFGetOrder_C", NULL)); 4283ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 429211a84d6SLisandro Dalcin } 430211a84d6SLisandro Dalcin 431d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetUp_BDF(TS ts) 432d71ae5a4SJacob Faibussowitsch { 433211a84d6SLisandro Dalcin TS_BDF *bdf = (TS_BDF *)ts->data; 434400f6f02SBarry Smith size_t n = PETSC_STATIC_ARRAY_LENGTH(bdf->work); 4352ffb9264SLisandro Dalcin PetscReal low, high, two = 2; 43636de76e7SStefano Zampini PetscInt cnt = 0; 437211a84d6SLisandro Dalcin 438211a84d6SLisandro Dalcin PetscFunctionBegin; 4399566063dSJacob Faibussowitsch PetscCall(TSHasTransientVariable(ts, &bdf->transientvar)); 44036de76e7SStefano Zampini for (size_t i = 0; i < n; i++) { 44136de76e7SStefano Zampini if (!bdf->work[i]) PetscCall(VecDuplicate(ts->vec_sol, &bdf->work[i])); 44236de76e7SStefano Zampini else cnt++; 44336de76e7SStefano Zampini if (i && bdf->transientvar && !bdf->tvwork[i]) PetscCall(VecDuplicate(ts->vec_sol, &bdf->tvwork[i])); 444e3c11fc1SJed Brown } 44536de76e7SStefano Zampini if (!cnt) bdf->k = bdf->n = 0; 4469566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &bdf->vec_dot)); 4479566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &bdf->vec_wrk)); 4489566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &bdf->vec_lte)); 4499566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &ts->dm)); 4509566063dSJacob Faibussowitsch PetscCall(DMCoarsenHookAdd(ts->dm, DMCoarsenHook_TSBDF, DMRestrictHook_TSBDF, ts)); 451211a84d6SLisandro Dalcin 4529566063dSJacob Faibussowitsch PetscCall(TSGetAdapt(ts, &ts->adapt)); 4539566063dSJacob Faibussowitsch PetscCall(TSAdaptCandidatesClear(ts->adapt)); 4549566063dSJacob Faibussowitsch PetscCall(TSAdaptGetClip(ts->adapt, &low, &high)); 4559566063dSJacob Faibussowitsch PetscCall(TSAdaptSetClip(ts->adapt, low, PetscMin(high, two))); 456211a84d6SLisandro Dalcin 4579566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts, &ts->snes)); 4583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 459211a84d6SLisandro Dalcin } 460211a84d6SLisandro Dalcin 461d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetFromOptions_BDF(TS ts, PetscOptionItems *PetscOptionsObject) 462d71ae5a4SJacob Faibussowitsch { 463211a84d6SLisandro Dalcin PetscFunctionBegin; 464d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "BDF ODE solver options"); 465211a84d6SLisandro Dalcin { 466211a84d6SLisandro Dalcin PetscBool flg; 467e5b8ffdfSLisandro Dalcin PetscInt order; 4689566063dSJacob Faibussowitsch PetscCall(TSBDFGetOrder(ts, &order)); 4699566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-ts_bdf_order", "Order of the BDF method", "TSBDFSetOrder", order, &order, &flg)); 4709566063dSJacob Faibussowitsch if (flg) PetscCall(TSBDFSetOrder(ts, order)); 471211a84d6SLisandro Dalcin } 472d0609cedSBarry Smith PetscOptionsHeadEnd(); 4733ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 474211a84d6SLisandro Dalcin } 475211a84d6SLisandro Dalcin 476d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSView_BDF(TS ts, PetscViewer viewer) 477d71ae5a4SJacob Faibussowitsch { 478211a84d6SLisandro Dalcin TS_BDF *bdf = (TS_BDF *)ts->data; 479211a84d6SLisandro Dalcin PetscBool iascii; 480211a84d6SLisandro Dalcin 481211a84d6SLisandro Dalcin PetscFunctionBegin; 4829566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 48348a46eb9SPierre Jolivet if (iascii) PetscCall(PetscViewerASCIIPrintf(viewer, " Order=%" PetscInt_FMT "\n", bdf->order)); 4843ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 485211a84d6SLisandro Dalcin } 486211a84d6SLisandro Dalcin 487211a84d6SLisandro Dalcin /* ------------------------------------------------------------ */ 488211a84d6SLisandro Dalcin 489d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSBDFSetOrder_BDF(TS ts, PetscInt order) 490d71ae5a4SJacob Faibussowitsch { 491211a84d6SLisandro Dalcin TS_BDF *bdf = (TS_BDF *)ts->data; 492211a84d6SLisandro Dalcin 493211a84d6SLisandro Dalcin PetscFunctionBegin; 4943ba16761SJacob Faibussowitsch if (order == bdf->order) PetscFunctionReturn(PETSC_SUCCESS); 49563a3b9bcSJacob Faibussowitsch PetscCheck(order >= 1 && order <= 6, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_OUTOFRANGE, "BDF Order %" PetscInt_FMT " not implemented", order); 496211a84d6SLisandro Dalcin bdf->order = order; 4973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 498211a84d6SLisandro Dalcin } 499211a84d6SLisandro Dalcin 500d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSBDFGetOrder_BDF(TS ts, PetscInt *order) 501d71ae5a4SJacob Faibussowitsch { 502211a84d6SLisandro Dalcin TS_BDF *bdf = (TS_BDF *)ts->data; 503211a84d6SLisandro Dalcin 504211a84d6SLisandro Dalcin PetscFunctionBegin; 505211a84d6SLisandro Dalcin *order = bdf->order; 5063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 507211a84d6SLisandro Dalcin } 508211a84d6SLisandro Dalcin 509211a84d6SLisandro Dalcin /* ------------------------------------------------------------ */ 510211a84d6SLisandro Dalcin 511211a84d6SLisandro Dalcin /*MC 512211a84d6SLisandro Dalcin TSBDF - DAE solver using BDF methods 513211a84d6SLisandro Dalcin 514211a84d6SLisandro Dalcin Level: beginner 515211a84d6SLisandro Dalcin 5161cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSCreate()`, `TSSetType()`, `TSType` 517211a84d6SLisandro Dalcin M*/ 518d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSCreate_BDF(TS ts) 519d71ae5a4SJacob Faibussowitsch { 520211a84d6SLisandro Dalcin TS_BDF *bdf; 521211a84d6SLisandro Dalcin 522211a84d6SLisandro Dalcin PetscFunctionBegin; 523211a84d6SLisandro Dalcin ts->ops->reset = TSReset_BDF; 524211a84d6SLisandro Dalcin ts->ops->destroy = TSDestroy_BDF; 525211a84d6SLisandro Dalcin ts->ops->view = TSView_BDF; 526211a84d6SLisandro Dalcin ts->ops->setup = TSSetUp_BDF; 527211a84d6SLisandro Dalcin ts->ops->setfromoptions = TSSetFromOptions_BDF; 528211a84d6SLisandro Dalcin ts->ops->step = TSStep_BDF; 529211a84d6SLisandro Dalcin ts->ops->evaluatewlte = TSEvaluateWLTE_BDF; 530211a84d6SLisandro Dalcin ts->ops->interpolate = TSInterpolate_BDF; 53136de76e7SStefano Zampini ts->ops->resizeregister = TSResizeRegister_BDF; 532211a84d6SLisandro Dalcin ts->ops->snesfunction = SNESTSFormFunction_BDF; 533211a84d6SLisandro Dalcin ts->ops->snesjacobian = SNESTSFormJacobian_BDF; 5342ffb9264SLisandro Dalcin ts->default_adapt_type = TSADAPTBASIC; 535211a84d6SLisandro Dalcin 536efd4aadfSBarry Smith ts->usessnes = PETSC_TRUE; 537efd4aadfSBarry Smith 5384dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&bdf)); 539211a84d6SLisandro Dalcin ts->data = (void *)bdf; 540211a84d6SLisandro Dalcin 541211a84d6SLisandro Dalcin bdf->status = TS_STEP_COMPLETE; 542c61711c8SStefano Zampini for (size_t i = 0; i < PETSC_STATIC_ARRAY_LENGTH(bdf->work); i++) { bdf->work[i] = bdf->tvwork[i] = NULL; } 543211a84d6SLisandro Dalcin 5449566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSBDFSetOrder_C", TSBDFSetOrder_BDF)); 5459566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSBDFGetOrder_C", TSBDFGetOrder_BDF)); 5469566063dSJacob Faibussowitsch PetscCall(TSBDFSetOrder(ts, 2)); 5473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 548211a84d6SLisandro Dalcin } 549211a84d6SLisandro Dalcin 550211a84d6SLisandro Dalcin /* ------------------------------------------------------------ */ 551211a84d6SLisandro Dalcin 552211a84d6SLisandro Dalcin /*@ 553bcf0153eSBarry Smith TSBDFSetOrder - Set the order of the `TSBDF` method 554211a84d6SLisandro Dalcin 555c3339decSBarry Smith Logically Collective 556211a84d6SLisandro Dalcin 557d8d19677SJose E. Roman Input Parameters: 558211a84d6SLisandro Dalcin + ts - timestepping context 559211a84d6SLisandro Dalcin - order - order of the method 560211a84d6SLisandro Dalcin 561bcf0153eSBarry Smith Options Database Key: 562147403d9SBarry Smith . -ts_bdf_order <order> - select the order 563211a84d6SLisandro Dalcin 564211a84d6SLisandro Dalcin Level: intermediate 565211a84d6SLisandro Dalcin 566bcf0153eSBarry Smith .seealso: `TSBDFGetOrder()`, `TS`, `TSBDF` 567211a84d6SLisandro Dalcin @*/ 568d71ae5a4SJacob Faibussowitsch PetscErrorCode TSBDFSetOrder(TS ts, PetscInt order) 569d71ae5a4SJacob Faibussowitsch { 570211a84d6SLisandro Dalcin PetscFunctionBegin; 571211a84d6SLisandro Dalcin PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 572211a84d6SLisandro Dalcin PetscValidLogicalCollectiveInt(ts, order, 2); 573cac4c232SBarry Smith PetscTryMethod(ts, "TSBDFSetOrder_C", (TS, PetscInt), (ts, order)); 5743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 575211a84d6SLisandro Dalcin } 576211a84d6SLisandro Dalcin 577211a84d6SLisandro Dalcin /*@ 578bcf0153eSBarry Smith TSBDFGetOrder - Get the order of the `TSBDF` method 579211a84d6SLisandro Dalcin 580211a84d6SLisandro Dalcin Not Collective 581211a84d6SLisandro Dalcin 582211a84d6SLisandro Dalcin Input Parameter: 583211a84d6SLisandro Dalcin . ts - timestepping context 584211a84d6SLisandro Dalcin 585211a84d6SLisandro Dalcin Output Parameter: 586211a84d6SLisandro Dalcin . order - order of the method 587211a84d6SLisandro Dalcin 588211a84d6SLisandro Dalcin Level: intermediate 589211a84d6SLisandro Dalcin 590bcf0153eSBarry Smith .seealso: `TSBDFSetOrder()`, `TS`, `TSBDF` 591211a84d6SLisandro Dalcin @*/ 592d71ae5a4SJacob Faibussowitsch PetscErrorCode TSBDFGetOrder(TS ts, PetscInt *order) 593d71ae5a4SJacob Faibussowitsch { 594211a84d6SLisandro Dalcin PetscFunctionBegin; 595211a84d6SLisandro Dalcin PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 5964f572ea9SToby Isaac PetscAssertPointer(order, 2); 597cac4c232SBarry Smith PetscUseMethod(ts, "TSBDFGetOrder_C", (TS, PetscInt *), (ts, order)); 5983ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 599211a84d6SLisandro Dalcin } 600