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; 8211a84d6SLisandro Dalcin static const char citation[] = 9211a84d6SLisandro Dalcin "@book{Brenan1995,\n" 10211a84d6SLisandro Dalcin " title = {Numerical Solution of Initial-Value Problems in Differential-Algebraic Equations},\n" 11211a84d6SLisandro Dalcin " author = {Brenan, K. and Campbell, S. and Petzold, L.},\n" 12211a84d6SLisandro Dalcin " publisher = {Society for Industrial and Applied Mathematics},\n" 13211a84d6SLisandro Dalcin " year = {1995},\n" 14211a84d6SLisandro Dalcin " doi = {10.1137/1.9781611971224},\n}\n"; 15211a84d6SLisandro Dalcin 16211a84d6SLisandro Dalcin typedef struct { 17211a84d6SLisandro Dalcin PetscInt k,n; 18211a84d6SLisandro Dalcin PetscReal time[6+2]; 19211a84d6SLisandro Dalcin Vec work[6+2]; 20e3c11fc1SJed Brown Vec tvwork[6+2]; 21211a84d6SLisandro Dalcin PetscReal shift; 22e3c11fc1SJed Brown Vec vec_dot; /* Xdot when !transientvar, else Cdot where C(X) is the transient variable. */ 231117961dSLisandro Dalcin Vec vec_wrk; 24211a84d6SLisandro Dalcin Vec vec_lte; 25211a84d6SLisandro Dalcin 26e3c11fc1SJed Brown PetscBool transientvar; 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 */ 34211a84d6SLisandro Dalcin PETSC_STATIC_INLINE void LagrangeBasisVals(PetscInt n,PetscReal t,const PetscReal T[],PetscScalar L[]) 35211a84d6SLisandro Dalcin { 36211a84d6SLisandro Dalcin PetscInt k,j; 37211a84d6SLisandro Dalcin for (k=0; k<n; k++) 38211a84d6SLisandro Dalcin for (L[k]=1, j=0; j<n; j++) 39211a84d6SLisandro Dalcin if (j != k) 40211a84d6SLisandro Dalcin L[k] *= (t - T[j])/(T[k] - T[j]); 41211a84d6SLisandro Dalcin } 42211a84d6SLisandro Dalcin 43211a84d6SLisandro Dalcin PETSC_STATIC_INLINE void LagrangeBasisDers(PetscInt n,PetscReal t,const PetscReal T[],PetscScalar dL[]) 44211a84d6SLisandro Dalcin { 45211a84d6SLisandro Dalcin PetscInt k,j,i; 46211a84d6SLisandro Dalcin for (k=0; k<n; k++) 47211a84d6SLisandro Dalcin for (dL[k]=0, j=0; j<n; j++) 48211a84d6SLisandro Dalcin if (j != k) { 49211a84d6SLisandro Dalcin PetscReal L = 1/(T[k] - T[j]); 50211a84d6SLisandro Dalcin for (i=0; i<n; i++) 51211a84d6SLisandro Dalcin if (i != j && i != k) 52211a84d6SLisandro Dalcin L *= (t - T[i])/(T[k] - T[i]); 53211a84d6SLisandro Dalcin dL[k] += L; 54211a84d6SLisandro Dalcin } 55211a84d6SLisandro Dalcin } 56211a84d6SLisandro Dalcin 571117961dSLisandro Dalcin static PetscErrorCode TSBDF_GetVecs(TS ts,DM dm,Vec *Xdot,Vec *Ydot) 581117961dSLisandro Dalcin { 591117961dSLisandro Dalcin TS_BDF *bdf = (TS_BDF*)ts->data; 601117961dSLisandro Dalcin PetscErrorCode ierr; 611117961dSLisandro Dalcin 621117961dSLisandro Dalcin PetscFunctionBegin; 631117961dSLisandro Dalcin if (dm && dm != ts->dm) { 641117961dSLisandro Dalcin ierr = DMGetNamedGlobalVector(dm,"TSBDF_Vec_Xdot",Xdot);CHKERRQ(ierr); 651117961dSLisandro Dalcin ierr = DMGetNamedGlobalVector(dm,"TSBDF_Vec_Ydot",Ydot);CHKERRQ(ierr); 661117961dSLisandro Dalcin } else { 671117961dSLisandro Dalcin *Xdot = bdf->vec_dot; 681117961dSLisandro Dalcin *Ydot = bdf->vec_wrk; 691117961dSLisandro Dalcin } 701117961dSLisandro Dalcin PetscFunctionReturn(0); 711117961dSLisandro Dalcin } 721117961dSLisandro Dalcin 731117961dSLisandro Dalcin static PetscErrorCode TSBDF_RestoreVecs(TS ts,DM dm,Vec *Xdot,Vec *Ydot) 741117961dSLisandro Dalcin { 751117961dSLisandro Dalcin TS_BDF *bdf = (TS_BDF*)ts->data; 761117961dSLisandro Dalcin PetscErrorCode ierr; 771117961dSLisandro Dalcin 781117961dSLisandro Dalcin PetscFunctionBegin; 791117961dSLisandro Dalcin if (dm && dm != ts->dm) { 801117961dSLisandro Dalcin ierr = DMRestoreNamedGlobalVector(dm,"TSBDF_Vec_Xdot",Xdot);CHKERRQ(ierr); 811117961dSLisandro Dalcin ierr = DMRestoreNamedGlobalVector(dm,"TSBDF_Vec_Ydot",Ydot);CHKERRQ(ierr); 821117961dSLisandro Dalcin } else { 831117961dSLisandro Dalcin if (*Xdot != bdf->vec_dot) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_INCOMP,"Vec does not match the cache"); 841117961dSLisandro Dalcin if (*Ydot != bdf->vec_wrk) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_INCOMP,"Vec does not match the cache"); 851117961dSLisandro Dalcin *Xdot = NULL; 861117961dSLisandro Dalcin *Ydot = NULL; 871117961dSLisandro Dalcin } 881117961dSLisandro Dalcin PetscFunctionReturn(0); 891117961dSLisandro Dalcin } 901117961dSLisandro Dalcin 911117961dSLisandro Dalcin static PetscErrorCode DMCoarsenHook_TSBDF(DM fine,DM coarse,void *ctx) 921117961dSLisandro Dalcin { 931117961dSLisandro Dalcin PetscFunctionBegin; 941117961dSLisandro Dalcin PetscFunctionReturn(0); 951117961dSLisandro Dalcin } 961117961dSLisandro Dalcin 971117961dSLisandro Dalcin static PetscErrorCode DMRestrictHook_TSBDF(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx) 981117961dSLisandro Dalcin { 991117961dSLisandro Dalcin TS ts = (TS)ctx; 1001117961dSLisandro Dalcin Vec Ydot,Ydot_c; 1011117961dSLisandro Dalcin Vec Xdot,Xdot_c; 1021117961dSLisandro Dalcin PetscErrorCode ierr; 1031117961dSLisandro Dalcin 1041117961dSLisandro Dalcin PetscFunctionBegin; 1051117961dSLisandro Dalcin ierr = TSBDF_GetVecs(ts,fine,&Xdot,&Ydot);CHKERRQ(ierr); 1061117961dSLisandro Dalcin ierr = TSBDF_GetVecs(ts,coarse,&Xdot_c,&Ydot_c);CHKERRQ(ierr); 1071117961dSLisandro Dalcin 1081117961dSLisandro Dalcin ierr = MatRestrict(restrct,Ydot,Ydot_c);CHKERRQ(ierr); 1091117961dSLisandro Dalcin ierr = VecPointwiseMult(Ydot_c,rscale,Ydot_c);CHKERRQ(ierr); 1101117961dSLisandro Dalcin 1111117961dSLisandro Dalcin ierr = TSBDF_RestoreVecs(ts,fine,&Xdot,&Ydot);CHKERRQ(ierr); 1121117961dSLisandro Dalcin ierr = TSBDF_RestoreVecs(ts,coarse,&Xdot_c,&Ydot_c);CHKERRQ(ierr); 1131117961dSLisandro Dalcin PetscFunctionReturn(0); 1141117961dSLisandro Dalcin } 1151117961dSLisandro Dalcin 116211a84d6SLisandro Dalcin static PetscErrorCode TSBDF_Advance(TS ts,PetscReal t,Vec X) 117211a84d6SLisandro Dalcin { 118211a84d6SLisandro Dalcin TS_BDF *bdf = (TS_BDF*)ts->data; 119211a84d6SLisandro Dalcin PetscInt i,n = (PetscInt)(sizeof(bdf->work)/sizeof(Vec)); 120e3c11fc1SJed Brown Vec tail = bdf->work[n-1],tvtail = bdf->tvwork[n-1]; 121211a84d6SLisandro Dalcin PetscErrorCode ierr; 122211a84d6SLisandro Dalcin 123211a84d6SLisandro Dalcin PetscFunctionBegin; 124211a84d6SLisandro Dalcin for (i=n-1; i>=2; i--) { 125211a84d6SLisandro Dalcin bdf->time[i] = bdf->time[i-1]; 126211a84d6SLisandro Dalcin bdf->work[i] = bdf->work[i-1]; 127e3c11fc1SJed Brown bdf->tvwork[i] = bdf->tvwork[i-1]; 128211a84d6SLisandro Dalcin } 129211a84d6SLisandro Dalcin bdf->n = PetscMin(bdf->n+1,n-1); 130211a84d6SLisandro Dalcin bdf->time[1] = t; 131211a84d6SLisandro Dalcin bdf->work[1] = tail; 132e3c11fc1SJed Brown bdf->tvwork[1] = tvtail; 133211a84d6SLisandro Dalcin ierr = VecCopy(X,tail);CHKERRQ(ierr); 134e3c11fc1SJed Brown ierr = TSComputeTransientVariable(ts,tail,tvtail);CHKERRQ(ierr); 135211a84d6SLisandro Dalcin PetscFunctionReturn(0); 136211a84d6SLisandro Dalcin } 137211a84d6SLisandro Dalcin 138211a84d6SLisandro Dalcin static PetscErrorCode TSBDF_VecLTE(TS ts,PetscInt order,Vec lte) 139211a84d6SLisandro Dalcin { 140211a84d6SLisandro Dalcin TS_BDF *bdf = (TS_BDF*)ts->data; 141211a84d6SLisandro Dalcin PetscInt i,n = order+1; 142211a84d6SLisandro Dalcin PetscReal *time = bdf->time; 143211a84d6SLisandro Dalcin Vec *vecs = bdf->work; 144211a84d6SLisandro Dalcin PetscScalar a[8],b[8],alpha[8]; 145211a84d6SLisandro Dalcin PetscErrorCode ierr; 146211a84d6SLisandro Dalcin 147211a84d6SLisandro Dalcin PetscFunctionBegin; 148211a84d6SLisandro Dalcin LagrangeBasisDers(n+0,time[0],time,a); a[n] =0; 149211a84d6SLisandro Dalcin LagrangeBasisDers(n+1,time[0],time,b); 150211a84d6SLisandro Dalcin for (i=0; i<n+1; i++) alpha[i] = (a[i]-b[i])/a[0]; 151211a84d6SLisandro Dalcin ierr = VecZeroEntries(lte);CHKERRQ(ierr); 152211a84d6SLisandro Dalcin ierr = VecMAXPY(lte,n+1,alpha,vecs);CHKERRQ(ierr); 153211a84d6SLisandro Dalcin PetscFunctionReturn(0); 154211a84d6SLisandro Dalcin } 155211a84d6SLisandro Dalcin 156211a84d6SLisandro Dalcin static PetscErrorCode TSBDF_Extrapolate(TS ts,PetscInt order,PetscReal t,Vec X) 157211a84d6SLisandro Dalcin { 158211a84d6SLisandro Dalcin TS_BDF *bdf = (TS_BDF*)ts->data; 159211a84d6SLisandro Dalcin PetscInt n = order+1; 160211a84d6SLisandro Dalcin PetscReal *time = bdf->time+1; 161211a84d6SLisandro Dalcin Vec *vecs = bdf->work+1; 162211a84d6SLisandro Dalcin PetscScalar alpha[7]; 163211a84d6SLisandro Dalcin PetscErrorCode ierr; 164211a84d6SLisandro Dalcin 165211a84d6SLisandro Dalcin PetscFunctionBegin; 166211a84d6SLisandro Dalcin n = PetscMin(n,bdf->n); 167211a84d6SLisandro Dalcin LagrangeBasisVals(n,t,time,alpha); 168211a84d6SLisandro Dalcin ierr = VecZeroEntries(X);CHKERRQ(ierr); 169211a84d6SLisandro Dalcin ierr = VecMAXPY(X,n,alpha,vecs);CHKERRQ(ierr); 170211a84d6SLisandro Dalcin PetscFunctionReturn(0); 171211a84d6SLisandro Dalcin } 172211a84d6SLisandro Dalcin 173211a84d6SLisandro Dalcin static PetscErrorCode TSBDF_Interpolate(TS ts,PetscInt order,PetscReal t,Vec X) 174211a84d6SLisandro Dalcin { 175211a84d6SLisandro Dalcin TS_BDF *bdf = (TS_BDF*)ts->data; 176211a84d6SLisandro Dalcin PetscInt n = order+1; 177211a84d6SLisandro Dalcin PetscReal *time = bdf->time; 178211a84d6SLisandro Dalcin Vec *vecs = bdf->work; 179211a84d6SLisandro Dalcin PetscScalar alpha[7]; 180211a84d6SLisandro Dalcin PetscErrorCode ierr; 181211a84d6SLisandro Dalcin 182211a84d6SLisandro Dalcin PetscFunctionBegin; 183211a84d6SLisandro Dalcin LagrangeBasisVals(n,t,time,alpha); 184211a84d6SLisandro Dalcin ierr = VecZeroEntries(X);CHKERRQ(ierr); 185211a84d6SLisandro Dalcin ierr = VecMAXPY(X,n,alpha,vecs);CHKERRQ(ierr); 186211a84d6SLisandro Dalcin PetscFunctionReturn(0); 187211a84d6SLisandro Dalcin } 188211a84d6SLisandro Dalcin 189e3c11fc1SJed Brown /* Compute the affine term V0 such that Xdot = shift*X + V0. 190e3c11fc1SJed Brown * 191e3c11fc1SJed Brown * When using transient variables, we're computing Cdot = shift*C(X) + V0, and thus choose a linear combination of tvwork. 192e3c11fc1SJed Brown */ 1931117961dSLisandro Dalcin static PetscErrorCode TSBDF_PreSolve(TS ts) 1941117961dSLisandro Dalcin { 1951117961dSLisandro Dalcin TS_BDF *bdf = (TS_BDF*)ts->data; 1961117961dSLisandro Dalcin PetscInt i,n = PetscMax(bdf->k,1) + 1; 1971117961dSLisandro Dalcin Vec V,V0; 1981117961dSLisandro Dalcin Vec vecs[7]; 1991117961dSLisandro Dalcin PetscScalar alpha[7]; 2001117961dSLisandro Dalcin PetscErrorCode ierr; 2011117961dSLisandro Dalcin 2021117961dSLisandro Dalcin PetscFunctionBegin; 2031117961dSLisandro Dalcin ierr = TSBDF_GetVecs(ts,NULL,&V,&V0);CHKERRQ(ierr); 2041117961dSLisandro Dalcin LagrangeBasisDers(n,bdf->time[0],bdf->time,alpha); 205e3c11fc1SJed Brown for (i=1; i<n; i++) { 206e3c11fc1SJed Brown vecs[i] = bdf->transientvar ? bdf->tvwork[i] : bdf->work[i]; 207e3c11fc1SJed Brown } 2081117961dSLisandro Dalcin ierr = VecZeroEntries(V0);CHKERRQ(ierr); 2091117961dSLisandro Dalcin ierr = VecMAXPY(V0,n-1,alpha+1,vecs+1);CHKERRQ(ierr); 2101117961dSLisandro Dalcin bdf->shift = PetscRealPart(alpha[0]); 2111117961dSLisandro Dalcin ierr = TSBDF_RestoreVecs(ts,NULL,&V,&V0);CHKERRQ(ierr); 2121117961dSLisandro Dalcin PetscFunctionReturn(0); 2131117961dSLisandro Dalcin } 2141117961dSLisandro Dalcin 215470880abSPatrick Sanan static PetscErrorCode TSBDF_SNESSolve(TS ts,Vec b,Vec x) 216211a84d6SLisandro Dalcin { 217211a84d6SLisandro Dalcin PetscInt nits,lits; 218211a84d6SLisandro Dalcin PetscErrorCode ierr; 219211a84d6SLisandro Dalcin 220211a84d6SLisandro Dalcin PetscFunctionBegin; 2211117961dSLisandro Dalcin ierr = TSBDF_PreSolve(ts);CHKERRQ(ierr); 222211a84d6SLisandro Dalcin ierr = SNESSolve(ts->snes,b,x);CHKERRQ(ierr); 223211a84d6SLisandro Dalcin ierr = SNESGetIterationNumber(ts->snes,&nits);CHKERRQ(ierr); 224211a84d6SLisandro Dalcin ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr); 225211a84d6SLisandro Dalcin ts->snes_its += nits; ts->ksp_its += lits; 226211a84d6SLisandro Dalcin PetscFunctionReturn(0); 227211a84d6SLisandro Dalcin } 228211a84d6SLisandro Dalcin 229211a84d6SLisandro Dalcin static PetscErrorCode TSBDF_Restart(TS ts,PetscBool *accept) 230211a84d6SLisandro Dalcin { 231211a84d6SLisandro Dalcin TS_BDF *bdf = (TS_BDF*)ts->data; 232211a84d6SLisandro Dalcin PetscErrorCode ierr; 233211a84d6SLisandro Dalcin 234211a84d6SLisandro Dalcin PetscFunctionBegin; 235211a84d6SLisandro Dalcin bdf->k = 1; bdf->n = 0; 236211a84d6SLisandro Dalcin ierr = TSBDF_Advance(ts,ts->ptime,ts->vec_sol);CHKERRQ(ierr); 237211a84d6SLisandro Dalcin 238211a84d6SLisandro Dalcin bdf->time[0] = ts->ptime + ts->time_step/2; 239211a84d6SLisandro Dalcin ierr = VecCopy(bdf->work[1],bdf->work[0]);CHKERRQ(ierr); 240211a84d6SLisandro Dalcin ierr = TSPreStage(ts,bdf->time[0]);CHKERRQ(ierr); 241470880abSPatrick Sanan ierr = TSBDF_SNESSolve(ts,NULL,bdf->work[0]);CHKERRQ(ierr); 242211a84d6SLisandro Dalcin ierr = TSPostStage(ts,bdf->time[0],0,&bdf->work[0]);CHKERRQ(ierr); 243211a84d6SLisandro Dalcin ierr = TSAdaptCheckStage(ts->adapt,ts,bdf->time[0],bdf->work[0],accept);CHKERRQ(ierr); 244211a84d6SLisandro Dalcin if (!*accept) PetscFunctionReturn(0); 245211a84d6SLisandro Dalcin 246211a84d6SLisandro Dalcin bdf->k = PetscMin(2,bdf->order); bdf->n++; 247211a84d6SLisandro Dalcin ierr = VecCopy(bdf->work[0],bdf->work[2]);CHKERRQ(ierr); 248211a84d6SLisandro Dalcin bdf->time[2] = bdf->time[0]; 249e3c11fc1SJed Brown ierr = TSComputeTransientVariable(ts,bdf->work[2],bdf->tvwork[2]);CHKERRQ(ierr); 250211a84d6SLisandro Dalcin PetscFunctionReturn(0); 251211a84d6SLisandro Dalcin } 252211a84d6SLisandro Dalcin 253211a84d6SLisandro Dalcin static const char *const BDF_SchemeName[] = {"", "1", "2", "3", "4", "5", "6"}; 254211a84d6SLisandro Dalcin 255211a84d6SLisandro Dalcin static PetscErrorCode TSStep_BDF(TS ts) 256211a84d6SLisandro Dalcin { 257211a84d6SLisandro Dalcin TS_BDF *bdf = (TS_BDF*)ts->data; 258211a84d6SLisandro Dalcin PetscInt rejections = 0; 259211a84d6SLisandro Dalcin PetscBool stageok,accept = PETSC_TRUE; 260211a84d6SLisandro Dalcin PetscReal next_time_step = ts->time_step; 261211a84d6SLisandro Dalcin PetscErrorCode ierr; 262211a84d6SLisandro Dalcin 263211a84d6SLisandro Dalcin PetscFunctionBegin; 264211a84d6SLisandro Dalcin ierr = PetscCitationsRegister(citation,&cited);CHKERRQ(ierr); 265211a84d6SLisandro Dalcin 266211a84d6SLisandro Dalcin if (!ts->steprollback && !ts->steprestart) { 267211a84d6SLisandro Dalcin bdf->k = PetscMin(bdf->k+1,bdf->order); 268211a84d6SLisandro Dalcin ierr = TSBDF_Advance(ts,ts->ptime,ts->vec_sol);CHKERRQ(ierr); 269211a84d6SLisandro Dalcin } 270211a84d6SLisandro Dalcin 271211a84d6SLisandro Dalcin bdf->status = TS_STEP_INCOMPLETE; 272211a84d6SLisandro Dalcin while (!ts->reason && bdf->status != TS_STEP_COMPLETE) { 273211a84d6SLisandro Dalcin 274211a84d6SLisandro Dalcin if (ts->steprestart) { 275211a84d6SLisandro Dalcin ierr = TSBDF_Restart(ts,&stageok);CHKERRQ(ierr); 276211a84d6SLisandro Dalcin if (!stageok) goto reject_step; 277211a84d6SLisandro Dalcin } 278211a84d6SLisandro Dalcin 279211a84d6SLisandro Dalcin bdf->time[0] = ts->ptime + ts->time_step; 280211a84d6SLisandro Dalcin ierr = TSBDF_Extrapolate(ts,bdf->k-(accept?0:1),bdf->time[0],bdf->work[0]);CHKERRQ(ierr); 281211a84d6SLisandro Dalcin ierr = TSPreStage(ts,bdf->time[0]);CHKERRQ(ierr); 282470880abSPatrick Sanan ierr = TSBDF_SNESSolve(ts,NULL,bdf->work[0]);CHKERRQ(ierr); 283211a84d6SLisandro Dalcin ierr = TSPostStage(ts,bdf->time[0],0,&bdf->work[0]);CHKERRQ(ierr); 284211a84d6SLisandro Dalcin ierr = TSAdaptCheckStage(ts->adapt,ts,bdf->time[0],bdf->work[0],&stageok);CHKERRQ(ierr); 285211a84d6SLisandro Dalcin if (!stageok) goto reject_step; 286211a84d6SLisandro Dalcin 287211a84d6SLisandro Dalcin bdf->status = TS_STEP_PENDING; 288211a84d6SLisandro Dalcin ierr = TSAdaptCandidatesClear(ts->adapt);CHKERRQ(ierr); 289211a84d6SLisandro Dalcin ierr = TSAdaptCandidateAdd(ts->adapt,BDF_SchemeName[bdf->k],bdf->k,1,1.0,1.0,PETSC_TRUE);CHKERRQ(ierr); 290211a84d6SLisandro Dalcin ierr = TSAdaptChoose(ts->adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr); 291211a84d6SLisandro Dalcin bdf->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE; 292211a84d6SLisandro Dalcin if (!accept) { ts->time_step = next_time_step; goto reject_step; } 293211a84d6SLisandro Dalcin 294211a84d6SLisandro Dalcin ierr = VecCopy(bdf->work[0],ts->vec_sol);CHKERRQ(ierr); 295211a84d6SLisandro Dalcin ts->ptime += ts->time_step; 296211a84d6SLisandro Dalcin ts->time_step = next_time_step; 297211a84d6SLisandro Dalcin break; 298211a84d6SLisandro Dalcin 299211a84d6SLisandro Dalcin reject_step: 300211a84d6SLisandro Dalcin ts->reject++; accept = PETSC_FALSE; 301211a84d6SLisandro Dalcin if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) { 302211a84d6SLisandro Dalcin ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);CHKERRQ(ierr); 303211a84d6SLisandro Dalcin ts->reason = TS_DIVERGED_STEP_REJECTED; 304211a84d6SLisandro Dalcin } 305211a84d6SLisandro Dalcin } 306211a84d6SLisandro Dalcin PetscFunctionReturn(0); 307211a84d6SLisandro Dalcin } 308211a84d6SLisandro Dalcin 309211a84d6SLisandro Dalcin static PetscErrorCode TSInterpolate_BDF(TS ts,PetscReal t,Vec X) 310211a84d6SLisandro Dalcin { 311211a84d6SLisandro Dalcin TS_BDF *bdf = (TS_BDF*)ts->data; 312211a84d6SLisandro Dalcin PetscErrorCode ierr; 313211a84d6SLisandro Dalcin 314211a84d6SLisandro Dalcin PetscFunctionBegin; 315211a84d6SLisandro Dalcin ierr = TSBDF_Interpolate(ts,bdf->k,t,X);CHKERRQ(ierr); 316211a84d6SLisandro Dalcin PetscFunctionReturn(0); 317211a84d6SLisandro Dalcin } 318211a84d6SLisandro Dalcin 319211a84d6SLisandro Dalcin static PetscErrorCode TSEvaluateWLTE_BDF(TS ts,NormType wnormtype,PetscInt *order,PetscReal *wlte) 320211a84d6SLisandro Dalcin { 321211a84d6SLisandro Dalcin TS_BDF *bdf = (TS_BDF*)ts->data; 322211a84d6SLisandro Dalcin PetscInt k = bdf->k; 3237453f775SEmil Constantinescu PetscReal wltea,wlter; 324211a84d6SLisandro Dalcin Vec X = bdf->work[0], Y = bdf->vec_lte; 325211a84d6SLisandro Dalcin PetscErrorCode ierr; 326211a84d6SLisandro Dalcin 327211a84d6SLisandro Dalcin PetscFunctionBegin; 328211a84d6SLisandro Dalcin k = PetscMin(k,bdf->n-1); 329211a84d6SLisandro Dalcin ierr = TSBDF_VecLTE(ts,k,Y);CHKERRQ(ierr); 330211a84d6SLisandro Dalcin ierr = VecAXPY(Y,1,X);CHKERRQ(ierr); 3317453f775SEmil Constantinescu ierr = TSErrorWeightedNorm(ts,X,Y,wnormtype,wlte,&wltea,&wlter);CHKERRQ(ierr); 332211a84d6SLisandro Dalcin if (order) *order = k + 1; 333211a84d6SLisandro Dalcin PetscFunctionReturn(0); 334211a84d6SLisandro Dalcin } 335211a84d6SLisandro Dalcin 336211a84d6SLisandro Dalcin static PetscErrorCode TSRollBack_BDF(TS ts) 337211a84d6SLisandro Dalcin { 338211a84d6SLisandro Dalcin TS_BDF *bdf = (TS_BDF*)ts->data; 339211a84d6SLisandro Dalcin PetscErrorCode ierr; 340211a84d6SLisandro Dalcin 341211a84d6SLisandro Dalcin PetscFunctionBegin; 342211a84d6SLisandro Dalcin ierr = VecCopy(bdf->work[1],ts->vec_sol);CHKERRQ(ierr); 343211a84d6SLisandro Dalcin PetscFunctionReturn(0); 344211a84d6SLisandro Dalcin } 345211a84d6SLisandro Dalcin 3461117961dSLisandro Dalcin static PetscErrorCode SNESTSFormFunction_BDF(SNES snes,Vec X,Vec F,TS ts) 347211a84d6SLisandro Dalcin { 348211a84d6SLisandro Dalcin TS_BDF *bdf = (TS_BDF*)ts->data; 3491117961dSLisandro Dalcin DM dm, dmsave = ts->dm; 350211a84d6SLisandro Dalcin PetscReal t = bdf->time[0]; 3511117961dSLisandro Dalcin PetscReal shift = bdf->shift; 3521117961dSLisandro Dalcin Vec V,V0; 353211a84d6SLisandro Dalcin PetscErrorCode ierr; 354211a84d6SLisandro Dalcin 355211a84d6SLisandro Dalcin PetscFunctionBegin; 3561117961dSLisandro Dalcin ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 3571117961dSLisandro Dalcin ierr = TSBDF_GetVecs(ts,dm,&V,&V0);CHKERRQ(ierr); 358e3c11fc1SJed Brown if (bdf->transientvar) { /* shift*C(X) + V0 */ 359e3c11fc1SJed Brown ierr = TSComputeTransientVariable(ts,X,V);CHKERRQ(ierr); 360e3c11fc1SJed Brown ierr = VecAYPX(V,shift,V0);CHKERRQ(ierr); 361e3c11fc1SJed Brown } else { /* shift*X + V0 */ 3621117961dSLisandro Dalcin ierr = VecWAXPY(V,shift,X,V0);CHKERRQ(ierr); 363e3c11fc1SJed Brown } 3641117961dSLisandro Dalcin 365211a84d6SLisandro Dalcin /* F = Function(t,X,V) */ 3661117961dSLisandro Dalcin ts->dm = dm; 367211a84d6SLisandro Dalcin ierr = TSComputeIFunction(ts,t,X,V,F,PETSC_FALSE);CHKERRQ(ierr); 3681117961dSLisandro Dalcin ts->dm = dmsave; 3691117961dSLisandro Dalcin 3701117961dSLisandro Dalcin ierr = TSBDF_RestoreVecs(ts,dm,&V,&V0);CHKERRQ(ierr); 371211a84d6SLisandro Dalcin PetscFunctionReturn(0); 372211a84d6SLisandro Dalcin } 373211a84d6SLisandro Dalcin 3741117961dSLisandro Dalcin static PetscErrorCode SNESTSFormJacobian_BDF(SNES snes,Vec X,Mat J,Mat P,TS ts) 375211a84d6SLisandro Dalcin { 376211a84d6SLisandro Dalcin TS_BDF *bdf = (TS_BDF*)ts->data; 3771117961dSLisandro Dalcin DM dm, dmsave = ts->dm; 378211a84d6SLisandro Dalcin PetscReal t = bdf->time[0]; 3791117961dSLisandro Dalcin PetscReal shift = bdf->shift; 3801117961dSLisandro Dalcin Vec V,V0; 381211a84d6SLisandro Dalcin PetscErrorCode ierr; 382211a84d6SLisandro Dalcin 383211a84d6SLisandro Dalcin PetscFunctionBegin; 3841117961dSLisandro Dalcin ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 3851117961dSLisandro Dalcin ierr = TSBDF_GetVecs(ts,dm,&V,&V0);CHKERRQ(ierr); 3861117961dSLisandro Dalcin 387211a84d6SLisandro Dalcin /* J,P = Jacobian(t,X,V) */ 3881117961dSLisandro Dalcin ts->dm = dm; 3891117961dSLisandro Dalcin ierr = TSComputeIJacobian(ts,t,X,V,shift,J,P,PETSC_FALSE);CHKERRQ(ierr); 3901117961dSLisandro Dalcin ts->dm = dmsave; 3911117961dSLisandro Dalcin 3921117961dSLisandro Dalcin ierr = TSBDF_RestoreVecs(ts,dm,&V,&V0);CHKERRQ(ierr); 393211a84d6SLisandro Dalcin PetscFunctionReturn(0); 394211a84d6SLisandro Dalcin } 395211a84d6SLisandro Dalcin 396211a84d6SLisandro Dalcin static PetscErrorCode TSReset_BDF(TS ts) 397211a84d6SLisandro Dalcin { 398211a84d6SLisandro Dalcin TS_BDF *bdf = (TS_BDF*)ts->data; 399211a84d6SLisandro Dalcin size_t i,n = sizeof(bdf->work)/sizeof(Vec); 400211a84d6SLisandro Dalcin PetscErrorCode ierr; 401211a84d6SLisandro Dalcin 402211a84d6SLisandro Dalcin PetscFunctionBegin; 4031117961dSLisandro Dalcin bdf->k = bdf->n = 0; 404e3c11fc1SJed Brown for (i=0; i<n; i++) { 405e3c11fc1SJed Brown ierr = VecDestroy(&bdf->work[i]);CHKERRQ(ierr); 406e3c11fc1SJed Brown ierr = VecDestroy(&bdf->tvwork[i]);CHKERRQ(ierr); 407e3c11fc1SJed Brown } 408211a84d6SLisandro Dalcin ierr = VecDestroy(&bdf->vec_dot);CHKERRQ(ierr); 4091117961dSLisandro Dalcin ierr = VecDestroy(&bdf->vec_wrk);CHKERRQ(ierr); 410211a84d6SLisandro Dalcin ierr = VecDestroy(&bdf->vec_lte);CHKERRQ(ierr); 4111117961dSLisandro Dalcin if (ts->dm) {ierr = DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSBDF,DMRestrictHook_TSBDF,ts);CHKERRQ(ierr);} 412211a84d6SLisandro Dalcin PetscFunctionReturn(0); 413211a84d6SLisandro Dalcin } 414211a84d6SLisandro Dalcin 415211a84d6SLisandro Dalcin static PetscErrorCode TSDestroy_BDF(TS ts) 416211a84d6SLisandro Dalcin { 417211a84d6SLisandro Dalcin PetscErrorCode ierr; 418211a84d6SLisandro Dalcin 419211a84d6SLisandro Dalcin PetscFunctionBegin; 420211a84d6SLisandro Dalcin ierr = TSReset_BDF(ts);CHKERRQ(ierr); 421211a84d6SLisandro Dalcin ierr = PetscFree(ts->data);CHKERRQ(ierr); 422211a84d6SLisandro Dalcin ierr = PetscObjectComposeFunction((PetscObject)ts,"TSBDFSetOrder_C",NULL);CHKERRQ(ierr); 423211a84d6SLisandro Dalcin ierr = PetscObjectComposeFunction((PetscObject)ts,"TSBDFGetOrder_C",NULL);CHKERRQ(ierr); 424211a84d6SLisandro Dalcin PetscFunctionReturn(0); 425211a84d6SLisandro Dalcin } 426211a84d6SLisandro Dalcin 427211a84d6SLisandro Dalcin static PetscErrorCode TSSetUp_BDF(TS ts) 428211a84d6SLisandro Dalcin { 429211a84d6SLisandro Dalcin TS_BDF *bdf = (TS_BDF*)ts->data; 430211a84d6SLisandro Dalcin size_t i,n = sizeof(bdf->work)/sizeof(Vec); 4312ffb9264SLisandro Dalcin PetscReal low,high,two = 2; 432211a84d6SLisandro Dalcin PetscErrorCode ierr; 433211a84d6SLisandro Dalcin 434211a84d6SLisandro Dalcin PetscFunctionBegin; 435e3c11fc1SJed Brown ierr = TSHasTransientVariable(ts,&bdf->transientvar);CHKERRQ(ierr); 436211a84d6SLisandro Dalcin bdf->k = bdf->n = 0; 437e3c11fc1SJed Brown for (i=0; i<n; i++) { 438e3c11fc1SJed Brown ierr = VecDuplicate(ts->vec_sol,&bdf->work[i]);CHKERRQ(ierr); 439e3c11fc1SJed Brown if (i && bdf->transientvar) { 440e3c11fc1SJed Brown ierr = VecDuplicate(ts->vec_sol,&bdf->tvwork[i]);CHKERRQ(ierr); 441e3c11fc1SJed Brown } 442e3c11fc1SJed Brown } 443211a84d6SLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&bdf->vec_dot);CHKERRQ(ierr); 4441117961dSLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&bdf->vec_wrk);CHKERRQ(ierr); 445211a84d6SLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&bdf->vec_lte);CHKERRQ(ierr); 4461117961dSLisandro Dalcin ierr = TSGetDM(ts,&ts->dm);CHKERRQ(ierr); 4471117961dSLisandro Dalcin ierr = DMCoarsenHookAdd(ts->dm,DMCoarsenHook_TSBDF,DMRestrictHook_TSBDF,ts);CHKERRQ(ierr); 448211a84d6SLisandro Dalcin 449211a84d6SLisandro Dalcin ierr = TSGetAdapt(ts,&ts->adapt);CHKERRQ(ierr); 450211a84d6SLisandro Dalcin ierr = TSAdaptCandidatesClear(ts->adapt);CHKERRQ(ierr); 4511917a363SLisandro Dalcin ierr = TSAdaptGetClip(ts->adapt,&low,&high);CHKERRQ(ierr); 4522ffb9264SLisandro Dalcin ierr = TSAdaptSetClip(ts->adapt,low,PetscMin(high,two));CHKERRQ(ierr); 453211a84d6SLisandro Dalcin 454211a84d6SLisandro Dalcin ierr = TSGetSNES(ts,&ts->snes);CHKERRQ(ierr); 455211a84d6SLisandro Dalcin PetscFunctionReturn(0); 456211a84d6SLisandro Dalcin } 457211a84d6SLisandro Dalcin 458211a84d6SLisandro Dalcin static PetscErrorCode TSSetFromOptions_BDF(PetscOptionItems *PetscOptionsObject,TS ts) 459211a84d6SLisandro Dalcin { 460211a84d6SLisandro Dalcin PetscErrorCode ierr; 461211a84d6SLisandro Dalcin 462211a84d6SLisandro Dalcin PetscFunctionBegin; 463211a84d6SLisandro Dalcin ierr = PetscOptionsHead(PetscOptionsObject,"BDF ODE solver options");CHKERRQ(ierr); 464211a84d6SLisandro Dalcin { 465211a84d6SLisandro Dalcin PetscBool flg; 466e5b8ffdfSLisandro Dalcin PetscInt order; 467e5b8ffdfSLisandro Dalcin ierr = TSBDFGetOrder(ts,&order);CHKERRQ(ierr); 468211a84d6SLisandro Dalcin ierr = PetscOptionsInt("-ts_bdf_order","Order of the BDF method","TSBDFSetOrder",order,&order,&flg);CHKERRQ(ierr); 469211a84d6SLisandro Dalcin if (flg) {ierr = TSBDFSetOrder(ts,order);CHKERRQ(ierr);} 470211a84d6SLisandro Dalcin } 471211a84d6SLisandro Dalcin ierr = PetscOptionsTail();CHKERRQ(ierr); 472211a84d6SLisandro Dalcin PetscFunctionReturn(0); 473211a84d6SLisandro Dalcin } 474211a84d6SLisandro Dalcin 475211a84d6SLisandro Dalcin static PetscErrorCode TSView_BDF(TS ts,PetscViewer viewer) 476211a84d6SLisandro Dalcin { 477211a84d6SLisandro Dalcin TS_BDF *bdf = (TS_BDF*)ts->data; 478211a84d6SLisandro Dalcin PetscBool iascii; 479211a84d6SLisandro Dalcin PetscErrorCode ierr; 480211a84d6SLisandro Dalcin 481211a84d6SLisandro Dalcin PetscFunctionBegin; 482211a84d6SLisandro Dalcin ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 483211a84d6SLisandro Dalcin if (iascii) { 484211a84d6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," Order=%D\n",bdf->order);CHKERRQ(ierr); 485211a84d6SLisandro Dalcin } 486211a84d6SLisandro Dalcin PetscFunctionReturn(0); 487211a84d6SLisandro Dalcin } 488211a84d6SLisandro Dalcin 489211a84d6SLisandro Dalcin /* ------------------------------------------------------------ */ 490211a84d6SLisandro Dalcin 491211a84d6SLisandro Dalcin static PetscErrorCode TSBDFSetOrder_BDF(TS ts,PetscInt order) 492211a84d6SLisandro Dalcin { 493211a84d6SLisandro Dalcin TS_BDF *bdf = (TS_BDF*)ts->data; 494211a84d6SLisandro Dalcin 495211a84d6SLisandro Dalcin PetscFunctionBegin; 496211a84d6SLisandro Dalcin if (order == bdf->order) PetscFunctionReturn(0); 497211a84d6SLisandro Dalcin if (order < 1 || order > 6) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"BDF Order %D not implemented",order); 498211a84d6SLisandro Dalcin bdf->order = order; 499211a84d6SLisandro Dalcin PetscFunctionReturn(0); 500211a84d6SLisandro Dalcin } 501211a84d6SLisandro Dalcin 502211a84d6SLisandro Dalcin static PetscErrorCode TSBDFGetOrder_BDF(TS ts,PetscInt *order) 503211a84d6SLisandro Dalcin { 504211a84d6SLisandro Dalcin TS_BDF *bdf = (TS_BDF*)ts->data; 505211a84d6SLisandro Dalcin 506211a84d6SLisandro Dalcin PetscFunctionBegin; 507211a84d6SLisandro Dalcin *order = bdf->order; 508211a84d6SLisandro Dalcin PetscFunctionReturn(0); 509211a84d6SLisandro Dalcin } 510211a84d6SLisandro Dalcin 511211a84d6SLisandro Dalcin /* ------------------------------------------------------------ */ 512211a84d6SLisandro Dalcin 513211a84d6SLisandro Dalcin /*MC 514211a84d6SLisandro Dalcin TSBDF - DAE solver using BDF methods 515211a84d6SLisandro Dalcin 516211a84d6SLisandro Dalcin Level: beginner 517211a84d6SLisandro Dalcin 518211a84d6SLisandro Dalcin .seealso: TS, TSCreate(), TSSetType() 519211a84d6SLisandro Dalcin M*/ 520211a84d6SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSCreate_BDF(TS ts) 521211a84d6SLisandro Dalcin { 522211a84d6SLisandro Dalcin TS_BDF *bdf; 523211a84d6SLisandro Dalcin PetscErrorCode ierr; 524211a84d6SLisandro Dalcin 525211a84d6SLisandro Dalcin PetscFunctionBegin; 526211a84d6SLisandro Dalcin ts->ops->reset = TSReset_BDF; 527211a84d6SLisandro Dalcin ts->ops->destroy = TSDestroy_BDF; 528211a84d6SLisandro Dalcin ts->ops->view = TSView_BDF; 529211a84d6SLisandro Dalcin ts->ops->setup = TSSetUp_BDF; 530211a84d6SLisandro Dalcin ts->ops->setfromoptions = TSSetFromOptions_BDF; 531211a84d6SLisandro Dalcin ts->ops->step = TSStep_BDF; 532211a84d6SLisandro Dalcin ts->ops->evaluatewlte = TSEvaluateWLTE_BDF; 533211a84d6SLisandro Dalcin ts->ops->rollback = TSRollBack_BDF; 534211a84d6SLisandro Dalcin ts->ops->interpolate = TSInterpolate_BDF; 535211a84d6SLisandro Dalcin ts->ops->snesfunction = SNESTSFormFunction_BDF; 536211a84d6SLisandro Dalcin ts->ops->snesjacobian = SNESTSFormJacobian_BDF; 5372ffb9264SLisandro Dalcin ts->default_adapt_type = TSADAPTBASIC; 538211a84d6SLisandro Dalcin 539efd4aadfSBarry Smith ts->usessnes = PETSC_TRUE; 540efd4aadfSBarry Smith 541211a84d6SLisandro Dalcin ierr = PetscNewLog(ts,&bdf);CHKERRQ(ierr); 542211a84d6SLisandro Dalcin ts->data = (void*)bdf; 543211a84d6SLisandro Dalcin 544211a84d6SLisandro Dalcin bdf->status = TS_STEP_COMPLETE; 545211a84d6SLisandro Dalcin 546211a84d6SLisandro Dalcin ierr = PetscObjectComposeFunction((PetscObject)ts,"TSBDFSetOrder_C",TSBDFSetOrder_BDF);CHKERRQ(ierr); 547211a84d6SLisandro Dalcin ierr = PetscObjectComposeFunction((PetscObject)ts,"TSBDFGetOrder_C",TSBDFGetOrder_BDF);CHKERRQ(ierr); 548e5b8ffdfSLisandro Dalcin ierr = TSBDFSetOrder(ts,2);CHKERRQ(ierr); 549211a84d6SLisandro Dalcin PetscFunctionReturn(0); 550211a84d6SLisandro Dalcin } 551211a84d6SLisandro Dalcin 552211a84d6SLisandro Dalcin /* ------------------------------------------------------------ */ 553211a84d6SLisandro Dalcin 554211a84d6SLisandro Dalcin /*@ 555211a84d6SLisandro Dalcin TSBDFSetOrder - Set the order of the BDF method 556211a84d6SLisandro Dalcin 557211a84d6SLisandro Dalcin Logically Collective on TS 558211a84d6SLisandro Dalcin 559*d8d19677SJose E. Roman Input Parameters: 560211a84d6SLisandro Dalcin + ts - timestepping context 561211a84d6SLisandro Dalcin - order - order of the method 562211a84d6SLisandro Dalcin 563211a84d6SLisandro Dalcin Options Database: 564211a84d6SLisandro Dalcin . -ts_bdf_order <order> 565211a84d6SLisandro Dalcin 566211a84d6SLisandro Dalcin Level: intermediate 567211a84d6SLisandro Dalcin 568211a84d6SLisandro Dalcin @*/ 569211a84d6SLisandro Dalcin PetscErrorCode TSBDFSetOrder(TS ts,PetscInt order) 570211a84d6SLisandro Dalcin { 571211a84d6SLisandro Dalcin PetscErrorCode ierr; 572211a84d6SLisandro Dalcin 573211a84d6SLisandro Dalcin PetscFunctionBegin; 574211a84d6SLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,1); 575211a84d6SLisandro Dalcin PetscValidLogicalCollectiveInt(ts,order,2); 576211a84d6SLisandro Dalcin ierr = PetscTryMethod(ts,"TSBDFSetOrder_C",(TS,PetscInt),(ts,order));CHKERRQ(ierr); 577211a84d6SLisandro Dalcin PetscFunctionReturn(0); 578211a84d6SLisandro Dalcin } 579211a84d6SLisandro Dalcin 580211a84d6SLisandro Dalcin /*@ 581211a84d6SLisandro Dalcin TSBDFGetOrder - Get the order of the BDF method 582211a84d6SLisandro Dalcin 583211a84d6SLisandro Dalcin Not Collective 584211a84d6SLisandro Dalcin 585211a84d6SLisandro Dalcin Input Parameter: 586211a84d6SLisandro Dalcin . ts - timestepping context 587211a84d6SLisandro Dalcin 588211a84d6SLisandro Dalcin Output Parameter: 589211a84d6SLisandro Dalcin . order - order of the method 590211a84d6SLisandro Dalcin 591211a84d6SLisandro Dalcin Level: intermediate 592211a84d6SLisandro Dalcin 593211a84d6SLisandro Dalcin @*/ 594211a84d6SLisandro Dalcin PetscErrorCode TSBDFGetOrder(TS ts,PetscInt *order) 595211a84d6SLisandro Dalcin { 596211a84d6SLisandro Dalcin PetscErrorCode ierr; 597211a84d6SLisandro Dalcin 598211a84d6SLisandro Dalcin PetscFunctionBegin; 599211a84d6SLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,1); 600211a84d6SLisandro Dalcin PetscValidIntPointer(order,2); 601211a84d6SLisandro Dalcin ierr = PetscUseMethod(ts,"TSBDFGetOrder_C",(TS,PetscInt*),(ts,order));CHKERRQ(ierr); 602211a84d6SLisandro Dalcin PetscFunctionReturn(0); 603211a84d6SLisandro Dalcin } 604