xref: /petsc/src/ts/impls/explicit/rk/mrk.c (revision 0fe4d17ebc7eff27e4b161545be6c5b3fed71f51)
1474dd773SHong Zhang /*
2474dd773SHong Zhang   Code for time stepping with the multi-rate Runge-Kutta method
3474dd773SHong Zhang 
4474dd773SHong Zhang   Notes:
5474dd773SHong Zhang   1) The general system is written as
6474dd773SHong Zhang      Udot = F(t,U) for the nonsplit version of multi-rate RK,
7474dd773SHong Zhang      user should give the indexes for both slow and fast components;
8474dd773SHong Zhang   2) The general system is written as
9474dd773SHong Zhang      Usdot = Fs(t,Us,Uf)
10474dd773SHong Zhang      Ufdot = Ff(t,Us,Uf) for multi-rate RK with RHS splits,
11474dd773SHong Zhang      user should partioned RHS by themselves and also provide the indexes for both slow and fast components.
12474dd773SHong Zhang */
13474dd773SHong Zhang 
14474dd773SHong Zhang #include <petsc/private/tsimpl.h>
15474dd773SHong Zhang #include <petscdm.h>
16474dd773SHong Zhang #include <../src/ts/impls/explicit/rk/rk.h>
17474dd773SHong Zhang 
18e5bffa22SHong Zhang static PetscErrorCode TSReset_MRKNONSPLIT(TS ts)
19e5bffa22SHong Zhang {
20e5bffa22SHong Zhang   TS_RK          *rk = (TS_RK*)ts->data;
21e5bffa22SHong Zhang   RKTableau      tab = rk->tableau;
22e5bffa22SHong Zhang   PetscErrorCode ierr;
23e5bffa22SHong Zhang 
24e5bffa22SHong Zhang   PetscFunctionBegin;
25e5bffa22SHong Zhang   ierr = VecDestroy(&rk->X0);CHKERRQ(ierr);
26e5bffa22SHong Zhang   ierr = VecDestroyVecs(tab->s,&rk->YdotRHS_slow);CHKERRQ(ierr);
27e5bffa22SHong Zhang   PetscFunctionReturn(0);
28e5bffa22SHong Zhang }
29e5bffa22SHong Zhang 
30e5bffa22SHong Zhang static PetscErrorCode TSInterpolate_MRKNONSPLIT(TS ts,PetscReal itime,Vec X)
31e5bffa22SHong Zhang {
32e5bffa22SHong Zhang   TS_RK            *rk = (TS_RK*)ts->data;
33e5bffa22SHong Zhang   PetscInt         s = rk->tableau->s,p = rk->tableau->p,i,j;
3463a6f1b4SHong Zhang   PetscReal        h = ts->time_step;
35e5bffa22SHong Zhang   PetscReal        tt,t;
36e5bffa22SHong Zhang   PetscScalar      *b;
37e5bffa22SHong Zhang   const PetscReal  *B = rk->tableau->binterp;
38e5bffa22SHong Zhang   PetscErrorCode   ierr;
39e5bffa22SHong Zhang 
40e5bffa22SHong Zhang   PetscFunctionBegin;
41e5bffa22SHong Zhang   if (!B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSRK %s does not have an interpolation formula",rk->tableau->name);
4263a6f1b4SHong Zhang   t = (itime - rk->ptime)/h;
43e5bffa22SHong Zhang   ierr = PetscMalloc1(s,&b);CHKERRQ(ierr);
44e5bffa22SHong Zhang   for (i=0; i<s; i++) b[i] = 0;
45e5bffa22SHong Zhang   for (j=0,tt=t; j<p; j++,tt*=t) {
46e5bffa22SHong Zhang     for (i=0; i<s; i++) {
47e5bffa22SHong Zhang       b[i]  += h * B[i*p+j] * tt;
48e5bffa22SHong Zhang     }
49e5bffa22SHong Zhang   }
50e5bffa22SHong Zhang   ierr = VecCopy(rk->X0,X);CHKERRQ(ierr);
51e5bffa22SHong Zhang   ierr = VecMAXPY(X,s,b,rk->YdotRHS_slow);CHKERRQ(ierr);
52e5bffa22SHong Zhang   ierr = PetscFree(b);CHKERRQ(ierr);
53e5bffa22SHong Zhang   PetscFunctionReturn(0);
54e5bffa22SHong Zhang }
55e5bffa22SHong Zhang 
5663a6f1b4SHong Zhang static PetscErrorCode TSStepRefine_MRKNONSPLIT(TS ts)
5763a6f1b4SHong Zhang {
5863a6f1b4SHong Zhang   TS              previousts,subts;
5963a6f1b4SHong Zhang   TS_RK           *rk = (TS_RK*)ts->data;
6063a6f1b4SHong Zhang   RKTableau       tab = rk->tableau;
6163a6f1b4SHong Zhang   Vec             *Y = rk->Y,*YdotRHS = rk->YdotRHS;
6263a6f1b4SHong Zhang   Vec             vec_fast,subvec_fast;
6363a6f1b4SHong Zhang   const PetscInt  s = tab->s;
6463a6f1b4SHong Zhang   const PetscReal *A = tab->A,*c = tab->c;
6563a6f1b4SHong Zhang   PetscScalar     *w = rk->work;
6663a6f1b4SHong Zhang   PetscInt        i,j,k;
6763a6f1b4SHong Zhang   PetscReal       t = ts->ptime,h = ts->time_step;
6863a6f1b4SHong Zhang   PetscErrorCode  ierr;
6963a6f1b4SHong Zhang 
7063a6f1b4SHong Zhang   ierr = VecDuplicate(ts->vec_sol,&vec_fast);CHKERRQ(ierr);
7163a6f1b4SHong Zhang   previousts = rk->subts_current;
7263a6f1b4SHong Zhang   ierr = TSRHSSplitGetSubTS(rk->subts_current,"fast",&subts);CHKERRQ(ierr);
7363a6f1b4SHong Zhang   ierr = TSRHSSplitGetSubTS(subts,"fast",&subts);CHKERRQ(ierr);
7463a6f1b4SHong Zhang   for (k=0; k<rk->dtratio; k++) {
7563a6f1b4SHong Zhang     for (i=0; i<s; i++) {
7663a6f1b4SHong Zhang       ierr = TSInterpolate_MRKNONSPLIT(ts,t+k*h/rk->dtratio+h/rk->dtratio*c[i],Y[i]);CHKERRQ(ierr);
7763a6f1b4SHong Zhang       for (j=0; j<i; j++) w[j] = h/rk->dtratio*A[i*s+j];
7863a6f1b4SHong Zhang       /* update the fast components in the stage value, the slow components will be ignored, so we do not care the slow part in vec_fast */
7963a6f1b4SHong Zhang       ierr = VecCopy(ts->vec_sol,vec_fast);CHKERRQ(ierr);
8063a6f1b4SHong Zhang       ierr = VecMAXPY(vec_fast,i,w,YdotRHS);CHKERRQ(ierr);
8163a6f1b4SHong Zhang       /* update the fast components in the stage value */
8263a6f1b4SHong Zhang       ierr = VecGetSubVector(vec_fast,rk->is_fast,&subvec_fast);CHKERRQ(ierr);
8363a6f1b4SHong Zhang       ierr = VecISCopy(Y[i],rk->is_fast,SCATTER_FORWARD,subvec_fast);CHKERRQ(ierr);
8463a6f1b4SHong Zhang       ierr = VecRestoreSubVector(vec_fast,rk->is_fast,&subvec_fast);CHKERRQ(ierr);
8563a6f1b4SHong Zhang       /* compute the stage RHS */
8663a6f1b4SHong Zhang       ierr = TSComputeRHSFunction(ts,t+k*h/rk->dtratio+h/rk->dtratio*c[i],Y[i],YdotRHS[i]);CHKERRQ(ierr);
8763a6f1b4SHong Zhang     }
8863a6f1b4SHong Zhang     ierr = VecCopy(ts->vec_sol,vec_fast);CHKERRQ(ierr);
8963a6f1b4SHong Zhang     ierr = TSEvaluateStep(ts,tab->order,vec_fast,NULL);CHKERRQ(ierr);
9063a6f1b4SHong Zhang     /* update the fast components in the solution */
9163a6f1b4SHong Zhang     ierr = VecGetSubVector(vec_fast,rk->is_fast,&subvec_fast);CHKERRQ(ierr);
9263a6f1b4SHong Zhang     ierr = VecISCopy(ts->vec_sol,rk->is_fast,SCATTER_FORWARD,subvec_fast);CHKERRQ(ierr);
9363a6f1b4SHong Zhang     ierr = VecRestoreSubVector(vec_fast,rk->is_fast,&subvec_fast);CHKERRQ(ierr);
9463a6f1b4SHong Zhang 
9563a6f1b4SHong Zhang     if (subts) {
9663a6f1b4SHong Zhang       Vec *YdotRHS_copy;
9763a6f1b4SHong Zhang       ierr = VecDuplicateVecs(ts->vec_sol,s,&YdotRHS_copy);CHKERRQ(ierr);
9863a6f1b4SHong Zhang       rk->subts_current = rk->subts_fast;
9963a6f1b4SHong Zhang       ts->ptime = t+k*h/rk->dtratio;
10063a6f1b4SHong Zhang       ts->time_step = h/rk->dtratio;
10163a6f1b4SHong Zhang       ierr = TSRHSSplitGetIS(rk->subts_current,"fast",&rk->is_fast);CHKERRQ(ierr);
10263a6f1b4SHong Zhang       for (i=0; i<s; i++) {
10363a6f1b4SHong Zhang         ierr = VecCopy(rk->YdotRHS_slow[i],YdotRHS_copy[i]);CHKERRQ(ierr);
10463a6f1b4SHong Zhang         ierr = VecCopy(YdotRHS[i],rk->YdotRHS_slow[i]);CHKERRQ(ierr);
10563a6f1b4SHong Zhang       }
10663a6f1b4SHong Zhang 
10763a6f1b4SHong Zhang       ierr = TSStepRefine_MRKNONSPLIT(ts);CHKERRQ(ierr);
10863a6f1b4SHong Zhang 
10963a6f1b4SHong Zhang       rk->subts_current = previousts;
11063a6f1b4SHong Zhang       ts->ptime = t;
11163a6f1b4SHong Zhang       ts->time_step = h;
11263a6f1b4SHong Zhang       ierr = TSRHSSplitGetIS(previousts,"fast",&rk->is_fast);CHKERRQ(ierr);
11363a6f1b4SHong Zhang       for (i=0; i<s; i++) {
11463a6f1b4SHong Zhang         ierr = VecCopy(YdotRHS_copy[i],rk->YdotRHS_slow[i]);CHKERRQ(ierr);
11563a6f1b4SHong Zhang       }
11663a6f1b4SHong Zhang       ierr = VecDestroyVecs(s,&YdotRHS_copy);CHKERRQ(ierr);
11763a6f1b4SHong Zhang     }
11863a6f1b4SHong Zhang   }
11963a6f1b4SHong Zhang   ierr = VecDestroy(&vec_fast);CHKERRQ(ierr);
12063a6f1b4SHong Zhang   PetscFunctionReturn(0);
12163a6f1b4SHong Zhang }
12263a6f1b4SHong Zhang 
123e5bffa22SHong Zhang static PetscErrorCode TSStep_MRKNONSPLIT(TS ts)
124e5bffa22SHong Zhang {
125e5bffa22SHong Zhang   TS_RK           *rk = (TS_RK*)ts->data;
126e5bffa22SHong Zhang   RKTableau       tab = rk->tableau;
127e5bffa22SHong Zhang   Vec             *Y = rk->Y,*YdotRHS = rk->YdotRHS,*YdotRHS_slow = rk->YdotRHS_slow;
128e5bffa22SHong Zhang   Vec             stage_slow,sol_slow; /* vectors store the slow components */
129e5bffa22SHong Zhang   Vec             subvec_slow; /* sub vector to store the slow components */
130e5bffa22SHong Zhang   IS              is_slow = rk->is_slow;
131e5bffa22SHong Zhang   const PetscInt  s = tab->s;
132e5bffa22SHong Zhang   const PetscReal *A = tab->A,*c = tab->c;
133e5bffa22SHong Zhang   PetscScalar     *w = rk->work;
13463a6f1b4SHong Zhang   PetscInt        i,j,dtratio = rk->dtratio;
135e5bffa22SHong Zhang   PetscReal       next_time_step = ts->time_step,t = ts->ptime,h = ts->time_step;
136e5bffa22SHong Zhang   PetscErrorCode  ierr;
137e5bffa22SHong Zhang 
138e5bffa22SHong Zhang   PetscFunctionBegin;
139e5bffa22SHong Zhang   rk->status = TS_STEP_INCOMPLETE;
140e5bffa22SHong Zhang   ierr = VecDuplicate(ts->vec_sol,&stage_slow);CHKERRQ(ierr);
141e5bffa22SHong Zhang   ierr = VecDuplicate(ts->vec_sol,&sol_slow);CHKERRQ(ierr);
142e5bffa22SHong Zhang   ierr = VecCopy(ts->vec_sol,rk->X0);CHKERRQ(ierr);
143e5bffa22SHong Zhang   for (i=0; i<s; i++) {
144e5bffa22SHong Zhang     rk->stage_time = t + h*c[i];
145e5bffa22SHong Zhang     ierr = TSPreStage(ts,rk->stage_time);CHKERRQ(ierr);
146e5bffa22SHong Zhang     ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr);
147e5bffa22SHong Zhang     for (j=0; j<i; j++) w[j] = h*A[i*s+j];
148e5bffa22SHong Zhang     ierr = VecMAXPY(Y[i],i,w,YdotRHS_slow);CHKERRQ(ierr);
149e5bffa22SHong Zhang     ierr = TSPostStage(ts,rk->stage_time,i,Y); CHKERRQ(ierr);
150e5bffa22SHong Zhang     /* compute the stage RHS */
151e5bffa22SHong Zhang     ierr = TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS_slow[i]);CHKERRQ(ierr);
152e5bffa22SHong Zhang   }
153e5bffa22SHong Zhang   /* update the slow components in the solution */
154e5bffa22SHong Zhang   rk->YdotRHS = YdotRHS_slow;
155e5bffa22SHong Zhang   rk->dtratio = 1;
156e5bffa22SHong Zhang   ierr = TSEvaluateStep(ts,tab->order,sol_slow,NULL);CHKERRQ(ierr);
157e5bffa22SHong Zhang   rk->dtratio = dtratio;
158e5bffa22SHong Zhang   rk->YdotRHS = YdotRHS;
159e5bffa22SHong Zhang   /* update the slow components in the solution */
160e5bffa22SHong Zhang   ierr = VecGetSubVector(sol_slow,is_slow,&subvec_slow);CHKERRQ(ierr);
161e5bffa22SHong Zhang   ierr = VecISCopy(ts->vec_sol,is_slow,SCATTER_FORWARD,subvec_slow);CHKERRQ(ierr);
162e5bffa22SHong Zhang   ierr = VecRestoreSubVector(sol_slow,is_slow,&subvec_slow);CHKERRQ(ierr);
163e5bffa22SHong Zhang 
16463a6f1b4SHong Zhang   rk->subts_current = ts;
16563a6f1b4SHong Zhang   rk->ptime = t;
16663a6f1b4SHong Zhang   rk->time_step = h;
16763a6f1b4SHong Zhang   ierr = TSStepRefine_MRKNONSPLIT(ts);CHKERRQ(ierr);
16863a6f1b4SHong Zhang 
169e5bffa22SHong Zhang   ts->ptime = t + ts->time_step;
170e5bffa22SHong Zhang   ts->time_step = next_time_step;
171e5bffa22SHong Zhang   rk->status = TS_STEP_COMPLETE;
17263a6f1b4SHong Zhang 
173e5bffa22SHong Zhang   /* free memory of work vectors */
174e5bffa22SHong Zhang   ierr = VecDestroy(&stage_slow);CHKERRQ(ierr);
175e5bffa22SHong Zhang   ierr = VecDestroy(&sol_slow);CHKERRQ(ierr);
176e5bffa22SHong Zhang   PetscFunctionReturn(0);
177e5bffa22SHong Zhang }
178e5bffa22SHong Zhang 
179*0fe4d17eSHong Zhang static PetscErrorCode TSSetUp_MRKNONSPLIT(TS ts)
180*0fe4d17eSHong Zhang {
181*0fe4d17eSHong Zhang   TS_RK          *rk = (TS_RK*)ts->data;
182*0fe4d17eSHong Zhang   RKTableau      tab = rk->tableau;
183*0fe4d17eSHong Zhang   PetscErrorCode ierr;
184*0fe4d17eSHong Zhang 
185*0fe4d17eSHong Zhang   PetscFunctionBegin;
186*0fe4d17eSHong Zhang   ierr = TSRHSSplitGetIS(ts,"slow",&rk->is_slow);CHKERRQ(ierr);
187*0fe4d17eSHong Zhang   ierr = TSRHSSplitGetIS(ts,"fast",&rk->is_fast);CHKERRQ(ierr);
188*0fe4d17eSHong Zhang   if (!rk->is_slow || !rk->is_fast) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Must set up RHSSplits with TSRHSSplitSetIS() using split names 'slow' and 'fast' respectively in order to use multirate RK");
189*0fe4d17eSHong Zhang   ierr = TSRHSSplitGetSubTS(ts,"slow",&rk->subts_slow);CHKERRQ(ierr);
190*0fe4d17eSHong Zhang   ierr = TSRHSSplitGetSubTS(ts,"fast",&rk->subts_fast);CHKERRQ(ierr);
191*0fe4d17eSHong Zhang   if (!rk->subts_slow || !rk->subts_fast) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Must set up the RHSFunctions for 'slow' and 'fast' components using TSRHSSplitSetRHSFunction() or calling TSSetRHSFunction() for each sub-TS");
192*0fe4d17eSHong Zhang   ierr = VecDuplicate(ts->vec_sol,&rk->X0);CHKERRQ(ierr);
193*0fe4d17eSHong Zhang   ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&rk->YdotRHS_slow);CHKERRQ(ierr);
194*0fe4d17eSHong Zhang   rk->subts_current = rk->subts_fast;
195*0fe4d17eSHong Zhang 
196*0fe4d17eSHong Zhang   ts->ops->step        = TSStep_MRKNONSPLIT;
197*0fe4d17eSHong Zhang   ts->ops->interpolate = TSInterpolate_MRKNONSPLIT;
198*0fe4d17eSHong Zhang   PetscFunctionReturn(0);
199*0fe4d17eSHong Zhang }
200*0fe4d17eSHong Zhang 
20163a6f1b4SHong Zhang /*
20263a6f1b4SHong Zhang   Copy DM from tssrc to tsdest, while keeping the original DMTS and DMSNES in tsdest.
20363a6f1b4SHong Zhang */
20463a6f1b4SHong Zhang static PetscErrorCode TSCopyDM(TS tssrc,TS tsdest)
20563a6f1b4SHong Zhang {
20663a6f1b4SHong Zhang   DM             newdm,dmsrc,dmdest;
20763a6f1b4SHong Zhang   PetscErrorCode ierr;
20863a6f1b4SHong Zhang 
20963a6f1b4SHong Zhang   PetscFunctionBegin;
21063a6f1b4SHong Zhang   ierr = TSGetDM(tssrc,&dmsrc);CHKERRQ(ierr);
21163a6f1b4SHong Zhang   ierr = DMClone(dmsrc,&newdm);CHKERRQ(ierr);
21263a6f1b4SHong Zhang   ierr = TSGetDM(tsdest,&dmdest);CHKERRQ(ierr);
21363a6f1b4SHong Zhang   ierr = DMCopyDMTS(dmdest,newdm);CHKERRQ(ierr);
21463a6f1b4SHong Zhang   ierr = DMCopyDMSNES(dmdest,newdm);CHKERRQ(ierr);
21563a6f1b4SHong Zhang   ierr = TSSetDM(tsdest,newdm);CHKERRQ(ierr);
21663a6f1b4SHong Zhang   ierr = DMDestroy(&newdm);CHKERRQ(ierr);
21763a6f1b4SHong Zhang   PetscFunctionReturn(0);
21863a6f1b4SHong Zhang }
21963a6f1b4SHong Zhang 
220474dd773SHong Zhang static PetscErrorCode TSReset_MRKSPLIT(TS ts)
221474dd773SHong Zhang {
222474dd773SHong Zhang   TS_RK          *rk = (TS_RK*)ts->data;
223474dd773SHong Zhang   PetscErrorCode ierr;
224474dd773SHong Zhang 
225474dd773SHong Zhang   PetscFunctionBegin;
226bb42530cSHong Zhang   if (rk->subts_slow) {
227bb42530cSHong Zhang     ierr = PetscFree(rk->subts_slow->data);CHKERRQ(ierr);
228bb42530cSHong Zhang     rk->subts_slow = NULL;
229bb42530cSHong Zhang   }
230bb42530cSHong Zhang   if (rk->subts_fast) {
23163a6f1b4SHong Zhang     ierr = PetscFree(rk->YdotRHS_fast);CHKERRQ(ierr);
23263a6f1b4SHong Zhang     ierr = PetscFree(rk->YdotRHS_slow);CHKERRQ(ierr);
233bb42530cSHong Zhang     ierr = VecDestroy(&rk->X0);CHKERRQ(ierr);
234bb42530cSHong Zhang     ierr = TSReset_MRKSPLIT(rk->subts_fast);
235bb42530cSHong Zhang     ierr = PetscFree(rk->subts_fast->data);CHKERRQ(ierr);
236bb42530cSHong Zhang     rk->subts_fast = NULL;
237bb42530cSHong Zhang   }
238474dd773SHong Zhang   PetscFunctionReturn(0);
239474dd773SHong Zhang }
240474dd773SHong Zhang 
241474dd773SHong Zhang static PetscErrorCode TSInterpolate_MRKSPLIT(TS ts,PetscReal itime,Vec X)
242474dd773SHong Zhang {
243474dd773SHong Zhang   TS_RK           *rk = (TS_RK*)ts->data;
24463a6f1b4SHong Zhang   Vec             Xslow;
245474dd773SHong Zhang   PetscInt        s = rk->tableau->s,p = rk->tableau->p,i,j;
246474dd773SHong Zhang   PetscReal       h;
247474dd773SHong Zhang   PetscReal       tt,t;
248474dd773SHong Zhang   PetscScalar     *b;
249474dd773SHong Zhang   const PetscReal *B = rk->tableau->binterp;
250474dd773SHong Zhang   PetscErrorCode  ierr;
251474dd773SHong Zhang 
252474dd773SHong Zhang   PetscFunctionBegin;
253474dd773SHong Zhang   if (!B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSRK %s does not have an interpolation formula",rk->tableau->name);
254474dd773SHong Zhang 
255474dd773SHong Zhang   switch (rk->status) {
256474dd773SHong Zhang     case TS_STEP_INCOMPLETE:
257474dd773SHong Zhang     case TS_STEP_PENDING:
258474dd773SHong Zhang       h = ts->time_step;
259474dd773SHong Zhang       t = (itime - ts->ptime)/h;
260474dd773SHong Zhang       break;
261474dd773SHong Zhang     case TS_STEP_COMPLETE:
262474dd773SHong Zhang       h = ts->ptime - ts->ptime_prev;
263474dd773SHong Zhang       t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */
264474dd773SHong Zhang       break;
265474dd773SHong Zhang     default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
266474dd773SHong Zhang   }
267474dd773SHong Zhang   ierr = PetscMalloc1(s,&b);CHKERRQ(ierr);
268474dd773SHong Zhang   for (i=0; i<s; i++) b[i] = 0;
269474dd773SHong Zhang   for (j=0,tt=t; j<p; j++,tt*=t) {
270474dd773SHong Zhang     for (i=0; i<s; i++) {
271474dd773SHong Zhang       b[i]  += h * B[i*p+j] * tt;
272474dd773SHong Zhang     }
273474dd773SHong Zhang   }
27463a6f1b4SHong Zhang   for (i=0; i<s; i++) {
27563a6f1b4SHong Zhang     ierr = VecGetSubVector(rk->YdotRHS[i],rk->is_slow,&rk->YdotRHS_slow[i]);CHKERRQ(ierr);
27663a6f1b4SHong Zhang   }
27763a6f1b4SHong Zhang   ierr = VecGetSubVector(X,rk->is_slow,&Xslow);CHKERRQ(ierr);
27863a6f1b4SHong Zhang   ierr = VecISCopy(rk->X0,rk->is_slow,SCATTER_REVERSE,Xslow);CHKERRQ(ierr);
27963a6f1b4SHong Zhang   ierr = VecMAXPY(Xslow,s,b,rk->YdotRHS_slow);CHKERRQ(ierr);
28063a6f1b4SHong Zhang   ierr = VecRestoreSubVector(X,rk->is_slow,&Xslow);CHKERRQ(ierr);
28163a6f1b4SHong Zhang   for (i=0; i<s; i++) {
28263a6f1b4SHong Zhang     ierr = VecRestoreSubVector(rk->YdotRHS[i],rk->is_slow,&rk->YdotRHS_slow[i]);CHKERRQ(ierr);
28363a6f1b4SHong Zhang   }
284474dd773SHong Zhang   ierr = PetscFree(b);CHKERRQ(ierr);
285474dd773SHong Zhang   PetscFunctionReturn(0);
286474dd773SHong Zhang }
287474dd773SHong Zhang 
288474dd773SHong Zhang /*
289474dd773SHong Zhang  This is for partitioned RHS multirate RK method
290474dd773SHong Zhang  The step completion formula is
291474dd773SHong Zhang 
292474dd773SHong Zhang  x1 = x0 + h b^T YdotRHS
293474dd773SHong Zhang 
294474dd773SHong Zhang */
295474dd773SHong Zhang static PetscErrorCode TSEvaluateStep_MRKSPLIT(TS ts,PetscInt order,Vec X,PetscBool *done)
296474dd773SHong Zhang {
297474dd773SHong Zhang   TS_RK          *rk = (TS_RK*)ts->data;
298474dd773SHong Zhang   RKTableau      tab = rk->tableau;
299474dd773SHong Zhang   Vec            Xslow,Xfast;                  /* subvectors of X which store slow components and fast components respectively */
300474dd773SHong Zhang   PetscScalar    *w = rk->work;
301474dd773SHong Zhang   PetscReal      h = ts->time_step;
302474dd773SHong Zhang   PetscInt       s = tab->s,j;
303474dd773SHong Zhang   PetscErrorCode ierr;
304474dd773SHong Zhang 
305474dd773SHong Zhang   PetscFunctionBegin;
306474dd773SHong Zhang   ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
307474dd773SHong Zhang   if (rk->slow) {
308474dd773SHong Zhang     for (j=0; j<s; j++) w[j] = h*tab->b[j];
309474dd773SHong Zhang     ierr = VecGetSubVector(ts->vec_sol,rk->is_slow,&Xslow);CHKERRQ(ierr);
31063a6f1b4SHong Zhang     ierr = VecMAXPY(Xslow,s,w,rk->YdotRHS_slow);CHKERRQ(ierr);
311474dd773SHong Zhang     ierr = VecRestoreSubVector(ts->vec_sol,rk->is_slow,&Xslow);CHKERRQ(ierr);;
312474dd773SHong Zhang   } else {
313474dd773SHong Zhang     for (j=0; j<s; j++) w[j] = h/rk->dtratio*tab->b[j];
314474dd773SHong Zhang     ierr = VecGetSubVector(X,rk->is_fast,&Xfast);CHKERRQ(ierr);
31563a6f1b4SHong Zhang     ierr = VecMAXPY(Xfast,s,w,rk->YdotRHS_fast);CHKERRQ(ierr);
316474dd773SHong Zhang     ierr = VecRestoreSubVector(X,rk->is_fast,&Xfast);CHKERRQ(ierr);
317474dd773SHong Zhang   }
318474dd773SHong Zhang   PetscFunctionReturn(0);
319474dd773SHong Zhang }
320474dd773SHong Zhang 
32163a6f1b4SHong Zhang static PetscErrorCode TSStepRefine_MRKSPLIT(TS ts)
32263a6f1b4SHong Zhang {
32363a6f1b4SHong Zhang   TS_RK           *rk = (TS_RK*)ts->data;
32463a6f1b4SHong Zhang   TS              subts_fast = rk->subts_fast,currentlevelts;
32563a6f1b4SHong Zhang   TS_RK           *subrk_fast = (TS_RK*)subts_fast->data;
32663a6f1b4SHong Zhang   RKTableau       tab = rk->tableau;
32763a6f1b4SHong Zhang   Vec             *Y = rk->Y;
32863a6f1b4SHong Zhang   Vec             *YdotRHS = rk->YdotRHS,*YdotRHS_fast = rk->YdotRHS_fast;
32963a6f1b4SHong Zhang   Vec             Yfast,Xfast;
33063a6f1b4SHong Zhang   const PetscInt  s = tab->s;
33163a6f1b4SHong Zhang   const PetscReal *A = tab->A,*c = tab->c;
33263a6f1b4SHong Zhang   PetscScalar     *w = rk->work;
33363a6f1b4SHong Zhang   PetscInt        i,j,k;
33463a6f1b4SHong Zhang   PetscReal       t = ts->ptime,h = ts->time_step;
33563a6f1b4SHong Zhang   PetscErrorCode  ierr;
33663a6f1b4SHong Zhang 
33763a6f1b4SHong Zhang   for (k=0; k<rk->dtratio; k++) {
33863a6f1b4SHong Zhang     ierr = VecGetSubVector(ts->vec_sol,rk->is_fast,&Xfast);CHKERRQ(ierr);
33963a6f1b4SHong Zhang     for (i=0; i<s; i++) {
34063a6f1b4SHong Zhang       ierr = VecGetSubVector(YdotRHS[i],rk->is_fast,&YdotRHS_fast[i]);CHKERRQ(ierr);
34163a6f1b4SHong Zhang     }
34263a6f1b4SHong Zhang     /* propogate fast component using small time steps */
34363a6f1b4SHong Zhang     for (i=0; i<s; i++) {
34463a6f1b4SHong Zhang       /* stage value for slow components */
34563a6f1b4SHong Zhang       ierr = TSInterpolate_MRKSPLIT(rk->ts_root,t+k*h/rk->dtratio+h/rk->dtratio*c[i],Y[i]);CHKERRQ(ierr);
34663a6f1b4SHong Zhang       currentlevelts = rk->ts_root;
34763a6f1b4SHong Zhang       while (currentlevelts != ts) { /* all the slow parts need to be interpolated separated */
34863a6f1b4SHong Zhang         currentlevelts = ((TS_RK*)currentlevelts->data)->subts_fast;
34963a6f1b4SHong Zhang         ierr = TSInterpolate_MRKSPLIT(currentlevelts,t+k*h/rk->dtratio+h/rk->dtratio*c[i],Y[i]);CHKERRQ(ierr);
35063a6f1b4SHong Zhang       }
35163a6f1b4SHong Zhang       for (j=0; j<i; j++) w[j] = h/rk->dtratio*A[i*s+j];
35263a6f1b4SHong Zhang       subrk_fast->stage_time = t + h/rk->dtratio*c[i];
35363a6f1b4SHong Zhang       ierr = TSPreStage(subts_fast,subrk_fast->stage_time);CHKERRQ(ierr);
35463a6f1b4SHong Zhang       /* stage value for fast components */
35563a6f1b4SHong Zhang       ierr = VecGetSubVector(Y[i],rk->is_fast,&Yfast);CHKERRQ(ierr);
35663a6f1b4SHong Zhang       ierr = VecCopy(Xfast,Yfast);CHKERRQ(ierr);
35763a6f1b4SHong Zhang       ierr = VecMAXPY(Yfast,i,w,YdotRHS_fast);CHKERRQ(ierr);
35863a6f1b4SHong Zhang       ierr = VecRestoreSubVector(Y[i],rk->is_fast,&Yfast);CHKERRQ(ierr);
35963a6f1b4SHong Zhang       ierr = TSPostStage(subts_fast,subrk_fast->stage_time,i,Y);CHKERRQ(ierr);
36063a6f1b4SHong Zhang       /* compute the stage RHS for fast components */
36163a6f1b4SHong Zhang       ierr = TSComputeRHSFunction(subts_fast,t+k*h*rk->dtratio+h/rk->dtratio*c[i],Y[i],YdotRHS_fast[i]);CHKERRQ(ierr);
36263a6f1b4SHong Zhang     }
36363a6f1b4SHong Zhang     ierr = VecRestoreSubVector(ts->vec_sol,rk->is_fast,&Xfast);CHKERRQ(ierr);
36463a6f1b4SHong Zhang     /* update the value of fast components using fast time step */
36563a6f1b4SHong Zhang     rk->slow = PETSC_FALSE;
36663a6f1b4SHong Zhang     ierr = TSEvaluateStep_MRKSPLIT(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr);
36763a6f1b4SHong Zhang     for (i=0; i<s; i++) {
36863a6f1b4SHong Zhang       ierr = VecRestoreSubVector(YdotRHS[i],rk->is_fast,&YdotRHS_fast[i]);CHKERRQ(ierr);
36963a6f1b4SHong Zhang     }
37063a6f1b4SHong Zhang 
37163a6f1b4SHong Zhang     if (subrk_fast->subts_fast) {
37263a6f1b4SHong Zhang       subts_fast->ptime = t+k*h/rk->dtratio;
37363a6f1b4SHong Zhang       subts_fast->time_step = h/rk->dtratio;
37463a6f1b4SHong Zhang       ierr = TSStepRefine_MRKSPLIT(subts_fast);CHKERRQ(ierr);
37563a6f1b4SHong Zhang     }
37663a6f1b4SHong Zhang     /* update the fast components of the solution */
37763a6f1b4SHong Zhang     ierr = VecGetSubVector(ts->vec_sol,rk->is_fast,&Xfast);CHKERRQ(ierr);
37863a6f1b4SHong Zhang     ierr = VecISCopy(rk->X0,rk->is_fast,SCATTER_FORWARD,Xfast);CHKERRQ(ierr);
37963a6f1b4SHong Zhang     ierr = VecRestoreSubVector(ts->vec_sol,rk->is_fast,&Xfast);CHKERRQ(ierr);
38063a6f1b4SHong Zhang   }
38163a6f1b4SHong Zhang   PetscFunctionReturn(0);
38263a6f1b4SHong Zhang }
38363a6f1b4SHong Zhang 
384474dd773SHong Zhang static PetscErrorCode TSStep_MRKSPLIT(TS ts)
385474dd773SHong Zhang {
386474dd773SHong Zhang   TS_RK           *rk = (TS_RK*)ts->data;
387474dd773SHong Zhang   RKTableau       tab = rk->tableau;
388474dd773SHong Zhang   Vec             *Y = rk->Y,*YdotRHS = rk->YdotRHS;
38963a6f1b4SHong Zhang   Vec             *YdotRHS_fast = rk->YdotRHS_fast,*YdotRHS_slow = rk->YdotRHS_slow;
390474dd773SHong Zhang   Vec             Yslow,Yfast; /* subvectors store the stges of slow components and fast components respectively                           */
391474dd773SHong Zhang   const PetscInt  s = tab->s;
392474dd773SHong Zhang   const PetscReal *A = tab->A,*c = tab->c;
393474dd773SHong Zhang   PetscScalar     *w = rk->work;
39463a6f1b4SHong Zhang   PetscInt        i,j;
395474dd773SHong Zhang   PetscReal       next_time_step = ts->time_step,t = ts->ptime,h = ts->time_step;
396474dd773SHong Zhang   PetscErrorCode  ierr;
397474dd773SHong Zhang 
398474dd773SHong Zhang   PetscFunctionBegin;
399474dd773SHong Zhang   rk->status = TS_STEP_INCOMPLETE;
400474dd773SHong Zhang   for (i=0; i<s; i++) {
401474dd773SHong Zhang     ierr = VecGetSubVector(YdotRHS[i],rk->is_slow,&YdotRHS_slow[i]);CHKERRQ(ierr);
402474dd773SHong Zhang     ierr = VecGetSubVector(YdotRHS[i],rk->is_fast,&YdotRHS_fast[i]);CHKERRQ(ierr);
403474dd773SHong Zhang   }
404474dd773SHong Zhang   ierr = VecCopy(ts->vec_sol,rk->X0);CHKERRQ(ierr);
405474dd773SHong Zhang   /* propogate both slow and fast components using large time steps */
406474dd773SHong Zhang   for (i=0; i<s; i++) {
407474dd773SHong Zhang     rk->stage_time = t + h*c[i];
408474dd773SHong Zhang     ierr = TSPreStage(ts,rk->stage_time);CHKERRQ(ierr);
409474dd773SHong Zhang     ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr);
410474dd773SHong Zhang     ierr = VecGetSubVector(Y[i],rk->is_fast,&Yfast);CHKERRQ(ierr);
411bb42530cSHong Zhang     ierr = VecGetSubVector(Y[i],rk->is_slow,&Yslow);CHKERRQ(ierr);
412474dd773SHong Zhang     for (j=0; j<i; j++) w[j] = h*A[i*s+j];
413474dd773SHong Zhang     ierr = VecMAXPY(Yfast,i,w,YdotRHS_fast);CHKERRQ(ierr);
414bb42530cSHong Zhang     ierr = VecMAXPY(Yslow,i,w,YdotRHS_slow);CHKERRQ(ierr);
415474dd773SHong Zhang     ierr = VecRestoreSubVector(Y[i],rk->is_fast,&Yfast);CHKERRQ(ierr);
416bb42530cSHong Zhang     ierr = VecRestoreSubVector(Y[i],rk->is_slow,&Yslow);CHKERRQ(ierr);
417474dd773SHong Zhang     ierr = TSPostStage(ts,rk->stage_time,i,Y); CHKERRQ(ierr);
418474dd773SHong Zhang     ierr = TSComputeRHSFunction(rk->subts_slow,t+h*c[i],Y[i],YdotRHS_slow[i]);CHKERRQ(ierr);
419474dd773SHong Zhang     ierr = TSComputeRHSFunction(rk->subts_fast,t+h*c[i],Y[i],YdotRHS_fast[i]);CHKERRQ(ierr);
420474dd773SHong Zhang   }
421474dd773SHong Zhang   rk->slow = PETSC_TRUE;
42263a6f1b4SHong Zhang   /* update the slow components of the solution using slow time step */
423474dd773SHong Zhang   ierr = TSEvaluateStep_MRKSPLIT(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr);
42463a6f1b4SHong Zhang   /* YdotRHS will be used for interpolation during refinement */
425474dd773SHong Zhang   for (i=0; i<s; i++) {
426474dd773SHong Zhang     ierr = VecRestoreSubVector(YdotRHS[i],rk->is_slow,&YdotRHS_slow[i]);CHKERRQ(ierr);
427474dd773SHong Zhang     ierr = VecRestoreSubVector(YdotRHS[i],rk->is_fast,&YdotRHS_fast[i]);CHKERRQ(ierr);
428474dd773SHong Zhang   }
42963a6f1b4SHong Zhang 
43063a6f1b4SHong Zhang   ierr = TSStepRefine_MRKSPLIT(ts);CHKERRQ(ierr);
43163a6f1b4SHong Zhang 
432bb42530cSHong Zhang   ts->ptime = t + ts->time_step;
433474dd773SHong Zhang   ts->time_step = next_time_step;
434474dd773SHong Zhang   rk->status = TS_STEP_COMPLETE;
435474dd773SHong Zhang   PetscFunctionReturn(0);
436474dd773SHong Zhang }
437474dd773SHong Zhang 
438*0fe4d17eSHong Zhang static PetscErrorCode TSSetUp_MRKSPLIT(TS ts)
439*0fe4d17eSHong Zhang {
440*0fe4d17eSHong Zhang   TS_RK          *rk = (TS_RK*)ts->data,*nextlevelrk,*currentlevelrk;
441*0fe4d17eSHong Zhang   TS             nextlevelts;
442*0fe4d17eSHong Zhang   Vec            X0;
443*0fe4d17eSHong Zhang   PetscErrorCode ierr;
444*0fe4d17eSHong Zhang 
445*0fe4d17eSHong Zhang   PetscFunctionBegin;
446*0fe4d17eSHong Zhang   ierr = TSRHSSplitGetIS(ts,"slow",&rk->is_slow);CHKERRQ(ierr);
447*0fe4d17eSHong Zhang   ierr = TSRHSSplitGetIS(ts,"fast",&rk->is_fast);CHKERRQ(ierr);
448*0fe4d17eSHong Zhang   if (!rk->is_slow || !rk->is_fast) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Must set up RHSSplits with TSRHSSplitSetIS() using split names 'slow' and 'fast' respectively in order to use -ts_type bsi");
449*0fe4d17eSHong Zhang   ierr = TSRHSSplitGetSubTS(ts,"slow",&rk->subts_slow);CHKERRQ(ierr);
450*0fe4d17eSHong Zhang   ierr = TSRHSSplitGetSubTS(ts,"fast",&rk->subts_fast);CHKERRQ(ierr);
451*0fe4d17eSHong Zhang   if (!rk->subts_slow || !rk->subts_fast) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Must set up the RHSFunctions for 'slow' and 'fast' components using TSRHSSplitSetRHSFunction() or calling TSSetRHSFunction() for each sub-TS");
452*0fe4d17eSHong Zhang 
453*0fe4d17eSHong Zhang   ierr = VecDuplicate(ts->vec_sol,&X0);CHKERRQ(ierr);
454*0fe4d17eSHong Zhang   /* The TS at each level share the same tableau, work array, solution vector, stage values and stage derivatives */
455*0fe4d17eSHong Zhang   currentlevelrk = rk;
456*0fe4d17eSHong Zhang   while (currentlevelrk->subts_fast) {
457*0fe4d17eSHong Zhang     ierr = PetscMalloc1(rk->tableau->s,&currentlevelrk->YdotRHS_fast);CHKERRQ(ierr);
458*0fe4d17eSHong Zhang     ierr = PetscMalloc1(rk->tableau->s,&currentlevelrk->YdotRHS_slow);CHKERRQ(ierr);
459*0fe4d17eSHong Zhang     ierr = PetscObjectReference((PetscObject)X0);CHKERRQ(ierr);
460*0fe4d17eSHong Zhang     currentlevelrk->X0 = X0;
461*0fe4d17eSHong Zhang     currentlevelrk->ts_root = ts;
462*0fe4d17eSHong Zhang 
463*0fe4d17eSHong Zhang     /* set up the ts for the slow part */
464*0fe4d17eSHong Zhang     nextlevelts = currentlevelrk->subts_slow;
465*0fe4d17eSHong Zhang     ierr = PetscNewLog(nextlevelts,&nextlevelrk);CHKERRQ(ierr);
466*0fe4d17eSHong Zhang     nextlevelrk->tableau = rk->tableau;
467*0fe4d17eSHong Zhang     nextlevelrk->work = rk->work;
468*0fe4d17eSHong Zhang     nextlevelrk->Y = rk->Y;
469*0fe4d17eSHong Zhang     nextlevelrk->YdotRHS = rk->YdotRHS;
470*0fe4d17eSHong Zhang     nextlevelts->data = (void*)nextlevelrk;
471*0fe4d17eSHong Zhang     ierr = TSCopyDM(ts,nextlevelts);CHKERRQ(ierr);
472*0fe4d17eSHong Zhang     ierr = TSSetSolution(nextlevelts,ts->vec_sol);CHKERRQ(ierr);
473*0fe4d17eSHong Zhang 
474*0fe4d17eSHong Zhang     /* set up the ts for the fast part */
475*0fe4d17eSHong Zhang     nextlevelts = currentlevelrk->subts_fast;
476*0fe4d17eSHong Zhang     ierr = PetscNewLog(nextlevelts,&nextlevelrk);CHKERRQ(ierr);
477*0fe4d17eSHong Zhang     nextlevelrk->tableau = rk->tableau;
478*0fe4d17eSHong Zhang     nextlevelrk->work = rk->work;
479*0fe4d17eSHong Zhang     nextlevelrk->Y = rk->Y;
480*0fe4d17eSHong Zhang     nextlevelrk->YdotRHS = rk->YdotRHS;
481*0fe4d17eSHong Zhang     nextlevelrk->dtratio = rk->dtratio;
482*0fe4d17eSHong Zhang     ierr = TSRHSSplitGetIS(nextlevelts,"slow",&nextlevelrk->is_slow);CHKERRQ(ierr);
483*0fe4d17eSHong Zhang     ierr = TSRHSSplitGetSubTS(nextlevelts,"slow",&nextlevelrk->subts_slow);CHKERRQ(ierr);
484*0fe4d17eSHong Zhang     ierr = TSRHSSplitGetIS(nextlevelts,"fast",&nextlevelrk->is_fast);CHKERRQ(ierr);
485*0fe4d17eSHong Zhang     ierr = TSRHSSplitGetSubTS(nextlevelts,"fast",&nextlevelrk->subts_fast);CHKERRQ(ierr);
486*0fe4d17eSHong Zhang     nextlevelts->data = (void*)nextlevelrk;
487*0fe4d17eSHong Zhang     ierr = TSCopyDM(ts,nextlevelts);CHKERRQ(ierr);
488*0fe4d17eSHong Zhang     ierr = TSSetSolution(nextlevelts,ts->vec_sol);CHKERRQ(ierr);
489*0fe4d17eSHong Zhang 
490*0fe4d17eSHong Zhang     currentlevelrk = nextlevelrk;
491*0fe4d17eSHong Zhang   }
492*0fe4d17eSHong Zhang   ierr = VecDestroy(&X0);CHKERRQ(ierr);
493*0fe4d17eSHong Zhang 
494*0fe4d17eSHong Zhang   ts->ops->step         = TSStep_MRKSPLIT;
495*0fe4d17eSHong Zhang   ts->ops->evaluatestep = TSEvaluateStep_MRKSPLIT;
496*0fe4d17eSHong Zhang   ts->ops->interpolate  = TSInterpolate_MRKSPLIT;
497*0fe4d17eSHong Zhang   PetscFunctionReturn(0);
498*0fe4d17eSHong Zhang }
499*0fe4d17eSHong Zhang 
500474dd773SHong Zhang /*@C
501*0fe4d17eSHong Zhang   TSRKSetMultirate - Use the interpolation-based multirate RK method
502474dd773SHong Zhang 
503474dd773SHong Zhang   Logically collective
504474dd773SHong Zhang 
505474dd773SHong Zhang   Input Parameter:
506474dd773SHong Zhang +  ts - timestepping context
507*0fe4d17eSHong Zhang -  use_multirate - PETSC_TRUE enables the multirate RK method, sets the basic method to be RK2A and sets the ratio between slow stepsize and fast stepsize to be 2
508474dd773SHong Zhang 
509474dd773SHong Zhang   Options Database:
510*0fe4d17eSHong Zhang .   -ts_rk_multirate - <true,false>
511*0fe4d17eSHong Zhang 
512*0fe4d17eSHong Zhang   Notes:
513*0fe4d17eSHong Zhang   The multirate method requires interpolation. The default interpolation works for 1st- and 2nd- order RK, but not for high-order RKs except TSRK5DP which comes with the interpolation coeffcients (binterp).
514474dd773SHong Zhang 
515474dd773SHong Zhang   Level: intermediate
516*0fe4d17eSHong Zhang 
517*0fe4d17eSHong Zhang .seealso: TSRKGetMultirate()
518474dd773SHong Zhang @*/
519*0fe4d17eSHong Zhang PetscErrorCode TSRKSetMultirate(TS ts,PetscBool use_multirate)
520474dd773SHong Zhang {
521474dd773SHong Zhang   TS_RK          *rk = (TS_RK*)ts->data;
522474dd773SHong Zhang   PetscErrorCode ierr;
523474dd773SHong Zhang 
524474dd773SHong Zhang   PetscFunctionBegin;
525474dd773SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
526*0fe4d17eSHong Zhang   rk->use_multirate = use_multirate;
527*0fe4d17eSHong Zhang   if (use_multirate) {
528474dd773SHong Zhang     rk->dtratio = 2;
529474dd773SHong Zhang     ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSetUp_MRKSPLIT_C",TSSetUp_MRKSPLIT);CHKERRQ(ierr);
530474dd773SHong Zhang     ierr = PetscObjectComposeFunction((PetscObject)ts,"TSReset_MRKSPLIT_C",TSReset_MRKSPLIT);CHKERRQ(ierr);
531*0fe4d17eSHong Zhang     ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSetUp_MRKNONSPLIT_C",TSSetUp_MRKNONSPLIT);CHKERRQ(ierr);
532*0fe4d17eSHong Zhang     ierr = PetscObjectComposeFunction((PetscObject)ts,"TSReset_MRKNONSPLIT_C",TSReset_MRKNONSPLIT);CHKERRQ(ierr);
533*0fe4d17eSHong Zhang   } else {
534*0fe4d17eSHong Zhang     rk->dtratio = 0;
535*0fe4d17eSHong Zhang     ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSetUp_MRKSPLIT_C",NULL);CHKERRQ(ierr);
536*0fe4d17eSHong Zhang     ierr = PetscObjectComposeFunction((PetscObject)ts,"TSReset_MRKSPLIT_C",NULL);CHKERRQ(ierr);
537*0fe4d17eSHong Zhang     ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSetUp_MRKNONSPLIT_C",NULL);CHKERRQ(ierr);
538*0fe4d17eSHong Zhang     ierr = PetscObjectComposeFunction((PetscObject)ts,"TSReset_MRKNONSPLIT_C",NULL);CHKERRQ(ierr);
539474dd773SHong Zhang   }
540474dd773SHong Zhang   PetscFunctionReturn(0);
541474dd773SHong Zhang }
542474dd773SHong Zhang 
543*0fe4d17eSHong Zhang /*@C
544*0fe4d17eSHong Zhang   TSRKGetMultirate - Gets whether to Use the interpolation-based multirate RK method
545*0fe4d17eSHong Zhang 
546*0fe4d17eSHong Zhang   Not collective
547*0fe4d17eSHong Zhang 
548*0fe4d17eSHong Zhang   Input Parameter:
549*0fe4d17eSHong Zhang .  ts - timestepping context
550*0fe4d17eSHong Zhang 
551*0fe4d17eSHong Zhang   Output Parameter:
552*0fe4d17eSHong Zhang .  use_multirate - PETSC_TRUE if the multirate RK method is enabled, PETSC_FALSE otherwise
553*0fe4d17eSHong Zhang 
554*0fe4d17eSHong Zhang   Level: intermediate
555*0fe4d17eSHong Zhang 
556*0fe4d17eSHong Zhang .seealso: TSRKSetMultirate()
557*0fe4d17eSHong Zhang @*/
558*0fe4d17eSHong Zhang PetscErrorCode TSRKGetMultirate(TS ts,PetscBool *use_multirate)
559*0fe4d17eSHong Zhang {
560*0fe4d17eSHong Zhang   TS_RK          *rk = (TS_RK*)ts->data;
561*0fe4d17eSHong Zhang 
562*0fe4d17eSHong Zhang   PetscFunctionBegin;
563*0fe4d17eSHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
564*0fe4d17eSHong Zhang   *use_multirate = rk->use_multirate;
565*0fe4d17eSHong Zhang   PetscFunctionReturn(0);
566*0fe4d17eSHong Zhang }
567