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