1211a84d6SLisandro Dalcin /* 2211a84d6SLisandro Dalcin Code for timestepping with BDF methods 3211a84d6SLisandro Dalcin */ 4211a84d6SLisandro Dalcin #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 5*1117961dSLisandro 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]; 20211a84d6SLisandro Dalcin PetscReal shift; 21211a84d6SLisandro Dalcin Vec vec_dot; 22*1117961dSLisandro Dalcin Vec vec_wrk; 23211a84d6SLisandro Dalcin Vec vec_lte; 24211a84d6SLisandro Dalcin 25211a84d6SLisandro Dalcin PetscInt order; 26211a84d6SLisandro Dalcin TSStepStatus status; 27211a84d6SLisandro Dalcin } TS_BDF; 28211a84d6SLisandro Dalcin 29211a84d6SLisandro Dalcin 30211a84d6SLisandro Dalcin PETSC_STATIC_INLINE void LagrangeBasisVals(PetscInt n,PetscReal t,const PetscReal T[],PetscScalar L[]) 31211a84d6SLisandro Dalcin { 32211a84d6SLisandro Dalcin PetscInt k,j; 33211a84d6SLisandro Dalcin for (k=0; k<n; k++) 34211a84d6SLisandro Dalcin for (L[k]=1, j=0; j<n; j++) 35211a84d6SLisandro Dalcin if (j != k) 36211a84d6SLisandro Dalcin L[k] *= (t - T[j])/(T[k] - T[j]); 37211a84d6SLisandro Dalcin } 38211a84d6SLisandro Dalcin 39211a84d6SLisandro Dalcin PETSC_STATIC_INLINE void LagrangeBasisDers(PetscInt n,PetscReal t,const PetscReal T[],PetscScalar dL[]) 40211a84d6SLisandro Dalcin { 41211a84d6SLisandro Dalcin PetscInt k,j,i; 42211a84d6SLisandro Dalcin for (k=0; k<n; k++) 43211a84d6SLisandro Dalcin for (dL[k]=0, j=0; j<n; j++) 44211a84d6SLisandro Dalcin if (j != k) { 45211a84d6SLisandro Dalcin PetscReal L = 1/(T[k] - T[j]); 46211a84d6SLisandro Dalcin for (i=0; i<n; i++) 47211a84d6SLisandro Dalcin if (i != j && i != k) 48211a84d6SLisandro Dalcin L *= (t - T[i])/(T[k] - T[i]); 49211a84d6SLisandro Dalcin dL[k] += L; 50211a84d6SLisandro Dalcin } 51211a84d6SLisandro Dalcin } 52211a84d6SLisandro Dalcin 53*1117961dSLisandro Dalcin static PetscErrorCode TSBDF_GetVecs(TS ts,DM dm,Vec *Xdot,Vec *Ydot) 54*1117961dSLisandro Dalcin { 55*1117961dSLisandro Dalcin TS_BDF *bdf = (TS_BDF*)ts->data; 56*1117961dSLisandro Dalcin PetscErrorCode ierr; 57*1117961dSLisandro Dalcin 58*1117961dSLisandro Dalcin PetscFunctionBegin; 59*1117961dSLisandro Dalcin if (dm && dm != ts->dm) { 60*1117961dSLisandro Dalcin ierr = DMGetNamedGlobalVector(dm,"TSBDF_Vec_Xdot",Xdot);CHKERRQ(ierr); 61*1117961dSLisandro Dalcin ierr = DMGetNamedGlobalVector(dm,"TSBDF_Vec_Ydot",Ydot);CHKERRQ(ierr); 62*1117961dSLisandro Dalcin } else { 63*1117961dSLisandro Dalcin *Xdot = bdf->vec_dot; 64*1117961dSLisandro Dalcin *Ydot = bdf->vec_wrk; 65*1117961dSLisandro Dalcin } 66*1117961dSLisandro Dalcin PetscFunctionReturn(0); 67*1117961dSLisandro Dalcin } 68*1117961dSLisandro Dalcin 69*1117961dSLisandro Dalcin static PetscErrorCode TSBDF_RestoreVecs(TS ts,DM dm,Vec *Xdot,Vec *Ydot) 70*1117961dSLisandro Dalcin { 71*1117961dSLisandro Dalcin TS_BDF *bdf = (TS_BDF*)ts->data; 72*1117961dSLisandro Dalcin PetscErrorCode ierr; 73*1117961dSLisandro Dalcin 74*1117961dSLisandro Dalcin PetscFunctionBegin; 75*1117961dSLisandro Dalcin if (dm && dm != ts->dm) { 76*1117961dSLisandro Dalcin ierr = DMRestoreNamedGlobalVector(dm,"TSBDF_Vec_Xdot",Xdot);CHKERRQ(ierr); 77*1117961dSLisandro Dalcin ierr = DMRestoreNamedGlobalVector(dm,"TSBDF_Vec_Ydot",Ydot);CHKERRQ(ierr); 78*1117961dSLisandro Dalcin } else { 79*1117961dSLisandro Dalcin if (*Xdot != bdf->vec_dot) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_INCOMP,"Vec does not match the cache"); 80*1117961dSLisandro Dalcin if (*Ydot != bdf->vec_wrk) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_INCOMP,"Vec does not match the cache"); 81*1117961dSLisandro Dalcin *Xdot = NULL; 82*1117961dSLisandro Dalcin *Ydot = NULL; 83*1117961dSLisandro Dalcin } 84*1117961dSLisandro Dalcin PetscFunctionReturn(0); 85*1117961dSLisandro Dalcin } 86*1117961dSLisandro Dalcin 87*1117961dSLisandro Dalcin static PetscErrorCode DMCoarsenHook_TSBDF(DM fine,DM coarse,void *ctx) 88*1117961dSLisandro Dalcin { 89*1117961dSLisandro Dalcin PetscFunctionBegin; 90*1117961dSLisandro Dalcin PetscFunctionReturn(0); 91*1117961dSLisandro Dalcin } 92*1117961dSLisandro Dalcin 93*1117961dSLisandro Dalcin static PetscErrorCode DMRestrictHook_TSBDF(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx) 94*1117961dSLisandro Dalcin { 95*1117961dSLisandro Dalcin TS ts = (TS)ctx; 96*1117961dSLisandro Dalcin Vec Ydot,Ydot_c; 97*1117961dSLisandro Dalcin Vec Xdot,Xdot_c; 98*1117961dSLisandro Dalcin PetscErrorCode ierr; 99*1117961dSLisandro Dalcin 100*1117961dSLisandro Dalcin PetscFunctionBegin; 101*1117961dSLisandro Dalcin ierr = TSBDF_GetVecs(ts,fine,&Xdot,&Ydot);CHKERRQ(ierr); 102*1117961dSLisandro Dalcin ierr = TSBDF_GetVecs(ts,coarse,&Xdot_c,&Ydot_c);CHKERRQ(ierr); 103*1117961dSLisandro Dalcin 104*1117961dSLisandro Dalcin ierr = MatRestrict(restrct,Ydot,Ydot_c);CHKERRQ(ierr); 105*1117961dSLisandro Dalcin ierr = VecPointwiseMult(Ydot_c,rscale,Ydot_c);CHKERRQ(ierr); 106*1117961dSLisandro Dalcin 107*1117961dSLisandro Dalcin ierr = TSBDF_RestoreVecs(ts,fine,&Xdot,&Ydot);CHKERRQ(ierr); 108*1117961dSLisandro Dalcin ierr = TSBDF_RestoreVecs(ts,coarse,&Xdot_c,&Ydot_c);CHKERRQ(ierr); 109*1117961dSLisandro Dalcin PetscFunctionReturn(0); 110*1117961dSLisandro Dalcin } 111*1117961dSLisandro Dalcin 112211a84d6SLisandro Dalcin static PetscErrorCode TSBDF_Advance(TS ts,PetscReal t,Vec X) 113211a84d6SLisandro Dalcin { 114211a84d6SLisandro Dalcin TS_BDF *bdf = (TS_BDF*)ts->data; 115211a84d6SLisandro Dalcin PetscInt i,n = (PetscInt)(sizeof(bdf->work)/sizeof(Vec)); 116211a84d6SLisandro Dalcin Vec tail = bdf->work[n-1]; 117211a84d6SLisandro Dalcin PetscErrorCode ierr; 118211a84d6SLisandro Dalcin 119211a84d6SLisandro Dalcin PetscFunctionBegin; 120211a84d6SLisandro Dalcin for (i=n-1; i>=2; i--) { 121211a84d6SLisandro Dalcin bdf->time[i] = bdf->time[i-1]; 122211a84d6SLisandro Dalcin bdf->work[i] = bdf->work[i-1]; 123211a84d6SLisandro Dalcin } 124211a84d6SLisandro Dalcin bdf->n = PetscMin(bdf->n+1,n-1); 125211a84d6SLisandro Dalcin bdf->time[1] = t; 126211a84d6SLisandro Dalcin bdf->work[1] = tail; 127211a84d6SLisandro Dalcin ierr = VecCopy(X,tail);CHKERRQ(ierr); 128211a84d6SLisandro Dalcin PetscFunctionReturn(0); 129211a84d6SLisandro Dalcin } 130211a84d6SLisandro Dalcin 131211a84d6SLisandro Dalcin static PetscErrorCode TSBDF_VecLTE(TS ts,PetscInt order,Vec lte) 132211a84d6SLisandro Dalcin { 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 PetscErrorCode ierr; 139211a84d6SLisandro Dalcin 140211a84d6SLisandro Dalcin PetscFunctionBegin; 141211a84d6SLisandro Dalcin LagrangeBasisDers(n+0,time[0],time,a); 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]; 144211a84d6SLisandro Dalcin ierr = VecZeroEntries(lte);CHKERRQ(ierr); 145211a84d6SLisandro Dalcin ierr = VecMAXPY(lte,n+1,alpha,vecs);CHKERRQ(ierr); 146211a84d6SLisandro Dalcin PetscFunctionReturn(0); 147211a84d6SLisandro Dalcin } 148211a84d6SLisandro Dalcin 149211a84d6SLisandro Dalcin static PetscErrorCode TSBDF_Extrapolate(TS ts,PetscInt order,PetscReal t,Vec X) 150211a84d6SLisandro Dalcin { 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 PetscErrorCode ierr; 157211a84d6SLisandro Dalcin 158211a84d6SLisandro Dalcin PetscFunctionBegin; 159211a84d6SLisandro Dalcin n = PetscMin(n,bdf->n); 160211a84d6SLisandro Dalcin LagrangeBasisVals(n,t,time,alpha); 161211a84d6SLisandro Dalcin ierr = VecZeroEntries(X);CHKERRQ(ierr); 162211a84d6SLisandro Dalcin ierr = VecMAXPY(X,n,alpha,vecs);CHKERRQ(ierr); 163211a84d6SLisandro Dalcin PetscFunctionReturn(0); 164211a84d6SLisandro Dalcin } 165211a84d6SLisandro Dalcin 166211a84d6SLisandro Dalcin static PetscErrorCode TSBDF_Interpolate(TS ts,PetscInt order,PetscReal t,Vec X) 167211a84d6SLisandro Dalcin { 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 PetscErrorCode ierr; 174211a84d6SLisandro Dalcin 175211a84d6SLisandro Dalcin PetscFunctionBegin; 176211a84d6SLisandro Dalcin LagrangeBasisVals(n,t,time,alpha); 177211a84d6SLisandro Dalcin ierr = VecZeroEntries(X);CHKERRQ(ierr); 178211a84d6SLisandro Dalcin ierr = VecMAXPY(X,n,alpha,vecs);CHKERRQ(ierr); 179211a84d6SLisandro Dalcin PetscFunctionReturn(0); 180211a84d6SLisandro Dalcin } 181211a84d6SLisandro Dalcin 182*1117961dSLisandro Dalcin static PetscErrorCode TSBDF_PreSolve(TS ts) 183*1117961dSLisandro Dalcin { 184*1117961dSLisandro Dalcin TS_BDF *bdf = (TS_BDF*)ts->data; 185*1117961dSLisandro Dalcin PetscInt i,n = PetscMax(bdf->k,1) + 1; 186*1117961dSLisandro Dalcin Vec V,V0; 187*1117961dSLisandro Dalcin Vec vecs[7]; 188*1117961dSLisandro Dalcin PetscScalar alpha[7]; 189*1117961dSLisandro Dalcin PetscErrorCode ierr; 190*1117961dSLisandro Dalcin 191*1117961dSLisandro Dalcin PetscFunctionBegin; 192*1117961dSLisandro Dalcin ierr = TSBDF_GetVecs(ts,NULL,&V,&V0);CHKERRQ(ierr); 193*1117961dSLisandro Dalcin LagrangeBasisDers(n,bdf->time[0],bdf->time,alpha); 194*1117961dSLisandro Dalcin for (i=1; i<n; i++) vecs[i] = bdf->work[i]; 195*1117961dSLisandro Dalcin ierr = VecZeroEntries(V0);CHKERRQ(ierr); 196*1117961dSLisandro Dalcin ierr = VecMAXPY(V0,n-1,alpha+1,vecs+1);CHKERRQ(ierr); 197*1117961dSLisandro Dalcin bdf->shift = PetscRealPart(alpha[0]); 198*1117961dSLisandro Dalcin ierr = TSBDF_RestoreVecs(ts,NULL,&V,&V0);CHKERRQ(ierr); 199*1117961dSLisandro Dalcin PetscFunctionReturn(0); 200*1117961dSLisandro Dalcin } 201*1117961dSLisandro Dalcin 202470880abSPatrick Sanan static PetscErrorCode TSBDF_SNESSolve(TS ts,Vec b,Vec x) 203211a84d6SLisandro Dalcin { 204211a84d6SLisandro Dalcin PetscInt nits,lits; 205211a84d6SLisandro Dalcin PetscErrorCode ierr; 206211a84d6SLisandro Dalcin 207211a84d6SLisandro Dalcin PetscFunctionBegin; 208*1117961dSLisandro Dalcin ierr = TSBDF_PreSolve(ts);CHKERRQ(ierr); 209211a84d6SLisandro Dalcin ierr = SNESSolve(ts->snes,b,x);CHKERRQ(ierr); 210211a84d6SLisandro Dalcin ierr = SNESGetIterationNumber(ts->snes,&nits);CHKERRQ(ierr); 211211a84d6SLisandro Dalcin ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr); 212211a84d6SLisandro Dalcin ts->snes_its += nits; ts->ksp_its += lits; 213211a84d6SLisandro Dalcin PetscFunctionReturn(0); 214211a84d6SLisandro Dalcin } 215211a84d6SLisandro Dalcin 216211a84d6SLisandro Dalcin static PetscErrorCode TSBDF_Restart(TS ts,PetscBool *accept) 217211a84d6SLisandro Dalcin { 218211a84d6SLisandro Dalcin TS_BDF *bdf = (TS_BDF*)ts->data; 219211a84d6SLisandro Dalcin PetscErrorCode ierr; 220211a84d6SLisandro Dalcin 221211a84d6SLisandro Dalcin PetscFunctionBegin; 222211a84d6SLisandro Dalcin bdf->k = 1; bdf->n = 0; 223211a84d6SLisandro Dalcin ierr = TSBDF_Advance(ts,ts->ptime,ts->vec_sol);CHKERRQ(ierr); 224211a84d6SLisandro Dalcin 225211a84d6SLisandro Dalcin bdf->time[0] = ts->ptime + ts->time_step/2; 226211a84d6SLisandro Dalcin ierr = VecCopy(bdf->work[1],bdf->work[0]);CHKERRQ(ierr); 227211a84d6SLisandro Dalcin ierr = TSPreStage(ts,bdf->time[0]);CHKERRQ(ierr); 228470880abSPatrick Sanan ierr = TSBDF_SNESSolve(ts,NULL,bdf->work[0]);CHKERRQ(ierr); 229211a84d6SLisandro Dalcin ierr = TSPostStage(ts,bdf->time[0],0,&bdf->work[0]);CHKERRQ(ierr); 230211a84d6SLisandro Dalcin ierr = TSAdaptCheckStage(ts->adapt,ts,bdf->time[0],bdf->work[0],accept);CHKERRQ(ierr); 231211a84d6SLisandro Dalcin if (!*accept) PetscFunctionReturn(0); 232211a84d6SLisandro Dalcin 233211a84d6SLisandro Dalcin bdf->k = PetscMin(2,bdf->order); bdf->n++; 234211a84d6SLisandro Dalcin ierr = VecCopy(bdf->work[0],bdf->work[2]);CHKERRQ(ierr); 235211a84d6SLisandro Dalcin bdf->time[2] = bdf->time[0]; 236211a84d6SLisandro Dalcin PetscFunctionReturn(0); 237211a84d6SLisandro Dalcin } 238211a84d6SLisandro Dalcin 239211a84d6SLisandro Dalcin static const char *const BDF_SchemeName[] = {"", "1", "2", "3", "4", "5", "6"}; 240211a84d6SLisandro Dalcin 241211a84d6SLisandro Dalcin static PetscErrorCode TSStep_BDF(TS ts) 242211a84d6SLisandro Dalcin { 243211a84d6SLisandro Dalcin TS_BDF *bdf = (TS_BDF*)ts->data; 244211a84d6SLisandro Dalcin PetscInt rejections = 0; 245211a84d6SLisandro Dalcin PetscBool stageok,accept = PETSC_TRUE; 246211a84d6SLisandro Dalcin PetscReal next_time_step = ts->time_step; 247211a84d6SLisandro Dalcin PetscErrorCode ierr; 248211a84d6SLisandro Dalcin 249211a84d6SLisandro Dalcin PetscFunctionBegin; 250211a84d6SLisandro Dalcin ierr = PetscCitationsRegister(citation,&cited);CHKERRQ(ierr); 251211a84d6SLisandro Dalcin 252211a84d6SLisandro Dalcin if (!ts->steprollback && !ts->steprestart) { 253211a84d6SLisandro Dalcin bdf->k = PetscMin(bdf->k+1,bdf->order); 254211a84d6SLisandro Dalcin ierr = TSBDF_Advance(ts,ts->ptime,ts->vec_sol);CHKERRQ(ierr); 255211a84d6SLisandro Dalcin } 256211a84d6SLisandro Dalcin 257211a84d6SLisandro Dalcin bdf->status = TS_STEP_INCOMPLETE; 258211a84d6SLisandro Dalcin while (!ts->reason && bdf->status != TS_STEP_COMPLETE) { 259211a84d6SLisandro Dalcin 260211a84d6SLisandro Dalcin if (ts->steprestart) { 261211a84d6SLisandro Dalcin ierr = TSBDF_Restart(ts,&stageok);CHKERRQ(ierr); 262211a84d6SLisandro Dalcin if (!stageok) goto reject_step; 263211a84d6SLisandro Dalcin } 264211a84d6SLisandro Dalcin 265211a84d6SLisandro Dalcin bdf->time[0] = ts->ptime + ts->time_step; 266211a84d6SLisandro Dalcin ierr = TSBDF_Extrapolate(ts,bdf->k-(accept?0:1),bdf->time[0],bdf->work[0]);CHKERRQ(ierr); 267211a84d6SLisandro Dalcin ierr = TSPreStage(ts,bdf->time[0]);CHKERRQ(ierr); 268470880abSPatrick Sanan ierr = TSBDF_SNESSolve(ts,NULL,bdf->work[0]);CHKERRQ(ierr); 269211a84d6SLisandro Dalcin ierr = TSPostStage(ts,bdf->time[0],0,&bdf->work[0]);CHKERRQ(ierr); 270211a84d6SLisandro Dalcin ierr = TSAdaptCheckStage(ts->adapt,ts,bdf->time[0],bdf->work[0],&stageok);CHKERRQ(ierr); 271211a84d6SLisandro Dalcin if (!stageok) goto reject_step; 272211a84d6SLisandro Dalcin 273211a84d6SLisandro Dalcin bdf->status = TS_STEP_PENDING; 274211a84d6SLisandro Dalcin ierr = TSAdaptCandidatesClear(ts->adapt);CHKERRQ(ierr); 275211a84d6SLisandro Dalcin ierr = TSAdaptCandidateAdd(ts->adapt,BDF_SchemeName[bdf->k],bdf->k,1,1.0,1.0,PETSC_TRUE);CHKERRQ(ierr); 276211a84d6SLisandro Dalcin ierr = TSAdaptChoose(ts->adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr); 277211a84d6SLisandro Dalcin bdf->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE; 278211a84d6SLisandro Dalcin if (!accept) { ts->time_step = next_time_step; goto reject_step; } 279211a84d6SLisandro Dalcin 280211a84d6SLisandro Dalcin ierr = VecCopy(bdf->work[0],ts->vec_sol);CHKERRQ(ierr); 281211a84d6SLisandro Dalcin ts->ptime += ts->time_step; 282211a84d6SLisandro Dalcin ts->time_step = next_time_step; 283211a84d6SLisandro Dalcin break; 284211a84d6SLisandro Dalcin 285211a84d6SLisandro Dalcin reject_step: 286211a84d6SLisandro Dalcin ts->reject++; accept = PETSC_FALSE; 287211a84d6SLisandro Dalcin if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) { 288211a84d6SLisandro Dalcin ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);CHKERRQ(ierr); 289211a84d6SLisandro Dalcin ts->reason = TS_DIVERGED_STEP_REJECTED; 290211a84d6SLisandro Dalcin } 291211a84d6SLisandro Dalcin } 292211a84d6SLisandro Dalcin PetscFunctionReturn(0); 293211a84d6SLisandro Dalcin } 294211a84d6SLisandro Dalcin 295211a84d6SLisandro Dalcin static PetscErrorCode TSInterpolate_BDF(TS ts,PetscReal t,Vec X) 296211a84d6SLisandro Dalcin { 297211a84d6SLisandro Dalcin TS_BDF *bdf = (TS_BDF*)ts->data; 298211a84d6SLisandro Dalcin PetscErrorCode ierr; 299211a84d6SLisandro Dalcin 300211a84d6SLisandro Dalcin PetscFunctionBegin; 301211a84d6SLisandro Dalcin ierr = TSBDF_Interpolate(ts,bdf->k,t,X);CHKERRQ(ierr); 302211a84d6SLisandro Dalcin PetscFunctionReturn(0); 303211a84d6SLisandro Dalcin } 304211a84d6SLisandro Dalcin 305211a84d6SLisandro Dalcin static PetscErrorCode TSEvaluateWLTE_BDF(TS ts,NormType wnormtype,PetscInt *order,PetscReal *wlte) 306211a84d6SLisandro Dalcin { 307211a84d6SLisandro Dalcin TS_BDF *bdf = (TS_BDF*)ts->data; 308211a84d6SLisandro Dalcin PetscInt k = bdf->k; 3097453f775SEmil Constantinescu PetscReal wltea,wlter; 310211a84d6SLisandro Dalcin Vec X = bdf->work[0], Y = bdf->vec_lte; 311211a84d6SLisandro Dalcin PetscErrorCode ierr; 312211a84d6SLisandro Dalcin 313211a84d6SLisandro Dalcin PetscFunctionBegin; 314211a84d6SLisandro Dalcin k = PetscMin(k,bdf->n-1); 315211a84d6SLisandro Dalcin ierr = TSBDF_VecLTE(ts,k,Y);CHKERRQ(ierr); 316211a84d6SLisandro Dalcin ierr = VecAXPY(Y,1,X);CHKERRQ(ierr); 3177453f775SEmil Constantinescu ierr = TSErrorWeightedNorm(ts,X,Y,wnormtype,wlte,&wltea,&wlter);CHKERRQ(ierr); 318211a84d6SLisandro Dalcin if (order) *order = k + 1; 319211a84d6SLisandro Dalcin PetscFunctionReturn(0); 320211a84d6SLisandro Dalcin } 321211a84d6SLisandro Dalcin 322211a84d6SLisandro Dalcin static PetscErrorCode TSRollBack_BDF(TS ts) 323211a84d6SLisandro Dalcin { 324211a84d6SLisandro Dalcin TS_BDF *bdf = (TS_BDF*)ts->data; 325211a84d6SLisandro Dalcin PetscErrorCode ierr; 326211a84d6SLisandro Dalcin 327211a84d6SLisandro Dalcin PetscFunctionBegin; 328211a84d6SLisandro Dalcin ierr = VecCopy(bdf->work[1],ts->vec_sol);CHKERRQ(ierr); 329211a84d6SLisandro Dalcin PetscFunctionReturn(0); 330211a84d6SLisandro Dalcin } 331211a84d6SLisandro Dalcin 332*1117961dSLisandro Dalcin static PetscErrorCode SNESTSFormFunction_BDF(SNES snes,Vec X,Vec F,TS ts) 333211a84d6SLisandro Dalcin { 334211a84d6SLisandro Dalcin TS_BDF *bdf = (TS_BDF*)ts->data; 335*1117961dSLisandro Dalcin DM dm, dmsave = ts->dm; 336211a84d6SLisandro Dalcin PetscReal t = bdf->time[0]; 337*1117961dSLisandro Dalcin PetscReal shift = bdf->shift; 338*1117961dSLisandro Dalcin Vec V,V0; 339211a84d6SLisandro Dalcin PetscErrorCode ierr; 340211a84d6SLisandro Dalcin 341211a84d6SLisandro Dalcin PetscFunctionBegin; 342*1117961dSLisandro Dalcin ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 343*1117961dSLisandro Dalcin ierr = TSBDF_GetVecs(ts,dm,&V,&V0);CHKERRQ(ierr); 344*1117961dSLisandro Dalcin ierr = VecWAXPY(V,shift,X,V0);CHKERRQ(ierr); 345*1117961dSLisandro Dalcin 346211a84d6SLisandro Dalcin /* F = Function(t,X,V) */ 347*1117961dSLisandro Dalcin ts->dm = dm; 348211a84d6SLisandro Dalcin ierr = TSComputeIFunction(ts,t,X,V,F,PETSC_FALSE);CHKERRQ(ierr); 349*1117961dSLisandro Dalcin ts->dm = dmsave; 350*1117961dSLisandro Dalcin 351*1117961dSLisandro Dalcin ierr = TSBDF_RestoreVecs(ts,dm,&V,&V0);CHKERRQ(ierr); 352211a84d6SLisandro Dalcin PetscFunctionReturn(0); 353211a84d6SLisandro Dalcin } 354211a84d6SLisandro Dalcin 355*1117961dSLisandro Dalcin static PetscErrorCode SNESTSFormJacobian_BDF(SNES snes,Vec X,Mat J,Mat P,TS ts) 356211a84d6SLisandro Dalcin { 357211a84d6SLisandro Dalcin TS_BDF *bdf = (TS_BDF*)ts->data; 358*1117961dSLisandro Dalcin DM dm, dmsave = ts->dm; 359211a84d6SLisandro Dalcin PetscReal t = bdf->time[0]; 360*1117961dSLisandro Dalcin PetscReal shift = bdf->shift; 361*1117961dSLisandro Dalcin Vec V,V0; 362211a84d6SLisandro Dalcin PetscErrorCode ierr; 363211a84d6SLisandro Dalcin 364211a84d6SLisandro Dalcin PetscFunctionBegin; 365*1117961dSLisandro Dalcin ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 366*1117961dSLisandro Dalcin ierr = TSBDF_GetVecs(ts,dm,&V,&V0);CHKERRQ(ierr); 367*1117961dSLisandro Dalcin 368211a84d6SLisandro Dalcin /* J,P = Jacobian(t,X,V) */ 369*1117961dSLisandro Dalcin ts->dm = dm; 370*1117961dSLisandro Dalcin ierr = TSComputeIJacobian(ts,t,X,V,shift,J,P,PETSC_FALSE);CHKERRQ(ierr); 371*1117961dSLisandro Dalcin ts->dm = dmsave; 372*1117961dSLisandro Dalcin 373*1117961dSLisandro Dalcin ierr = TSBDF_RestoreVecs(ts,dm,&V,&V0);CHKERRQ(ierr); 374211a84d6SLisandro Dalcin PetscFunctionReturn(0); 375211a84d6SLisandro Dalcin } 376211a84d6SLisandro Dalcin 377211a84d6SLisandro Dalcin static PetscErrorCode TSReset_BDF(TS ts) 378211a84d6SLisandro Dalcin { 379211a84d6SLisandro Dalcin TS_BDF *bdf = (TS_BDF*)ts->data; 380211a84d6SLisandro Dalcin size_t i,n = sizeof(bdf->work)/sizeof(Vec); 381211a84d6SLisandro Dalcin PetscErrorCode ierr; 382211a84d6SLisandro Dalcin 383211a84d6SLisandro Dalcin PetscFunctionBegin; 384*1117961dSLisandro Dalcin bdf->k = bdf->n = 0; 385211a84d6SLisandro Dalcin for (i=0; i<n; i++) {ierr = VecDestroy(&bdf->work[i]);CHKERRQ(ierr);} 386211a84d6SLisandro Dalcin ierr = VecDestroy(&bdf->vec_dot);CHKERRQ(ierr); 387*1117961dSLisandro Dalcin ierr = VecDestroy(&bdf->vec_wrk);CHKERRQ(ierr); 388211a84d6SLisandro Dalcin ierr = VecDestroy(&bdf->vec_lte);CHKERRQ(ierr); 389*1117961dSLisandro Dalcin if (ts->dm) {ierr = DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSBDF,DMRestrictHook_TSBDF,ts);CHKERRQ(ierr);} 390211a84d6SLisandro Dalcin PetscFunctionReturn(0); 391211a84d6SLisandro Dalcin } 392211a84d6SLisandro Dalcin 393211a84d6SLisandro Dalcin static PetscErrorCode TSDestroy_BDF(TS ts) 394211a84d6SLisandro Dalcin { 395211a84d6SLisandro Dalcin PetscErrorCode ierr; 396211a84d6SLisandro Dalcin 397211a84d6SLisandro Dalcin PetscFunctionBegin; 398211a84d6SLisandro Dalcin ierr = TSReset_BDF(ts);CHKERRQ(ierr); 399211a84d6SLisandro Dalcin ierr = PetscFree(ts->data);CHKERRQ(ierr); 400211a84d6SLisandro Dalcin ierr = PetscObjectComposeFunction((PetscObject)ts,"TSBDFSetOrder_C",NULL);CHKERRQ(ierr); 401211a84d6SLisandro Dalcin ierr = PetscObjectComposeFunction((PetscObject)ts,"TSBDFGetOrder_C",NULL);CHKERRQ(ierr); 402211a84d6SLisandro Dalcin PetscFunctionReturn(0); 403211a84d6SLisandro Dalcin } 404211a84d6SLisandro Dalcin 405211a84d6SLisandro Dalcin static PetscErrorCode TSSetUp_BDF(TS ts) 406211a84d6SLisandro Dalcin { 407211a84d6SLisandro Dalcin TS_BDF *bdf = (TS_BDF*)ts->data; 408211a84d6SLisandro Dalcin size_t i,n = sizeof(bdf->work)/sizeof(Vec); 4092ffb9264SLisandro Dalcin PetscReal low,high,two = 2; 410211a84d6SLisandro Dalcin PetscErrorCode ierr; 411211a84d6SLisandro Dalcin 412211a84d6SLisandro Dalcin PetscFunctionBegin; 413211a84d6SLisandro Dalcin bdf->k = bdf->n = 0; 414211a84d6SLisandro Dalcin for (i=0; i<n; i++) {ierr = VecDuplicate(ts->vec_sol,&bdf->work[i]);CHKERRQ(ierr);} 415211a84d6SLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&bdf->vec_dot);CHKERRQ(ierr); 416*1117961dSLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&bdf->vec_wrk);CHKERRQ(ierr); 417211a84d6SLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&bdf->vec_lte);CHKERRQ(ierr); 418*1117961dSLisandro Dalcin ierr = TSGetDM(ts,&ts->dm);CHKERRQ(ierr); 419*1117961dSLisandro Dalcin ierr = DMCoarsenHookAdd(ts->dm,DMCoarsenHook_TSBDF,DMRestrictHook_TSBDF,ts);CHKERRQ(ierr); 420211a84d6SLisandro Dalcin 421211a84d6SLisandro Dalcin ierr = TSGetAdapt(ts,&ts->adapt);CHKERRQ(ierr); 422211a84d6SLisandro Dalcin ierr = TSAdaptCandidatesClear(ts->adapt);CHKERRQ(ierr); 4231917a363SLisandro Dalcin ierr = TSAdaptGetClip(ts->adapt,&low,&high);CHKERRQ(ierr); 4242ffb9264SLisandro Dalcin ierr = TSAdaptSetClip(ts->adapt,low,PetscMin(high,two));CHKERRQ(ierr); 425211a84d6SLisandro Dalcin 426211a84d6SLisandro Dalcin ierr = TSGetSNES(ts,&ts->snes);CHKERRQ(ierr); 427211a84d6SLisandro Dalcin PetscFunctionReturn(0); 428211a84d6SLisandro Dalcin } 429211a84d6SLisandro Dalcin 430211a84d6SLisandro Dalcin static PetscErrorCode TSSetFromOptions_BDF(PetscOptionItems *PetscOptionsObject,TS ts) 431211a84d6SLisandro Dalcin { 432211a84d6SLisandro Dalcin TS_BDF *bdf = (TS_BDF*)ts->data; 433211a84d6SLisandro Dalcin PetscErrorCode ierr; 434211a84d6SLisandro Dalcin 435211a84d6SLisandro Dalcin PetscFunctionBegin; 436211a84d6SLisandro Dalcin ierr = PetscOptionsHead(PetscOptionsObject,"BDF ODE solver options");CHKERRQ(ierr); 437211a84d6SLisandro Dalcin { 438211a84d6SLisandro Dalcin PetscBool flg; 439211a84d6SLisandro Dalcin PetscInt order = bdf->order; 440211a84d6SLisandro Dalcin ierr = PetscOptionsInt("-ts_bdf_order","Order of the BDF method","TSBDFSetOrder",order,&order,&flg);CHKERRQ(ierr); 441211a84d6SLisandro Dalcin if (flg) {ierr = TSBDFSetOrder(ts,order);CHKERRQ(ierr);} 442211a84d6SLisandro Dalcin } 443211a84d6SLisandro Dalcin ierr = PetscOptionsTail();CHKERRQ(ierr); 444211a84d6SLisandro Dalcin PetscFunctionReturn(0); 445211a84d6SLisandro Dalcin } 446211a84d6SLisandro Dalcin 447211a84d6SLisandro Dalcin static PetscErrorCode TSView_BDF(TS ts,PetscViewer viewer) 448211a84d6SLisandro Dalcin { 449211a84d6SLisandro Dalcin TS_BDF *bdf = (TS_BDF*)ts->data; 450211a84d6SLisandro Dalcin PetscBool iascii; 451211a84d6SLisandro Dalcin PetscErrorCode ierr; 452211a84d6SLisandro Dalcin 453211a84d6SLisandro Dalcin PetscFunctionBegin; 454211a84d6SLisandro Dalcin ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 455211a84d6SLisandro Dalcin if (iascii) { 456211a84d6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," Order=%D\n",bdf->order);CHKERRQ(ierr); 457211a84d6SLisandro Dalcin } 458211a84d6SLisandro Dalcin PetscFunctionReturn(0); 459211a84d6SLisandro Dalcin } 460211a84d6SLisandro Dalcin 461211a84d6SLisandro Dalcin /* ------------------------------------------------------------ */ 462211a84d6SLisandro Dalcin 463211a84d6SLisandro Dalcin static PetscErrorCode TSBDFSetOrder_BDF(TS ts,PetscInt order) 464211a84d6SLisandro Dalcin { 465211a84d6SLisandro Dalcin TS_BDF *bdf = (TS_BDF*)ts->data; 466211a84d6SLisandro Dalcin 467211a84d6SLisandro Dalcin PetscFunctionBegin; 468211a84d6SLisandro Dalcin if (order == bdf->order) PetscFunctionReturn(0); 469211a84d6SLisandro Dalcin if (order < 1 || order > 6) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"BDF Order %D not implemented",order); 470211a84d6SLisandro Dalcin bdf->order = order; 471211a84d6SLisandro Dalcin PetscFunctionReturn(0); 472211a84d6SLisandro Dalcin } 473211a84d6SLisandro Dalcin 474211a84d6SLisandro Dalcin static PetscErrorCode TSBDFGetOrder_BDF(TS ts,PetscInt *order) 475211a84d6SLisandro Dalcin { 476211a84d6SLisandro Dalcin TS_BDF *bdf = (TS_BDF*)ts->data; 477211a84d6SLisandro Dalcin 478211a84d6SLisandro Dalcin PetscFunctionBegin; 479211a84d6SLisandro Dalcin *order = bdf->order; 480211a84d6SLisandro Dalcin PetscFunctionReturn(0); 481211a84d6SLisandro Dalcin } 482211a84d6SLisandro Dalcin 483211a84d6SLisandro Dalcin /* ------------------------------------------------------------ */ 484211a84d6SLisandro Dalcin 485211a84d6SLisandro Dalcin /*MC 486211a84d6SLisandro Dalcin TSBDF - DAE solver using BDF methods 487211a84d6SLisandro Dalcin 488211a84d6SLisandro Dalcin Level: beginner 489211a84d6SLisandro Dalcin 490211a84d6SLisandro Dalcin .seealso: TS, TSCreate(), TSSetType() 491211a84d6SLisandro Dalcin M*/ 492211a84d6SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSCreate_BDF(TS ts) 493211a84d6SLisandro Dalcin { 494211a84d6SLisandro Dalcin TS_BDF *bdf; 495211a84d6SLisandro Dalcin PetscErrorCode ierr; 496211a84d6SLisandro Dalcin 497211a84d6SLisandro Dalcin PetscFunctionBegin; 498211a84d6SLisandro Dalcin ts->ops->reset = TSReset_BDF; 499211a84d6SLisandro Dalcin ts->ops->destroy = TSDestroy_BDF; 500211a84d6SLisandro Dalcin ts->ops->view = TSView_BDF; 501211a84d6SLisandro Dalcin ts->ops->setup = TSSetUp_BDF; 502211a84d6SLisandro Dalcin ts->ops->setfromoptions = TSSetFromOptions_BDF; 503211a84d6SLisandro Dalcin ts->ops->step = TSStep_BDF; 504211a84d6SLisandro Dalcin ts->ops->evaluatewlte = TSEvaluateWLTE_BDF; 505211a84d6SLisandro Dalcin ts->ops->rollback = TSRollBack_BDF; 506211a84d6SLisandro Dalcin ts->ops->interpolate = TSInterpolate_BDF; 507211a84d6SLisandro Dalcin ts->ops->snesfunction = SNESTSFormFunction_BDF; 508211a84d6SLisandro Dalcin ts->ops->snesjacobian = SNESTSFormJacobian_BDF; 5092ffb9264SLisandro Dalcin ts->default_adapt_type = TSADAPTBASIC; 510211a84d6SLisandro Dalcin 511efd4aadfSBarry Smith ts->usessnes = PETSC_TRUE; 512efd4aadfSBarry Smith 513211a84d6SLisandro Dalcin ierr = PetscNewLog(ts,&bdf);CHKERRQ(ierr); 514211a84d6SLisandro Dalcin ts->data = (void*)bdf; 515211a84d6SLisandro Dalcin 516211a84d6SLisandro Dalcin bdf->order = 2; 517211a84d6SLisandro Dalcin bdf->status = TS_STEP_COMPLETE; 518211a84d6SLisandro Dalcin 519211a84d6SLisandro Dalcin ierr = PetscObjectComposeFunction((PetscObject)ts,"TSBDFSetOrder_C",TSBDFSetOrder_BDF);CHKERRQ(ierr); 520211a84d6SLisandro Dalcin ierr = PetscObjectComposeFunction((PetscObject)ts,"TSBDFGetOrder_C",TSBDFGetOrder_BDF);CHKERRQ(ierr); 521211a84d6SLisandro Dalcin PetscFunctionReturn(0); 522211a84d6SLisandro Dalcin } 523211a84d6SLisandro Dalcin 524211a84d6SLisandro Dalcin /* ------------------------------------------------------------ */ 525211a84d6SLisandro Dalcin 526211a84d6SLisandro Dalcin /*@ 527211a84d6SLisandro Dalcin TSBDFSetOrder - Set the order of the BDF method 528211a84d6SLisandro Dalcin 529211a84d6SLisandro Dalcin Logically Collective on TS 530211a84d6SLisandro Dalcin 531211a84d6SLisandro Dalcin Input Parameter: 532211a84d6SLisandro Dalcin + ts - timestepping context 533211a84d6SLisandro Dalcin - order - order of the method 534211a84d6SLisandro Dalcin 535211a84d6SLisandro Dalcin Options Database: 536211a84d6SLisandro Dalcin . -ts_bdf_order <order> 537211a84d6SLisandro Dalcin 538211a84d6SLisandro Dalcin Level: intermediate 539211a84d6SLisandro Dalcin 540211a84d6SLisandro Dalcin @*/ 541211a84d6SLisandro Dalcin PetscErrorCode TSBDFSetOrder(TS ts,PetscInt order) 542211a84d6SLisandro Dalcin { 543211a84d6SLisandro Dalcin PetscErrorCode ierr; 544211a84d6SLisandro Dalcin 545211a84d6SLisandro Dalcin PetscFunctionBegin; 546211a84d6SLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,1); 547211a84d6SLisandro Dalcin PetscValidLogicalCollectiveInt(ts,order,2); 548211a84d6SLisandro Dalcin ierr = PetscTryMethod(ts,"TSBDFSetOrder_C",(TS,PetscInt),(ts,order));CHKERRQ(ierr); 549211a84d6SLisandro Dalcin PetscFunctionReturn(0); 550211a84d6SLisandro Dalcin } 551211a84d6SLisandro Dalcin 552211a84d6SLisandro Dalcin /*@ 553211a84d6SLisandro Dalcin TSBDFGetOrder - Get the order of the BDF method 554211a84d6SLisandro Dalcin 555211a84d6SLisandro Dalcin Not Collective 556211a84d6SLisandro Dalcin 557211a84d6SLisandro Dalcin Input Parameter: 558211a84d6SLisandro Dalcin . ts - timestepping context 559211a84d6SLisandro Dalcin 560211a84d6SLisandro Dalcin Output Parameter: 561211a84d6SLisandro Dalcin . order - order of the method 562211a84d6SLisandro Dalcin 563211a84d6SLisandro Dalcin Level: intermediate 564211a84d6SLisandro Dalcin 565211a84d6SLisandro Dalcin @*/ 566211a84d6SLisandro Dalcin PetscErrorCode TSBDFGetOrder(TS ts,PetscInt *order) 567211a84d6SLisandro Dalcin { 568211a84d6SLisandro Dalcin PetscErrorCode ierr; 569211a84d6SLisandro Dalcin 570211a84d6SLisandro Dalcin PetscFunctionBegin; 571211a84d6SLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,1); 572211a84d6SLisandro Dalcin PetscValidIntPointer(order,2); 573211a84d6SLisandro Dalcin ierr = PetscUseMethod(ts,"TSBDFGetOrder_C",(TS,PetscInt*),(ts,order));CHKERRQ(ierr); 574211a84d6SLisandro Dalcin PetscFunctionReturn(0); 575211a84d6SLisandro Dalcin } 576