xref: /petsc/src/ts/impls/explicit/rk/mrk.c (revision 030f984af8d8bb4c203755d35bded3c05b3d83ce)
1 /*
2   Code for time stepping with the multi-rate Runge-Kutta method
3 
4   Notes:
5   1) The general system is written as
6      Udot = F(t,U) for the nonsplit version of multi-rate RK,
7      user should give the indexes for both slow and fast components;
8   2) The general system is written as
9      Usdot = Fs(t,Us,Uf)
10      Ufdot = Ff(t,Us,Uf) for multi-rate RK with RHS splits,
11      user should partioned RHS by themselves and also provide the indexes for both slow and fast components.
12 */
13 
14 #include <petsc/private/tsimpl.h>
15 #include <petscdm.h>
16 #include <../src/ts/impls/explicit/rk/rk.h>
17 #include <../src/ts/impls/explicit/rk/mrk.h>
18 
19 static PetscErrorCode TSReset_RK_MultirateNonsplit(TS ts)
20 {
21   TS_RK          *rk = (TS_RK*)ts->data;
22   RKTableau      tab = rk->tableau;
23   PetscErrorCode ierr;
24 
25   PetscFunctionBegin;
26   ierr = VecDestroy(&rk->X0);CHKERRQ(ierr);
27   ierr = VecDestroyVecs(tab->s,&rk->YdotRHS_slow);CHKERRQ(ierr);
28   PetscFunctionReturn(0);
29 }
30 
31 static PetscErrorCode TSInterpolate_RK_MultirateNonsplit(TS ts,PetscReal itime,Vec X)
32 {
33   TS_RK            *rk = (TS_RK*)ts->data;
34   PetscInt         s = rk->tableau->s,p = rk->tableau->p,i,j;
35   PetscReal        h = ts->time_step;
36   PetscReal        tt,t;
37   PetscScalar      *b;
38   const PetscReal  *B = rk->tableau->binterp;
39   PetscErrorCode   ierr;
40 
41   PetscFunctionBegin;
42   if (!B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSRK %s does not have an interpolation formula",rk->tableau->name);
43   t = (itime - rk->ptime)/h;
44   ierr = PetscMalloc1(s,&b);CHKERRQ(ierr);
45   for (i=0; i<s; i++) b[i] = 0;
46   for (j=0,tt=t; j<p; j++,tt*=t) {
47     for (i=0; i<s; i++) {
48       b[i]  += h * B[i*p+j] * tt;
49     }
50   }
51   ierr = VecCopy(rk->X0,X);CHKERRQ(ierr);
52   ierr = VecMAXPY(X,s,b,rk->YdotRHS_slow);CHKERRQ(ierr);
53   ierr = PetscFree(b);CHKERRQ(ierr);
54   PetscFunctionReturn(0);
55 }
56 
57 static PetscErrorCode TSStepRefine_RK_MultirateNonsplit(TS ts)
58 {
59   TS              previousts,subts;
60   TS_RK           *rk = (TS_RK*)ts->data;
61   RKTableau       tab = rk->tableau;
62   Vec             *Y = rk->Y,*YdotRHS = rk->YdotRHS;
63   Vec             vec_fast,subvec_fast;
64   const PetscInt  s = tab->s;
65   const PetscReal *A = tab->A,*c = tab->c;
66   PetscScalar     *w = rk->work;
67   PetscInt        i,j,k;
68   PetscReal       t = ts->ptime,h = ts->time_step;
69   PetscErrorCode  ierr;
70 
71   PetscFunctionBegin;
72   ierr = VecDuplicate(ts->vec_sol,&vec_fast);CHKERRQ(ierr);
73   previousts = rk->subts_current;
74   ierr = TSRHSSplitGetSubTS(rk->subts_current,"fast",&subts);CHKERRQ(ierr);
75   ierr = TSRHSSplitGetSubTS(subts,"fast",&subts);CHKERRQ(ierr);
76   for (k=0; k<rk->dtratio; k++) {
77     for (i=0; i<s; i++) {
78       ierr = TSInterpolate_RK_MultirateNonsplit(ts,t+k*h/rk->dtratio+h/rk->dtratio*c[i],Y[i]);CHKERRQ(ierr);
79       for (j=0; j<i; j++) w[j] = h/rk->dtratio*A[i*s+j];
80       /* 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 */
81       ierr = VecCopy(ts->vec_sol,vec_fast);CHKERRQ(ierr);
82       ierr = VecMAXPY(vec_fast,i,w,YdotRHS);CHKERRQ(ierr);
83       /* update the fast components in the stage value */
84       ierr = VecGetSubVector(vec_fast,rk->is_fast,&subvec_fast);CHKERRQ(ierr);
85       ierr = VecISCopy(Y[i],rk->is_fast,SCATTER_FORWARD,subvec_fast);CHKERRQ(ierr);
86       ierr = VecRestoreSubVector(vec_fast,rk->is_fast,&subvec_fast);CHKERRQ(ierr);
87       /* compute the stage RHS */
88       ierr = TSComputeRHSFunction(ts,t+k*h/rk->dtratio+h/rk->dtratio*c[i],Y[i],YdotRHS[i]);CHKERRQ(ierr);
89     }
90     ierr = VecCopy(ts->vec_sol,vec_fast);CHKERRQ(ierr);
91     ierr = TSEvaluateStep(ts,tab->order,vec_fast,NULL);CHKERRQ(ierr);
92     /* update the fast components in the solution */
93     ierr = VecGetSubVector(vec_fast,rk->is_fast,&subvec_fast);CHKERRQ(ierr);
94     ierr = VecISCopy(ts->vec_sol,rk->is_fast,SCATTER_FORWARD,subvec_fast);CHKERRQ(ierr);
95     ierr = VecRestoreSubVector(vec_fast,rk->is_fast,&subvec_fast);CHKERRQ(ierr);
96 
97     if (subts) {
98       Vec *YdotRHS_copy;
99       ierr = VecDuplicateVecs(ts->vec_sol,s,&YdotRHS_copy);CHKERRQ(ierr);
100       rk->subts_current = rk->subts_fast;
101       ts->ptime = t+k*h/rk->dtratio;
102       ts->time_step = h/rk->dtratio;
103       ierr = TSRHSSplitGetIS(rk->subts_current,"fast",&rk->is_fast);CHKERRQ(ierr);
104       for (i=0; i<s; i++) {
105         ierr = VecCopy(rk->YdotRHS_slow[i],YdotRHS_copy[i]);CHKERRQ(ierr);
106         ierr = VecCopy(YdotRHS[i],rk->YdotRHS_slow[i]);CHKERRQ(ierr);
107       }
108 
109       ierr = TSStepRefine_RK_MultirateNonsplit(ts);CHKERRQ(ierr);
110 
111       rk->subts_current = previousts;
112       ts->ptime = t;
113       ts->time_step = h;
114       ierr = TSRHSSplitGetIS(previousts,"fast",&rk->is_fast);CHKERRQ(ierr);
115       for (i=0; i<s; i++) {
116         ierr = VecCopy(YdotRHS_copy[i],rk->YdotRHS_slow[i]);CHKERRQ(ierr);
117       }
118       ierr = VecDestroyVecs(s,&YdotRHS_copy);CHKERRQ(ierr);
119     }
120   }
121   ierr = VecDestroy(&vec_fast);CHKERRQ(ierr);
122   PetscFunctionReturn(0);
123 }
124 
125 static PetscErrorCode TSStep_RK_MultirateNonsplit(TS ts)
126 {
127   TS_RK           *rk = (TS_RK*)ts->data;
128   RKTableau       tab = rk->tableau;
129   Vec             *Y = rk->Y,*YdotRHS = rk->YdotRHS,*YdotRHS_slow = rk->YdotRHS_slow;
130   Vec             stage_slow,sol_slow; /* vectors store the slow components */
131   Vec             subvec_slow; /* sub vector to store the slow components */
132   IS              is_slow = rk->is_slow;
133   const PetscInt  s = tab->s;
134   const PetscReal *A = tab->A,*c = tab->c;
135   PetscScalar     *w = rk->work;
136   PetscInt        i,j,dtratio = rk->dtratio;
137   PetscReal       next_time_step = ts->time_step,t = ts->ptime,h = ts->time_step;
138   PetscErrorCode  ierr;
139 
140   PetscFunctionBegin;
141   rk->status = TS_STEP_INCOMPLETE;
142   ierr = VecDuplicate(ts->vec_sol,&stage_slow);CHKERRQ(ierr);
143   ierr = VecDuplicate(ts->vec_sol,&sol_slow);CHKERRQ(ierr);
144   ierr = VecCopy(ts->vec_sol,rk->X0);CHKERRQ(ierr);
145   for (i=0; i<s; i++) {
146     rk->stage_time = t + h*c[i];
147     ierr = TSPreStage(ts,rk->stage_time);CHKERRQ(ierr);
148     ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr);
149     for (j=0; j<i; j++) w[j] = h*A[i*s+j];
150     ierr = VecMAXPY(Y[i],i,w,YdotRHS_slow);CHKERRQ(ierr);
151     ierr = TSPostStage(ts,rk->stage_time,i,Y);CHKERRQ(ierr);
152     /* compute the stage RHS */
153     ierr = TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS_slow[i]);CHKERRQ(ierr);
154   }
155   /* update the slow components in the solution */
156   rk->YdotRHS = YdotRHS_slow;
157   rk->dtratio = 1;
158   ierr = TSEvaluateStep(ts,tab->order,sol_slow,NULL);CHKERRQ(ierr);
159   rk->dtratio = dtratio;
160   rk->YdotRHS = YdotRHS;
161   /* update the slow components in the solution */
162   ierr = VecGetSubVector(sol_slow,is_slow,&subvec_slow);CHKERRQ(ierr);
163   ierr = VecISCopy(ts->vec_sol,is_slow,SCATTER_FORWARD,subvec_slow);CHKERRQ(ierr);
164   ierr = VecRestoreSubVector(sol_slow,is_slow,&subvec_slow);CHKERRQ(ierr);
165 
166   rk->subts_current = ts;
167   rk->ptime = t;
168   rk->time_step = h;
169   ierr = TSStepRefine_RK_MultirateNonsplit(ts);CHKERRQ(ierr);
170 
171   ts->ptime = t + ts->time_step;
172   ts->time_step = next_time_step;
173   rk->status = TS_STEP_COMPLETE;
174 
175   /* free memory of work vectors */
176   ierr = VecDestroy(&stage_slow);CHKERRQ(ierr);
177   ierr = VecDestroy(&sol_slow);CHKERRQ(ierr);
178   PetscFunctionReturn(0);
179 }
180 
181 static PetscErrorCode TSSetUp_RK_MultirateNonsplit(TS ts)
182 {
183   TS_RK          *rk = (TS_RK*)ts->data;
184   RKTableau      tab = rk->tableau;
185   PetscErrorCode ierr;
186 
187   PetscFunctionBegin;
188   ierr = TSRHSSplitGetIS(ts,"slow",&rk->is_slow);CHKERRQ(ierr);
189   ierr = TSRHSSplitGetIS(ts,"fast",&rk->is_fast);CHKERRQ(ierr);
190   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");
191   ierr = TSRHSSplitGetSubTS(ts,"slow",&rk->subts_slow);CHKERRQ(ierr);
192   ierr = TSRHSSplitGetSubTS(ts,"fast",&rk->subts_fast);CHKERRQ(ierr);
193   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");
194   ierr = VecDuplicate(ts->vec_sol,&rk->X0);CHKERRQ(ierr);
195   ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&rk->YdotRHS_slow);CHKERRQ(ierr);
196   rk->subts_current = rk->subts_fast;
197 
198   ts->ops->step        = TSStep_RK_MultirateNonsplit;
199   ts->ops->interpolate = TSInterpolate_RK_MultirateNonsplit;
200   PetscFunctionReturn(0);
201 }
202 
203 /*
204   Copy DM from tssrc to tsdest, while keeping the original DMTS and DMSNES in tsdest.
205 */
206 static PetscErrorCode TSCopyDM(TS tssrc,TS tsdest)
207 {
208   DM             newdm,dmsrc,dmdest;
209   PetscErrorCode ierr;
210 
211   PetscFunctionBegin;
212   ierr = TSGetDM(tssrc,&dmsrc);CHKERRQ(ierr);
213   ierr = DMClone(dmsrc,&newdm);CHKERRQ(ierr);
214   ierr = TSGetDM(tsdest,&dmdest);CHKERRQ(ierr);
215   ierr = DMCopyDMTS(dmdest,newdm);CHKERRQ(ierr);
216   ierr = DMCopyDMSNES(dmdest,newdm);CHKERRQ(ierr);
217   ierr = TSSetDM(tsdest,newdm);CHKERRQ(ierr);
218   ierr = DMDestroy(&newdm);CHKERRQ(ierr);
219   PetscFunctionReturn(0);
220 }
221 
222 static PetscErrorCode TSReset_RK_MultirateSplit(TS ts)
223 {
224   TS_RK          *rk = (TS_RK*)ts->data;
225   PetscErrorCode ierr;
226 
227   PetscFunctionBegin;
228   if (rk->subts_slow) {
229     ierr = PetscFree(rk->subts_slow->data);CHKERRQ(ierr);
230     rk->subts_slow = NULL;
231   }
232   if (rk->subts_fast) {
233     ierr = PetscFree(rk->YdotRHS_fast);CHKERRQ(ierr);
234     ierr = PetscFree(rk->YdotRHS_slow);CHKERRQ(ierr);
235     ierr = VecDestroy(&rk->X0);CHKERRQ(ierr);
236     ierr = TSReset_RK_MultirateSplit(rk->subts_fast);CHKERRQ(ierr);
237     ierr = PetscFree(rk->subts_fast->data);CHKERRQ(ierr);
238     rk->subts_fast = NULL;
239   }
240   PetscFunctionReturn(0);
241 }
242 
243 static PetscErrorCode TSInterpolate_RK_MultirateSplit(TS ts,PetscReal itime,Vec X)
244 {
245   TS_RK           *rk = (TS_RK*)ts->data;
246   Vec             Xslow;
247   PetscInt        s = rk->tableau->s,p = rk->tableau->p,i,j;
248   PetscReal       h;
249   PetscReal       tt,t;
250   PetscScalar     *b;
251   const PetscReal *B = rk->tableau->binterp;
252   PetscErrorCode  ierr;
253 
254   PetscFunctionBegin;
255   if (!B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSRK %s does not have an interpolation formula",rk->tableau->name);
256 
257   switch (rk->status) {
258     case TS_STEP_INCOMPLETE:
259     case TS_STEP_PENDING:
260       h = ts->time_step;
261       t = (itime - ts->ptime)/h;
262       break;
263     case TS_STEP_COMPLETE:
264       h = ts->ptime - ts->ptime_prev;
265       t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */
266       break;
267     default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
268   }
269   ierr = PetscMalloc1(s,&b);CHKERRQ(ierr);
270   for (i=0; i<s; i++) b[i] = 0;
271   for (j=0,tt=t; j<p; j++,tt*=t) {
272     for (i=0; i<s; i++) {
273       b[i]  += h * B[i*p+j] * tt;
274     }
275   }
276   for (i=0; i<s; i++) {
277     ierr = VecGetSubVector(rk->YdotRHS[i],rk->is_slow,&rk->YdotRHS_slow[i]);CHKERRQ(ierr);
278   }
279   ierr = VecGetSubVector(X,rk->is_slow,&Xslow);CHKERRQ(ierr);
280   ierr = VecISCopy(rk->X0,rk->is_slow,SCATTER_REVERSE,Xslow);CHKERRQ(ierr);
281   ierr = VecMAXPY(Xslow,s,b,rk->YdotRHS_slow);CHKERRQ(ierr);
282   ierr = VecRestoreSubVector(X,rk->is_slow,&Xslow);CHKERRQ(ierr);
283   for (i=0; i<s; i++) {
284     ierr = VecRestoreSubVector(rk->YdotRHS[i],rk->is_slow,&rk->YdotRHS_slow[i]);CHKERRQ(ierr);
285   }
286   ierr = PetscFree(b);CHKERRQ(ierr);
287   PetscFunctionReturn(0);
288 }
289 
290 /*
291  This is for partitioned RHS multirate RK method
292  The step completion formula is
293 
294  x1 = x0 + h b^T YdotRHS
295 
296 */
297 static PetscErrorCode TSEvaluateStep_RK_MultirateSplit(TS ts,PetscInt order,Vec X,PetscBool *done)
298 {
299   TS_RK          *rk = (TS_RK*)ts->data;
300   RKTableau      tab = rk->tableau;
301   Vec            Xslow,Xfast;                  /* subvectors of X which store slow components and fast components respectively */
302   PetscScalar    *w = rk->work;
303   PetscReal      h = ts->time_step;
304   PetscInt       s = tab->s,j;
305   PetscErrorCode ierr;
306 
307   PetscFunctionBegin;
308   ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
309   if (rk->slow) {
310     for (j=0; j<s; j++) w[j] = h*tab->b[j];
311     ierr = VecGetSubVector(ts->vec_sol,rk->is_slow,&Xslow);CHKERRQ(ierr);
312     ierr = VecMAXPY(Xslow,s,w,rk->YdotRHS_slow);CHKERRQ(ierr);
313     ierr = VecRestoreSubVector(ts->vec_sol,rk->is_slow,&Xslow);CHKERRQ(ierr);
314   } else {
315     for (j=0; j<s; j++) w[j] = h/rk->dtratio*tab->b[j];
316     ierr = VecGetSubVector(X,rk->is_fast,&Xfast);CHKERRQ(ierr);
317     ierr = VecMAXPY(Xfast,s,w,rk->YdotRHS_fast);CHKERRQ(ierr);
318     ierr = VecRestoreSubVector(X,rk->is_fast,&Xfast);CHKERRQ(ierr);
319   }
320   PetscFunctionReturn(0);
321 }
322 
323 static PetscErrorCode TSStepRefine_RK_MultirateSplit(TS ts)
324 {
325   TS_RK           *rk = (TS_RK*)ts->data;
326   TS              subts_fast = rk->subts_fast,currentlevelts;
327   TS_RK           *subrk_fast = (TS_RK*)subts_fast->data;
328   RKTableau       tab = rk->tableau;
329   Vec             *Y = rk->Y;
330   Vec             *YdotRHS = rk->YdotRHS,*YdotRHS_fast = rk->YdotRHS_fast;
331   Vec             Yfast,Xfast;
332   const PetscInt  s = tab->s;
333   const PetscReal *A = tab->A,*c = tab->c;
334   PetscScalar     *w = rk->work;
335   PetscInt        i,j,k;
336   PetscReal       t = ts->ptime,h = ts->time_step;
337   PetscErrorCode  ierr;
338 
339   PetscFunctionBegin;
340   for (k=0; k<rk->dtratio; k++) {
341     ierr = VecGetSubVector(ts->vec_sol,rk->is_fast,&Xfast);CHKERRQ(ierr);
342     for (i=0; i<s; i++) {
343       ierr = VecGetSubVector(YdotRHS[i],rk->is_fast,&YdotRHS_fast[i]);CHKERRQ(ierr);
344     }
345     /* propogate fast component using small time steps */
346     for (i=0; i<s; i++) {
347       /* stage value for slow components */
348       ierr = TSInterpolate_RK_MultirateSplit(rk->ts_root,t+k*h/rk->dtratio+h/rk->dtratio*c[i],Y[i]);CHKERRQ(ierr);
349       currentlevelts = rk->ts_root;
350       while (currentlevelts != ts) { /* all the slow parts need to be interpolated separated */
351         currentlevelts = ((TS_RK*)currentlevelts->data)->subts_fast;
352         ierr = TSInterpolate_RK_MultirateSplit(currentlevelts,t+k*h/rk->dtratio+h/rk->dtratio*c[i],Y[i]);CHKERRQ(ierr);
353       }
354       for (j=0; j<i; j++) w[j] = h/rk->dtratio*A[i*s+j];
355       subrk_fast->stage_time = t + h/rk->dtratio*c[i];
356       ierr = TSPreStage(subts_fast,subrk_fast->stage_time);CHKERRQ(ierr);
357       /* stage value for fast components */
358       ierr = VecGetSubVector(Y[i],rk->is_fast,&Yfast);CHKERRQ(ierr);
359       ierr = VecCopy(Xfast,Yfast);CHKERRQ(ierr);
360       ierr = VecMAXPY(Yfast,i,w,YdotRHS_fast);CHKERRQ(ierr);
361       ierr = VecRestoreSubVector(Y[i],rk->is_fast,&Yfast);CHKERRQ(ierr);
362       ierr = TSPostStage(subts_fast,subrk_fast->stage_time,i,Y);CHKERRQ(ierr);
363       /* compute the stage RHS for fast components */
364       ierr = TSComputeRHSFunction(subts_fast,t+k*h*rk->dtratio+h/rk->dtratio*c[i],Y[i],YdotRHS_fast[i]);CHKERRQ(ierr);
365     }
366     ierr = VecRestoreSubVector(ts->vec_sol,rk->is_fast,&Xfast);CHKERRQ(ierr);
367     /* update the value of fast components using fast time step */
368     rk->slow = PETSC_FALSE;
369     ierr = TSEvaluateStep_RK_MultirateSplit(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr);
370     for (i=0; i<s; i++) {
371       ierr = VecRestoreSubVector(YdotRHS[i],rk->is_fast,&YdotRHS_fast[i]);CHKERRQ(ierr);
372     }
373 
374     if (subrk_fast->subts_fast) {
375       subts_fast->ptime = t+k*h/rk->dtratio;
376       subts_fast->time_step = h/rk->dtratio;
377       ierr = TSStepRefine_RK_MultirateSplit(subts_fast);CHKERRQ(ierr);
378     }
379     /* update the fast components of the solution */
380     ierr = VecGetSubVector(ts->vec_sol,rk->is_fast,&Xfast);CHKERRQ(ierr);
381     ierr = VecISCopy(rk->X0,rk->is_fast,SCATTER_FORWARD,Xfast);CHKERRQ(ierr);
382     ierr = VecRestoreSubVector(ts->vec_sol,rk->is_fast,&Xfast);CHKERRQ(ierr);
383   }
384   PetscFunctionReturn(0);
385 }
386 
387 static PetscErrorCode TSStep_RK_MultirateSplit(TS ts)
388 {
389   TS_RK           *rk = (TS_RK*)ts->data;
390   RKTableau       tab = rk->tableau;
391   Vec             *Y = rk->Y,*YdotRHS = rk->YdotRHS;
392   Vec             *YdotRHS_fast = rk->YdotRHS_fast,*YdotRHS_slow = rk->YdotRHS_slow;
393   Vec             Yslow,Yfast; /* subvectors store the stges of slow components and fast components respectively                           */
394   const PetscInt  s = tab->s;
395   const PetscReal *A = tab->A,*c = tab->c;
396   PetscScalar     *w = rk->work;
397   PetscInt        i,j;
398   PetscReal       next_time_step = ts->time_step,t = ts->ptime,h = ts->time_step;
399   PetscErrorCode  ierr;
400 
401   PetscFunctionBegin;
402   rk->status = TS_STEP_INCOMPLETE;
403   for (i=0; i<s; i++) {
404     ierr = VecGetSubVector(YdotRHS[i],rk->is_slow,&YdotRHS_slow[i]);CHKERRQ(ierr);
405     ierr = VecGetSubVector(YdotRHS[i],rk->is_fast,&YdotRHS_fast[i]);CHKERRQ(ierr);
406   }
407   ierr = VecCopy(ts->vec_sol,rk->X0);CHKERRQ(ierr);
408   /* propogate both slow and fast components using large time steps */
409   for (i=0; i<s; i++) {
410     rk->stage_time = t + h*c[i];
411     ierr = TSPreStage(ts,rk->stage_time);CHKERRQ(ierr);
412     ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr);
413     ierr = VecGetSubVector(Y[i],rk->is_fast,&Yfast);CHKERRQ(ierr);
414     ierr = VecGetSubVector(Y[i],rk->is_slow,&Yslow);CHKERRQ(ierr);
415     for (j=0; j<i; j++) w[j] = h*A[i*s+j];
416     ierr = VecMAXPY(Yfast,i,w,YdotRHS_fast);CHKERRQ(ierr);
417     ierr = VecMAXPY(Yslow,i,w,YdotRHS_slow);CHKERRQ(ierr);
418     ierr = VecRestoreSubVector(Y[i],rk->is_fast,&Yfast);CHKERRQ(ierr);
419     ierr = VecRestoreSubVector(Y[i],rk->is_slow,&Yslow);CHKERRQ(ierr);
420     ierr = TSPostStage(ts,rk->stage_time,i,Y);CHKERRQ(ierr);
421     ierr = TSComputeRHSFunction(rk->subts_slow,t+h*c[i],Y[i],YdotRHS_slow[i]);CHKERRQ(ierr);
422     ierr = TSComputeRHSFunction(rk->subts_fast,t+h*c[i],Y[i],YdotRHS_fast[i]);CHKERRQ(ierr);
423   }
424   rk->slow = PETSC_TRUE;
425   /* update the slow components of the solution using slow time step */
426   ierr = TSEvaluateStep_RK_MultirateSplit(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr);
427   /* YdotRHS will be used for interpolation during refinement */
428   for (i=0; i<s; i++) {
429     ierr = VecRestoreSubVector(YdotRHS[i],rk->is_slow,&YdotRHS_slow[i]);CHKERRQ(ierr);
430     ierr = VecRestoreSubVector(YdotRHS[i],rk->is_fast,&YdotRHS_fast[i]);CHKERRQ(ierr);
431   }
432 
433   ierr = TSStepRefine_RK_MultirateSplit(ts);CHKERRQ(ierr);
434 
435   ts->ptime = t + ts->time_step;
436   ts->time_step = next_time_step;
437   rk->status = TS_STEP_COMPLETE;
438   PetscFunctionReturn(0);
439 }
440 
441 static PetscErrorCode TSSetUp_RK_MultirateSplit(TS ts)
442 {
443   TS_RK          *rk = (TS_RK*)ts->data,*nextlevelrk,*currentlevelrk;
444   TS             nextlevelts;
445   Vec            X0;
446   PetscErrorCode ierr;
447 
448   PetscFunctionBegin;
449   ierr = TSRHSSplitGetIS(ts,"slow",&rk->is_slow);CHKERRQ(ierr);
450   ierr = TSRHSSplitGetIS(ts,"fast",&rk->is_fast);CHKERRQ(ierr);
451   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");
452   ierr = TSRHSSplitGetSubTS(ts,"slow",&rk->subts_slow);CHKERRQ(ierr);
453   ierr = TSRHSSplitGetSubTS(ts,"fast",&rk->subts_fast);CHKERRQ(ierr);
454   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");
455 
456   ierr = VecDuplicate(ts->vec_sol,&X0);CHKERRQ(ierr);
457   /* The TS at each level share the same tableau, work array, solution vector, stage values and stage derivatives */
458   currentlevelrk = rk;
459   while (currentlevelrk->subts_fast) {
460     ierr = PetscMalloc1(rk->tableau->s,&currentlevelrk->YdotRHS_fast);CHKERRQ(ierr);
461     ierr = PetscMalloc1(rk->tableau->s,&currentlevelrk->YdotRHS_slow);CHKERRQ(ierr);
462     ierr = PetscObjectReference((PetscObject)X0);CHKERRQ(ierr);
463     currentlevelrk->X0 = X0;
464     currentlevelrk->ts_root = ts;
465 
466     /* set up the ts for the slow part */
467     nextlevelts = currentlevelrk->subts_slow;
468     ierr = PetscNewLog(nextlevelts,&nextlevelrk);CHKERRQ(ierr);
469     nextlevelrk->tableau = rk->tableau;
470     nextlevelrk->work = rk->work;
471     nextlevelrk->Y = rk->Y;
472     nextlevelrk->YdotRHS = rk->YdotRHS;
473     nextlevelts->data = (void*)nextlevelrk;
474     ierr = TSCopyDM(ts,nextlevelts);CHKERRQ(ierr);
475     ierr = TSSetSolution(nextlevelts,ts->vec_sol);CHKERRQ(ierr);
476 
477     /* set up the ts for the fast part */
478     nextlevelts = currentlevelrk->subts_fast;
479     ierr = PetscNewLog(nextlevelts,&nextlevelrk);CHKERRQ(ierr);
480     nextlevelrk->tableau = rk->tableau;
481     nextlevelrk->work = rk->work;
482     nextlevelrk->Y = rk->Y;
483     nextlevelrk->YdotRHS = rk->YdotRHS;
484     nextlevelrk->dtratio = rk->dtratio;
485     ierr = TSRHSSplitGetIS(nextlevelts,"slow",&nextlevelrk->is_slow);CHKERRQ(ierr);
486     ierr = TSRHSSplitGetSubTS(nextlevelts,"slow",&nextlevelrk->subts_slow);CHKERRQ(ierr);
487     ierr = TSRHSSplitGetIS(nextlevelts,"fast",&nextlevelrk->is_fast);CHKERRQ(ierr);
488     ierr = TSRHSSplitGetSubTS(nextlevelts,"fast",&nextlevelrk->subts_fast);CHKERRQ(ierr);
489     nextlevelts->data = (void*)nextlevelrk;
490     ierr = TSCopyDM(ts,nextlevelts);CHKERRQ(ierr);
491     ierr = TSSetSolution(nextlevelts,ts->vec_sol);CHKERRQ(ierr);
492 
493     currentlevelrk = nextlevelrk;
494   }
495   ierr = VecDestroy(&X0);CHKERRQ(ierr);
496 
497   ts->ops->step         = TSStep_RK_MultirateSplit;
498   ts->ops->evaluatestep = TSEvaluateStep_RK_MultirateSplit;
499   ts->ops->interpolate  = TSInterpolate_RK_MultirateSplit;
500   PetscFunctionReturn(0);
501 }
502 
503 PetscErrorCode TSRKSetMultirate_RK(TS ts,PetscBool use_multirate)
504 {
505   TS_RK          *rk = (TS_RK*)ts->data;
506   PetscErrorCode ierr;
507 
508   PetscFunctionBegin;
509   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
510   rk->use_multirate = use_multirate;
511   if (use_multirate) {
512     rk->dtratio = 2;
513     ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSetUp_RK_MultirateSplit_C",TSSetUp_RK_MultirateSplit);CHKERRQ(ierr);
514     ierr = PetscObjectComposeFunction((PetscObject)ts,"TSReset_RK_MultirateSplit_C",TSReset_RK_MultirateSplit);CHKERRQ(ierr);
515     ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSetUp_RK_MultirateNonsplit_C",TSSetUp_RK_MultirateNonsplit);CHKERRQ(ierr);
516     ierr = PetscObjectComposeFunction((PetscObject)ts,"TSReset_RK_MultirateNonsplit_C",TSReset_RK_MultirateNonsplit);CHKERRQ(ierr);
517   } else {
518     rk->dtratio = 0;
519     ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSetUp_RK_MultirateSplit_C",NULL);CHKERRQ(ierr);
520     ierr = PetscObjectComposeFunction((PetscObject)ts,"TSReset_RK_MultirateSplit_C",NULL);CHKERRQ(ierr);
521     ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSetUp_RK_MultirateNonsplit_C",NULL);CHKERRQ(ierr);
522     ierr = PetscObjectComposeFunction((PetscObject)ts,"TSReset_RK_MultirateNonsplit_C",NULL);CHKERRQ(ierr);
523   }
524   PetscFunctionReturn(0);
525 }
526 
527 PetscErrorCode TSRKGetMultirate_RK(TS ts,PetscBool *use_multirate)
528 {
529   TS_RK *rk = (TS_RK*)ts->data;
530 
531   PetscFunctionBegin;
532   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
533   *use_multirate = rk->use_multirate;
534   PetscFunctionReturn(0);
535 }
536